![]() |
Gowdy solver
|
Time evolution integrator for the Method of Lines (MoL). More...
#include <integrator.h>
Public Member Functions | |
MoLIntegrator (Parameters ¶ms, UniformGrid &ug, GridOutputWriter &out) | |
Constructs the MoL integrator as specified in the parameter file. | |
void | quit (Int m=-1) |
Moves the integration final time to the current time effectively stopping the integration. | |
void | addToEvolution (IntegFace *eomInterface) |
Adds the EoM to the list of all evolved subsystems. | |
void | keepConstant (const std::vector< Int > &gfs) |
The given fields will be kept constant by the integrator. | |
void | keepEvolved (const std::vector< fld::EvolvedBy > &gfs) |
Add the given functions to the evolution list. | |
Real | dt () const |
Returns the time step. | |
bool | evolveEquations () |
Integrates the evolution equations in the given grid realm. | |
Private Member Functions | |
void | reportIntegrationTime (Int m) |
Reports the integration time and the estimated real-time to finish. | |
void | integStep_Begin (Int m) |
Executed at the beginning of each integration step. | |
void | integStep_End (Int m1, Int m) |
Executed at the end of each integration step. | |
void | integStep_Checkpoint (Int m) |
Executed at each checkpoint. | |
void | integStep_Euler (Int mNext, Int m) |
Evolves the state variables using the FT scheme. | |
void | integStep_AvgICN (Int mNext, Int mMid, Int m) |
Evolves the sate variables using the ICN averaging. | |
void | integrate_AvgICN (int ICN_steps) |
Time evolution using the iterated Crank-Nicolson (ICN). | |
void | integStep_MoLInit (Int m1, Int m, Real beta_delta_t) |
Prepares the state variables for the intermediate MoL steps. | |
void | integStep_MoL (Int m1, Int m, Real alpha) |
The intermediate method of lines (MoL) steps. | |
void | integrate_MoL (const MoLDescriptor &MoL) |
Time evolution using the method of lines (MoL). | |
Static Private Member Functions | |
static void | signalHandler (int signum) |
Forward the signal to all integrators (signaling them to quit). | |
Private Attributes | |
std::vector< IntegFace * > | eomList |
Equations of motion that we evolve. | |
std::vector< Int > | constantGF |
Grid functions that are kept constant. | |
std::vector< fld::EvolvedBy > | evolvedGF |
Grid functions that are integrated in time. | |
GridOutputWriter * | output |
The integration output is directed here. | |
time_hr | rt_start |
The realtime when the integration started. | |
Int | method |
The integration method (Euler, ICN, MoL, ...) | |
Real | t_0 |
The initial t . | |
Real | t_1 |
Integrate up to this t . | |
Real | delta_t |
The integration step. | |
Real | dissip |
Kreiss-Oliger dissipation coefficient. | |
Real | dissip_delta_r |
dissip / delta_r | |
Real | cur_t |
Current time. | |
Static Private Attributes | |
static std::map< std::string, int > | knownMethods |
Implemented integration methods. | |
static std::vector < MoLIntegrator * > | knownIntegrators |
Keep the list of all created integrators (used by signalHandler). |
Time evolution integrator for the Method of Lines (MoL).
The implemented methods are declared in MoLIntegrator::knownMethods.
The integrator can evolve equations from several modules at the same time. All such modules should implement IntegFace ans sign in to the integrator. The modules should also specify which grid functions are either evolved or kept constant by the integrator. (Those grid functions that are unknown to the integrator have to be maintained by the the modules themselves.)
The system provides the coordinates t
and r
as global grid functions (they are defined in fld::sysIndex). The integrator keeps time t
evolved and the spatial coordinate r
constant.
The integrator handles process signals (e.g., Ctrl-C), which will cause a premature end of the integration.
Definition at line 31 of file integrator.h.
MoLIntegrator::MoLIntegrator | ( | Parameters & | params, |
UniformGrid & | ug, | ||
GridOutputWriter & | out | ||
) | [inline] |
Constructs the MoL integrator as specified in the parameter file.
Definition at line 191 of file integrator.h.
References abs(), constantGF, cur_t, GridUser::delta_r, delta_t, dissip, dissip_delta_r, evolvedGF, Parameters::get(), GridOutputWriter::get_mSkip(), knownIntegrators, knownMethods, method, output, fld::r, GridOutputWriter::recordSizeInBytes(), signalHandler(), slog, t_0, and t_1.
: GridUser( ug ), output( &out ) { constantGF.reserve( 100 ); evolvedGF.reserve( 100 ); constantGF.push_back( fld::r ); // Radial coordinate is kept constant by default params.get( "integ.t_0", t_0, 0.0 ); params.get( "integ.t_1", t_1, 1.0 ); params.get( "integ.delta_t", delta_t, 0.1 ); params.get( "integ.dissip", dissip, 0.0 ); std::string name = params.get( "integ.method", method, 0, knownMethods ); // Fix the sign of delta_t // delta_t = t_1 < t_0 ? -abs(delta_t) : abs(delta_t); cur_t = t_0; // The initial `t` // Calculate Kreiss-Oliger dissipation corrected delta r // dissip_delta_r = dissip / delta_r; // Expected output file data size (in bytes) // const Real sz = ( ( t_1 - t_0 ) / delta_t / output->get_mSkip() + 1 ) * output->recordSizeInBytes(); slog << "Integrator:" << std::endl << std::endl << " t_0 = " << t_0 << ", t_1 = " << t_1 << ", delta_t = " << delta_t << ", method = " << name << " (#" << method << ")" << std::endl << " Kreiss-Oliger dissipation = " << dissip << std::endl << " Expected output data size = " << round( ( sz < 1e3 ? sz : sz < 1e6 ? sz/1e3 : sz/1e6 ) * 10 ) / 10 << ( sz < 1e3 ? "" : sz < 1e6 ? " kB" : " MB" ) << std::endl << std::endl; // Register to receive the signals from the system // knownIntegrators.push_back( this ); if ( knownIntegrators.size () == 1 ) { signal( SIGINT, signalHandler ); signal( SIGTERM, signalHandler ); } }
void MoLIntegrator::addToEvolution | ( | IntegFace * | eomInterface | ) | [inline] |
Adds the EoM to the list of all evolved subsystems.
Definition at line 249 of file integrator.h.
References eomList.
Referenced by GowdyEvolve::GowdyEvolve().
{ eomList.push_back( eomInterface ); }
Real MoLIntegrator::dt | ( | ) | const [inline] |
Returns the time step.
Definition at line 270 of file integrator.h.
References delta_t.
{ return delta_t; }
bool MoLIntegrator::evolveEquations | ( | ) | [inline] |
Integrates the evolution equations in the given grid realm.
Definition at line 276 of file integrator.h.
References ICN2, ICN3, integrate_AvgICN(), integrate_MoL(), method, RK1, RK2, RK3, RK4, and rt_start.
Referenced by main().
{ // Start timing the real-time rt_start = std::chrono::high_resolution_clock::now (); switch( method ) { case 0: integrate_AvgICN ( 0 ); break; case 1: integrate_AvgICN ( 1 ); break; case 2: integrate_AvgICN ( 2 ); break; case 3: integrate_AvgICN ( 3 ); break; case 4: integrate_MoL ( ICN2 ); break; case 5: integrate_MoL ( ICN3 ); break; case 6: integrate_MoL ( RK1 ); break; case 7: integrate_MoL ( RK2 ); break; case 8: integrate_MoL ( RK3 ); break; case 9: integrate_MoL ( RK4 ); break; } return true; }
void MoLIntegrator::integStep_Begin | ( | Int | m | ) | [inline, private] |
Executed at the beginning of each integration step.
Prepares the intermediate variables, exchanges or calculates boundaries, and finally, calculates the RHS of time evolution equations.
Definition at line 105 of file integrator.h.
References eomList, UniformGrid::exchangeBoundaries(), GridUser::gridDriver, MPIWorld::isFirstInRank(), and MPIWorld::isLastInRank().
Referenced by integrate_AvgICN(), and integrate_MoL().
{ for( auto eom: eomList ) { eom->integStep_Prepare( m ); } // Exchange boundaries with the neighboring MPI ranks (if any) // gridDriver->exchangeBoundaries( m ); if ( gridDriver->isFirstInRank () ) { for( auto eom: eomList ) { eom->applyLeftBoundaryCondition( m ); } } if ( gridDriver->isLastInRank () ) { for( auto eom: eomList ) { eom->applyRightBoundaryCondition( m ); } } for( auto eom: eomList ) { eom->integStep_CalcEvolutionRHS( m ); } }
void MoLIntegrator::integStep_Checkpoint | ( | Int | m | ) | [inline, private] |
Executed at each checkpoint.
Reports the integration time, writes the output data and checks for eventual NaNs.
Definition at line 146 of file integrator.h.
References cur_t, delta_t, eomList, GridOutputWriter::get_nOut(), GridUser::nGhost, GridUser::nLen, output, quit(), reportIntegrationTime(), and GridOutputWriter::write().
Referenced by integrate_AvgICN(), and integrate_MoL().
{ reportIntegrationTime( m ); output->write( m, cur_t, delta_t ); // Check for NaNs in the central/lower region of the grid // Int nFrom = nGhost + nLen/10; Int nTo = nGhost + std::min( output->get_nOut(), nLen/2 ); for( auto eom: eomList ) { if ( eom->integStep_CheckForNaNs( m, nFrom, nTo ) ) { quit( m ); } } }
void MoLIntegrator::integStep_End | ( | Int | m1, |
Int | m | ||
) | [inline, private] |
Executed at the end of each integration step.
Finalizes the step (e.g., post processing of EoM) and waits exchange of the boundaries.
Definition at line 135 of file integrator.h.
References eomList, GridUser::gridDriver, and MPIWorld::waitExchange().
Referenced by integrate_AvgICN(), and integrate_MoL().
{ for( auto eom: eomList ) { eom->integStep_Finalize( m1, m ); } gridDriver->waitExchange (); }
void MoLIntegrator::keepConstant | ( | const std::vector< Int > & | gfs | ) | [inline] |
The given fields will be kept constant by the integrator.
Definition at line 256 of file integrator.h.
References constantGF.
{ constantGF.insert( constantGF.end(), gfs.begin(), gfs.end() ); }
void MoLIntegrator::keepEvolved | ( | const std::vector< fld::EvolvedBy > & | gfs | ) | [inline] |
Add the given functions to the evolution list.
Definition at line 263 of file integrator.h.
References evolvedGF.
Referenced by GowdyEvolve::GowdyEvolve().
void MoLIntegrator::quit | ( | Int | m = -1 | ) | [inline] |
Moves the integration final time to the current time effectively stopping the integration.
Definition at line 241 of file integrator.h.
References cur_t, GF, GridUser::nGhost, slog, fld::t, and t_1.
Referenced by integStep_Checkpoint().
void MoLIntegrator::reportIntegrationTime | ( | Int | m | ) | [inline, private] |
Reports the integration time and the estimated real-time to finish.
Definition at line 72 of file integrator.h.
References delta_t, GridOutputWriter::get_mSkip(), GF, GridUser::nGhost, output, rt_start, fld::t, t_0, and t_1.
Referenced by integStep_Checkpoint().
{ Real tnow = GF( fld::t, m, nGhost ); // current integration time auto elapsed = 1e-3 * std::chrono::duration_cast<std::chrono::milliseconds> ( std::chrono::high_resolution_clock::now () - rt_start ).count (); auto rt_factor = elapsed / ( tnow - t_0 ); // real time/integration time auto rt_total = round( ( t_1 - t_0 ) * rt_factor ); // projected total time auto rt_left = round( ( t_1 - tnow ) * rt_factor ); // estimated time to cover #if 0 std::ios oldState( NULL ); oldState.copyfmt( std::cout ); std::cout.setf( std::ios::fixed, std::ios::floatfield ); std::cout << " t = " << std::setw( 1 + round(log10( t_1 / output->get_mSkip() / delta_t )) ) << std::setprecision( - round(log10( output->get_mSkip() * delta_t )) ) << tnow; std::cout.copyfmt( oldState ); #endif std::cout << " t = " << std::setw(6) << std::setprecision(2) << std::fixed << tnow << std::setw(-1) << std::setprecision(-1) << std::defaultfloat << ", real = " << round( elapsed ) << " s (" << round( 100 * ( tnow - t_0 ) / ( t_1 - t_0 ) ) << "% of " << rt_total << " s), ETA " << rt_left << " s \r"; std::flush( std::cout ); }
static void MoLIntegrator::signalHandler | ( | int | signum | ) | [inline, static, private] |
Forward the signal to all integrators (signaling them to quit).
Definition at line 58 of file integrator.h.
References knownIntegrators, and slog.
Referenced by MoLIntegrator().
{ slog << "*** Signal " << signum << " received ***"; if( knownIntegrators.size () == 0 ) { exit( -1 ); } for( auto i: knownIntegrators ) { i->quit (); } }
std::vector<Int> MoLIntegrator::constantGF [private] |
Grid functions that are kept constant.
Definition at line 34 of file integrator.h.
Referenced by integStep_Euler(), integStep_MoLInit(), keepConstant(), and MoLIntegrator().
Real MoLIntegrator::cur_t [private] |
Current time.
Definition at line 46 of file integrator.h.
Referenced by integrate_AvgICN(), integrate_MoL(), integStep_Checkpoint(), MoLIntegrator(), and quit().
Real MoLIntegrator::delta_t [private] |
The integration step.
Definition at line 43 of file integrator.h.
Referenced by dt(), integrate_MoL(), integStep_Checkpoint(), integStep_Euler(), MoLIntegrator(), and reportIntegrationTime().
Real MoLIntegrator::dissip [private] |
Kreiss-Oliger dissipation coefficient.
Definition at line 44 of file integrator.h.
Referenced by MoLIntegrator().
Real MoLIntegrator::dissip_delta_r [private] |
std::vector<IntegFace*> MoLIntegrator::eomList [private] |
Equations of motion that we evolve.
Definition at line 33 of file integrator.h.
Referenced by addToEvolution(), integStep_Begin(), integStep_Checkpoint(), and integStep_End().
std::vector<fld::EvolvedBy> MoLIntegrator::evolvedGF [private] |
Grid functions that are integrated in time.
Definition at line 35 of file integrator.h.
Referenced by integStep_Euler(), integStep_MoL(), integStep_MoLInit(), keepEvolved(), and MoLIntegrator().
Int MoLIntegrator::method [private] |
The integration method (Euler, ICN, MoL, ...)
Definition at line 40 of file integrator.h.
Referenced by evolveEquations(), and MoLIntegrator().
GridOutputWriter* MoLIntegrator::output [private] |
The integration output is directed here.
Definition at line 37 of file integrator.h.
Referenced by integrate_AvgICN(), integrate_MoL(), integStep_Checkpoint(), MoLIntegrator(), and reportIntegrationTime().
time_hr MoLIntegrator::rt_start [private] |
The realtime when the integration started.
Definition at line 39 of file integrator.h.
Referenced by evolveEquations(), and reportIntegrationTime().
Real MoLIntegrator::t_0 [private] |
The initial t
.
Definition at line 41 of file integrator.h.
Referenced by integrate_AvgICN(), integrate_MoL(), MoLIntegrator(), and reportIntegrationTime().
Real MoLIntegrator::t_1 [private] |
Integrate up to this t
.
Definition at line 42 of file integrator.h.
Referenced by integrate_AvgICN(), integrate_MoL(), MoLIntegrator(), quit(), and reportIntegrationTime().