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