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), 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 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 // @TODO: Refactor to take coeffs (FwdTrans) if boolean flag (in constructor
248 // function) says to.
249 for (i = 0; i < nvariables; ++i)
250 {
251 fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
252 m_fields[m_intVariables[i]]->SetPhysState(false);
253 }
254
255 // @TODO: Virtual function that allows to transform the field space, embed
256 // the MultiplyMassMatrix in here.
257 // @TODO: Specify what the fields variables are physical or coefficient,
258 // boolean in UnsteadySystem class...
259
260 v_ALEPreMultiplyMass(fields);
261
262 // Initialise time integration scheme.
263 m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
264
265 // Initialise filters.
266 for (auto &x : m_filters)
267 {
268 x.second->Initialise(m_fields, m_time);
269 }
270
272 bool doCheckTime = false;
273 int step = m_initialStep;
274 int stepCounter = 0;
275 NekDouble intTime = 0.0;
276 NekDouble cpuTime = 0.0;
277 NekDouble cpuPrevious = 0.0;
278 NekDouble elapsed = 0.0;
279 NekDouble totFilterTime = 0.0;
280
281 m_lastCheckTime = 0.0;
282
283 Array<OneD, int> abortFlags(2, 0);
284 string abortFile = "abort";
285 if (m_session->DefinesSolverInfo("CheckAbortFile"))
286 {
287 abortFile = m_session->GetSolverInfo("CheckAbortFile");
288 }
289
290 NekDouble tmp_cflSafetyFactor = m_cflSafetyFactor;
291 while ((step < m_steps || m_time < m_fintime - NekConstants::kNekZeroTol) &&
292 abortFlags[1] == 0)
293 {
295 {
296 tmp_cflSafetyFactor =
297 min(m_CFLEnd, m_CFLGrowth * tmp_cflSafetyFactor);
298 }
299
300 // Frozen preconditioner checks.
301 if (!m_comm->IsParallelInTime())
302 {
304 {
305 m_cflSafetyFactor = tmp_cflSafetyFactor;
306
308 {
309 m_timestep = GetTimeStep(fields);
310 }
311
312 // Ensure that the final timestep finishes at the final
313 // time, or at a prescribed IO_CheckTime.
314 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
315 {
317 }
318 else if (m_checktime &&
320 {
323 doCheckTime = true;
324 }
325 }
326 }
327
328 if (m_TimeIncrementFactor > 1.0)
329 {
330 NekDouble timeincrementFactor = m_TimeIncrementFactor;
331 m_timestep *= timeincrementFactor;
332
333 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
334 {
336 }
337 }
338
339 // Perform any solver-specific pre-integration steps.
340 timer.Start();
341 if (v_PreIntegrate(
342 step)) // Could be possible to put a preintegrate step in the
343 // ALEHelper class, put in the Unsteady Advection class
344 {
345 break;
346 }
347
348 ASSERTL0(m_timestep > 0, "m_timestep < 0");
349
350 fields = m_intScheme->TimeIntegrate(stepCounter, m_timestep);
351 timer.Stop();
352
354 elapsed = timer.TimePerTest(1);
355 intTime += elapsed;
356 cpuTime += elapsed;
357
358 // Write out status information.
359 v_PrintStatusInformation(step, cpuTime);
360 if (m_infosteps &&
361 m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
362 !((step + 1) % m_infosteps))
363 {
364 cpuPrevious = cpuTime;
365 cpuTime = 0.0;
366 }
367
368 // @TODO: Another virtual function with this in it based on if ALE or
369 // not.
370 if (m_ALESolver) // Change to advect coeffs, change flag to physical vs
371 // coefficent space
372 {
375 }
376 else
377 {
378 // Transform data into coefficient space
379 for (i = 0; i < nvariables; ++i)
380 {
381 // copy fields into ExpList::m_phys and assign the new
382 // array to fields
383 m_fields[m_intVariables[i]]->SetPhys(fields[i]);
384 fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
385 if (v_RequireFwdTrans())
386 {
387 if (m_comm->IsParallelInTime())
388 {
389 m_fields[m_intVariables[i]]->FwdTrans(
390 m_fields[m_intVariables[i]]->GetPhys(),
391 m_fields[m_intVariables[i]]->UpdateCoeffs());
392 }
393 else
394 {
395 m_fields[m_intVariables[i]]->FwdTransLocalElmt(
396 m_fields[m_intVariables[i]]->GetPhys(),
397 m_fields[m_intVariables[i]]->UpdateCoeffs());
398 }
399 }
400 m_fields[m_intVariables[i]]->SetPhysState(false);
401 }
402 }
403
404 // Perform any solver-specific post-integration steps.
405 if (v_PostIntegrate(step))
406 {
407 break;
408 }
409
410 // Check for steady-state.
412 (!((step + 1) % m_steadyStateSteps)))
413 {
414 if (CheckSteadyState(step, intTime))
415 {
416 if (m_comm->GetRank() == 0)
417 {
418 cout << "Reached Steady State to tolerance "
419 << m_steadyStateTol << endl;
420 }
421 break;
422 }
423 }
424
425 // Test for abort conditions (nan, or abort file).
426 if (m_abortSteps && !((step + 1) % m_abortSteps))
427 {
428 abortFlags[0] = 0;
429 for (i = 0; i < nvariables; ++i)
430 {
431 if (Vmath::Nnan(m_fields[m_intVariables[i]]->GetPhys().size(),
432 m_fields[m_intVariables[i]]->GetPhys(), 1) > 0)
433 {
434 abortFlags[0] = 1;
435 }
436 }
437
438 // Rank zero looks for abort file and deletes it
439 // if it exists. Communicates the abort.
440 if (m_session->GetComm()->GetSpaceComm()->GetRank() == 0)
441 {
442 if (fs::exists(abortFile))
443 {
444 fs::remove(abortFile);
445 abortFlags[1] = 1;
446 }
447 }
448
449 m_session->GetComm()->GetSpaceComm()->AllReduce(
450 abortFlags, LibUtilities::ReduceMax);
451
452 ASSERTL0(!abortFlags[0], "NaN found during time integration.");
453 }
454
455 // Update filters.
456 for (auto &x : m_filters)
457 {
458 timer.Start();
459 x.second->Update(m_fields, m_time);
460 timer.Stop();
461 elapsed = timer.TimePerTest(1);
462 totFilterTime += elapsed;
463
464 // Write out individual filter status information.
465 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
466 !((step + 1) % m_filtersInfosteps) && !m_filters.empty() &&
467 m_session->DefinesCmdLineArgument("verbose"))
468 {
469 stringstream s0;
470 s0 << x.first << ":";
471 stringstream s1;
472 s1 << elapsed << "s";
473 stringstream s2;
474 s2 << elapsed / cpuPrevious * 100 << "%";
475 cout << "CPU time for filter " << setw(25) << left << s0.str()
476 << setw(12) << left << s1.str() << endl
477 << "\t Percentage of time integration: " << setw(10)
478 << left << s2.str() << endl;
479 }
480 }
481
482 // Write out overall filter status information.
483 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
484 !((step + 1) % m_filtersInfosteps) && !m_filters.empty())
485 {
486 stringstream ss;
487 ss << totFilterTime << "s";
488 cout << "Total filters CPU Time:\t\t\t " << setw(10) << left
489 << ss.str() << endl;
490 }
491 totFilterTime = 0.0;
492
493 // Write out checkpoint files.
494 if ((m_checksteps && !((step + 1) % m_checksteps)) || doCheckTime)
495 {
497 {
498 // Transform to physical space for output.
499 vector<bool> transformed(nfields, false);
500 for (i = 0; i < nfields; i++)
501 {
502 if (m_fields[i]->GetWaveSpace())
503 {
504 m_fields[i]->SetWaveSpace(false);
505 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
506 m_fields[i]->UpdatePhys());
507 transformed[i] = true;
508 }
509 }
511 m_nchk++;
512
513 // Transform back to wave-space after output.
514 for (i = 0; i < nfields; i++)
515 {
516 if (transformed[i])
517 {
518 m_fields[i]->SetWaveSpace(true);
519 m_fields[i]->HomogeneousFwdTrans(
520 m_fields[i]->GetTotPoints(), m_fields[i]->GetPhys(),
521 m_fields[i]->UpdatePhys());
522 m_fields[i]->SetPhysState(false);
523 }
524 }
525 }
526 else
527 {
529 m_nchk++;
530 }
531 doCheckTime = false;
532 }
533
534 // Step advance.
535 ++step;
536 ++stepCounter;
537 }
538
539 // Print out summary statistics.
541
542 // If homogeneous, transform back into physical space if necessary.
543 if (!m_ALESolver)
544 {
546 {
547 for (i = 0; i < nfields; i++)
548 {
549 if (m_fields[i]->GetWaveSpace())
550 {
551 m_fields[i]->SetWaveSpace(false);
552 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
553 m_fields[i]->UpdatePhys());
554 }
555 }
556 }
557 else
558 {
559 for (i = 0; i < nvariables; ++i)
560 {
561 m_fields[m_intVariables[i]]->SetPhysState(true);
562 }
563 }
564 }
565 // Finalise filters.
566 for (auto &x : m_filters)
567 {
568 x.second->Finalise(m_fields, m_time);
569 }
570
571 // Print for 1D problems.
572 if (m_spacedim == 1)
573 {
575 }
576}
577
579 const NekDouble cpuTime)
580{
581 if (m_infosteps && m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
582 !((step + 1) % m_infosteps))
583 {
584 if (m_comm->IsParallelInTime())
585 {
586 cout << "RANK " << m_session->GetComm()->GetTimeComm()->GetRank()
587 << " Steps: " << setw(8) << left << step + 1 << " "
588 << "Time: " << setw(12) << left << m_time;
589 }
590 else
591 {
592 cout << "Steps: " << setw(8) << left << step + 1 << " "
593 << "Time: " << setw(12) << left << m_time;
594 }
595
597 {
598 cout << " Time-step: " << setw(12) << left << m_timestep;
599 }
600
601 stringstream ss;
602 ss << cpuTime << "s";
603 cout << " CPU Time: " << setw(8) << left << ss.str() << endl;
604 }
605}
606
608{
609 if (m_session->GetComm()->GetRank() == 0)
610 {
611 if (m_cflSafetyFactor > 0.0)
612 {
613 cout << "CFL safety factor : " << m_cflSafetyFactor << endl
614 << "CFL time-step : " << m_timestep << endl;
615 }
616
617 if (m_session->GetSolverInfo("Driver") != "SteadyState" &&
618 m_session->GetSolverInfo("Driver") != "Parareal" &&
619 m_session->GetSolverInfo("Driver") != "PFASST")
620 {
621 cout << "Time-integration : " << intTime << "s" << endl;
622 }
623 }
624}
625
626/**
627 * @brief Sets the initial conditions.
628 */
629void UnsteadySystem::v_DoInitialise(bool dumpInitialConditions)
630{
633 SetInitialConditions(m_time, dumpInitialConditions);
634
636
638}
639
640/**
641 * @brief Prints a summary with some information regards the
642 * time-stepping.
643 */
645{
647 AddSummaryItem(s, "Advect. advancement",
648 m_explicitAdvection ? "explicit" : "implicit");
649
650 AddSummaryItem(s, "Diffuse. advancement",
651 m_explicitDiffusion ? "explicit" : "implicit");
652
653 if (m_session->GetSolverInfo("EQTYPE") ==
654 "SteadyAdvectionDiffusionReaction")
655 {
656 AddSummaryItem(s, "React. advancement",
657 m_explicitReaction ? "explicit" : "implicit");
658 }
659
660 AddSummaryItem(s, "Time Step", m_timestep);
661 AddSummaryItem(s, "No. of Steps", m_steps);
662 AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
663 if (m_intScheme)
664 {
665 AddSummaryItem(s, "Integration Type", m_intScheme->GetName());
666 }
667}
668
669/**
670 * Stores the solution in a file for 1D problems only. This method has
671 * been implemented to facilitate the post-processing for 1D problems.
672 */
674{
675 // Coordinates of the quadrature points in the real physical space.
679 m_fields[0]->GetCoords(x, y, z);
680
681 // Print out the solution in a txt file.
682 ofstream outfile;
683 outfile.open("solution1D.txt");
684 for (int i = 0; i < GetNpoints(); i++)
685 {
686 outfile << scientific << setw(17) << setprecision(16) << x[i] << " "
687 << m_fields[m_intVariables[0]]->GetPhys()[i] << endl;
688 }
689 outfile << endl << endl;
690 outfile.close();
691}
692
693/**
694 *
695 */
697{
698 if (m_session->DefinesFunction("InitialConditions"))
699 {
700 for (int i = 0; i < m_fields.size(); ++i)
701 {
703
704 vType = m_session->GetFunctionType("InitialConditions",
705 m_session->GetVariable(i));
706
708 {
709 std::string filename = m_session->GetFunctionFilename(
710 "InitialConditions", m_session->GetVariable(i));
711
712 fs::path pfilename(filename);
713
714 // Redefine path for parallel file which is in directory.
715 if (fs::is_directory(pfilename))
716 {
717 fs::path metafile("Info.xml");
718 fs::path fullpath = pfilename / metafile;
719 filename = LibUtilities::PortablePath(fullpath);
720 }
723 fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
724
725 // Check to see if time defined.
727 {
728 auto iter = m_fieldMetaDataMap.find("Time");
729 if (iter != m_fieldMetaDataMap.end())
730 {
731 time = std::stod(iter->second);
732 }
733
734 iter = m_fieldMetaDataMap.find("ChkFileNum");
735 if (iter != m_fieldMetaDataMap.end())
736 {
737 nchk = std::stod(iter->second);
738 }
739 }
740
741 break;
742 }
743 }
744 }
745 if (m_session->DefinesCmdLineArgument("set-start-time"))
746 {
747 time = std::stod(
748 m_session->GetCmdLineArgument<std::string>("set-start-time")
749 .c_str());
750 }
751 if (m_session->DefinesCmdLineArgument("set-start-chknumber"))
752 {
753 nchk = boost::lexical_cast<int>(
754 m_session->GetCmdLineArgument<std::string>("set-start-chknumber"));
755 }
756 ASSERTL0(time >= 0 && nchk >= 0,
757 "Starting time and checkpoint number should be >= 0");
758}
759
760/**
761 * @brief Return the timestep to be used for the next step in the
762 * time-marching loop.
763 *
764 * @see UnsteadySystem::GetTimeStep
765 */
767 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &inarray)
768{
769 NEKERROR(ErrorUtil::efatal, "Not defined for this class");
770 return 0.0;
771}
772
773/**
774 *
775 */
776bool UnsteadySystem::v_PreIntegrate([[maybe_unused]] int step)
777{
778 return false;
779}
780
781/**
782 *
783 */
784bool UnsteadySystem::v_PostIntegrate([[maybe_unused]] int step)
785{
786 return false;
787}
788
789/**
790 *
791 */
793{
794 return true;
795}
796
797/**
798 *
799 */
801{
802 return true;
803}
804
805/**
806 *
807 */
810 StdRegions::VarCoeffMap &varCoeffMap)
811{
812 int phystot = m_fields[0]->GetTotPoints();
813 int nvel = vel.size();
814
815 Array<OneD, NekDouble> varcoeff(phystot), tmp;
816
817 // Calculate magnitude of v.
818 Vmath::Vmul(phystot, vel[0], 1, vel[0], 1, varcoeff, 1);
819 for (int n = 1; n < nvel; ++n)
820 {
821 Vmath::Vvtvp(phystot, vel[n], 1, vel[n], 1, varcoeff, 1, varcoeff, 1);
822 }
823 Vmath::Vsqrt(phystot, varcoeff, 1, varcoeff, 1);
824
825 for (int i = 0; i < m_fields[0]->GetNumElmts(); ++i)
826 {
827 int offset = m_fields[0]->GetPhys_Offset(i);
828 int nq = m_fields[0]->GetExp(i)->GetTotPoints();
829 Array<OneD, NekDouble> unit(nq, 1.0);
830
831 int nmodes = 0;
832
833 for (int n = 0; n < m_fields[0]->GetExp(i)->GetNumBases(); ++n)
834 {
835 nmodes = max(nmodes, m_fields[0]->GetExp(i)->GetBasisNumModes(n));
836 }
837
838 NekDouble h = m_fields[0]->GetExp(i)->Integral(unit);
839 h = pow(h, (NekDouble)(1.0 / nvel)) / ((NekDouble)nmodes);
840
841 Vmath::Smul(nq, h, varcoeff + offset, 1, tmp = varcoeff + offset, 1);
842 }
843
844 // Set up map with eVarCoffLaplacian key.
845 varCoeffMap[StdRegions::eVarCoeffLaplacian] = varcoeff;
846}
847
848/**
849 *
850 */
852{
853 if (m_steadyStateTol > 0.0)
854 {
855 const int nPoints = m_fields[0]->GetTotPoints();
858
859 for (int i = 0; i < m_fields.size(); ++i)
860 {
862 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1,
863 m_previousSolution[i], 1);
864 }
865
866 if (m_comm->GetRank() == 0)
867 {
868 std::string fName =
869 m_session->GetSessionName() + std::string(".resd");
870 m_errFile.open(fName.c_str());
871 m_errFile << setw(26) << left << "# Time";
872
873 m_errFile << setw(26) << left << "CPU_Time";
874
875 m_errFile << setw(26) << left << "Step";
876
877 for (int i = 0; i < m_fields.size(); ++i)
878 {
879 m_errFile << setw(26) << m_session->GetVariables()[i];
880 }
881
882 m_errFile << endl;
883 }
884 }
885}
886
887/**
888 * @brief Calculate whether the system has reached a steady state by
889 * observing residuals to a user-defined tolerance.
890 */
891bool UnsteadySystem::CheckSteadyState(int step, const NekDouble &totCPUTime)
892{
893 const int nPoints = GetTotPoints();
894 const int nFields = m_fields.size();
895
896 // Holds L2 errors.
897 Array<OneD, NekDouble> L2(nFields);
898
899 SteadyStateResidual(step, L2);
900
901 if (m_infosteps && m_comm->GetRank() == 0 &&
902 (((step + 1) % m_infosteps == 0) || ((step == m_initialStep))))
903 {
904 // Output time.
905 m_errFile << boost::format("%25.19e") % m_time;
906
907 m_errFile << " " << boost::format("%25.19e") % totCPUTime;
908
909 int stepp = step + 1;
910
911 m_errFile << " " << boost::format("%25.19e") % stepp;
912
913 // Output residuals.
914 for (int i = 0; i < nFields; ++i)
915 {
916 m_errFile << " " << boost::format("%25.19e") % L2[i];
917 }
918
919 m_errFile << endl;
920 }
921
922 // Calculate maximum L2 error.
923 NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
924
925 if (m_infosteps && m_session->DefinesCmdLineArgument("verbose") &&
926 m_comm->GetRank() == 0 && ((step + 1) % m_infosteps == 0))
927 {
928 cout << "-- Maximum L^2 residual: " << maxL2 << endl;
929 }
930
931 if (maxL2 <= m_steadyStateTol)
932 {
933 return true;
934 }
935
936 for (int i = 0; i < m_fields.size(); ++i)
937 {
938 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
939 1);
940 }
941
942 return false;
943}
944
945/**
946 *
947 */
948void UnsteadySystem::v_SteadyStateResidual([[maybe_unused]] int step,
950{
951 const int nPoints = GetTotPoints();
952 const int nFields = m_fields.size();
953
954 // Holds L2 errors.
955 Array<OneD, NekDouble> RHSL2(nFields);
956 Array<OneD, NekDouble> residual(nFields);
957 Array<OneD, NekDouble> reference(nFields);
958
959 for (int i = 0; i < nFields; ++i)
960 {
961 Array<OneD, NekDouble> tmp(nPoints);
962
963 Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
964 1, tmp, 1);
965 Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
966 residual[i] = Vmath::Vsum(nPoints, tmp, 1);
967
969 tmp, 1);
970 reference[i] = Vmath::Vsum(nPoints, tmp, 1);
971 }
972
973 m_comm->GetSpaceComm()->AllReduce(residual, LibUtilities::ReduceSum);
974 m_comm->GetSpaceComm()->AllReduce(reference, LibUtilities::ReduceSum);
975
976 // L2 error.
977 for (int i = 0; i < nFields; ++i)
978 {
979 reference[i] = (reference[i] == 0) ? 1 : reference[i];
980 L2[i] = sqrt(residual[i] / reference[i]);
981 }
982}
983
984/**
985 *
986 */
988 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
990 [[maybe_unused]] const NekDouble time)
991{
992
993 if (&inarray != &outarray)
994 {
995 for (int i = 0; i < inarray.size(); ++i)
996 {
997 Vmath::Vcopy(GetNpoints(), inarray[i], 1, outarray[i], 1);
998 }
999 }
1000}
1001
1002} // 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.
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)
Definition: ALEHelper.cpp:108
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....
Definition: ALEHelper.cpp:131
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity(const NekDouble &time)
Definition: ALEHelper.cpp:90
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.
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.
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()
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:375
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
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
std::vector< NekDouble > freeParams
Definition: SessionReader.h:87