39 #include <boost/core/ignore_unused.hpp>
40 #include <boost/format.hpp>
68 UnsteadySystem::UnsteadySystem(
87 m_session->MatchSolverInfo(
"DIFFUSIONADVANCEMENT",
"Explicit",
89 m_session->MatchSolverInfo(
"ADVECTIONADVANCEMENT",
"Explicit",
91 m_session->MatchSolverInfo(
"REACTIONADVANCEMENT",
"Explicit",
98 m_session->LoadParameter(
"SteadyStateSteps",
102 if (
m_session->DefinesSolverInfo(
"TIMEINTEGRATIONMETHOD"))
106 "TIMEINTEGRATIONMETHOD"));
110 m_session->LoadParameter(
"IO_FiltersInfoSteps",
122 "Final condition not unique: "
123 "fintime > 0.0 and Nsteps > 0");
126 "Timestep not unique: timestep > 0.0 & CFL > 0.0");
131 "Need to set either TimeStep or CFL");
136 boost::lexical_cast<std::string>(
m_time);
146 x.first,
m_session, shared_from_this(), x.second)));
169 TimeStability = 2.784;
191 "No CFL control implementation for this time"
192 "integration scheme");
195 return TimeStability;
208 int nfields =
m_fields.num_elements();
212 for (i = 0; i < nfields; ++i)
216 nvariables = nfields;
226 for(i = 0; i < nfields; ++i)
239 for(i = 0; i < nvariables; ++i)
256 bool doCheckTime =
false;
267 string abortFile =
"abort";
268 if (
m_session->DefinesSolverInfo(
"CheckAbortFile"))
270 abortFile =
m_session->GetSolverInfo(
"CheckAbortFile");
313 if (
m_session->GetComm()->GetRank() == 0 &&
316 cout <<
"Steps: " << setw(8) << left << step+1 <<
" "
317 <<
"Time: " << setw(12) << left <<
m_time;
321 cout <<
" Time-step: " << setw(12)
326 ss << cpuTime <<
"s";
327 cout <<
" CPU Time: " << setw(8) << left
329 cpuPrevious = cpuTime;
334 for (i = 0; i < nvariables; ++i)
360 if (
m_comm->GetRank() == 0)
362 cout <<
"Reached Steady State to tolerance "
373 for (i = 0; i < nvariables; ++i)
386 if(boost::filesystem::exists(abortFile))
388 boost::filesystem::remove(abortFile);
393 m_session->GetComm()->AllReduce(abortFlags,
397 "NaN found during time integration.");
407 totFilterTime += elapsed;
410 if(
m_session->GetComm()->GetRank() == 0 &&
412 m_session->DefinesCmdLineArgument(
"verbose"))
415 s0 << x.first <<
":";
417 s1 << elapsed <<
"s";
419 s2 << elapsed / cpuPrevious * 100 <<
"%";
420 cout <<
"CPU time for filter " << setw(25) << left
421 << s0.str() << setw(12) << left << s1.str() <<
422 endl <<
"\t Percentage of time integration: "
423 << setw(10) << left << s2.str() << endl;
428 if (
m_session->GetComm()->GetRank() == 0 &&
432 ss << totFilterTime <<
"s";
433 cout <<
"Total filters CPU Time:\t\t\t " << setw(10)
434 << left << ss.str() << endl;
444 vector<bool> transformed(nfields,
false);
445 for(i = 0; i < nfields; i++)
453 transformed[i] =
true;
458 for(i = 0; i < nfields; i++)
484 if (
m_session->GetComm()->GetRank() == 0)
492 if (
m_session->GetSolverInfo(
"Driver") !=
"SteadyState")
494 cout <<
"Time-integration : " << intTime <<
"s" << endl;
501 for(i = 0; i < nfields; i++)
514 for(i = 0; i < nvariables; ++i)
555 if(
m_session->DefinesSolverInfo(
"AdvectionType"))
558 m_session->GetSolverInfo(
"AdvectionType"));
565 ==
"SteadyAdvectionDiffusionReaction")
594 outfile.open(
"solution1D.txt");
597 outfile << scientific << setw (17) << setprecision(16) << x[i]
598 <<
" " << solution1D[0][i] << endl;
600 outfile << endl << endl;
606 if (
m_session->DefinesFunction(
"InitialConditions"))
608 for (
int i = 0; i <
m_fields.num_elements(); ++i)
613 "InitialConditions",
m_session->GetVariable(i));
619 "InitialConditions",
m_session->GetVariable(i));
621 fs::path pfilename(filename);
624 if(fs::is_directory(pfilename))
626 fs::path metafile(
"Info.xml");
627 fs::path fullpath = pfilename / metafile;
642 time = boost::lexical_cast<NekDouble>(
649 nchk = boost::lexical_cast<NekDouble>(
683 boost::ignore_unused(inarray);
690 boost::ignore_unused(step);
696 boost::ignore_unused(step);
704 int phystot =
m_fields[0]->GetTotPoints();
705 int nvel = vel.num_elements();
711 for(
int n = 1; n < nvel; ++n)
713 Vmath::Vvtvp(phystot,vel[n],1,vel[n],1,varcoeff,1,varcoeff,1);
717 for(
int i = 0; i <
m_fields[0]->GetNumElmts(); ++i)
719 int offset =
m_fields[0]->GetPhys_Offset(i);
720 int nq =
m_fields[0]->GetExp(i)->GetTotPoints();
725 for(
int n = 0; n <
m_fields[0]->GetExp(i)->GetNumBases(); ++n)
728 m_fields[0]->GetExp(i)->GetBasisNumModes(n));
734 Vmath::Smul(nq,h,varcoeff+offset,1,tmp = varcoeff+offset,1);
745 const int nPoints =
m_fields[0]->GetTotPoints();
749 for (
int i = 0; i <
m_fields.num_elements(); ++i)
756 if (
m_comm->GetRank() == 0)
758 std::string fName =
m_session->GetSessionName() +
761 m_errFile << setw(26) << left <<
"# Time";
763 for (
int i = 0; i <
m_fields.num_elements(); ++i)
780 const int nFields =
m_fields.num_elements();
787 for (
int i = 0; i < nFields; ++i)
805 for (
int i = 0; i < nFields; ++i)
807 reference[i] = (reference[i] == 0) ? 1 : reference[i];
808 L2[i] = sqrt(residual[i] / reference[i]);
817 for (
int i = 0; i < nFields; ++i)
828 if (
m_session->DefinesCmdLineArgument(
"verbose") &&
831 cout <<
"-- Maximum L^2 residual: " << maxL2 << endl;
839 for (
int i = 0; i <
m_fields.num_elements(); ++i)
#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.
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
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.
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()
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.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Initialisation object for EquationSystem.
SOLVER_UTILS_EXPORT int GetTotPoints()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitReaction
Indicates if explicit or implicit treatment of reaction is used.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Sets up initial conditions.
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
Storage for previous solution for steady-state check.
virtual SOLVER_UTILS_EXPORT void v_DoSolve()
Solves an unsteady problem.
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem()
Destructor.
NekDouble m_filterTimeWarning
Number of time steps between outputting status information.
NekDouble m_steadyStateTol
Tolerance to which steady state should be evaluated at.
int m_abortSteps
Number of steps between checks for abort conditions.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator()
Get the maximum timestep estimator for cfl control.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble >> &inarray)
Return the timestep to be used for the next step in the time-marching loop.
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans()
bool CheckSteadyState(int step)
Calculate whether the system has reached a steady state by observing residuals to a user-defined tole...
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate(int step)
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff(const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h t...
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
void InitializeSteadyState()
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime(NekDouble &time, int &nchk)
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
std::vector< int > m_intVariables
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
int m_steadyStateSteps
Check for steady state at step interval.
int m_filtersInfosteps
Number of time steps between outputting filters information.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D(Array< OneD, Array< OneD, NekDouble >> &solution1D)
Print the solution at each solution point in a txt file.
int m_infosteps
Number of time steps between outputting status information.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
@ eRungeKutta4
Classical RungeKutta4 method (new name for eClassicalRungeKutta4)
@ eRungeKutta3_SSP
Nonlinear SSP RungeKutta3 explicit.
@ eAdamsBashforthOrder2
Adams-Bashforth Forward multi-step scheme of order 2.
@ eMidpoint
midpoint method (old name)
@ eForwardEuler
Forward Euler scheme.
@ eRungeKutta2_SSP
Nonlinear SSP RungeKutta2 explicit (surrogate for eRungeKutta2_ImprovedEuler)
@ eRungeKutta2
Classical RungeKutta2 method (new name for eMidpoint)
@ eRungeKutta2_ImprovedEuler
Improved RungeKutta2 explicit (old name meaning Heun's method)
@ eClassicalRungeKutta4
Runge-Kutta multi-stage scheme 4th order explicit (old name)
@ eAdamsBashforthOrder1
Adams-Bashforth Forward multi-step scheme of order 1.
const char *const TimeIntegrationMethodMap[]
std::shared_ptr< FieldIO > FieldIOSharedPtr
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
TimeIntegrationWrapperFactory & GetTimeIntegrationWrapperFactory()
static const NekDouble kNekZeroTol
std::vector< std::pair< std::string, std::string > > SummaryList
FilterFactory & GetFilterFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
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.
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
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.