Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Incompressible Navier Stokes solver
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
38 
39 namespace Nektar
40 {
41 namespace SolverUtils
42 {
43 
46  "SteadyState", DriverSteadyState::create);
49  "Driver","SteadyState",0);
50 
51 /**
52  *
53  */
56  : DriverModifiedArnoldi(pSession)
57 {
58 }
59 
60 
61 /**
62  *
63  */
65 {
66 }
67 
68 
69 /**
70  *
71  */
73 {
75 }
76 
77 
78 void DriverSteadyState::v_Execute(ostream &out)
79 
80 {
81  // With a loop over "DoSolve", this Driver implements the
82  // "encaplulated" Selective Frequency Damping method
83  // to find the steady state of a flow above the critical Reynolds
84  // number.
85  m_equ[m_nequ - 1]->PrintSummary(out);
86  m_equ[m_nequ - 1]->DoInitialise();
87 
88  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 1000);
89  m_session->LoadParameter("IO_CheckSteps", m_checksteps, 100000);
90  m_session->LoadParameter("ControlCoeff",m_X, 1);
91  m_session->LoadParameter("FilterWidth", m_Delta, 1);
92 
93  // To evaluate optimum SFD parameters if growth rate provided in the
94  // xml file
95  m_session->LoadParameter("GrowthRateEV", GrowthRateEV, 0.0);
96 
97  // To evaluate optimum SFD parameters if frequency provided in the xml
98  // file
99  m_session->LoadParameter("FrequencyEV", FrequencyEV, 0.0);
100 
101  // Determine when SFD method is converged
102  m_session->LoadParameter("TOL", TOL, 1.0e-08);
103 
104  // Used only for the Adaptive SFD method
105  m_session->LoadParameter("AdaptiveTOL", AdaptiveTOL, 1.0e-02);
106 
107  // Used only for the Adaptive SFD method
108  m_session->LoadParameter("AdaptiveTime", AdaptiveTime, 25.0);
109 
110  if (m_comm->GetRank() == 0)
111  {
112  PrintSummarySFD();
113  }
114 
115  timer.Start();
116 
117  // Definition of shared pointer (used only for the Adaptive SFD method)
119 
120  // Condition necessary to run SFD for the compressible case
121  NumVar_SFD = m_equ[m_nequ - 1]->UpdateFields()[0]->GetCoordim(0);
122  if (m_session->GetSolverInfo("EqType") == "EulerCFE" ||
123  m_session->GetSolverInfo("EqType") == "NavierStokesCFE")
124  {
125  // Number of variables for the compressible equations
126  NumVar_SFD += 2;
127  }
128  if(m_session->DefinesSolverInfo("HOMOGENEOUS"))
129  {
130  if (m_session->GetSolverInfo("HOMOGENEOUS") == "1D")
131  {
132  NumVar_SFD += 1;
133  }
134  }
135 
136  // We store the time step
137  m_dt = m_equ[m_nequ - 1]->GetTimeStep();
138 
139  // Evaluate optimum SFD parameters if dominent EV given by xml file
140  if (GrowthRateEV != 0.0 && FrequencyEV != 0.0)
141  {
142  cout << "Besed on the dominant EV given in the xml file,"
143  << "a 1D model is used to evaluate the optumum parameters"
144  << "of the SFD method:" << endl;
145  complex<NekDouble> EV = polar(exp(GrowthRateEV), FrequencyEV);
147  }
148 
149 
150  // We set up the elements of the operator of the encapsulated
151  // formulation of the selective frequencive damping method
153 
154  // m_steps is set to 1. Then "m_equ[m_nequ - 1]->DoSolve()" will run
155  // for only one time step
156  m_equ[m_nequ - 1]->SetStepsToOne();
157  ofstream m_file("ConvergenceHistory.txt", ios::out | ios::trunc);
158 
159  Array<OneD, Array<OneD, NekDouble> > q0(NumVar_SFD);
160  Array<OneD, Array<OneD, NekDouble> > q1(NumVar_SFD);
161  Array<OneD, Array<OneD, NekDouble> > qBar0(NumVar_SFD);
162  Array<OneD, Array<OneD, NekDouble> > qBar1(NumVar_SFD);
163 
164  for(int i = 0; i < NumVar_SFD; ++i)
165  {
166  q0[i] = Array<OneD, NekDouble> (m_equ[m_nequ-1]->GetTotPoints(),
167  0.0); //q0 is initialised
168  qBar0[i] = Array<OneD, NekDouble> (m_equ[m_nequ-1]->GetTotPoints(),
169  0.0); //qBar0 initially set to zero
170  m_equ[m_nequ - 1]->CopyFromPhysField(i, qBar0[i]);
171  }
172 
173  ///Definition of variables used in this algorithm
174  m_stepCounter = 0;
175  m_Check = 0;
176  m_Check_BaseFlow = 1;
178  Diff_q_qBar = 1.0;
179  Diff_q1_q0 = 1.0;
180  cpuTime = 0.0;
181  elapsed = 0.0;
182  totalTime = 0.0;
183  FlowPartiallyConverged = false;
184 
185  while (max(Diff_q_qBar, Diff_q1_q0) > TOL)
186  {
187  ///Call the Navier-Stokes solver for one time step
188  m_equ[m_nequ - 1]->DoSolve();
189 
190  for(int i = 0; i < NumVar_SFD; ++i)
191  {
192  ///Copy the current flow field into q0
193  m_equ[m_nequ - 1]->CopyFromPhysField(i, q0[i]);
194 
195  ///Apply the linear operator to the outcome of the solver
196  ComputeSFD(i, q0, qBar0, q1, qBar1);
197 
198  ///Update qBar
199  qBar0[i] = qBar1[i];
200 
201  ///Copy the output of the SFD method into the current flow field
202  m_equ[m_nequ - 1]->CopyToPhysField(i, q1[i]);
203  }
204 
206  {
208 
209  ///Loop for the adaptive SFD method
211  FlowPartiallyConverged == false)
212  {
213 
215 
216  if (Diff_q_qBar < AdaptiveTOL)
217  {
218  if (m_comm->GetRank() == 0)
219  {
220  cout << "\n\t The SFD method is converging: we compute "
221  << "stability analysis using the 'partially "
222  << "converged' steady state as base flow:\n" << endl;
223  }
224 
225  m_equ[m_nequ - 1]->Checkpoint_BaseFlow(m_Check_BaseFlow);
227 
228  A->GetAdvObject()->SetBaseFlow(q0);
230 
231  if (m_comm->GetRank() == 0)
232  {
233  // Compute the update of the SFD parameters only on
234  // one processor
236  }
237  else
238  {
239  // On all the other processors, the parameters are set
240  // to 0
241  m_X = 0;
242  m_Delta = 0;
243  }
244  // The we give to all the processors the value of X and
245  // Delta of the first processor
248 
249  FlowPartiallyConverged = true;
250  }
252  >= AdaptiveTime)
253  {
254  if (m_comm->GetRank() == 0)
255  {
256  cout << "\n\t We compute stability analysis using"
257  << " the current flow field as base flow:\n"
258  << endl;
259  }
260 
261  m_equ[m_nequ - 1]->Checkpoint_BaseFlow(m_Check_BaseFlow);
263 
264  A->GetAdvObject()->SetBaseFlow(q0);
266 
267  if (m_comm->GetRank() == 0)
268  {
269  // Compute the update of the SFD parameters only on
270  // one processor
272  }
273  else
274  {
275  // On all the other processors, the parameters are set
276  // to 0
277  m_X = 0;
278  m_Delta = 0;
279  }
280  // The we give to all the processors the value of X and
281  // Delta of the first processor
284 
286  }
287  }
288  }
289 
291  {
292  m_Check++;
293  m_equ[m_nequ - 1]->Checkpoint_Output(m_Check);
294  }
295  m_stepCounter++;
296  }
297 
298  m_file.close();
299 
300  ///We save the final solution into a .fld file
301  m_equ[m_nequ - 1]->Output();
302 }
303 
304 
305 /**
306  * This routine defines the encapsulated SFD operator with first-order
307  * splitting and exact resolution of the second subproblem.
308  *(See http://scitation.aip.org/content/aip/journal/pof2/26/3/10.1063/1.4867482 for details)
309  */
311  const NekDouble Delta_input)
312 {
313  NekDouble X = X_input*m_dt;
314  NekDouble Delta = Delta_input/m_dt;
315  NekDouble coeff = 1.0/(1.0 + X*Delta);
316  M11 = coeff*(1.0 + X*Delta*exp(-(X + 1.0/Delta)));
317  M12 = coeff*(X*Delta*(1.0 - exp(-(X + 1.0/Delta))));
318  M21 = coeff*(1.0 - exp(-(X + 1.0/Delta)));
319  M22 = coeff*(X*Delta + exp(-(X + 1.0/Delta)));
320 }
321 
322 
324  const Array<OneD, const Array<OneD, NekDouble> > &q0,
325  const Array<OneD, const Array<OneD, NekDouble> > &qBar0,
326  Array<OneD, Array<OneD, NekDouble> > &q1,
327  Array<OneD, Array<OneD, NekDouble> > &qBar1)
328 {
329  q1[i] = Array<OneD, NekDouble> (m_equ[m_nequ - 1]->GetTotPoints(),0.0);
330  qBar1[i] = Array<OneD, NekDouble> (m_equ[m_nequ - 1]->GetTotPoints(),0.0);
331 
332  ///Encapsulated SFD method
333  Vmath::Svtvp(q1[i].num_elements(), M11, q0[i], 1, q1[i], 1, q1[i], 1 );
334  Vmath::Svtvp(q1[i].num_elements(), M12, qBar0[i], 1, q1[i], 1, q1[i], 1 );
335 
336  Vmath::Svtvp(qBar1[i].num_elements(), M21, q0[i], 1, qBar1[i], 1,
337  qBar1[i], 1 );
338  Vmath::Svtvp(qBar1[i].num_elements(), M22, qBar0[i], 1, qBar1[i], 1,
339  qBar1[i], 1 );
340 }
341 
342 
344 {
345  NekDouble growthEV(0.0);
346  NekDouble frequencyEV(0.0);
347 
348  // m_kdim is the dimension of Krylov subspace (defined in the xml file and
349  // used in DriverArnoldi.cpp)
350  ReadEVfile(m_kdim, growthEV, frequencyEV);
351 
352  cout << "\n\tgrowthEV = " << growthEV << endl;
353  cout << "\tfrequencyEV = " << frequencyEV << endl;
354 
355  complex<NekDouble> ApproxEV = polar(exp(growthEV), frequencyEV);
356 
357  NekDouble X_new = m_X;
358  NekDouble Delta_new = m_Delta;
359 
360  GradientDescentMethod(ApproxEV, X_new, Delta_new);
361 
362  m_X = X_new;
363  m_Delta = Delta_new;
364 
366 }
367 
368 /**
369  * This routine implements a gradient descent method to find the parameters X
370  * end Delta which give the minimum eigenlavue of the SFD problem applied to
371  * the scalar case u(n+1) = \alpha*u(n).
372  */
374  const complex<NekDouble> &alpha,
375  NekDouble &X_output,
376  NekDouble &Delta_output)
377 {
378  cout << "\n\tWe enter the Gradient Descent Method [...]" << endl;
379  bool OptParmFound = false;
380  bool Descending = true;
381  NekDouble X_input = X_output;
382  NekDouble Delta_input = Delta_output;
383 
384  NekDouble X_init = X_output;
385  NekDouble Delta_init = Delta_output;
386  int stepCounter(0);
387 
388  NekDouble F0(0.0);
389  NekDouble F0x(0.0);
390  NekDouble F0y(0.0);
391  NekDouble F1(0.0);
392  NekDouble dx = 0.00000001;
393  NekDouble dirx(0.0);
394  NekDouble diry(0.0);
395  NekDouble s(0.0);
396  NekDouble CurrentMin = 1.0;
397 
398  while (OptParmFound == false)
399  {
400  Descending = true;
401  EvalEV_ScalarSFD(X_input, Delta_input, alpha, F0);
402  EvalEV_ScalarSFD(X_input + dx, Delta_input, alpha, F0x);
403  EvalEV_ScalarSFD(X_input, Delta_input + dx, alpha, F0y);
404  dirx = (F0 - F0x);
405  diry = (F0 - F0y);
406  s = abs(0.000001/dirx);
407  X_output = X_input + s*dirx;
408  Delta_output = Delta_input + s*diry;
409  F1 = F0;
410 
411  while (Descending == true)
412  {
413  CurrentMin = F1;
414  X_input = X_output;
415  Delta_input = Delta_output;
416  EvalEV_ScalarSFD(X_output, Delta_output, alpha, F1);
417 
418  if (F1 > CurrentMin)
419  {
420  Descending = false;
421  }
422  else
423  {
424  s = s + s*0.01;
425  X_output = X_input + s*dirx;
426  Delta_output = Delta_input + s*diry;
427  }
428 
429  if(stepCounter > 9999999)
430  {
431  //We are stuck in this loop..
432  //Then we restart it with different initail conditions
433  Descending = false;
434  X_input = X_init;
435  Delta_init = Delta_init + Delta_init*0.1;
436  Delta_input = Delta_init;
437  stepCounter = 0;
438  }
439  stepCounter++;
440  }
441 
442  if (abs(F0-F1) < dx)
443  {
444  cout << "\tThe Gradient Descent Method has converged!" << endl;
445  EvalEV_ScalarSFD(X_output, Delta_output, alpha, F1);
446  cout << "\n\tThe updated parameters are: X_tilde = " << X_output
447  << " and Delta_tilde = " << Delta_output << endl;
448  OptParmFound = true;
449  }
450 
451  }
452 }
453 
454 
455 /**
456  * This routine evaluates the maximum eigenvalue of the SFD system when applied
457  * to the 1D model u(n+1) = alpha*u(n)
458  */
460  const NekDouble &X_input,
461  const NekDouble &Delta_input,
462  const complex<NekDouble> &alpha,
463  NekDouble &MaxEV)
464 {
465  NekDouble A11 = ( 1.0 + X_input * Delta_input *
466  exp(-(X_input + 1.0/Delta_input)) )/(1.0 + X_input*Delta_input);
467  NekDouble A12 = ( X_input*Delta_input - X_input * Delta_input *
468  exp(-(X_input + 1.0/Delta_input)) )/(1.0 + X_input*Delta_input);
469  NekDouble A21 = ( 1.0 - 1.0 *
470  exp(-(X_input + 1.0/Delta_input)) )/(1 + X_input*Delta_input);
471  NekDouble A22 = ( X_input*Delta_input + 1.0 *
472  exp(-(X_input + 1.0/Delta_input)) )/(1.0 + X_input*Delta_input);
473 
474  complex<NekDouble> B11 = alpha;
475  NekDouble B12 = 0.0;
476  NekDouble B21 = 0.0;
477  NekDouble B22 = 1.0;
478 
479  complex<NekDouble> a = A11*B11 + A12*B21;
480  NekDouble b = A11*B12 + A12*B22;
481  complex<NekDouble> c = A21*B11 + A22*B21;
482  NekDouble d = A21*B12 + A22*B22;
483 
484  complex<NekDouble> delt = sqrt((a-d)*(a-d) + 4.0*b*c);
485  complex<NekDouble> lambda_1 = 0.5*(a+d + delt);
486  complex<NekDouble> lambda_2 = 0.5*(a+d - delt);
487 
488  NekDouble NORM_1 = abs(lambda_1);
489  NekDouble NORM_2 = abs(lambda_2);
490 
491  MaxEV = max(NORM_1, NORM_2);
492 }
493 
494 
496  int &KrylovSubspaceDim,
497  NekDouble &growthEV,
498  NekDouble &frequencyEV)
499 {
500  // This routine reads the .evl file written by the Arnoldi algorithm
501  // (written in September 2014)
502  std::string EVfileName = m_session->GetSessionName() + ".evl";
503  std::ifstream EVfile(EVfileName.c_str());
504 
505  int NumLinesInFile(0);
506  NekDouble NonReleventNumber(0.0);
507 
508  if(EVfile)
509  {
510  std::string line;
511 
512  // This block counts the total number of lines of the .evl file
513  // We keep going util we reach the end of the file
514  while(getline(EVfile, line))
515  {
516  NumLinesInFile += 1;
517  }
518  EVfile.close();
519 
520  // It may happen that the Stability method that have produced the .elv
521  // file converges in less than m_kdim iterations. In this case,
522  // KrylovSubspaceDim has to be changed here
523  if(NumLinesInFile < KrylovSubspaceDim*2.0
524  + KrylovSubspaceDim*(KrylovSubspaceDim+1.0)/2.0)
525  {
526  for(int i = 1; i <= KrylovSubspaceDim; ++i)
527  {
528  if(NumLinesInFile == i*2.0 + i*(i+1.0)/2.0)
529  {
530  KrylovSubspaceDim = i;
531  }
532  }
533  }
534 
535  // go back to the beginning of the file
536  std::ifstream EVfile(EVfileName.c_str());
537 
538  // We now want to go to the line where the most unstable eigenlavue was
539  // written
540  for(int i = 0; i < (NumLinesInFile - KrylovSubspaceDim); ++i)
541  {
542  std::getline(EVfile, line);
543  }
544 
545  // Then we read this line by skipping the first three values written
546  EVfile >> NonReleventNumber;
547  EVfile >> NonReleventNumber;
548  EVfile >> NonReleventNumber;
549 
550  // The growth rate and the frequency of the EV are at the 4th and 5th
551  // colums of the .evl file
552  EVfile >> growthEV;
553  EVfile >> frequencyEV;
554  }
555  else
556  {
557  cout << "An error occured when openning the .evl file" << endl;
558  }
559  EVfile.close();
560 }
561 
562 
563 /**
564  * This routine evaluates |q-qBar|_inf (and |q1-q0|_inf) and writes the values
565  * in "ConvergenceHistory.txt"
566  */
568  const Array<OneD, const Array<OneD, NekDouble> > &qBar1,
569  const Array<OneD, const Array<OneD, NekDouble> > &q0,
570  NekDouble &MaxNormDiff_q_qBar,
571  NekDouble &MaxNormDiff_q1_q0)
572 {
573  Array<OneD, NekDouble > NormDiff_q_qBar(NumVar_SFD, 1.0);
574  Array<OneD, NekDouble > NormDiff_q1_q0(NumVar_SFD, 1.0);
575  MaxNormDiff_q_qBar=0.0;
576  MaxNormDiff_q1_q0=0.0;
577 
578  for(int i = 0; i < NumVar_SFD; ++i)
579  {
580  NormDiff_q_qBar[i] = m_equ[m_nequ - 1]->LinfError(i, qBar1[i]);
581  NormDiff_q1_q0[i] = m_equ[m_nequ - 1]->LinfError(i, q0[i]);
582 
583  if (MaxNormDiff_q_qBar < NormDiff_q_qBar[i])
584  {
585  MaxNormDiff_q_qBar = NormDiff_q_qBar[i];
586  }
587  if (MaxNormDiff_q1_q0 < NormDiff_q1_q0[i])
588  {
589  MaxNormDiff_q1_q0 = NormDiff_q1_q0[i];
590  }
591  }
592 
593  timer.Stop();
595  cpuTime += elapsed;
596  totalTime += elapsed;
597 
598  if (m_comm->GetRank() == 0)
599  {
600  cout << "SFD - Step: " << left << m_stepCounter+1
601  << ";\tTime: " << left << m_equ[m_nequ - 1]->GetFinalTime()
602  << ";\tCPU time = " << left << cpuTime << " s"
603  << ";\tTot time = " << left << totalTime << " s"
604  << ";\tX = " << left << m_X
605  << ";\tDelta = " << left << m_Delta
606  << ";\t|q-qBar|inf = " << left << MaxNormDiff_q_qBar << endl;
607  std::ofstream m_file( "ConvergenceHistory.txt", std::ios_base::app);
608  m_file << m_stepCounter+1 << "\t"
609  << m_equ[m_nequ - 1]->GetFinalTime() << "\t"
610  << totalTime << "\t"
611  << MaxNormDiff_q_qBar << "\t"
612  << MaxNormDiff_q1_q0 << endl;
613  m_file.close();
614  }
615 
616  cpuTime = 0.0;
617  timer.Start();
618 }
619 
620 
622 {
623  cout << "\n====================================="
624  "=================================="
625  << endl;
626  cout << "Parameters for the SFD method:" << endl;
627  cout << "\tControl Coefficient: X = " << m_X << endl;
628  cout << "\tFilter Width: Delta = " << m_Delta << endl;
629  cout << "The simulation is stopped when |q-qBar|inf < " << TOL << endl;
631  {
632  cout << "\nWe use the adaptive SFD method:" << endl;
633  cout << " The parameters are updated every " << AdaptiveTime
634  << " time units;" << endl;
635  cout << " until |q-qBar|inf becomes smaller than " << AdaptiveTOL
636  << endl;
637  }
638  cout << "====================================="
639  "==================================\n" << endl;
640 }
641 
642 }
643 }
644