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::SolverUtils
45{
46
51 LibUtilities::SessionReader::RegisterEnumValue("Driver", "SteadyState", 0);
52
53/**
54 *
55 */
59 : DriverModifiedArnoldi(pSession, pGraph)
60{
61}
62
63/**
64 *
65 */
67{
69}
70
71/**
72 *
73 */
75
76{
77 // With a loop over "DoSolve", this Driver implements the
78 // "encapsulated" Selective Frequency Damping method(SFD)
79 // to find the steady state of a flow above the critical Reynolds
80 // number.
81 m_equ[m_nequ - 1]->PrintSummary(out);
82
83 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 1000);
84 m_session->LoadParameter("IO_CheckSteps", m_checksteps, 100000);
85 m_session->LoadParameter("ControlCoeff", m_X, 1);
86 m_session->LoadParameter("FilterWidth", m_Delta, 2);
87
88 // To evaluate optimum SFD parameters if growth rate provided in the
89 // xml file
90 m_session->LoadParameter("GrowthRateEV", GrowthRateEV, 0.0);
91
92 // To evaluate optimum SFD parameters if frequency provided in the xml
93 // file
94 m_session->LoadParameter("FrequencyEV", FrequencyEV, 0.0);
95
96 // Determine when SFD method is converged
97 m_session->LoadParameter("TOL", TOL, 1.0e-08);
98
99 // Used only for the Adaptive SFD method
100 m_session->LoadParameter("AdaptiveTOL", AdaptiveTOL, 1.0e-02);
101
102 // Used only for the Adaptive SFD method
103 m_session->LoadParameter("AdaptiveTime", AdaptiveTime, 50.0 * m_Delta);
104
105 if (m_comm->GetRank() == 0)
106 {
108 }
109
110 timer.Start();
111
112 // Definition of shared pointer (used only for the Adaptive SFD method)
114
115 NumVar_SFD = m_equ[m_nequ - 1]->UpdateFields()[0]->GetCoordim(0);
116 // SFD to run for incompressible case with scalar field
117 if (m_session->GetSolverInfo("EqType") == "UnsteadyNavierStokes" ||
118 m_session->GetSolverInfo("EqType") == "SteadyNavierStokes")
119 {
120 int nConvectiveFields = m_session->GetVariables().size();
121 if (boost::iequals(m_session->GetVariable(nConvectiveFields - 1), "p"))
122 {
123 nConvectiveFields -= 1;
124 }
125 NumVar_SFD = nConvectiveFields;
126 }
127 // Condition necessary to run SFD for the compressible case
128 if (m_session->GetSolverInfo("EqType") == "EulerCFE" ||
129 m_session->GetSolverInfo("EqType") == "NavierStokesCFE" ||
130 m_session->GetSolverInfo("EqType") == "NavierStokesImplicitCFE")
131 {
132 // Number of variables for the compressible equations
133 NumVar_SFD += 2;
134 }
135 if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
136 {
137 if (m_session->GetSolverInfo("HOMOGENEOUS") == "1D")
138 {
139 NumVar_SFD += 1;
140 }
141 }
142
143 // We store the time step
144 m_dt = m_equ[m_nequ - 1]->GetTimeStep();
145
146 // Evaluate optimum SFD parameters if dominent EV given by xml file
147 if (GrowthRateEV != 0.0 && FrequencyEV != 0.0)
148 {
149 cout << "Besed on the dominant EV given in the xml file,"
150 << "a 1D model is used to evaluate the optumum parameters"
151 << "of the SFD method:" << endl;
152 complex<NekDouble> EV = polar(exp(GrowthRateEV), FrequencyEV);
154 }
155
156 // We set up the elements of the operator of the encapsulated
157 // formulation of the selective frequencive damping method
159
160 // m_steps is set to 1. Then "m_equ[m_nequ - 1]->DoSolve()" will run
161 // for only one time step
162 m_equ[m_nequ - 1]->SetSteps(1);
163 ofstream m_file("ConvergenceHistory.txt", ios::out | ios::trunc);
164
169
170 for (int i = 0; i < NumVar_SFD; ++i)
171 {
172 q0[i] = Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(),
173 0.0); // q0 is initialised
174 qBar0[i] =
175 Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(), 0.0);
176 m_equ[m_nequ - 1]->CopyFromPhysField(i, qBar0[i]);
177 }
178
179 /// Definition of variables used in this algorithm
180 m_stepCounter = 0;
181 m_Check = 0;
184 Diff_q_qBar = 1.0;
185 Diff_q1_q0 = 1.0;
186 cpuTime = 0.0;
187 elapsed = 0.0;
188 totalTime = 0.0;
190
191 while (max(Diff_q_qBar, Diff_q1_q0) > TOL)
192 {
193 /// Call the Navier-Stokes solver for one time step
194 m_equ[m_nequ - 1]->DoSolve();
195
196 for (int i = 0; i < NumVar_SFD; ++i)
197 {
198 /// Copy the current flow field into q0
199 m_equ[m_nequ - 1]->CopyFromPhysField(i, q0[i]);
200
201 /// Apply the linear operator to the outcome of the solver
202 ComputeSFD(i, q0, qBar0, q1, qBar1);
203
204 /// Update qBar
205 qBar0[i] = qBar1[i];
206
207 /// Copy the output of the SFD method into the current flow field
208 m_equ[m_nequ - 1]->CopyToPhysField(i, q1[i]);
209 }
210
211 if (m_infosteps && !((m_stepCounter + 1) % m_infosteps))
212 {
214
215 /// Loop for the adaptive SFD method
217 FlowPartiallyConverged == false)
218 {
219
221
223 {
224 if (m_comm->GetRank() == 0)
225 {
226 cout << "\n\t The SFD method is converging: we compute "
227 << "stability analysis using the 'partially "
228 << "converged' steady state as base flow:\n"
229 << endl;
230 }
231
232 m_equ[m_nequ - 1]->Checkpoint_BaseFlow(m_Check_BaseFlow);
234
235 A->GetAdvObject()->SetBaseFlow(q0,
236 m_equ[0]->UpdateFields());
238
239 if (m_comm->GetRank() == 0)
240 {
241 // Compute the update of the SFD parameters only on
242 // one processor
244 }
245 else
246 {
247 // On all the other processors, the parameters are set
248 // to 0
249 m_X = 0;
250 m_Delta = 0;
251 }
252 // The we give to all the processors the value of X and
253 // Delta of the first processor
256
258 }
261 {
262 if (m_comm->GetRank() == 0)
263 {
264 cout << "\n\t We compute stability analysis using"
265 << " the current flow field as base flow:\n"
266 << endl;
267 }
268
269 m_equ[m_nequ - 1]->Checkpoint_BaseFlow(m_Check_BaseFlow);
271
272 A->GetAdvObject()->SetBaseFlow(q0,
273 m_equ[0]->UpdateFields());
275
276 if (m_comm->GetRank() == 0)
277 {
278 // Compute the update of the SFD parameters only on
279 // one processor
281 }
282 else
283 {
284 // On all the other processors, the parameters are set
285 // to 0
286 m_X = 0;
287 m_Delta = 0;
288 }
289 // The we give to all the processors the value of X and
290 // Delta of the first processor
293
295 }
296 }
297 }
298
300 (!((m_stepCounter + 1) % m_checksteps)))
301 {
302 m_Check++;
303 m_equ[m_nequ - 1]->Checkpoint_Output(m_Check);
304 }
306 }
307
308 m_file.close();
309
310 /// We save the final solution into a .fld file
311 m_equ[m_nequ - 1]->Output();
312
313 for (int j = 0; j < m_equ[m_nequ - 1]->GetNvariables(); ++j)
314 {
315 NekDouble vL2Error = m_equ[m_nequ - 1]->L2Error(j, false);
316 NekDouble vLinfError = m_equ[m_nequ - 1]->LinfError(j);
317 if (m_comm->GetRank() == 0)
318 {
319 out << "L 2 error (variable " << m_equ[m_nequ - 1]->GetVariable(j)
320 << ") : " << vL2Error << endl;
321 out << "L inf error (variable " << m_equ[m_nequ - 1]->GetVariable(j)
322 << ") : " << vLinfError << endl;
323 }
324 }
325}
326
327/**
328 * This routine defines the encapsulated SFD operator with first-order
329 * splitting and exact resolution of the second subproblem.
330 *(See http://scitation.aip.org/content/aip/journal/pof2/26/3/10.1063/1.4867482
331 *for details)
332 */
334 const NekDouble Delta_input)
335{
336 NekDouble X = X_input * m_dt;
337 NekDouble Delta = Delta_input / m_dt;
338 NekDouble coeff = 1.0 / (1.0 + X * Delta);
339 M11 = coeff * (1.0 + X * Delta * exp(-(X + 1.0 / Delta)));
340 M12 = coeff * (X * Delta * (1.0 - exp(-(X + 1.0 / Delta))));
341 M21 = coeff * (1.0 - exp(-(X + 1.0 / Delta)));
342 M22 = coeff * (X * Delta + exp(-(X + 1.0 / Delta)));
343}
344
345/**
346 *
347 */
349 const int i, const Array<OneD, const Array<OneD, NekDouble>> &q0,
350 const Array<OneD, const Array<OneD, NekDouble>> &qBar0,
353{
354 q1[i] = Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(), 0.0);
355 qBar1[i] = Array<OneD, NekDouble>(m_equ[m_nequ - 1]->GetTotPoints(), 0.0);
356
357 /// Encapsulated SFD method
358 Vmath::Svtvp(q1[i].size(), M11, q0[i], 1, q1[i], 1, q1[i], 1);
359 Vmath::Svtvp(q1[i].size(), M12, qBar0[i], 1, q1[i], 1, q1[i], 1);
360
361 Vmath::Svtvp(qBar1[i].size(), M21, q0[i], 1, qBar1[i], 1, qBar1[i], 1);
362 Vmath::Svtvp(qBar1[i].size(), M22, qBar0[i], 1, qBar1[i], 1, qBar1[i], 1);
363}
364
365/**
366 *
367 */
369{
370 NekDouble growthEV(0.0);
371 NekDouble frequencyEV(0.0);
372
373 // m_kdim is the dimension of Krylov subspace (defined in the xml file and
374 // used in DriverArnoldi.cpp)
375 ReadEVfile(m_kdim, growthEV, frequencyEV);
376
377 cout << "\n\tgrowthEV = " << growthEV << endl;
378 cout << "\tfrequencyEV = " << frequencyEV << endl;
379
380 complex<NekDouble> ApproxEV = polar(exp(growthEV), frequencyEV);
381
382 NekDouble X_new = m_X;
383 NekDouble Delta_new = m_Delta;
384
385 GradientDescentMethod(ApproxEV, X_new, Delta_new);
386
387 m_X = X_new;
388 m_Delta = Delta_new;
389
391}
392
393/**
394 * This routine implements a gradient descent method to find the parameters X
395 * end Delta which give the minimum eigenlavue of the SFD problem applied to
396 * the scalar case u(n+1) = \alpha*u(n).
397 */
398void DriverSteadyState::GradientDescentMethod(const complex<NekDouble> &alpha,
399 NekDouble &X_output,
400 NekDouble &Delta_output)
401{
402 cout << "\n\tWe enter the Gradient Descent Method [...]" << endl;
403 bool OptParmFound = false;
404 bool Descending = true;
405 NekDouble X_input = X_output;
406 NekDouble Delta_input = Delta_output;
407
408 NekDouble X_init = X_output;
409 NekDouble Delta_init = Delta_output;
410 int stepCounter(0);
411
412 NekDouble F0(0.0);
413 NekDouble F0x(0.0);
414 NekDouble F0y(0.0);
415 NekDouble F1(0.0);
416 NekDouble dx = 0.00000001;
417 NekDouble dirx(0.0);
418 NekDouble diry(0.0);
419 NekDouble s(0.0);
420 NekDouble CurrentMin = 1.0;
421
422 while (OptParmFound == false)
423 {
424 Descending = true;
425 EvalEV_ScalarSFD(X_input, Delta_input, alpha, F0);
426 EvalEV_ScalarSFD(X_input + dx, Delta_input, alpha, F0x);
427 EvalEV_ScalarSFD(X_input, Delta_input + dx, alpha, F0y);
428 dirx = (F0 - F0x);
429 diry = (F0 - F0y);
430 s = abs(0.000001 / dirx);
431 X_output = X_input + s * dirx;
432 Delta_output = Delta_input + s * diry;
433 F1 = F0;
434
435 while (Descending == true)
436 {
437 CurrentMin = F1;
438 X_input = X_output;
439 Delta_input = Delta_output;
440 EvalEV_ScalarSFD(X_output, Delta_output, alpha, F1);
441
442 if (F1 > CurrentMin)
443 {
444 Descending = false;
445 }
446 else
447 {
448 s = s + s * 0.01;
449 X_output = X_input + s * dirx;
450 Delta_output = Delta_input + s * diry;
451 }
452
453 if (stepCounter > 9999999)
454 {
455 // We are stuck in this loop..
456 // Then we restart it with different initail conditions
457 Descending = false;
458 X_input = X_init;
459 Delta_init = Delta_init + Delta_init * 0.1;
460 Delta_input = Delta_init;
461 stepCounter = 0;
462 }
463 stepCounter++;
464 }
465
466 if (abs(F0 - F1) < dx)
467 {
468 cout << "\tThe Gradient Descent Method has converged!" << endl;
469 EvalEV_ScalarSFD(X_output, Delta_output, alpha, F1);
470 cout << "\n\tThe updated parameters are: X_tilde = " << X_output
471 << " and Delta_tilde = " << Delta_output << endl;
472 OptParmFound = true;
473 }
474 }
475}
476
477/**
478 * This routine evaluates the maximum eigenvalue of the SFD system when applied
479 * to the 1D model u(n+1) = alpha*u(n)
480 */
482 const NekDouble &Delta_input,
483 const complex<NekDouble> &alpha,
484 NekDouble &MaxEV)
485{
486 NekDouble A11 =
487 (1.0 + X_input * Delta_input * exp(-(X_input + 1.0 / Delta_input))) /
488 (1.0 + X_input * Delta_input);
489 NekDouble A12 =
490 (X_input * Delta_input -
491 X_input * Delta_input * exp(-(X_input + 1.0 / Delta_input))) /
492 (1.0 + X_input * Delta_input);
493 NekDouble A21 = (1.0 - 1.0 * exp(-(X_input + 1.0 / Delta_input))) /
494 (1 + X_input * Delta_input);
495 NekDouble A22 =
496 (X_input * Delta_input + 1.0 * exp(-(X_input + 1.0 / Delta_input))) /
497 (1.0 + X_input * Delta_input);
498
499 complex<NekDouble> B11 = alpha;
500 NekDouble B12 = 0.0;
501 NekDouble B21 = 0.0;
502 NekDouble B22 = 1.0;
503
504 complex<NekDouble> a = A11 * B11 + A12 * B21;
505 NekDouble b = A11 * B12 + A12 * B22;
506 complex<NekDouble> c = A21 * B11 + A22 * B21;
507 NekDouble d = A21 * B12 + A22 * B22;
508
509 complex<NekDouble> delt = sqrt((a - d) * (a - d) + 4.0 * b * c);
510 complex<NekDouble> lambda_1 = 0.5 * (a + d + delt);
511 complex<NekDouble> lambda_2 = 0.5 * (a + d - delt);
512
513 NekDouble NORM_1 = abs(lambda_1);
514 NekDouble NORM_2 = abs(lambda_2);
515
516 MaxEV = max(NORM_1, NORM_2);
517}
518
519/**
520 *
521 */
522void DriverSteadyState::ReadEVfile(int &KrylovSubspaceDim, NekDouble &growthEV,
523 NekDouble &frequencyEV)
524{
525 // This routine reads the .evl file written by the Arnoldi algorithm
526 // (written in September 2014)
527 std::string line;
528 int NumLinesInFile = 0;
529 std::string EVfileName = m_session->GetSessionName() + ".evl";
530 std::ifstream EVfile(EVfileName.c_str());
531 ASSERTL0(EVfile.good(), "Cannot open .evl file.");
532
533 if (EVfile)
534 {
535 // This block counts the total number of lines of the .evl file
536 // We keep going util we reach the end of the file
537 while (getline(EVfile, line))
538 {
539 ++NumLinesInFile;
540 }
541
542 // It may happen that the Stability method that have produced the .elv
543 // file converges in less than m_kdim iterations. In this case,
544 // KrylovSubspaceDim has to be changed here
545 if (NumLinesInFile <
546 KrylovSubspaceDim * 2.0 +
547 KrylovSubspaceDim * (KrylovSubspaceDim + 1.0) / 2.0)
548 {
549 for (int i = 1; i <= KrylovSubspaceDim; ++i)
550 {
551 if (NumLinesInFile == i * 2.0 + i * (i + 1.0) / 2.0)
552 {
553 KrylovSubspaceDim = i;
554 }
555 }
556 }
557
558 // go back to the beginning of the file
559 EVfile.clear();
560 EVfile.seekg(0, ios::beg);
561
562 // We now want to go to the line where the most unstable eigenlavue was
563 // written
564 for (int i = 0; i < (NumLinesInFile - KrylovSubspaceDim); ++i)
565 {
566 std::getline(EVfile, line);
567 cout << "Discard line: " << line << endl;
568 }
569
570 std::vector<std::string> tokens;
571 std::getline(EVfile, line);
572 boost::algorithm::split(tokens, line, boost::is_any_of("\t "),
573 boost::token_compress_on);
574
575 ASSERTL0(tokens.size() >= 5,
576 "Unexpected formatting of .evl file while reading line:\n" +
577 line);
578 growthEV = boost::lexical_cast<NekDouble>(tokens[4]);
579 frequencyEV = boost::lexical_cast<NekDouble>(tokens[5]);
580 }
581 else
582 {
583 cout << "An error occurred when opening the .evl file" << endl;
584 }
585 EVfile.close();
586}
587
588/**
589 * This routine evaluates |q-qBar|_inf (and |q1-q0|_inf) and writes the values
590 * in "ConvergenceHistory.txt"
591 */
593 const Array<OneD, const Array<OneD, NekDouble>> &qBar1,
594 const Array<OneD, const Array<OneD, NekDouble>> &q0,
595 NekDouble &MaxNormDiff_q_qBar, NekDouble &MaxNormDiff_q1_q0)
596{
597 Array<OneD, NekDouble> NormDiff_q_qBar(NumVar_SFD, 1.0);
598 Array<OneD, NekDouble> NormDiff_q1_q0(NumVar_SFD, 1.0);
599 MaxNormDiff_q_qBar = 0.0;
600 MaxNormDiff_q1_q0 = 0.0;
601
602 for (int i = 0; i < NumVar_SFD; ++i)
603 {
604 NormDiff_q_qBar[i] = m_equ[m_nequ - 1]->LinfError(i, qBar1[i]);
605 NormDiff_q1_q0[i] = m_equ[m_nequ - 1]->LinfError(i, q0[i]);
606
607 if (MaxNormDiff_q_qBar < NormDiff_q_qBar[i])
608 {
609 MaxNormDiff_q_qBar = NormDiff_q_qBar[i];
610 }
611 if (MaxNormDiff_q1_q0 < NormDiff_q1_q0[i])
612 {
613 MaxNormDiff_q1_q0 = NormDiff_q1_q0[i];
614 }
615 }
616
617 timer.Stop();
619 cpuTime += elapsed;
621
622 if (m_comm->GetRank() == 0)
623 {
624 cout << "SFD - Step: " << left << m_stepCounter + 1
625 << ";\tTime: " << left << m_equ[m_nequ - 1]->GetTime()
626 << ";\tCPU time = " << left << cpuTime << " s"
627 << ";\tTot time = " << left << totalTime << " s"
628 << ";\tX = " << left << m_X << ";\tDelta = " << left << m_Delta
629 << ";\t|q-qBar|inf = " << left << MaxNormDiff_q_qBar << endl;
630 std::ofstream m_file("ConvergenceHistory.txt", std::ios_base::app);
631 m_file << m_stepCounter + 1 << "\t" << m_equ[m_nequ - 1]->GetTime()
632 << "\t" << totalTime << "\t" << MaxNormDiff_q_qBar << "\t"
633 << MaxNormDiff_q1_q0 << endl;
634 m_file.close();
635 }
636
637 cpuTime = 0.0;
638 timer.Start();
639}
640
641/**
642 *
643 */
645{
646 cout << "\n====================================="
647 "=================================="
648 << endl;
649 cout << "Parameters for the SFD method:" << endl;
650 cout << "\tControl Coefficient: X = " << m_X << endl;
651 cout << "\tFilter Width: Delta = " << m_Delta << endl;
652 cout << "The simulation is stopped when |q-qBar|inf < " << TOL << endl;
654 {
655 cout << "\nWe use the adaptive SFD method:" << endl;
656 cout << " The parameters are updated every " << AdaptiveTime
657 << " time units;" << endl;
658 cout << " until |q-qBar|inf becomes smaller than " << AdaptiveTOL
659 << endl;
660 }
661 cout << "====================================="
662 "==================================\n"
663 << endl;
664}
665
666} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
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:65
A base class for PDEs which include an advection component.
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:83
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:80
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:98
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:92
int m_nequ
number of equations
Definition: Driver.h:95
void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
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)
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.
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:65
std::shared_ptr< AdvectionSystem > AdvectionSystemSharedPtr
Shared pointer to an AdvectionSystem class.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::vector< double > d(NPUPPER *NPUPPER)
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.hpp:396
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294