40 #include <boost/core/ignore_unused.hpp>
41 #include <boost/format.hpp>
69 UnsteadySystem::UnsteadySystem(
88 m_session->MatchSolverInfo(
"DIFFUSIONADVANCEMENT",
"Explicit",
90 m_session->MatchSolverInfo(
"ADVECTIONADVANCEMENT",
"Explicit",
92 m_session->MatchSolverInfo(
"REACTIONADVANCEMENT",
"Explicit",
95 m_session->MatchSolverInfo(
"FLAGIMPLICITITSSTATISTICS",
"True",
102 m_session->LoadParameter(
"SteadyStateSteps",
106 if (
m_session->DefinesSolverInfo(
"TimeIntegrationMethod") ||
117 "TimeIntegrationMethod");
128 m_session->LoadParameter(
"IO_FiltersInfoSteps",
137 m_session->LoadParameter(
"FilterTimeWarning",
145 "Final condition not unique: "
146 "fintime > 0.0 and Nsteps > 0");
149 "Timestep not unique: timestep > 0.0 & CFL > 0.0");
154 "Need to set either TimeStep or CFL");
162 "m_CFLEnd >= m_cflSafetyFactor required");
167 boost::lexical_cast<std::string>(
m_time);
177 x.first,
m_session, shared_from_this(), x.second)));
210 for (i = 0; i < nfields; ++i)
214 nvariables = nfields;
224 for(i = 0; i < nfields; ++i)
237 for(i = 0; i < nvariables; ++i)
253 bool doCheckTime =
false;
256 int restartStep = -1;
270 string abortFile =
"abort";
271 if (
m_session->DefinesSolverInfo(
"CheckAbortFile"))
273 abortFile =
m_session->GetSolverInfo(
"CheckAbortFile");
287 tmp_cflSafetyFactor =
353 if (
m_session->GetComm()->GetRank() == 0 &&
356 cout <<
"Steps: " << setw(8) << left << step+1 <<
" "
357 <<
"Time: " << setw(12) << left <<
m_time;
361 cout <<
" Time-step: " << setw(12)
366 ss << cpuTime <<
"s";
367 cout <<
" CPU Time: " << setw(8) << left
369 cpuPrevious = cpuTime;
384 for (i = 0; i < nvariables; ++i)
410 if (
m_comm->GetRank() == 0)
412 cout <<
"Reached Steady State to tolerance "
423 for (i = 0; i < nvariables; ++i)
436 if(boost::filesystem::exists(abortFile))
438 boost::filesystem::remove(abortFile);
443 m_session->GetComm()->AllReduce(abortFlags,
447 "NaN found during time integration.");
457 totFilterTime += elapsed;
460 if(
m_session->GetComm()->GetRank() == 0 &&
462 m_session->DefinesCmdLineArgument(
"verbose"))
465 s0 << x.first <<
":";
467 s1 << elapsed <<
"s";
469 s2 << elapsed / cpuPrevious * 100 <<
"%";
470 cout <<
"CPU time for filter " << setw(25) << left
471 << s0.str() << setw(12) << left << s1.str() <<
472 endl <<
"\t Percentage of time integration: "
473 << setw(10) << left << s2.str() << endl;
478 if (
m_session->GetComm()->GetRank() == 0 &&
482 ss << totFilterTime <<
"s";
483 cout <<
"Total filters CPU Time:\t\t\t " << setw(10)
484 << left << ss.str() << endl;
494 vector<bool> transformed(nfields,
false);
495 for(i = 0; i < nfields; i++)
503 transformed[i] =
true;
508 for(i = 0; i < nfields; i++)
534 if (
m_session->GetComm()->GetRank() == 0)
542 if (
m_session->GetSolverInfo(
"Driver") !=
"SteadyState")
544 cout <<
"Time-integration : " << intTime <<
"s" << endl;
550 <<
"-------------------------------------------" << endl
554 <<
"-------------------------------------------" << endl;
561 for(i = 0; i < nfields; i++)
574 for(i = 0; i < nvariables; ++i)
615 if(
m_session->DefinesSolverInfo(
"AdvectionType"))
618 m_session->GetSolverInfo(
"AdvectionType"));
625 ==
"SteadyAdvectionDiffusionReaction")
652 outfile.open(
"solution1D.txt");
655 outfile << scientific << setw (17) << setprecision(16) << x[i]
656 <<
" " << solution1D[0][i] << endl;
658 outfile << endl << endl;
664 if (
m_session->DefinesFunction(
"InitialConditions"))
666 for (
int i = 0; i <
m_fields.size(); ++i)
671 "InitialConditions",
m_session->GetVariable(i));
677 "InitialConditions",
m_session->GetVariable(i));
679 fs::path pfilename(filename);
682 if(fs::is_directory(pfilename))
684 fs::path metafile(
"Info.xml");
685 fs::path fullpath = pfilename / metafile;
700 time = boost::lexical_cast<NekDouble>(
707 nchk = boost::lexical_cast<NekDouble>(
741 boost::ignore_unused(inarray);
748 boost::ignore_unused(step);
754 boost::ignore_unused(step);
762 int phystot =
m_fields[0]->GetTotPoints();
763 int nvel = vel.size();
769 for(
int n = 1; n < nvel; ++n)
771 Vmath::Vvtvp(phystot,vel[n],1,vel[n],1,varcoeff,1,varcoeff,1);
775 for(
int i = 0; i <
m_fields[0]->GetNumElmts(); ++i)
777 int offset =
m_fields[0]->GetPhys_Offset(i);
778 int nq =
m_fields[0]->GetExp(i)->GetTotPoints();
783 for(
int n = 0; n <
m_fields[0]->GetExp(i)->GetNumBases(); ++n)
786 m_fields[0]->GetExp(i)->GetBasisNumModes(n));
792 Vmath::Smul(nq,h,varcoeff+offset,1,tmp = varcoeff+offset,1);
803 const int nPoints =
m_fields[0]->GetTotPoints();
807 for (
int i = 0; i <
m_fields.size(); ++i)
814 if (
m_comm->GetRank() == 0)
816 std::string fName =
m_session->GetSessionName() +
817 std::string(
".resd");
819 m_errFile << setw(26) << left <<
"# Time";
821 m_errFile << setw(26) << left <<
"CPU_Time";
825 for (
int i = 0; i <
m_fields.size(); ++i)
847 const int nFields =
m_fields.size();
866 for (
int i = 0; i < nFields; ++i)
877 if (
m_session->DefinesCmdLineArgument(
"verbose") &&
880 cout <<
"-- Maximum L^2 residual: " << maxL2 << endl;
888 for (
int i = 0; i <
m_fields.size(); ++i)
904 boost::ignore_unused(step);
906 const int nFields =
m_fields.size();
913 for (
int i = 0; i < nFields; ++i)
931 for (
int i = 0; i < nFields; ++i)
933 reference[i] = (reference[i] == 0) ? 1 : reference[i];
934 L2[i] =
sqrt(residual[i] / reference[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.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
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.
NekDouble m_lastCheckTime
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()
NekDouble m_timestepMax
Time step size.
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.
NekDouble m_TimeIncrementFactor
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.
NekDouble m_CFLGrowth
CFL growth rate.
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.
bool m_CalcPhysicalAV
flag to update artificial viscosity
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.
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)
NekDouble m_CFLEnd
maximun cfl in cfl growth
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 void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)
bool m_flagImplicitSolver
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime(NekDouble &time, int &nchk)
NekDouble m_steadyStateRes
SOLVER_UTILS_EXPORT void SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)
NekDouble m_steadyStateRes0
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
std::vector< int > m_intVariables
bool m_flagUpdatePreconMat
int m_steadyStateSteps
Check for steady state at step interval.
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck()
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.
bool m_flagImplicitItsStatistics
TimeIntegrationSchemeFactory & GetTimeIntegrationSchemeFactory()
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
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
The above copyright notice and this permission notice shall be included.
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*x.
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.
scalarT< T > sqrt(scalarT< T > in)
std::vector< NekDouble > freeParams