Nektar++
UnsteadySystem.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadySystem.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: Generic timestepping for Unsteady solvers
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 
36 #include <iostream>
37 #include <iomanip>
38 using namespace std;
39 
40 #include <boost/core/ignore_unused.hpp>
41 #include <boost/format.hpp>
42 
47 
48 namespace Nektar
49 {
50  namespace SolverUtils
51  {
52  /**
53  * @class UnsteadySystem
54  *
55  * Provides the underlying timestepping framework for unsteady solvers
56  * including the general timestepping routines. This class is not
57  * intended to be directly instantiated, but rather is a base class
58  * on which to define unsteady solvers.
59  *
60  * For details on implementing unsteady solvers see
61  * \ref sectionADRSolverModuleImplementation here
62  */
63 
64  /**
65  * Processes SolverInfo parameters from the session file and sets up
66  * timestepping-specific code.
67  * @param pSession Session object to read parameters from.
68  */
69  UnsteadySystem::UnsteadySystem(
72  : EquationSystem(pSession, pGraph),
73  m_infosteps(10)
74 
75  {
76  }
77 
78  /**
79  * Initialization object for UnsteadySystem class.
80  */
82  {
84 
85  m_initialStep = 0;
86 
87  // Load SolverInfo parameters
88  m_session->MatchSolverInfo("DIFFUSIONADVANCEMENT","Explicit",
89  m_explicitDiffusion,true);
90  m_session->MatchSolverInfo("ADVECTIONADVANCEMENT","Explicit",
91  m_explicitAdvection,true);
92  m_session->MatchSolverInfo("REACTIONADVANCEMENT", "Explicit",
93  m_explicitReaction, true);
94 
95  m_session->MatchSolverInfo("FLAGIMPLICITITSSTATISTICS", "True",
97 
98  m_session->LoadParameter("CheckAbortSteps", m_abortSteps, 1);
99  // Steady state tolerance
100  m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 0.0);
101  // Frequency for checking steady state
102  m_session->LoadParameter("SteadyStateSteps",
103  m_steadyStateSteps, 1);
104 
105  // For steady problems, we do not initialise the time integration
106  if (m_session->DefinesSolverInfo("TimeIntegrationMethod") ||
107  m_session->DefinesTimeIntScheme())
108  {
110  if (m_session->DefinesTimeIntScheme())
111  {
112  timeInt = m_session->GetTimeIntScheme();
113  }
114  else
115  {
116  timeInt.method = m_session->GetSolverInfo(
117  "TimeIntegrationMethod");
118  }
119 
121  .CreateInstance(timeInt.method,
122  timeInt.variant,
123  timeInt.order,
124  timeInt.freeParams);
125 
126  // Load generic input parameters
127  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
128  m_session->LoadParameter("IO_FiltersInfoSteps",
130  m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
131  m_session->LoadParameter("CFLEnd", m_CFLEnd, 0.0);
132  m_session->LoadParameter("CFLGrowth", m_CFLGrowth, 1.0);
133  m_session->LoadParameter("CFLGrowth", m_CFLGrowth, 1.0);
134 
135 
136  // Time tolerance between filter update time and time integration
137  m_session->LoadParameter("FilterTimeWarning",
139 
140  // Ensure that there is no conflict of parameters
141  if(m_cflSafetyFactor > 0.0)
142  {
143  // Check final condition
144  ASSERTL0(m_fintime == 0.0 || m_steps == 0,
145  "Final condition not unique: "
146  "fintime > 0.0 and Nsteps > 0");
147  // Check timestep condition
148  ASSERTL0(m_timestep == 0.0,
149  "Timestep not unique: timestep > 0.0 & CFL > 0.0");
150  }
151  else
152  {
153  ASSERTL0(m_timestep != 0.0,
154  "Need to set either TimeStep or CFL");
155  }
156 
157  // Ensure that there is no conflict of parameters
158  if (m_CFLGrowth > 1.0)
159  {
160  // Check final condition
162  "m_CFLEnd >= m_cflSafetyFactor required");
163  }
164 
165  // Set up time to be dumped in field information
166  m_fieldMetaDataMap["Time"] =
167  boost::lexical_cast<std::string>(m_time);
168  }
169 
170  // By default attempt to forward transform initial condition.
171  m_homoInitialFwd = true;
172 
173  // Set up filters
174  for (auto &x : m_session->GetFilters())
175  {
176  m_filters.push_back(make_pair(x.first, GetFilterFactory().CreateInstance(
177  x.first, m_session, shared_from_this(), x.second)));
178  }
179  }
180 
181  /**
182  * Destructor for the class UnsteadyAdvection.
183  */
185  {
186  }
187 
188  /**
189  * @brief Returns the maximum time estimator for CFL control.
190  */
192  {
193  return m_intScheme->GetTimeStability();
194  }
195 
196  /**
197  * @brief Initialises the time integration scheme (as specified in the
198  * session file), and perform the time integration.
199  */
201  {
202  ASSERTL0(m_intScheme != 0, "No time integration scheme.");
203 
204  int i = 1;
205  int nvariables = 0;
206  int nfields = m_fields.size();
207 
208  if (m_intVariables.empty())
209  {
210  for (i = 0; i < nfields; ++i)
211  {
212  m_intVariables.push_back(i);
213  }
214  nvariables = nfields;
215  }
216  else
217  {
218  nvariables = m_intVariables.size();
219  }
220 
221  // Integrate in wave-space if using homogeneous1D
223  {
224  for(i = 0; i < nfields; ++i)
225  {
226  m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetPhys(),
227  m_fields[i]->UpdatePhys());
228  m_fields[i]->SetWaveSpace(true);
229  m_fields[i]->SetPhysState(false);
230  }
231  }
232 
233  // Set up wrapper to fields data storage.
234  Array<OneD, Array<OneD, NekDouble> > fields(nvariables);
235 
236  // Order storage to list time-integrated fields first.
237  for(i = 0; i < nvariables; ++i)
238  {
239  fields[i] = m_fields[m_intVariables[i]]->GetPhys();
240  m_fields[m_intVariables[i]]->SetPhysState(false);
241  }
242 
243  // Initialise time integration scheme
244  m_intScheme->InitializeScheme( m_timestep, fields, m_time, m_ode );
245 
246  // Initialise filters
247  for( auto &x : m_filters )
248  {
249  x.second->Initialise(m_fields, m_time);
250  }
251 
252  LibUtilities::Timer timer;
253  bool doCheckTime = false;
254  int step = m_initialStep;
255  int stepCounter = 0;
256  int restartStep = -1;
257  NekDouble intTime = 0.0;
258  NekDouble cpuTime = 0.0;
259  NekDouble cpuPrevious = 0.0;
260  NekDouble elapsed = 0.0;
261  NekDouble totFilterTime = 0.0;
262 
263  m_lastCheckTime = 0.0;
264 
265  m_TotNewtonIts = 0;
266  m_TotLinIts = 0;
267  m_TotImpStages = 0;
268 
269  Array<OneD, int> abortFlags(2, 0);
270  string abortFile = "abort";
271  if (m_session->DefinesSolverInfo("CheckAbortFile"))
272  {
273  abortFile = m_session->GetSolverInfo("CheckAbortFile");
274  }
275 
276  NekDouble tmp_cflSafetyFactor = m_cflSafetyFactor;
277 
279  while ((step < m_steps ||
281  abortFlags[1] == 0)
282  {
283  restartStep++;
284 
286  {
287  tmp_cflSafetyFactor =
288  min(m_CFLEnd,m_CFLGrowth*tmp_cflSafetyFactor);
289  }
290 
291  m_flagUpdatePreconMat = true;
292 
293  // Flag to update AV
294  m_CalcPhysicalAV = true;
295  // Frozen preconditioner checks
296  if (UpdateTimeStepCheck())
297  {
298  m_cflSafetyFactor = tmp_cflSafetyFactor;
299 
300  if (m_cflSafetyFactor)
301  {
302  m_timestep = GetTimeStep(fields);
303  }
304 
305  // Ensure that the final timestep finishes at the final
306  // time, or at a prescribed IO_CheckTime.
307  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
308  {
310  }
311  else if (m_checktime &&
313  {
316  doCheckTime = true;
317  }
318  }
319 
320  if (m_TimeIncrementFactor > 1.0)
321  {
322  NekDouble timeincrementFactor = m_TimeIncrementFactor;
323  m_timestep *= timeincrementFactor;
324 
325  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
326  {
328  }
329  }
330 
331  // Perform any solver-specific pre-integration steps
332  timer.Start();
333  if (v_PreIntegrate(step))
334  {
335  break;
336  }
337 
338  m_StagesPerStep = 0;
339  m_TotLinItePerStep = 0;
340 
341  ASSERTL0(m_timestep > 0, "m_timestep < 0");
342 
343  fields =
344  m_intScheme->TimeIntegrate( stepCounter, m_timestep, m_ode);
345  timer.Stop();
346 
347  m_time += m_timestep;
348  elapsed = timer.TimePerTest(1);
349  intTime += elapsed;
350  cpuTime += elapsed;
351 
352  // Write out status information
353  if (m_session->GetComm()->GetRank() == 0 &&
354  !((step+1) % m_infosteps))
355  {
356  cout << "Steps: " << setw(8) << left << step+1 << " "
357  << "Time: " << setw(12) << left << m_time;
358 
359  if (m_cflSafetyFactor)
360  {
361  cout << " Time-step: " << setw(12)
362  << left << m_timestep;
363  }
364 
365  stringstream ss;
366  ss << cpuTime << "s";
367  cout << " CPU Time: " << setw(8) << left
368  << ss.str() << endl;
369  cpuPrevious = cpuTime;
370  cpuTime = 0.0;
371 
373  {
374  cout
375  << " &&"
376  << " TotImpStages= " << m_TotImpStages
377  << " TotNewtonIts= " << m_TotNewtonIts
378  << " TotLinearIts = " << m_TotLinIts
379  << endl;
380  }
381  }
382 
383  // Transform data into coefficient space
384  for (i = 0; i < nvariables; ++i)
385  {
386  // copy fields into ExpList::m_phys and assign the new
387  // array to fields
388  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
389  fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
390  if( v_RequireFwdTrans() )
391  {
392  m_fields[m_intVariables[i]]->FwdTrans_IterPerExp(
393  fields[i],
394  m_fields[m_intVariables[i]]->UpdateCoeffs());
395  }
396  m_fields[m_intVariables[i]]->SetPhysState(false);
397  }
398 
399  // Perform any solver-specific post-integration steps
400  if (v_PostIntegrate(step))
401  {
402  break;
403  }
404 
405  // Check for steady-state
406  if (m_steadyStateTol > 0.0 && (!((step+1)%m_steadyStateSteps)) )
407  {
408  if (CheckSteadyState(step,intTime))
409  {
410  if (m_comm->GetRank() == 0)
411  {
412  cout << "Reached Steady State to tolerance "
413  << m_steadyStateTol << endl;
414  }
415  break;
416  }
417  }
418 
419  // test for abort conditions (nan, or abort file)
420  if (m_abortSteps && !((step+1) % m_abortSteps) )
421  {
422  abortFlags[0] = 0;
423  for (i = 0; i < nvariables; ++i)
424  {
425  if (Vmath::Nnan(fields[i].size(),
426  fields[i], 1) > 0)
427  {
428  abortFlags[0] = 1;
429  }
430  }
431 
432  //rank zero looks for abort file and deltes it
433  //if it exists. The communicates the abort
434  if(m_session->GetComm()->GetRank() == 0)
435  {
436  if(boost::filesystem::exists(abortFile))
437  {
438  boost::filesystem::remove(abortFile);
439  abortFlags[1] = 1;
440  }
441  }
442 
443  m_session->GetComm()->AllReduce(abortFlags,
445 
446  ASSERTL0 (!abortFlags[0],
447  "NaN found during time integration.");
448  }
449 
450  // Update filters
451  for (auto &x : m_filters)
452  {
453  timer.Start();
454  x.second->Update(m_fields, m_time);
455  timer.Stop();
456  elapsed = timer.TimePerTest(1);
457  totFilterTime += elapsed;
458 
459  // Write out individual filter status information
460  if(m_session->GetComm()->GetRank() == 0 &&
461  !((step+1) % m_filtersInfosteps) && !m_filters.empty() &&
462  m_session->DefinesCmdLineArgument("verbose"))
463  {
464  stringstream s0;
465  s0 << x.first << ":";
466  stringstream s1;
467  s1 << elapsed << "s";
468  stringstream s2;
469  s2 << elapsed / cpuPrevious * 100 << "%";
470  cout << "CPU time for filter " << setw(25) << left
471  << s0.str() << setw(12) << left << s1.str() <<
472  endl << "\t Percentage of time integration: "
473  << setw(10) << left << s2.str() << endl;
474  }
475  }
476 
477  // Write out overall filter status information
478  if (m_session->GetComm()->GetRank() == 0 &&
479  !((step+1) % m_filtersInfosteps) && !m_filters.empty())
480  {
481  stringstream ss;
482  ss << totFilterTime << "s";
483  cout << "Total filters CPU Time:\t\t\t " << setw(10)
484  << left << ss.str() << endl;
485  }
486  totFilterTime = 0.0;
487 
488  // Write out checkpoint files
489  if ((m_checksteps && !((step + 1) % m_checksteps)) ||
490  doCheckTime)
491  {
493  {
494  vector<bool> transformed(nfields, false);
495  for(i = 0; i < nfields; i++)
496  {
497  if (m_fields[i]->GetWaveSpace())
498  {
499  m_fields[i]->SetWaveSpace(false);
500  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
501  m_fields[i]->UpdatePhys());
502  m_fields[i]->SetPhysState(true);
503  transformed[i] = true;
504  }
505  }
507  m_nchk++;
508  for(i = 0; i < nfields; i++)
509  {
510  if (transformed[i])
511  {
512  m_fields[i]->SetWaveSpace(true);
513  m_fields[i]->HomogeneousFwdTrans(
514  m_fields[i]->GetPhys(),
515  m_fields[i]->UpdatePhys());
516  m_fields[i]->SetPhysState(false);
517  }
518  }
519  }
520  else
521  {
523  m_nchk++;
524  }
525  doCheckTime = false;
526  }
527 
528  // Step advance
529  ++step;
530  ++stepCounter;
531  }
532 
533  // Print out summary statistics
534  if (m_session->GetComm()->GetRank() == 0)
535  {
536  if (m_cflSafetyFactor > 0.0)
537  {
538  cout << "CFL safety factor : " << m_cflSafetyFactor << endl
539  << "CFL time-step : " << m_timestep << endl;
540  }
541 
542  if (m_session->GetSolverInfo("Driver") != "SteadyState")
543  {
544  cout << "Time-integration : " << intTime << "s" << endl;
545  }
546 
548  {
549  cout
550  << "-------------------------------------------" << endl
551  << "Total Implicit Stages: " << m_TotImpStages << endl
552  << "Total Newton Its : " << m_TotNewtonIts << endl
553  << "Total Linear Its : " << m_TotLinIts << endl
554  << "-------------------------------------------" << endl;
555  }
556  }
557 
558  // If homogeneous, transform back into physical space if necessary.
560  {
561  for(i = 0; i < nfields; i++)
562  {
563  if (m_fields[i]->GetWaveSpace())
564  {
565  m_fields[i]->SetWaveSpace(false);
566  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
567  m_fields[i]->UpdatePhys());
568  m_fields[i]->SetPhysState(true);
569  }
570  }
571  }
572  else
573  {
574  for(i = 0; i < nvariables; ++i)
575  {
576  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
577  m_fields[m_intVariables[i]]->SetPhysState(true);
578  }
579  }
580 
581  // Finalise filters
582  for (auto &x : m_filters)
583  {
584  x.second->Finalise(m_fields, m_time);
585  }
586 
587  // Print for 1D problems
588  if(m_spacedim == 1)
589  {
590  v_AppendOutput1D(fields);
591  }
592  }
593 
594  /**
595  * @brief Sets the initial conditions.
596  */
598  {
603  }
604 
605  /**
606  * @brief Prints a summary with some information regards the
607  * time-stepping.
608  */
610  {
612  AddSummaryItem(s, "Advection",
613  m_explicitAdvection ? "explicit" : "implicit");
614 
615  if(m_session->DefinesSolverInfo("AdvectionType"))
616  {
617  AddSummaryItem(s, "AdvectionType",
618  m_session->GetSolverInfo("AdvectionType"));
619  }
620 
621  AddSummaryItem(s, "Diffusion",
622  m_explicitDiffusion ? "explicit" : "implicit");
623 
624  if (m_session->GetSolverInfo("EQTYPE")
625  == "SteadyAdvectionDiffusionReaction")
626  {
627  AddSummaryItem(s, "Reaction",
628  m_explicitReaction ? "explicit" : "implicit");
629  }
630 
631  AddSummaryItem( s, "Time Step", m_timestep );
632  AddSummaryItem( s, "No. of Steps", m_steps );
633  AddSummaryItem( s, "Checkpoints (steps)", m_checksteps );
634  AddSummaryItem( s, "Integration Type", m_intScheme->GetName() );
635  }
636 
637  /**
638  * Stores the solution in a file for 1D problems only. This method has
639  * been implemented to facilitate the post-processing for 1D problems.
640  */
642  Array<OneD, Array<OneD, NekDouble> > &solution1D)
643  {
644  // Coordinates of the quadrature points in the real physical space
648  m_fields[0]->GetCoords(x, y, z);
649 
650  // Print out the solution in a txt file
651  ofstream outfile;
652  outfile.open("solution1D.txt");
653  for(int i = 0; i < GetNpoints(); i++)
654  {
655  outfile << scientific << setw (17) << setprecision(16) << x[i]
656  << " " << solution1D[0][i] << endl;
657  }
658  outfile << endl << endl;
659  outfile.close();
660  }
661 
663  {
664  if (m_session->DefinesFunction("InitialConditions"))
665  {
666  for (int i = 0; i < m_fields.size(); ++i)
667  {
669 
670  vType = m_session->GetFunctionType(
671  "InitialConditions", m_session->GetVariable(i));
672 
673  if (vType == LibUtilities::eFunctionTypeFile)
674  {
675  std::string filename
676  = m_session->GetFunctionFilename(
677  "InitialConditions", m_session->GetVariable(i));
678 
679  fs::path pfilename(filename);
680 
681  // redefine path for parallel file which is in directory
682  if(fs::is_directory(pfilename))
683  {
684  fs::path metafile("Info.xml");
685  fs::path fullpath = pfilename / metafile;
686  filename = LibUtilities::PortablePath(fullpath);
687  }
690  m_session, filename);
691  fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
692 
693  // check to see if time defined
694  if (m_fieldMetaDataMap !=
696  {
697  auto iter = m_fieldMetaDataMap.find("Time");
698  if (iter != m_fieldMetaDataMap.end())
699  {
700  time = boost::lexical_cast<NekDouble>(
701  iter->second);
702  }
703 
704  iter = m_fieldMetaDataMap.find("ChkFileNum");
705  if (iter != m_fieldMetaDataMap.end())
706  {
707  nchk = boost::lexical_cast<NekDouble>(
708  iter->second);
709  }
710  }
711 
712  break;
713  }
714  }
715  }
716  }
717 
718  /**
719  * @brief Return the timestep to be used for the next step in the
720  * time-marching loop.
721  *
722  * This function can be overloaded to facilitate solver which utilise a
723  * CFL (or other) parameter to determine a maximum timestep under which
724  * the problem remains stable.
725  */
727  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
728  {
729  return v_GetTimeStep(inarray);
730  }
731 
732  /**
733  * @brief Return the timestep to be used for the next step in the
734  * time-marching loop.
735  *
736  * @see UnsteadySystem::GetTimeStep
737  */
739  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
740  {
741  boost::ignore_unused(inarray);
742  NEKERROR(ErrorUtil::efatal, "Not defined for this class");
743  return 0.0;
744  }
745 
747  {
748  boost::ignore_unused(step);
749  return false;
750  }
751 
753  {
754  boost::ignore_unused(step);
755  return false;
756  }
757 
759  const Array<OneD, Array<OneD, NekDouble> > vel,
760  StdRegions::VarCoeffMap &varCoeffMap)
761  {
762  int phystot = m_fields[0]->GetTotPoints();
763  int nvel = vel.size();
764 
765  Array<OneD, NekDouble> varcoeff(phystot),tmp;
766 
767  // calculate magnitude of v
768  Vmath::Vmul(phystot,vel[0],1,vel[0],1,varcoeff,1);
769  for(int n = 1; n < nvel; ++n)
770  {
771  Vmath::Vvtvp(phystot,vel[n],1,vel[n],1,varcoeff,1,varcoeff,1);
772  }
773  Vmath::Vsqrt(phystot,varcoeff,1,varcoeff,1);
774 
775  for(int i = 0; i < m_fields[0]->GetNumElmts(); ++i)
776  {
777  int offset = m_fields[0]->GetPhys_Offset(i);
778  int nq = m_fields[0]->GetExp(i)->GetTotPoints();
779  Array<OneD, NekDouble> unit(nq,1.0);
780 
781  int nmodes = 0;
782 
783  for(int n = 0; n < m_fields[0]->GetExp(i)->GetNumBases(); ++n)
784  {
785  nmodes = max(nmodes,
786  m_fields[0]->GetExp(i)->GetBasisNumModes(n));
787  }
788 
789  NekDouble h = m_fields[0]->GetExp(i)->Integral(unit);
790  h = pow(h,(NekDouble) (1.0/nvel))/((NekDouble) nmodes);
791 
792  Vmath::Smul(nq,h,varcoeff+offset,1,tmp = varcoeff+offset,1);
793  }
794 
795  // set up map with eVarCoffLaplacian key
796  varCoeffMap[StdRegions::eVarCoeffLaplacian] = varcoeff;
797  }
798 
800  {
801  if (m_steadyStateTol > 0.0)
802  {
803  const int nPoints = m_fields[0]->GetTotPoints();
805  m_fields.size());
806 
807  for (int i = 0; i < m_fields.size(); ++i)
808  {
810  Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1,
811  m_previousSolution[i], 1);
812  }
813 
814  if (m_comm->GetRank() == 0)
815  {
816  std::string fName = m_session->GetSessionName() +
817  std::string(".resd");
818  m_errFile.open(fName.c_str());
819  m_errFile << setw(26) << left << "# Time";
820 
821  m_errFile << setw(26) << left << "CPU_Time";
822 
823  m_errFile << setw(26) << left << "Step";
824 
825  for (int i = 0; i < m_fields.size(); ++i)
826  {
827  m_errFile << setw(26) << m_session->GetVariables()[i];
828  }
829 
830  m_errFile << endl;
831  }
832  }
833  }
834 
835  /**
836  * @brief Calculate whether the system has reached a steady state by
837  * observing residuals to a user-defined tolerance.
838  */
840  {
841  return CheckSteadyState(step,0.0);
842  }
843 
844  bool UnsteadySystem::CheckSteadyState(int step, NekDouble totCPUTime)
845  {
846  const int nPoints = GetTotPoints();
847  const int nFields = m_fields.size();
848 
849  // Holds L2 errors.
850  Array<OneD, NekDouble> L2 (nFields);
851 
852  SteadyStateResidual(step, L2);
853 
854  if (m_comm->GetRank() == 0 && ( ((step+1) % m_infosteps == 0)||((step== m_initialStep)) ))
855  {
856  // Output time
857  m_errFile << boost::format("%25.19e") % m_time;
858 
859  m_errFile << " "<< boost::format("%25.19e") % totCPUTime;
860 
861  int stepp = step +1;
862 
863  m_errFile << " "<< boost::format("%25.19e") % stepp;
864 
865  // Output residuals
866  for (int i = 0; i < nFields; ++i)
867  {
868  m_errFile << " " << boost::format("%25.19e") % L2[i];
869  }
870 
871  m_errFile << endl;
872  }
873 
874  // Calculate maximum L2 error
875  NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
876 
877  if (m_session->DefinesCmdLineArgument("verbose") &&
878  m_comm->GetRank() == 0 && ((step+1) % m_infosteps == 0))
879  {
880  cout << "-- Maximum L^2 residual: " << maxL2 << endl;
881  }
882 
883  if (maxL2 <= m_steadyStateTol)
884  {
885  return true;
886  }
887 
888  for (int i = 0; i < m_fields.size(); ++i)
889  {
890  Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1,
891  m_previousSolution[i], 1);
892  }
893 
895  m_steadyStateRes = maxL2;
896 
897  return false;
898  }
899 
901  int step,
903  {
904  boost::ignore_unused(step);
905  const int nPoints = GetTotPoints();
906  const int nFields = m_fields.size();
907 
908  // Holds L2 errors.
909  Array<OneD, NekDouble> RHSL2 (nFields);
910  Array<OneD, NekDouble> residual (nFields);
911  Array<OneD, NekDouble> reference(nFields);
912 
913  for (int i = 0; i < nFields; ++i)
914  {
915  Array<OneD, NekDouble> tmp(nPoints);
916 
917  Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1,
918  m_previousSolution[i], 1, tmp, 1);
919  Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
920  residual[i] = Vmath::Vsum(nPoints, tmp, 1);
921 
922  Vmath::Vmul(nPoints, m_previousSolution[i], 1,
923  m_previousSolution[i], 1, tmp, 1);
924  reference[i] = Vmath::Vsum(nPoints, tmp, 1);
925  }
926 
927  m_comm->AllReduce(residual , LibUtilities::ReduceSum);
928  m_comm->AllReduce(reference, LibUtilities::ReduceSum);
929 
930  // L2 error
931  for (int i = 0; i < nFields; ++i)
932  {
933  reference[i] = (reference[i] == 0) ? 1 : reference[i];
934  L2[i] = sqrt(residual[i] / reference[i]);
935  }
936  }
937  }
938 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:226
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition: Timer.cpp:62
A base class for describing how to solve specific equations.
int m_spacedim
Spatial dimension (>= expansion dim).
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
NekDouble m_fintime
Finish time of the simulation.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
NekDouble m_checktime
Time between checkpoints.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()
NekDouble m_timestepMax
Time step size.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
int m_initialStep
Number of the step where the simulation should begin.
enum HomogeneousType m_HomogeneousType
SOLVER_UTILS_EXPORT int GetNpoints()
int m_nchk
Number of checkpoints written so far.
int m_steps
Number of steps to take.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &l)
Virtual function for generating summary information.
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Initialisation object for EquationSystem.
SOLVER_UTILS_EXPORT int GetTotPoints()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_CFLGrowth
CFL growth rate.
bool m_explicitReaction
Indicates if explicit or implicit treatment of reaction is used.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Sets up initial conditions.
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
Storage for previous solution for steady-state check.
virtual SOLVER_UTILS_EXPORT void v_DoSolve()
Solves an unsteady problem.
bool m_CalcPhysicalAV
flag to update artificial viscosity
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem()
Destructor.
NekDouble m_filterTimeWarning
Number of time steps between outputting status information.
NekDouble m_steadyStateTol
Tolerance to which steady state should be evaluated at.
int m_abortSteps
Number of steps between checks for abort conditions.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator()
Get the maximum timestep estimator for cfl control.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble >> &inarray)
Return the timestep to be used for the next step in the time-marching loop.
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans()
bool CheckSteadyState(int step)
Calculate whether the system has reached a steady state by observing residuals to a user-defined tole...
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate(int step)
NekDouble m_CFLEnd
maximun cfl in cfl growth
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff(const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h t...
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime(NekDouble &time, int &nchk)
SOLVER_UTILS_EXPORT void SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
int m_steadyStateSteps
Check for steady state at step interval.
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck()
int m_filtersInfosteps
Number of time steps between outputting filters information.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D(Array< OneD, Array< OneD, NekDouble >> &solution1D)
Print the solution at each solution point in a txt file.
int m_infosteps
Number of time steps between outputting status information.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
TimeIntegrationSchemeFactory & GetTimeIntegrationSchemeFactory()
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:306
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
static const NekDouble kNekZeroTol
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:475
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:192
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:513
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:846
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
Definition: Vmath.cpp:1015
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.cpp:892
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:372
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267
std::vector< NekDouble > freeParams
Definition: SessionReader.h:90