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