Nektar++
DriverSteadyState.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: DriverSteadyState.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Driver class for the steady-state solver
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
38#include <boost/algorithm/string/classification.hpp>
39#include <boost/algorithm/string/predicate.hpp>
40#include <boost/algorithm/string/split.hpp>
41
42using namespace std;
43
44namespace Nektar
45{
46namespace SolverUtils
47{
48
53 LibUtilities::SessionReader::RegisterEnumValue("Driver", "SteadyState", 0);
54
55/**
56 *
57 */
61 : DriverModifiedArnoldi(pSession, pGraph)
62{
63}
64
65/**
66 *
67 */
69{
71}
72
73/**
74 *
75 */
77
78{
79 // With a loop over "DoSolve", this Driver implements the
80 // "encapsulated" Selective Frequency Damping method(SFD)
81 // to find the steady state of a flow above the critical Reynolds
82 // number.
83 m_equ[m_nequ - 1]->PrintSummary(out);
84
85 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 1000);
86 m_session->LoadParameter("IO_CheckSteps", m_checksteps, 100000);
87 m_session->LoadParameter("ControlCoeff", m_X, 1);
88 m_session->LoadParameter("FilterWidth", m_Delta, 2);
89
90 // To evaluate optimum SFD parameters if growth rate provided in the
91 // xml file
92 m_session->LoadParameter("GrowthRateEV", GrowthRateEV, 0.0);
93
94 // To evaluate optimum SFD parameters if frequency provided in the xml
95 // file
96 m_session->LoadParameter("FrequencyEV", FrequencyEV, 0.0);
97
98 // Determine when SFD method is converged
99 m_session->LoadParameter("TOL", TOL, 1.0e-08);
100
101 // Used only for the Adaptive SFD method
102 m_session->LoadParameter("AdaptiveTOL", AdaptiveTOL, 1.0e-02);
103
104 // Used only for the Adaptive SFD method
105 m_session->LoadParameter("AdaptiveTime", AdaptiveTime, 50.0 * m_Delta);
106
107 if (m_comm->GetRank() == 0)
108 {
110 }
111
112 timer.Start();
113
114 // Definition of shared pointer (used only for the Adaptive SFD method)
116
117 NumVar_SFD = m_equ[m_nequ - 1]->UpdateFields()[0]->GetCoordim(0);
118 // SFD to run for incompressible case with scalar field
119 if (m_session->GetSolverInfo("EqType") == "UnsteadyNavierStokes" ||
120 m_session->GetSolverInfo("EqType") == "SteadyNavierStokes")
121 {
122 int nConvectiveFields = m_session->GetVariables().size();
123 if (boost::iequals(m_session->GetVariable(nConvectiveFields - 1), "p"))
124 {
125 nConvectiveFields -= 1;
126 }
127 NumVar_SFD = nConvectiveFields;
128 }
129 // Condition necessary to run SFD for the compressible case
130 if (m_session->GetSolverInfo("EqType") == "EulerCFE" ||
131 m_session->GetSolverInfo("EqType") == "NavierStokesCFE" ||
132 m_session->GetSolverInfo("EqType") == "NavierStokesImplicitCFE")
133 {
134 // Number of variables for the compressible equations
135 NumVar_SFD += 2;
136 }
137 if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
138 {
139 if (m_session->GetSolverInfo("HOMOGENEOUS") == "1D")
140 {
141 NumVar_SFD += 1;
142 }
143 }
144
145 // We store the time step
146 m_dt = m_equ[m_nequ - 1]->GetTimeStep();
147
148 // Evaluate optimum SFD parameters if dominent EV given by xml file
149 if (GrowthRateEV != 0.0 && FrequencyEV != 0.0)
150 {
151 cout << "Besed on the dominant EV given in the xml file,"
152 << "a 1D model is used to evaluate the optumum parameters"
153 << "of the SFD method:" << endl;
154 complex<NekDouble> EV = polar(exp(GrowthRateEV), FrequencyEV);
156 }
157
158 // We set up the elements of the operator of the encapsulated
159 // formulation of the selective frequencive damping method
161
162 // m_steps is set to 1. Then "m_equ[m_nequ - 1]->DoSolve()" will run
163 // for only one time step
164 m_equ[m_nequ - 1]->SetSteps(1);
165 ofstream m_file("ConvergenceHistory.txt", ios::out | ios::trunc);
166
171
172 for (int i = 0; i < NumVar_SFD; ++i)
173 {
174 q0[i] = Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(),
175 0.0); // q0 is initialised
176 qBar0[i] =
177 Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(), 0.0);
178 m_equ[m_nequ - 1]->CopyFromPhysField(i, qBar0[i]);
179 }
180
181 /// Definition of variables used in this algorithm
182 m_stepCounter = 0;
183 m_Check = 0;
186 Diff_q_qBar = 1.0;
187 Diff_q1_q0 = 1.0;
188 cpuTime = 0.0;
189 elapsed = 0.0;
190 totalTime = 0.0;
192
193 while (max(Diff_q_qBar, Diff_q1_q0) > TOL)
194 {
195 /// Call the Navier-Stokes solver for one time step
196 m_equ[m_nequ - 1]->DoSolve();
197
198 for (int i = 0; i < NumVar_SFD; ++i)
199 {
200 /// Copy the current flow field into q0
201 m_equ[m_nequ - 1]->CopyFromPhysField(i, q0[i]);
202
203 /// Apply the linear operator to the outcome of the solver
204 ComputeSFD(i, q0, qBar0, q1, qBar1);
205
206 /// Update qBar
207 qBar0[i] = qBar1[i];
208
209 /// Copy the output of the SFD method into the current flow field
210 m_equ[m_nequ - 1]->CopyToPhysField(i, q1[i]);
211 }
212
213 if (m_infosteps && !((m_stepCounter + 1) % m_infosteps))
214 {
216
217 /// Loop for the adaptive SFD method
219 FlowPartiallyConverged == false)
220 {
221
223
225 {
226 if (m_comm->GetRank() == 0)
227 {
228 cout << "\n\t The SFD method is converging: we compute "
229 << "stability analysis using the 'partially "
230 << "converged' steady state as base flow:\n"
231 << endl;
232 }
233
234 m_equ[m_nequ - 1]->Checkpoint_BaseFlow(m_Check_BaseFlow);
236
237 A->GetAdvObject()->SetBaseFlow(q0,
238 m_equ[0]->UpdateFields());
240
241 if (m_comm->GetRank() == 0)
242 {
243 // Compute the update of the SFD parameters only on
244 // one processor
246 }
247 else
248 {
249 // On all the other processors, the parameters are set
250 // to 0
251 m_X = 0;
252 m_Delta = 0;
253 }
254 // The we give to all the processors the value of X and
255 // Delta of the first processor
258
260 }
263 {
264 if (m_comm->GetRank() == 0)
265 {
266 cout << "\n\t We compute stability analysis using"
267 << " the current flow field as base flow:\n"
268 << endl;
269 }
270
271 m_equ[m_nequ - 1]->Checkpoint_BaseFlow(m_Check_BaseFlow);
273
274 A->GetAdvObject()->SetBaseFlow(q0,
275 m_equ[0]->UpdateFields());
277
278 if (m_comm->GetRank() == 0)
279 {
280 // Compute the update of the SFD parameters only on
281 // one processor
283 }
284 else
285 {
286 // On all the other processors, the parameters are set
287 // to 0
288 m_X = 0;
289 m_Delta = 0;
290 }
291 // The we give to all the processors the value of X and
292 // Delta of the first processor
295
297 }
298 }
299 }
300
302 (!((m_stepCounter + 1) % m_checksteps)))
303 {
304 m_Check++;
305 m_equ[m_nequ - 1]->Checkpoint_Output(m_Check);
306 }
308 }
309
310 m_file.close();
311
312 /// We save the final solution into a .fld file
313 m_equ[m_nequ - 1]->Output();
314
315 for (int j = 0; j < m_equ[m_nequ - 1]->GetNvariables(); ++j)
316 {
317 NekDouble vL2Error = m_equ[m_nequ - 1]->L2Error(j, false);
318 NekDouble vLinfError = m_equ[m_nequ - 1]->LinfError(j);
319 if (m_comm->GetRank() == 0)
320 {
321 out << "L 2 error (variable " << m_equ[m_nequ - 1]->GetVariable(j)
322 << ") : " << vL2Error << endl;
323 out << "L inf error (variable " << m_equ[m_nequ - 1]->GetVariable(j)
324 << ") : " << vLinfError << endl;
325 }
326 }
327}
328
329/**
330 * This routine defines the encapsulated SFD operator with first-order
331 * splitting and exact resolution of the second subproblem.
332 *(See http://scitation.aip.org/content/aip/journal/pof2/26/3/10.1063/1.4867482
333 *for details)
334 */
336 const NekDouble Delta_input)
337{
338 NekDouble X = X_input * m_dt;
339 NekDouble Delta = Delta_input / m_dt;
340 NekDouble coeff = 1.0 / (1.0 + X * Delta);
341 M11 = coeff * (1.0 + X * Delta * exp(-(X + 1.0 / Delta)));
342 M12 = coeff * (X * Delta * (1.0 - exp(-(X + 1.0 / Delta))));
343 M21 = coeff * (1.0 - exp(-(X + 1.0 / Delta)));
344 M22 = coeff * (X * Delta + exp(-(X + 1.0 / Delta)));
345}
346
347/**
348 *
349 */
351 const int i, const Array<OneD, const Array<OneD, NekDouble>> &q0,
352 const Array<OneD, const Array<OneD, NekDouble>> &qBar0,
355{
356 q1[i] = Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(), 0.0);
357 qBar1[i] = Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(), 0.0);
358
359 /// Encapsulated SFD method
360 Vmath::Svtvp(q1[i].size(), M11, q0[i], 1, q1[i], 1, q1[i], 1);
361 Vmath::Svtvp(q1[i].size(), M12, qBar0[i], 1, q1[i], 1, q1[i], 1);
362
363 Vmath::Svtvp(qBar1[i].size(), M21, q0[i], 1, qBar1[i], 1, qBar1[i], 1);
364 Vmath::Svtvp(qBar1[i].size(), M22, qBar0[i], 1, qBar1[i], 1, qBar1[i], 1);
365}
366
367/**
368 *
369 */
371{
372 NekDouble growthEV(0.0);
373 NekDouble frequencyEV(0.0);
374
375 // m_kdim is the dimension of Krylov subspace (defined in the xml file and
376 // used in DriverArnoldi.cpp)
377 ReadEVfile(m_kdim, growthEV, frequencyEV);
378
379 cout << "\n\tgrowthEV = " << growthEV << endl;
380 cout << "\tfrequencyEV = " << frequencyEV << endl;
381
382 complex<NekDouble> ApproxEV = polar(exp(growthEV), frequencyEV);
383
384 NekDouble X_new = m_X;
385 NekDouble Delta_new = m_Delta;
386
387 GradientDescentMethod(ApproxEV, X_new, Delta_new);
388
389 m_X = X_new;
390 m_Delta = Delta_new;
391
393}
394
395/**
396 * This routine implements a gradient descent method to find the parameters X
397 * end Delta which give the minimum eigenlavue of the SFD problem applied to
398 * the scalar case u(n+1) = \alpha*u(n).
399 */
400void DriverSteadyState::GradientDescentMethod(const complex<NekDouble> &alpha,
401 NekDouble &X_output,
402 NekDouble &Delta_output)
403{
404 cout << "\n\tWe enter the Gradient Descent Method [...]" << endl;
405 bool OptParmFound = false;
406 bool Descending = true;
407 NekDouble X_input = X_output;
408 NekDouble Delta_input = Delta_output;
409
410 NekDouble X_init = X_output;
411 NekDouble Delta_init = Delta_output;
412 int stepCounter(0);
413
414 NekDouble F0(0.0);
415 NekDouble F0x(0.0);
416 NekDouble F0y(0.0);
417 NekDouble F1(0.0);
418 NekDouble dx = 0.00000001;
419 NekDouble dirx(0.0);
420 NekDouble diry(0.0);
421 NekDouble s(0.0);
422 NekDouble CurrentMin = 1.0;
423
424 while (OptParmFound == false)
425 {
426 Descending = true;
427 EvalEV_ScalarSFD(X_input, Delta_input, alpha, F0);
428 EvalEV_ScalarSFD(X_input + dx, Delta_input, alpha, F0x);
429 EvalEV_ScalarSFD(X_input, Delta_input + dx, alpha, F0y);
430 dirx = (F0 - F0x);
431 diry = (F0 - F0y);
432 s = abs(0.000001 / dirx);
433 X_output = X_input + s * dirx;
434 Delta_output = Delta_input + s * diry;
435 F1 = F0;
436
437 while (Descending == true)
438 {
439 CurrentMin = F1;
440 X_input = X_output;
441 Delta_input = Delta_output;
442 EvalEV_ScalarSFD(X_output, Delta_output, alpha, F1);
443
444 if (F1 > CurrentMin)
445 {
446 Descending = false;
447 }
448 else
449 {
450 s = s + s * 0.01;
451 X_output = X_input + s * dirx;
452 Delta_output = Delta_input + s * diry;
453 }
454
455 if (stepCounter > 9999999)
456 {
457 // We are stuck in this loop..
458 // Then we restart it with different initail conditions
459 Descending = false;
460 X_input = X_init;
461 Delta_init = Delta_init + Delta_init * 0.1;
462 Delta_input = Delta_init;
463 stepCounter = 0;
464 }
465 stepCounter++;
466 }
467
468 if (abs(F0 - F1) < dx)
469 {
470 cout << "\tThe Gradient Descent Method has converged!" << endl;
471 EvalEV_ScalarSFD(X_output, Delta_output, alpha, F1);
472 cout << "\n\tThe updated parameters are: X_tilde = " << X_output
473 << " and Delta_tilde = " << Delta_output << endl;
474 OptParmFound = true;
475 }
476 }
477}
478
479/**
480 * This routine evaluates the maximum eigenvalue of the SFD system when applied
481 * to the 1D model u(n+1) = alpha*u(n)
482 */
484 const NekDouble &Delta_input,
485 const complex<NekDouble> &alpha,
486 NekDouble &MaxEV)
487{
488 NekDouble A11 =
489 (1.0 + X_input * Delta_input * exp(-(X_input + 1.0 / Delta_input))) /
490 (1.0 + X_input * Delta_input);
491 NekDouble A12 =
492 (X_input * Delta_input -
493 X_input * Delta_input * exp(-(X_input + 1.0 / Delta_input))) /
494 (1.0 + X_input * Delta_input);
495 NekDouble A21 = (1.0 - 1.0 * exp(-(X_input + 1.0 / Delta_input))) /
496 (1 + X_input * Delta_input);
497 NekDouble A22 =
498 (X_input * Delta_input + 1.0 * exp(-(X_input + 1.0 / Delta_input))) /
499 (1.0 + X_input * Delta_input);
500
501 complex<NekDouble> B11 = alpha;
502 NekDouble B12 = 0.0;
503 NekDouble B21 = 0.0;
504 NekDouble B22 = 1.0;
505
506 complex<NekDouble> a = A11 * B11 + A12 * B21;
507 NekDouble b = A11 * B12 + A12 * B22;
508 complex<NekDouble> c = A21 * B11 + A22 * B21;
509 NekDouble d = A21 * B12 + A22 * B22;
510
511 complex<NekDouble> delt = sqrt((a - d) * (a - d) + 4.0 * b * c);
512 complex<NekDouble> lambda_1 = 0.5 * (a + d + delt);
513 complex<NekDouble> lambda_2 = 0.5 * (a + d - delt);
514
515 NekDouble NORM_1 = abs(lambda_1);
516 NekDouble NORM_2 = abs(lambda_2);
517
518 MaxEV = max(NORM_1, NORM_2);
519}
520
521/**
522 *
523 */
524void DriverSteadyState::ReadEVfile(int &KrylovSubspaceDim, NekDouble &growthEV,
525 NekDouble &frequencyEV)
526{
527 // This routine reads the .evl file written by the Arnoldi algorithm
528 // (written in September 2014)
529 std::string line;
530 int NumLinesInFile = 0;
531 std::string EVfileName = m_session->GetSessionName() + ".evl";
532 std::ifstream EVfile(EVfileName.c_str());
533 ASSERTL0(EVfile.good(), "Cannot open .evl file.");
534
535 if (EVfile)
536 {
537 // This block counts the total number of lines of the .evl file
538 // We keep going util we reach the end of the file
539 while (getline(EVfile, line))
540 {
541 ++NumLinesInFile;
542 }
543
544 // It may happen that the Stability method that have produced the .elv
545 // file converges in less than m_kdim iterations. In this case,
546 // KrylovSubspaceDim has to be changed here
547 if (NumLinesInFile <
548 KrylovSubspaceDim * 2.0 +
549 KrylovSubspaceDim * (KrylovSubspaceDim + 1.0) / 2.0)
550 {
551 for (int i = 1; i <= KrylovSubspaceDim; ++i)
552 {
553 if (NumLinesInFile == i * 2.0 + i * (i + 1.0) / 2.0)
554 {
555 KrylovSubspaceDim = i;
556 }
557 }
558 }
559
560 // go back to the beginning of the file
561 EVfile.clear();
562 EVfile.seekg(0, ios::beg);
563
564 // We now want to go to the line where the most unstable eigenlavue was
565 // written
566 for (int i = 0; i < (NumLinesInFile - KrylovSubspaceDim); ++i)
567 {
568 std::getline(EVfile, line);
569 cout << "Discard line: " << line << endl;
570 }
571
572 std::vector<std::string> tokens;
573 std::getline(EVfile, line);
574 boost::algorithm::split(tokens, line, boost::is_any_of("\t "),
575 boost::token_compress_on);
576
577 ASSERTL0(tokens.size() >= 5,
578 "Unexpected formatting of .evl file while reading line:\n" +
579 line);
580 growthEV = boost::lexical_cast<NekDouble>(tokens[4]);
581 frequencyEV = boost::lexical_cast<NekDouble>(tokens[5]);
582 }
583 else
584 {
585 cout << "An error occurred when opening the .evl file" << endl;
586 }
587 EVfile.close();
588}
589
590/**
591 * This routine evaluates |q-qBar|_inf (and |q1-q0|_inf) and writes the values
592 * in "ConvergenceHistory.txt"
593 */
595 const Array<OneD, const Array<OneD, NekDouble>> &qBar1,
596 const Array<OneD, const Array<OneD, NekDouble>> &q0,
597 NekDouble &MaxNormDiff_q_qBar, NekDouble &MaxNormDiff_q1_q0)
598{
599 Array<OneD, NekDouble> NormDiff_q_qBar(NumVar_SFD, 1.0);
600 Array<OneD, NekDouble> NormDiff_q1_q0(NumVar_SFD, 1.0);
601 MaxNormDiff_q_qBar = 0.0;
602 MaxNormDiff_q1_q0 = 0.0;
603
604 for (int i = 0; i < NumVar_SFD; ++i)
605 {
606 NormDiff_q_qBar[i] = m_equ[m_nequ - 1]->LinfError(i, qBar1[i]);
607 NormDiff_q1_q0[i] = m_equ[m_nequ - 1]->LinfError(i, q0[i]);
608
609 if (MaxNormDiff_q_qBar < NormDiff_q_qBar[i])
610 {
611 MaxNormDiff_q_qBar = NormDiff_q_qBar[i];
612 }
613 if (MaxNormDiff_q1_q0 < NormDiff_q1_q0[i])
614 {
615 MaxNormDiff_q1_q0 = NormDiff_q1_q0[i];
616 }
617 }
618
619 timer.Stop();
621 cpuTime += elapsed;
623
624 if (m_comm->GetRank() == 0)
625 {
626 cout << "SFD - Step: " << left << m_stepCounter + 1
627 << ";\tTime: " << left << m_equ[m_nequ - 1]->GetFinalTime()
628 << ";\tCPU time = " << left << cpuTime << " s"
629 << ";\tTot time = " << left << totalTime << " s"
630 << ";\tX = " << left << m_X << ";\tDelta = " << left << m_Delta
631 << ";\t|q-qBar|inf = " << left << MaxNormDiff_q_qBar << endl;
632 std::ofstream m_file("ConvergenceHistory.txt", std::ios_base::app);
633 m_file << m_stepCounter + 1 << "\t" << m_equ[m_nequ - 1]->GetFinalTime()
634 << "\t" << totalTime << "\t" << MaxNormDiff_q_qBar << "\t"
635 << MaxNormDiff_q1_q0 << endl;
636 m_file.close();
637 }
638
639 cpuTime = 0.0;
640 timer.Start();
641}
642
643/**
644 *
645 */
647{
648 cout << "\n====================================="
649 "=================================="
650 << endl;
651 cout << "Parameters for the SFD method:" << endl;
652 cout << "\tControl Coefficient: X = " << m_X << endl;
653 cout << "\tFilter Width: Delta = " << m_Delta << endl;
654 cout << "The simulation is stopped when |q-qBar|inf < " << TOL << endl;
656 {
657 cout << "\nWe use the adaptive SFD method:" << endl;
658 cout << " The parameters are updated every " << AdaptiveTime
659 << " time units;" << endl;
660 cout << " until |q-qBar|inf becomes smaller than " << AdaptiveTOL
661 << endl;
662 }
663 cout << "====================================="
664 "==================================\n"
665 << endl;
666}
667
668} // namespace SolverUtils
669} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition: Timer.cpp:67
A base class for PDEs which include an advection component.
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:85
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:82
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:100
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
int m_nequ
number of equations
Definition: Driver.h:97
virtual void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
virtual void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
void ComputeSFD(const int i, const Array< OneD, const Array< OneD, NekDouble > > &q0, const Array< OneD, const Array< OneD, NekDouble > > &qBar0, Array< OneD, Array< OneD, NekDouble > > &q1, Array< OneD, Array< OneD, NekDouble > > &qBar1)
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout) override
Initialises EquationSystem class members.
void ConvergenceHistory(const Array< OneD, const Array< OneD, NekDouble > > &qBar1, const Array< OneD, const Array< OneD, NekDouble > > &q0, NekDouble &MaxNormDiff_q_qBar, NekDouble &MaxNormDiff_q1_q0)
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
SOLVER_UTILS_EXPORT DriverSteadyState(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.
NekDouble m_Delta
Definition of the SFD parameters:
void SetSFDOperator(const NekDouble X_input, const NekDouble Delta_input)
NekDouble M11
Definition of the SFD operator.
virtual SOLVER_UTILS_EXPORT void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
void GradientDescentMethod(const std::complex< NekDouble > &alpha, NekDouble &X_output, NekDouble &Delta_output)
static std::string className
Name of the class.
void EvalEV_ScalarSFD(const NekDouble &X_input, const NekDouble &Delta_input, const std::complex< NekDouble > &alpha, NekDouble &MaxEV)
bool FlowPartiallyConverged
For adaptive SFD method.
void ReadEVfile(int &KrylovSubspaceDim, NekDouble &growthEV, NekDouble &frequencyEV)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:67
std::shared_ptr< AdvectionSystem > AdvectionSystemSharedPtr
Shared pointer to an AdvectionSystem class.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::vector< double > d(NPUPPER *NPUPPER)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:617
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294