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>
37using namespace std;
38
39#include <boost/format.hpp>
40
45
46namespace Nektar::SolverUtils
47{
50 "set-start-time", "", "Set the starting time of the simulation.");
51
54 "set-start-chknumber", "",
55 "Set the starting number of the checkpoint file.");
56
57/**
58 * @class UnsteadySystem
59 *
60 * Provides the underlying timestepping framework for unsteady solvers
61 * including the general timestepping routines. This class is not
62 * intended to be directly instantiated, but rather is a base class
63 * on which to define unsteady solvers.
64 *
65 * For details on implementing unsteady solvers see
66 * \ref sectionADRSolverModuleImplementation here
67 */
68
69/**
70 * Processes SolverInfo parameters from the session file and sets up
71 * timestepping-specific code.
72 * @param pSession Session object to read parameters from.
73 */
77 : EquationSystem(pSession, pGraph)
78
79{
80}
81
82/**
83 * Initialization object for UnsteadySystem class.
84 */
85void UnsteadySystem::v_InitObject(bool DeclareField)
86{
87 EquationSystem::v_InitObject(DeclareField);
88
89 m_initialStep = 0;
90
91 // Load SolverInfo parameters.
92 m_session->MatchSolverInfo("DIFFUSIONADVANCEMENT", "Explicit",
94 m_session->MatchSolverInfo("ADVECTIONADVANCEMENT", "Explicit",
96 m_session->MatchSolverInfo("REACTIONADVANCEMENT", "Explicit",
97 m_explicitReaction, true);
98 m_session->LoadParameter("CheckAbortSteps", m_abortSteps, 1);
99
100 // Steady state tolerance.
101 m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 0.0);
102
103 // Frequency for checking steady state.
104 m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 1);
105
106 // For steady problems, we do not initialise the time integration.
107 if (m_session->DefinesSolverInfo("TimeIntegrationMethod") ||
108 m_session->DefinesTimeIntScheme())
109 {
111 if (m_session->DefinesTimeIntScheme())
112 {
113 timeInt = m_session->GetTimeIntScheme();
114 }
115 else
116 {
117 timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
118 }
119
122 timeInt.method, timeInt.variant, timeInt.order,
123 timeInt.freeParams);
124
125 // Load generic input parameters.
126 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
127 m_session->LoadParameter("IO_FiltersInfoSteps", m_filtersInfosteps,
128 10 * m_infosteps);
129 m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
130 m_session->LoadParameter("CFLEnd", m_CFLEnd, 0.0);
131 m_session->LoadParameter("CFLGrowth", m_CFLGrowth, 1.0);
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 Returns the time integration scheme.
190 */
193{
194 return m_intScheme;
195}
196
197/**
198 * @brief Returns the time integration scheme operators.
199 */
202{
203 return m_ode;
204}
205
206/**
207 * @brief Initialises the time integration scheme (as specified in the
208 * session file), and perform the time integration.
209 */
211{
212 ASSERTL0(m_intScheme != nullptr, "No time integration scheme.");
213
214 int i = 1;
215 int nvariables = 0;
216 int nfields = m_fields.size();
217 if (m_intVariables.empty())
218 {
219 for (i = 0; i < nfields; ++i)
220 {
221 m_intVariables.push_back(i);
222 }
223 nvariables = nfields;
224 }
225 else
226 {
227 nvariables = m_intVariables.size();
228 }
229
230 // Integrate in wave-space if using homogeneous1D.
232 {
233 for (i = 0; i < nfields; ++i)
234 {
235 m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetTotPoints(),
236 m_fields[i]->GetPhys(),
237 m_fields[i]->UpdatePhys());
238 m_fields[i]->SetWaveSpace(true);
239 m_fields[i]->SetPhysState(false);
240 }
241 }
242
243 // Set up wrapper to fields data storage.
244 Array<OneD, Array<OneD, NekDouble>> fields(nvariables);
245
246 // Order storage to list time-integrated fields first.
247 for (i = 0; i < nvariables; ++i)
248 {
249 fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
250 m_fields[m_intVariables[i]]->SetPhysState(false);
251 }
252
253 // Initialise time integration scheme.
254 m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
255
256 // Initialise filters.
257 for (auto &x : m_filters)
258 {
259 x.second->Initialise(m_fields, m_time);
260 }
261
263 bool doCheckTime = false;
264 int step = m_initialStep;
265 int stepCounter = 0;
266 NekDouble intTime = 0.0;
267 NekDouble cpuTime = 0.0;
268 NekDouble cpuPrevious = 0.0;
269 NekDouble elapsed = 0.0;
270 NekDouble totFilterTime = 0.0;
271
272 m_lastCheckTime = 0.0;
273
274 Array<OneD, int> abortFlags(2, 0);
275 string abortFile = "abort";
276 if (m_session->DefinesSolverInfo("CheckAbortFile"))
277 {
278 abortFile = m_session->GetSolverInfo("CheckAbortFile");
279 }
280
281 NekDouble tmp_cflSafetyFactor = m_cflSafetyFactor;
282 while ((step < m_steps || m_time < m_fintime - NekConstants::kNekZeroTol) &&
283 abortFlags[1] == 0)
284 {
286 {
287 tmp_cflSafetyFactor =
288 min(m_CFLEnd, m_CFLGrowth * tmp_cflSafetyFactor);
289 }
290
291 // Frozen preconditioner checks.
292 if (!m_comm->IsParallelInTime())
293 {
295 {
296 m_cflSafetyFactor = tmp_cflSafetyFactor;
297
299 {
300 m_timestep = GetTimeStep(fields);
301 }
302
303 // Ensure that the final timestep finishes at the final
304 // time, or at a prescribed IO_CheckTime.
305 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
306 {
308 }
309 else if (m_checktime &&
311 {
314 doCheckTime = true;
315 }
316 }
317 }
318
319 if (m_TimeIncrementFactor > 1.0)
320 {
321 NekDouble timeincrementFactor = m_TimeIncrementFactor;
322 m_timestep *= timeincrementFactor;
323
324 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
325 {
327 }
328 }
329
330 // Perform any solver-specific pre-integration steps.
331 timer.Start();
332 if (v_PreIntegrate(step))
333 {
334 break;
335 }
336
337 ASSERTL0(m_timestep > 0, "m_timestep < 0");
338
339 fields = m_intScheme->TimeIntegrate(stepCounter, m_timestep);
340 timer.Stop();
341
343 elapsed = timer.TimePerTest(1);
344 intTime += elapsed;
345 cpuTime += elapsed;
346
347 // Write out status information.
348 v_PrintStatusInformation(step, cpuTime);
349 if (m_infosteps &&
350 m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
351 !((step + 1) % m_infosteps))
352 {
353 cpuPrevious = cpuTime;
354 cpuTime = 0.0;
355 }
356
357 // Transform data into coefficient space.
358 for (i = 0; i < nvariables; ++i)
359 {
360 // Copy fields into ExpList::m_phys.
361 m_fields[m_intVariables[i]]->SetPhys(fields[i]);
362 if (v_RequireFwdTrans())
363 {
364 m_fields[m_intVariables[i]]->FwdTransLocalElmt(
365 m_fields[m_intVariables[i]]->GetPhys(),
366 m_fields[m_intVariables[i]]->UpdateCoeffs());
367 }
368 m_fields[m_intVariables[i]]->SetPhysState(false);
369 }
370
371 // Perform any solver-specific post-integration steps.
372 if (v_PostIntegrate(step))
373 {
374 break;
375 }
376
377 // Check for steady-state.
379 (!((step + 1) % m_steadyStateSteps)))
380 {
381 if (CheckSteadyState(step, intTime))
382 {
383 if (m_comm->GetRank() == 0)
384 {
385 cout << "Reached Steady State to tolerance "
386 << m_steadyStateTol << endl;
387 }
388 break;
389 }
390 }
391
392 // Test for abort conditions (nan, or abort file).
393 if (m_abortSteps && !((step + 1) % m_abortSteps))
394 {
395 abortFlags[0] = 0;
396 for (i = 0; i < nvariables; ++i)
397 {
398 if (Vmath::Nnan(m_fields[m_intVariables[i]]->GetPhys().size(),
399 m_fields[m_intVariables[i]]->GetPhys(), 1) > 0)
400 {
401 abortFlags[0] = 1;
402 }
403 }
404
405 // Rank zero looks for abort file and deletes it
406 // if it exists. Communicates the abort.
407 if (m_session->GetComm()->GetSpaceComm()->GetRank() == 0)
408 {
409 if (fs::exists(abortFile))
410 {
411 fs::remove(abortFile);
412 abortFlags[1] = 1;
413 }
414 }
415
416 m_session->GetComm()->GetSpaceComm()->AllReduce(
417 abortFlags, LibUtilities::ReduceMax);
418
419 ASSERTL0(!abortFlags[0], "NaN found during time integration.");
420 }
421
422 // Update filters.
423 for (auto &x : m_filters)
424 {
425 timer.Start();
426 x.second->Update(m_fields, m_time);
427 timer.Stop();
428 elapsed = timer.TimePerTest(1);
429 totFilterTime += elapsed;
430
431 // Write out individual filter status information.
432 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
433 !((step + 1) % m_filtersInfosteps) && !m_filters.empty() &&
434 m_session->DefinesCmdLineArgument("verbose"))
435 {
436 stringstream s0;
437 s0 << x.first << ":";
438 stringstream s1;
439 s1 << elapsed << "s";
440 stringstream s2;
441 s2 << elapsed / cpuPrevious * 100 << "%";
442 cout << "CPU time for filter " << setw(25) << left << s0.str()
443 << setw(12) << left << s1.str() << endl
444 << "\t Percentage of time integration: " << setw(10)
445 << left << s2.str() << endl;
446 }
447 }
448
449 // Write out overall filter status information.
450 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
451 !((step + 1) % m_filtersInfosteps) && !m_filters.empty())
452 {
453 stringstream ss;
454 ss << totFilterTime << "s";
455 cout << "Total filters CPU Time:\t\t\t " << setw(10) << left
456 << ss.str() << endl;
457 }
458 totFilterTime = 0.0;
459
460 // Write out checkpoint files.
461 if ((m_checksteps && !((step + 1) % m_checksteps)) || doCheckTime)
462 {
464 {
465 // Transform to physical space for output.
466 vector<bool> transformed(nfields, false);
467 for (i = 0; i < nfields; i++)
468 {
469 if (m_fields[i]->GetWaveSpace())
470 {
471 m_fields[i]->SetWaveSpace(false);
472 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
473 m_fields[i]->UpdatePhys());
474 transformed[i] = true;
475 }
476 }
478 m_nchk++;
479
480 // Transform back to wave-space after output.
481 for (i = 0; i < nfields; i++)
482 {
483 if (transformed[i])
484 {
485 m_fields[i]->SetWaveSpace(true);
486 m_fields[i]->HomogeneousFwdTrans(
487 m_fields[i]->GetTotPoints(), m_fields[i]->GetPhys(),
488 m_fields[i]->UpdatePhys());
489 m_fields[i]->SetPhysState(false);
490 }
491 }
492 }
493 else
494 {
496 m_nchk++;
497 }
498 doCheckTime = false;
499 }
500
501 // Step advance.
502 ++step;
503 ++stepCounter;
504 }
505
506 // Print out summary statistics.
508
509 // If homogeneous, transform back into physical space if necessary.
511 {
512 for (i = 0; i < nfields; i++)
513 {
514 if (m_fields[i]->GetWaveSpace())
515 {
516 m_fields[i]->SetWaveSpace(false);
517 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
518 m_fields[i]->UpdatePhys());
519 }
520 }
521 }
522 else
523 {
524 for (i = 0; i < nvariables; ++i)
525 {
526 m_fields[m_intVariables[i]]->SetPhysState(true);
527 }
528 }
529
530 // Finalise filters.
531 for (auto &x : m_filters)
532 {
533 x.second->Finalise(m_fields, m_time);
534 }
535
536 // Print for 1D problems.
537 if (m_spacedim == 1)
538 {
540 }
541}
542
544 const NekDouble cpuTime)
545{
546 if (m_infosteps && m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
547 !((step + 1) % m_infosteps))
548 {
549 if (m_comm->IsParallelInTime())
550 {
551 cout << "RANK " << m_session->GetComm()->GetTimeComm()->GetRank()
552 << " Steps: " << setw(8) << left << step + 1 << " "
553 << "Time: " << setw(12) << left << m_time;
554 }
555 else
556 {
557 cout << "Steps: " << setw(8) << left << step + 1 << " "
558 << "Time: " << setw(12) << left << m_time;
559 }
560
562 {
563 cout << " Time-step: " << setw(12) << left << m_timestep;
564 }
565
566 stringstream ss;
567 ss << cpuTime << "s";
568 cout << " CPU Time: " << setw(8) << left << ss.str() << endl;
569 }
570}
571
573{
574 if (m_session->GetComm()->GetRank() == 0)
575 {
576 if (m_cflSafetyFactor > 0.0)
577 {
578 cout << "CFL safety factor : " << m_cflSafetyFactor << endl
579 << "CFL time-step : " << m_timestep << endl;
580 }
581
582 if (m_session->GetSolverInfo("Driver") != "SteadyState" &&
583 m_session->GetSolverInfo("Driver") != "Parareal" &&
584 m_session->GetSolverInfo("Driver") != "PFASST")
585 {
586 cout << "Time-integration : " << intTime << "s" << endl;
587 }
588 }
589}
590
591/**
592 * @brief Sets the initial conditions.
593 */
594void UnsteadySystem::v_DoInitialise(bool dumpInitialConditions)
595{
598 SetInitialConditions(m_time, dumpInitialConditions);
600}
601
602/**
603 * @brief Prints a summary with some information regards the
604 * time-stepping.
605 */
607{
609 AddSummaryItem(s, "Advect. advancement",
610 m_explicitAdvection ? "explicit" : "implicit");
611
612 AddSummaryItem(s, "Diffuse. advancement",
613 m_explicitDiffusion ? "explicit" : "implicit");
614
615 if (m_session->GetSolverInfo("EQTYPE") ==
616 "SteadyAdvectionDiffusionReaction")
617 {
618 AddSummaryItem(s, "React. advancement",
619 m_explicitReaction ? "explicit" : "implicit");
620 }
621
622 AddSummaryItem(s, "Time Step", m_timestep);
623 AddSummaryItem(s, "No. of Steps", m_steps);
624 AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
625 if (m_intScheme)
626 {
627 AddSummaryItem(s, "Integration Type", m_intScheme->GetName());
628 }
629}
630
631/**
632 * Stores the solution in a file for 1D problems only. This method has
633 * been implemented to facilitate the post-processing for 1D problems.
634 */
636{
637 // Coordinates of the quadrature points in the real physical space.
641 m_fields[0]->GetCoords(x, y, z);
642
643 // Print out the solution in a txt file.
644 ofstream outfile;
645 outfile.open("solution1D.txt");
646 for (int i = 0; i < GetNpoints(); i++)
647 {
648 outfile << scientific << setw(17) << setprecision(16) << x[i] << " "
649 << m_fields[m_intVariables[0]]->GetPhys()[i] << endl;
650 }
651 outfile << endl << endl;
652 outfile.close();
653}
654
655/**
656 *
657 */
659{
660 if (m_session->DefinesFunction("InitialConditions"))
661 {
662 for (int i = 0; i < m_fields.size(); ++i)
663 {
665
666 vType = m_session->GetFunctionType("InitialConditions",
667 m_session->GetVariable(i));
668
670 {
671 std::string filename = m_session->GetFunctionFilename(
672 "InitialConditions", m_session->GetVariable(i));
673
674 fs::path pfilename(filename);
675
676 // Redefine path for parallel file which is in directory.
677 if (fs::is_directory(pfilename))
678 {
679 fs::path metafile("Info.xml");
680 fs::path fullpath = pfilename / metafile;
681 filename = LibUtilities::PortablePath(fullpath);
682 }
685 fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
686
687 // Check to see if time defined.
689 {
690 auto iter = m_fieldMetaDataMap.find("Time");
691 if (iter != m_fieldMetaDataMap.end())
692 {
693 time = boost::lexical_cast<NekDouble>(iter->second);
694 }
695
696 iter = m_fieldMetaDataMap.find("ChkFileNum");
697 if (iter != m_fieldMetaDataMap.end())
698 {
699 nchk = boost::lexical_cast<NekDouble>(iter->second);
700 }
701 }
702
703 break;
704 }
705 }
706 }
707 if (m_session->DefinesCmdLineArgument("set-start-time"))
708 {
709 time = boost::lexical_cast<NekDouble>(
710 m_session->GetCmdLineArgument<std::string>("set-start-time")
711 .c_str());
712 }
713 if (m_session->DefinesCmdLineArgument("set-start-chknumber"))
714 {
715 nchk = boost::lexical_cast<int>(
716 m_session->GetCmdLineArgument<std::string>("set-start-chknumber"));
717 }
718 ASSERTL0(time >= 0 && nchk >= 0,
719 "Starting time and checkpoint number should be >= 0");
720}
721
722/**
723 * @brief Return the timestep to be used for the next step in the
724 * time-marching loop.
725 *
726 * @see UnsteadySystem::GetTimeStep
727 */
729 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &inarray)
730{
731 NEKERROR(ErrorUtil::efatal, "Not defined for this class");
732 return 0.0;
733}
734
735/**
736 *
737 */
738bool UnsteadySystem::v_PreIntegrate([[maybe_unused]] int step)
739{
740 return false;
741}
742
743/**
744 *
745 */
746bool UnsteadySystem::v_PostIntegrate([[maybe_unused]] int step)
747{
748 return false;
749}
750
751/**
752 *
753 */
755{
756 return true;
757}
758
759/**
760 *
761 */
763{
764 return true;
765}
766
767/**
768 *
769 */
772 StdRegions::VarCoeffMap &varCoeffMap)
773{
774 int phystot = m_fields[0]->GetTotPoints();
775 int nvel = vel.size();
776
777 Array<OneD, NekDouble> varcoeff(phystot), tmp;
778
779 // Calculate magnitude of v.
780 Vmath::Vmul(phystot, vel[0], 1, vel[0], 1, varcoeff, 1);
781 for (int n = 1; n < nvel; ++n)
782 {
783 Vmath::Vvtvp(phystot, vel[n], 1, vel[n], 1, varcoeff, 1, varcoeff, 1);
784 }
785 Vmath::Vsqrt(phystot, varcoeff, 1, varcoeff, 1);
786
787 for (int i = 0; i < m_fields[0]->GetNumElmts(); ++i)
788 {
789 int offset = m_fields[0]->GetPhys_Offset(i);
790 int nq = m_fields[0]->GetExp(i)->GetTotPoints();
791 Array<OneD, NekDouble> unit(nq, 1.0);
792
793 int nmodes = 0;
794
795 for (int n = 0; n < m_fields[0]->GetExp(i)->GetNumBases(); ++n)
796 {
797 nmodes = max(nmodes, m_fields[0]->GetExp(i)->GetBasisNumModes(n));
798 }
799
800 NekDouble h = m_fields[0]->GetExp(i)->Integral(unit);
801 h = pow(h, (NekDouble)(1.0 / nvel)) / ((NekDouble)nmodes);
802
803 Vmath::Smul(nq, h, varcoeff + offset, 1, tmp = varcoeff + offset, 1);
804 }
805
806 // Set up map with eVarCoffLaplacian key.
807 varCoeffMap[StdRegions::eVarCoeffLaplacian] = varcoeff;
808}
809
810/**
811 *
812 */
814{
815 if (m_steadyStateTol > 0.0)
816 {
817 const int nPoints = m_fields[0]->GetTotPoints();
820
821 for (int i = 0; i < m_fields.size(); ++i)
822 {
824 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1,
825 m_previousSolution[i], 1);
826 }
827
828 if (m_comm->GetRank() == 0)
829 {
830 std::string fName =
831 m_session->GetSessionName() + std::string(".resd");
832 m_errFile.open(fName.c_str());
833 m_errFile << setw(26) << left << "# Time";
834
835 m_errFile << setw(26) << left << "CPU_Time";
836
837 m_errFile << setw(26) << left << "Step";
838
839 for (int i = 0; i < m_fields.size(); ++i)
840 {
841 m_errFile << setw(26) << m_session->GetVariables()[i];
842 }
843
844 m_errFile << endl;
845 }
846 }
847}
848
849/**
850 * @brief Calculate whether the system has reached a steady state by
851 * observing residuals to a user-defined tolerance.
852 */
853bool UnsteadySystem::CheckSteadyState(int step, const NekDouble &totCPUTime)
854{
855 const int nPoints = GetTotPoints();
856 const int nFields = m_fields.size();
857
858 // Holds L2 errors.
859 Array<OneD, NekDouble> L2(nFields);
860
861 SteadyStateResidual(step, L2);
862
863 if (m_infosteps && m_comm->GetRank() == 0 &&
864 (((step + 1) % m_infosteps == 0) || ((step == m_initialStep))))
865 {
866 // Output time.
867 m_errFile << boost::format("%25.19e") % m_time;
868
869 m_errFile << " " << boost::format("%25.19e") % totCPUTime;
870
871 int stepp = step + 1;
872
873 m_errFile << " " << boost::format("%25.19e") % stepp;
874
875 // Output residuals.
876 for (int i = 0; i < nFields; ++i)
877 {
878 m_errFile << " " << boost::format("%25.19e") % L2[i];
879 }
880
881 m_errFile << endl;
882 }
883
884 // Calculate maximum L2 error.
885 NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
886
887 if (m_infosteps && m_session->DefinesCmdLineArgument("verbose") &&
888 m_comm->GetRank() == 0 && ((step + 1) % m_infosteps == 0))
889 {
890 cout << "-- Maximum L^2 residual: " << maxL2 << endl;
891 }
892
893 if (maxL2 <= m_steadyStateTol)
894 {
895 return true;
896 }
897
898 for (int i = 0; i < m_fields.size(); ++i)
899 {
900 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
901 1);
902 }
903
904 return false;
905}
906
907/**
908 *
909 */
910void UnsteadySystem::v_SteadyStateResidual([[maybe_unused]] int step,
912{
913 const int nPoints = GetTotPoints();
914 const int nFields = m_fields.size();
915
916 // Holds L2 errors.
917 Array<OneD, NekDouble> RHSL2(nFields);
918 Array<OneD, NekDouble> residual(nFields);
919 Array<OneD, NekDouble> reference(nFields);
920
921 for (int i = 0; i < nFields; ++i)
922 {
923 Array<OneD, NekDouble> tmp(nPoints);
924
925 Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
926 1, tmp, 1);
927 Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
928 residual[i] = Vmath::Vsum(nPoints, tmp, 1);
929
931 tmp, 1);
932 reference[i] = Vmath::Vsum(nPoints, tmp, 1);
933 }
934
935 m_comm->GetSpaceComm()->AllReduce(residual, LibUtilities::ReduceSum);
936 m_comm->GetSpaceComm()->AllReduce(reference, LibUtilities::ReduceSum);
937
938 // L2 error.
939 for (int i = 0; i < nFields; ++i)
940 {
941 reference[i] = (reference[i] == 0) ? 1 : reference[i];
942 L2[i] = sqrt(residual[i] / reference[i]);
943 }
944}
945
946/**
947 *
948 */
950 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
952 [[maybe_unused]] const NekDouble time)
953{
954
955 if (&inarray != &outarray)
956 {
957 for (int i = 0; i < inarray.size(); ++i)
958 {
959 Vmath::Vcopy(GetNpoints(), inarray[i], 1, outarray[i], 1);
960 }
961 }
962}
963
964} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
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:224
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
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.
Binds a set of functions for use by time integration schemes.
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 describing how to solve specific equations.
int m_spacedim
Spatial dimension (>= expansion dim).
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.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
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.
SOLVER_UTILS_EXPORT void v_DoInitialise(bool dumpInitialConditions=true) override
Sets up initial conditions.
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck()
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.
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.
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 bool v_PostIntegrate(int step)
NekDouble m_CFLEnd
Maximun cfl in cfl growth.
SOLVER_UTILS_EXPORT ~UnsteadySystem() override
Destructor.
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics(const NekDouble intTime)
Print Summary Statistics.
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation(const int step, const NekDouble cpuTime)
Print Status Information.
bool CheckSteadyState(int step, const NekDouble &totCPUTime=0.0)
Calculate whether the system has reached a steady state by observing residuals to a user-defined tole...
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.
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)
SOLVER_UTILS_EXPORT void AppendOutput1D(void)
Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime(NekDouble &time, int &nchk)
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans()
SOLVER_UTILS_EXPORT void SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperators & GetTimeIntegrationSchemeOperators()
Returns the time integration scheme operators.
int m_steadyStateSteps
Check for steady state at step interval.
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtr & GetTimeIntegrationScheme()
Returns the time integration scheme.
int m_filtersInfosteps
Number of time steps between outputting filters information.
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...
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()
SOLVER_UTILS_EXPORT void v_DoSolve() override
Solves an unsteady problem.
TimeIntegrationSchemeFactory & GetTimeIntegrationSchemeFactory()
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
static const NekDouble kNekZeroTol
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39
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, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:346
std::vector< double > z(NPUPPER)
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.hpp:340
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.hpp:72
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.hpp:366
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608
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.hpp:100
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
Definition: Vmath.hpp:743
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.hpp:644
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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.hpp:220
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
std::vector< NekDouble > freeParams
Definition: SessionReader.h:87