Nektar++
Loading...
Searching...
No Matches
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), SolverUtils::ALEHelper()
78
79{
80}
81
82/**
83 * Initialization object for UnsteadySystem class.
84 */
85void UnsteadySystem::v_InitObject(bool DeclareField)
86{
87 EquationSystem::v_InitObject(DeclareField);
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 std::cout << "TimeIntegrationMethod is deprecated, please use "
118 "TIMEINTEGRATIONSCHEME";
119 timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
120 }
121
124 timeInt.method, timeInt.variant, timeInt.order,
125 timeInt.freeParams);
126
127 // Load generic input parameters.
128 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
129 m_session->LoadParameter("IO_FiltersInfoSteps", m_filtersInfosteps,
130 10 * m_infosteps);
131 m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
132 m_session->LoadParameter("CFLEnd", m_CFLEnd, 0.0);
133 m_session->LoadParameter("CFLGrowth", m_CFLGrowth, 1.0);
134
135 // Ensure that there is no conflict of parameters.
136 if (m_cflSafetyFactor > 0.0)
137 {
138 // Check final condition.
139 ASSERTL0(m_fintime == 0.0 || m_steps == 0,
140 "Final condition not unique: "
141 "fintime > 0.0 and Nsteps > 0");
142 // Check timestep condition.
143 ASSERTL0(m_timestep == 0.0,
144 "Timestep not unique: timestep > 0.0 & CFL > 0.0");
145 }
146 else
147 {
148 ASSERTL0(m_timestep != 0.0, "Need to set either TimeStep or CFL");
149 }
150
151 // Ensure that there is no conflict of parameters.
152 if (m_CFLGrowth > 1.0)
153 {
154 // Check final condition.
156 "m_CFLEnd >= m_cflSafetyFactor required");
157 }
158
159 // Set up time to be dumped in field information.
160 m_fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(m_time);
161 }
162
163 // By default attempt to forward transform initial condition.
164 m_homoInitialFwd = true;
165
166 // Set up filters.
167 for (auto &x : m_session->GetFilters())
168 {
169 m_filters.push_back(make_pair(
170 x.name, GetFilterFactory().CreateInstance(
171 x.name, m_session, shared_from_this(), x.params)));
172 }
173}
174
175/**
176 * @brief Returns the maximum time estimator for CFL control.
177 */
179{
180 return m_intScheme->GetTimeStability();
181}
182
183/**
184 * @brief Returns the time integration scheme.
185 */
191
192/**
193 * @brief Returns the time integration scheme operators.
194 */
200
201/**
202 * @brief Initialises the time integration scheme (as specified in the
203 * session file), and perform the time integration.
204 */
206{
207 ASSERTL0(m_intScheme != nullptr, "No time integration scheme.");
208
209 int i = 1;
210 int nvariables = 0;
211 int nfields = m_fields.size();
212 if (m_intVariables.empty())
213 {
214 for (i = 0; i < nfields; ++i)
215 {
216 m_intVariables.push_back(i);
217 }
218 nvariables = nfields;
219 }
220 else
221 {
222 nvariables = m_intVariables.size();
223 }
224
225 // Integrate in wave-space if using homogeneous1D.
227 {
228 for (i = 0; i < nfields; ++i)
229 {
230 m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetTotPoints(),
231 m_fields[i]->GetPhys(),
232 m_fields[i]->UpdatePhys());
233 m_fields[i]->SetWaveSpace(true);
234 m_fields[i]->SetPhysState(false);
235 }
236 }
237
238 // Set up wrapper to fields data storage.
239 Array<OneD, Array<OneD, NekDouble>> fields(nvariables);
240
241 // Order storage to list time-integrated fields first.
242 // @TODO: Refactor to take coeffs (FwdTrans) if boolean flag (in constructor
243 // function) says to.
244 for (i = 0; i < nvariables; ++i)
245 {
246 fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
247 m_fields[m_intVariables[i]]->SetPhysState(false);
248 }
249
250 // @TODO: Virtual function that allows to transform the field space, embed
251 // the MultiplyMassMatrix in here.
252 // @TODO: Specify what the fields variables are physical or coefficient,
253 // boolean in UnsteadySystem class...
254
255 if (m_meshDistorted)
256 {
257 v_ALEPreMultiplyMass(fields);
258 }
259
260 // Initialise time integration scheme.
261 m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
262
263 // Initialise filters.
264 for (auto &x : m_filters)
265 {
266 x.second->Initialise(m_fields, m_time);
267 }
268
270 bool doCheckTime = false;
271 int step = m_initialStep;
272 int stepCounter = 0;
273 NekDouble intTime = 0.0;
274 NekDouble cpuTime = 0.0;
275 NekDouble cpuPrevious = 0.0;
276 NekDouble elapsed = 0.0;
277 NekDouble totFilterTime = 0.0;
278
279 m_lastCheckTime = 0.0;
280
281 Array<OneD, int> abortFlags(2, 0);
282 string abortFile = "abort";
283 if (m_session->DefinesSolverInfo("CheckAbortFile"))
284 {
285 abortFile = m_session->GetSolverInfo("CheckAbortFile");
286 }
287
288 NekDouble tmp_cflSafetyFactor = m_cflSafetyFactor;
289 while ((step < m_steps || m_time < m_fintime - NekConstants::kNekZeroTol) &&
290 abortFlags[1] == 0)
291 {
293 {
294 tmp_cflSafetyFactor =
295 min(m_CFLEnd, m_CFLGrowth * tmp_cflSafetyFactor);
296 }
297
298 // Frozen preconditioner checks.
299 if (!m_comm->IsParallelInTime())
300 {
302 {
303 m_cflSafetyFactor = tmp_cflSafetyFactor;
304
306 {
307 m_timestep = GetTimeStep(fields);
308 }
309
310 // Ensure that the final timestep finishes at the final
311 // time, or at a prescribed IO_CheckTime.
312 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
313 {
315 }
316 else if (m_checktime &&
318 {
321 doCheckTime = true;
322 }
323 }
324 }
325
326 if (m_TimeIncrementFactor > 1.0)
327 {
328 NekDouble timeincrementFactor = m_TimeIncrementFactor;
329 m_timestep *= timeincrementFactor;
330
331 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
332 {
334 }
335 }
336
337 // Perform any solver-specific pre-integration steps.
338 timer.Start();
339 if (v_PreIntegrate(
340 step)) // Could be possible to put a preintegrate step in the
341 // ALEHelper class, put in the Unsteady Advection class
342 {
343 break;
344 }
345
346 ASSERTL0(m_timestep > 0, "m_timestep < 0");
347
348 fields = m_intScheme->TimeIntegrate(stepCounter, m_timestep);
349 timer.Stop();
350
352 elapsed = timer.TimePerTest(1);
353 intTime += elapsed;
354 cpuTime += elapsed;
355
356 // Write out status information.
357 v_PrintStatusInformation(step, cpuTime);
358 if (m_infosteps &&
359 m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
360 !((step + 1) % m_infosteps))
361 {
362 cpuPrevious = cpuTime;
363 cpuTime = 0.0;
364 }
365
366 // Change to advect coeffs
367 if (m_meshDistorted)
368 {
371 }
372 else
373 {
374 // Transform data into coefficient space
375 for (i = 0; i < nvariables; ++i)
376 {
377 // copy fields into ExpList::m_phys and assign the new
378 // array to fields
379 m_fields[m_intVariables[i]]->SetPhys(fields[i]);
380 fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
381 if (v_RequireFwdTrans())
382 {
383 if (m_comm->IsParallelInTime())
384 {
385 m_fields[m_intVariables[i]]->FwdTrans(
386 m_fields[m_intVariables[i]]->GetPhys(),
387 m_fields[m_intVariables[i]]->UpdateCoeffs());
388 }
389 else
390 {
391 m_fields[m_intVariables[i]]->FwdTransLocalElmt(
392 m_fields[m_intVariables[i]]->GetPhys(),
393 m_fields[m_intVariables[i]]->UpdateCoeffs());
394 }
395 }
396 m_fields[m_intVariables[i]]->SetPhysState(false);
397 }
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.
408 (!((step + 1) % m_steadyStateSteps)))
409 {
410 if (CheckSteadyState(step, intTime))
411 {
412 if (m_comm->GetRank() == 0)
413 {
414 cout << "Reached Steady State to tolerance "
415 << m_steadyStateTol << endl;
416 }
417 break;
418 }
419 }
420
421 // Test for abort conditions (nan, or abort file).
422 if (m_abortSteps && !((step + 1) % m_abortSteps))
423 {
424 abortFlags[0] = 0;
425 for (i = 0; i < nvariables; ++i)
426 {
427 if (Vmath::Nnan(m_fields[m_intVariables[i]]->GetPhys().size(),
428 m_fields[m_intVariables[i]]->GetPhys(), 1) > 0)
429 {
430 abortFlags[0] = 1;
431 }
432 }
433
434 // Rank zero looks for abort file and deletes it
435 // if it exists. Communicates the abort.
436 if (m_session->GetComm()->GetSpaceComm()->GetRank() == 0)
437 {
438 if (fs::exists(abortFile))
439 {
440 fs::remove(abortFile);
441 abortFlags[1] = 1;
442 }
443 }
444
445 m_session->GetComm()->GetSpaceComm()->AllReduce(
446 abortFlags, LibUtilities::ReduceMax);
447
448 ASSERTL0(!abortFlags[0], "NaN found during time integration.");
449 }
450
451 // Update filters.
452 for (auto &x : m_filters)
453 {
454 timer.Start();
455 x.second->Update(m_fields, m_time);
456 timer.Stop();
457 elapsed = timer.TimePerTest(1);
458 totFilterTime += elapsed;
459
460 // Write out individual filter status information.
461 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
462 !((step + 1) % m_filtersInfosteps) && !m_filters.empty() &&
463 m_session->DefinesCmdLineArgument("verbose"))
464 {
465 stringstream s0;
466 s0 << x.first << ":";
467 stringstream s1;
468 s1 << elapsed << "s";
469 stringstream s2;
470 s2 << elapsed / cpuPrevious * 100 << "%";
471 cout << "CPU time for filter " << setw(25) << left << s0.str()
472 << setw(12) << left << s1.str() << endl
473 << "\t Percentage of time integration: " << setw(10)
474 << left << s2.str() << endl;
475 }
476 }
477
478 // Write out overall filter status information.
479 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
480 !((step + 1) % m_filtersInfosteps) && !m_filters.empty())
481 {
482 stringstream ss;
483 ss << totFilterTime << "s";
484 cout << "Total filters CPU Time:\t\t\t " << setw(10) << left
485 << ss.str() << endl;
486 }
487 totFilterTime = 0.0;
488
489 // Write out checkpoint files.
490 if ((m_checksteps && !((step + 1) % m_checksteps)) || doCheckTime)
491 {
493 {
494 // Transform to physical space for output.
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 transformed[i] = true;
504 }
505 }
507 m_nchk++;
508
509 // Transform back to wave-space after output.
510 for (i = 0; i < nfields; i++)
511 {
512 if (transformed[i])
513 {
514 m_fields[i]->SetWaveSpace(true);
515 m_fields[i]->HomogeneousFwdTrans(
516 m_fields[i]->GetTotPoints(), m_fields[i]->GetPhys(),
517 m_fields[i]->UpdatePhys());
518 m_fields[i]->SetPhysState(false);
519 }
520 }
521 }
522 else
523 {
525 m_nchk++;
526 }
527 doCheckTime = false;
528 }
529
530 // Step advance.
531 ++step;
532 ++stepCounter;
533 }
534
535 // Print out summary statistics.
537
538 // If homogeneous, transform back into physical space if necessary.
539 if (!m_meshDistorted)
540 {
542 {
543 for (i = 0; i < nfields; i++)
544 {
545 if (m_fields[i]->GetWaveSpace())
546 {
547 m_fields[i]->SetWaveSpace(false);
548 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
549 m_fields[i]->UpdatePhys());
550 }
551 }
552 }
553 else
554 {
555 for (i = 0; i < nvariables; ++i)
556 {
557 m_fields[m_intVariables[i]]->SetPhysState(true);
558 }
559 }
560 }
561 // Finalise filters.
562 for (auto &x : m_filters)
563 {
564 x.second->Finalise(m_fields, m_time);
565 }
566
567 // Print for 1D problems.
568 if (m_spacedim == 1)
569 {
571 }
572}
573
575 const NekDouble cpuTime)
576{
577 if (m_infosteps && m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
578 !((step + 1) % m_infosteps))
579 {
580 if (m_comm->IsParallelInTime())
581 {
582 cout << "RANK " << m_session->GetComm()->GetTimeComm()->GetRank()
583 << " Steps: " << setw(8) << left << step + 1 << " "
584 << "Time: " << setw(12) << left << m_time;
585 }
586 else
587 {
588 cout << "Steps: " << setw(8) << left << step + 1 << " "
589 << "Time: " << setw(12) << left << m_time;
590 }
591
593 {
594 cout << " Time-step: " << setw(12) << left << m_timestep;
595 }
596
597 stringstream ss;
598 ss << cpuTime << "s";
599 cout << " CPU Time: " << setw(8) << left << ss.str() << endl;
600 }
601}
602
604{
605 if (m_session->GetComm()->GetRank() == 0)
606 {
607 if (m_cflSafetyFactor > 0.0)
608 {
609 cout << "CFL safety factor : " << m_cflSafetyFactor << endl
610 << "CFL time-step : " << m_timestep << endl;
611 }
612
613 if (m_session->GetSolverInfo("Driver") != "SteadyState" &&
614 m_session->GetSolverInfo("Driver") != "Parareal" &&
615 m_session->GetSolverInfo("Driver") != "PFASST")
616 {
617 cout << "Time-integration : " << intTime << "s" << endl;
618 }
619 }
620}
621
622/**
623 * @brief Sets the initial conditions.
624 */
625void UnsteadySystem::v_DoInitialise(bool dumpInitialConditions)
626{
629 SetInitialConditions(m_time, dumpInitialConditions);
630
632
634}
635
636/**
637 * @brief Prints a summary with some information regards the
638 * time-stepping.
639 */
641{
643 AddSummaryItem(s, "Advect. advancement",
644 m_explicitAdvection ? "explicit" : "implicit");
645
646 AddSummaryItem(s, "Diffuse. advancement",
647 m_explicitDiffusion ? "explicit" : "implicit");
648
649 if (m_session->GetSolverInfo("EQTYPE") ==
650 "SteadyAdvectionDiffusionReaction")
651 {
652 AddSummaryItem(s, "React. advancement",
653 m_explicitReaction ? "explicit" : "implicit");
654 }
655
656 AddSummaryItem(s, "Time Step", m_timestep);
657 AddSummaryItem(s, "No. of Steps", m_steps);
658 AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
659 if (m_intScheme)
660 {
661 AddSummaryItem(s, "Integration Type", m_intScheme->GetName());
662 }
663}
664
665/**
666 * Stores the solution in a file for 1D problems only. This method has
667 * been implemented to facilitate the post-processing for 1D problems.
668 */
670{
671 // Coordinates of the quadrature points in the real physical space.
675 m_fields[0]->GetCoords(x, y, z);
676
677 // Print out the solution in a txt file.
678 ofstream outfile;
679 outfile.open("solution1D.txt");
680 for (int i = 0; i < GetNpoints(); i++)
681 {
682 outfile << scientific << setw(17) << setprecision(16) << x[i] << " "
683 << m_fields[m_intVariables[0]]->GetPhys()[i] << endl;
684 }
685 outfile << endl << endl;
686 outfile.close();
687}
688
689/**
690 *
691 */
693{
694 if (m_session->DefinesFunction("InitialConditions"))
695 {
696 for (int i = 0; i < m_fields.size(); ++i)
697 {
699
700 vType = m_session->GetFunctionType("InitialConditions",
701 m_session->GetVariable(i));
702
704 {
705 std::string filename = m_session->GetFunctionFilename(
706 "InitialConditions", m_session->GetVariable(i));
707
708 fs::path pfilename(filename);
709
710 // Redefine path for parallel file which is in directory.
711 if (fs::is_directory(pfilename))
712 {
713 fs::path metafile("Info.xml");
714 fs::path fullpath = pfilename / metafile;
715 filename = LibUtilities::PortablePath(fullpath);
716 }
719 fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
720
721 // Check to see if time defined.
723 {
724 auto iter = m_fieldMetaDataMap.find("Time");
725 if (iter != m_fieldMetaDataMap.end())
726 {
727 time = std::stod(iter->second);
728 }
729
730 iter = m_fieldMetaDataMap.find("ChkFileNum");
731 if (iter != m_fieldMetaDataMap.end())
732 {
733 nchk = std::stod(iter->second);
734 }
735 }
736
737 break;
738 }
739 }
740 }
741 if (m_session->DefinesCmdLineArgument("set-start-time"))
742 {
743 time = std::stod(
744 m_session->GetCmdLineArgument<std::string>("set-start-time")
745 .c_str());
746 }
747 if (m_session->DefinesCmdLineArgument("set-start-chknumber"))
748 {
749 nchk = std::stoi(
750 m_session->GetCmdLineArgument<std::string>("set-start-chknumber"));
751 }
752 ASSERTL0(time >= 0 && nchk >= 0,
753 "Starting time and checkpoint number should be >= 0");
754}
755
756/**
757 * @brief Return the timestep to be used for the next step in the
758 * time-marching loop.
759 *
760 * @see UnsteadySystem::GetTimeStep
761 */
763 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &inarray)
764{
765 NEKERROR(ErrorUtil::efatal, "Not defined for this class");
766 return 0.0;
767}
768
769/**
770 *
771 */
772bool UnsteadySystem::v_PreIntegrate([[maybe_unused]] int step)
773{
774 return false;
775}
776
777/**
778 *
779 */
780bool UnsteadySystem::v_PostIntegrate([[maybe_unused]] int step)
781{
782 return false;
783}
784
785/**
786 *
787 */
789{
790 return true;
791}
792
793/**
794 *
795 */
797{
798 return true;
799}
800
801/**
802 *
803 */
806 StdRegions::VarCoeffMap &varCoeffMap)
807{
808 int phystot = m_fields[0]->GetTotPoints();
809 int nvel = vel.size();
810
811 Array<OneD, NekDouble> varcoeff(phystot), tmp;
812
813 // Calculate magnitude of v.
814 Vmath::Vmul(phystot, vel[0], 1, vel[0], 1, varcoeff, 1);
815 for (int n = 1; n < nvel; ++n)
816 {
817 Vmath::Vvtvp(phystot, vel[n], 1, vel[n], 1, varcoeff, 1, varcoeff, 1);
818 }
819 Vmath::Vsqrt(phystot, varcoeff, 1, varcoeff, 1);
820
821 for (int i = 0; i < m_fields[0]->GetNumElmts(); ++i)
822 {
823 int offset = m_fields[0]->GetPhys_Offset(i);
824 int nq = m_fields[0]->GetExp(i)->GetTotPoints();
825 Array<OneD, NekDouble> unit(nq, 1.0);
826
827 int nmodes = 0;
828
829 for (int n = 0; n < m_fields[0]->GetExp(i)->GetNumBases(); ++n)
830 {
831 nmodes = max(nmodes, m_fields[0]->GetExp(i)->GetBasisNumModes(n));
832 }
833
834 NekDouble h = m_fields[0]->GetExp(i)->Integral(unit);
835 h = pow(h, (NekDouble)(1.0 / nvel)) / ((NekDouble)nmodes);
836
837 Vmath::Smul(nq, h, varcoeff + offset, 1, tmp = varcoeff + offset, 1);
838 }
839
840 // Set up map with eVarCoffLaplacian key.
841 varCoeffMap[StdRegions::eVarCoeffLaplacian] = varcoeff;
842}
843
844/**
845 *
846 */
848{
849 if (m_steadyStateTol > 0.0)
850 {
851 const int nPoints = m_fields[0]->GetTotPoints();
854
855 for (int i = 0; i < m_fields.size(); ++i)
856 {
858 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1,
859 m_previousSolution[i], 1);
860 }
861
862 if (m_comm->GetRank() == 0)
863 {
864 std::string fName =
865 m_session->GetSessionName() + std::string(".resd");
866 m_errFile.open(fName.c_str());
867 m_errFile << setw(26) << left << "# Time";
868
869 m_errFile << setw(26) << left << "CPU_Time";
870
871 m_errFile << setw(26) << left << "Step";
872
873 for (int i = 0; i < m_fields.size(); ++i)
874 {
875 m_errFile << setw(26) << m_session->GetVariables()[i];
876 }
877
878 m_errFile << endl;
879 }
880 }
881}
882
883/**
884 * @brief Calculate whether the system has reached a steady state by
885 * observing residuals to a user-defined tolerance.
886 */
887bool UnsteadySystem::CheckSteadyState(int step, const NekDouble &totCPUTime)
888{
889 const int nPoints = GetTotPoints();
890 const int nFields = m_fields.size();
891
892 // Holds L2 errors.
893 Array<OneD, NekDouble> L2(nFields);
894
895 SteadyStateResidual(step, L2);
896
897 if (m_infosteps && m_comm->GetRank() == 0 &&
898 (((step + 1) % m_infosteps == 0) || ((step == m_initialStep))))
899 {
900 // Output time.
901 m_errFile << boost::format("%25.19e") % m_time;
902
903 m_errFile << " " << boost::format("%25.19e") % totCPUTime;
904
905 int stepp = step + 1;
906
907 m_errFile << " " << boost::format("%25.19e") % stepp;
908
909 // Output residuals.
910 for (int i = 0; i < nFields; ++i)
911 {
912 m_errFile << " " << boost::format("%25.19e") % L2[i];
913 }
914
915 m_errFile << endl;
916 }
917
918 // Calculate maximum L2 error.
919 NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
920
921 if (m_infosteps && m_session->DefinesCmdLineArgument("verbose") &&
922 m_comm->GetRank() == 0 && ((step + 1) % m_infosteps == 0))
923 {
924 cout << "-- Maximum L^2 residual: " << maxL2 << endl;
925 }
926
927 if (maxL2 <= m_steadyStateTol)
928 {
929 return true;
930 }
931
932 for (int i = 0; i < m_fields.size(); ++i)
933 {
934 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
935 1);
936 }
937
938 return false;
939}
940
941/**
942 *
943 */
944void UnsteadySystem::v_SteadyStateResidual([[maybe_unused]] int step,
946{
947 const int nPoints = GetTotPoints();
948 const int nFields = m_fields.size();
949
950 // Holds L2 errors.
951 Array<OneD, NekDouble> RHSL2(nFields);
952 Array<OneD, NekDouble> residual(nFields);
953 Array<OneD, NekDouble> reference(nFields);
954
955 for (int i = 0; i < nFields; ++i)
956 {
957 Array<OneD, NekDouble> tmp(nPoints);
958
959 Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
960 1, tmp, 1);
961 Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
962 residual[i] = Vmath::Vsum(nPoints, tmp, 1);
963
965 tmp, 1);
966 reference[i] = Vmath::Vsum(nPoints, tmp, 1);
967 }
968
969 m_comm->GetSpaceComm()->AllReduce(residual, LibUtilities::ReduceSum);
970 m_comm->GetSpaceComm()->AllReduce(reference, LibUtilities::ReduceSum);
971
972 // L2 error.
973 for (int i = 0; i < nFields; ++i)
974 {
975 reference[i] = (reference[i] == 0) ? 1 : reference[i];
976 L2[i] = sqrt(residual[i] / reference[i]);
977 }
978}
979
980/**
981 *
982 */
984 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
986 [[maybe_unused]] const NekDouble time)
987{
988
989 if (&inarray != &outarray)
990 {
991 for (int i = 0; i < inarray.size(); ++i)
992 {
993 Vmath::Vcopy(GetNpoints(), inarray[i], 1, outarray[i], 1);
994 }
995 }
996}
997
998} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
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:223
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
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
virtual SOLVER_UTILS_EXPORT void v_ALEPreMultiplyMass(Array< OneD, Array< OneD, NekDouble > > &fields)
virtual SOLVER_UTILS_EXPORT void v_ALEInitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Definition ALEHelper.cpp:42
SOLVER_UTILS_EXPORT void ALEDoElmtInvMass(Array< OneD, Array< OneD, NekDouble > > &traceNormals, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble time)
Update m_fields with u^n by multiplying by inverse mass matrix. That's then used in e....
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity(const NekDouble &time)
Definition ALEHelper.cpp:97
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.
SOLVER_UTILS_EXPORT int GetNpoints()
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.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
int m_initialStep
Number of the step where the simulation should begin.
SOLVER_UTILS_EXPORT int GetTotPoints()
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.
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.
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.
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()
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:217
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
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
STL namespace.
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:305
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:300
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290
std::vector< NekDouble > freeParams