![]() |
Gowdy solver
|
00001 /** 00002 * @file integrator.h 00003 * @brief Time evolution integrator for the Method of Lines (MoL). 00004 * @authors Mikica Kocic, Wilkinson & Reinsch, Press et al 00005 * @copyright GNU General Public License (GPLv3). 00006 */ 00007 00008 #include <csignal> 00009 00010 ///////////////////////////////////////////////////////////////////////////////////////// 00011 /** @defgroup g11 Numerical Methods */ 00012 ///////////////////////////////////////////////////////////////////////////////////////// 00013 /*@{*/ 00014 /** Time evolution integrator for the Method of Lines (MoL). 00015 * 00016 * The implemented methods are declared in MoLIntegrator::knownMethods. 00017 * 00018 * The integrator can evolve equations from several modules at the same time. 00019 * All such modules should implement @ref IntegFace ans sign in to the integrator. 00020 * The modules should also specify which grid functions are either evolved or 00021 * kept constant by the integrator. (Those grid functions that are unknown to 00022 * the integrator have to be maintained by the the modules themselves.) 00023 * 00024 * The system provides the coordinates `t` and `r` as global grid functions 00025 * (they are defined in @ref fld::sysIndex). The integrator keeps time `t` evolved and 00026 * the spatial coordinate `r` constant. 00027 * 00028 * The integrator handles process signals (e.g., Ctrl-C), which will cause 00029 * a premature end of the integration. 00030 */ 00031 class MoLIntegrator : GridUser 00032 { 00033 std::vector<IntegFace*> eomList; //!< Equations of motion that we evolve 00034 std::vector<Int> constantGF; //!< Grid functions that are kept constant 00035 std::vector<fld::EvolvedBy> evolvedGF; //!< Grid functions that are integrated in time 00036 00037 GridOutputWriter* output; //!< The integration output is directed here 00038 00039 time_hr rt_start; //!< The realtime when the integration started 00040 Int method; //!< The integration method (Euler, ICN, MoL, ...) 00041 Real t_0; //!< The initial `t`. 00042 Real t_1; //!< Integrate up to this `t`. 00043 Real delta_t; //!< The integration step 00044 Real dissip; //!< Kreiss-Oliger dissipation coefficient. 00045 Real dissip_delta_r; //!< `dissip / delta_r` 00046 Real cur_t; //!< Current time 00047 00048 /** Implemented integration methods 00049 */ 00050 static std::map<std::string,int> knownMethods; 00051 00052 /** Keep the list of all created integrators (used by signalHandler). 00053 */ 00054 static std::vector<MoLIntegrator*> knownIntegrators; 00055 00056 /** Forward the signal to all integrators (signaling them to quit). 00057 */ 00058 static void signalHandler( int signum ) 00059 { 00060 slog << "*** Signal " << signum << " received ***"; 00061 00062 if( knownIntegrators.size () == 0 ) { 00063 exit( -1 ); 00064 } 00065 for( auto i: knownIntegrators ) { 00066 i->quit (); 00067 } 00068 } 00069 00070 /** Reports the integration time and the estimated real-time to finish. 00071 */ 00072 void reportIntegrationTime( Int m ) 00073 { 00074 Real tnow = GF( fld::t, m, nGhost ); // current integration time 00075 00076 auto elapsed = 1e-3 * std::chrono::duration_cast<std::chrono::milliseconds> 00077 ( std::chrono::high_resolution_clock::now () - rt_start ).count (); 00078 00079 auto rt_factor = elapsed / ( tnow - t_0 ); // real time/integration time 00080 auto rt_total = round( ( t_1 - t_0 ) * rt_factor ); // projected total time 00081 auto rt_left = round( ( t_1 - tnow ) * rt_factor ); // estimated time to cover 00082 #if 0 00083 std::ios oldState( NULL ); 00084 oldState.copyfmt( std::cout ); 00085 std::cout.setf( std::ios::fixed, std::ios::floatfield ); 00086 std::cout << " t = " 00087 << std::setw( 1 + round(log10( t_1 / output->get_mSkip() / delta_t )) ) 00088 << std::setprecision( - round(log10( output->get_mSkip() * delta_t )) ) 00089 << tnow; 00090 std::cout.copyfmt( oldState ); 00091 #endif 00092 std::cout << " t = " << std::setw(6) << std::setprecision(2) << std::fixed 00093 << tnow << std::setw(-1) << std::setprecision(-1) << std::defaultfloat 00094 << ", real = " << round( elapsed ) 00095 << " s (" << round( 100 * ( tnow - t_0 ) / ( t_1 - t_0 ) ) 00096 << "% of " << rt_total 00097 << " s), ETA " << rt_left << " s \r"; 00098 std::flush( std::cout ); 00099 } 00100 00101 /** Executed at the beginning of each integration step. 00102 * Prepares the intermediate variables, exchanges or calculates boundaries, 00103 * and finally, calculates the RHS of time evolution equations. 00104 */ 00105 void integStep_Begin( Int m ) 00106 { 00107 for( auto eom: eomList ) { 00108 eom->integStep_Prepare( m ); 00109 } 00110 00111 // Exchange boundaries with the neighboring MPI ranks (if any) 00112 // 00113 gridDriver->exchangeBoundaries( m ); 00114 00115 if ( gridDriver->isFirstInRank () ) { 00116 for( auto eom: eomList ) { 00117 eom->applyLeftBoundaryCondition( m ); 00118 } 00119 } 00120 if ( gridDriver->isLastInRank () ) { 00121 for( auto eom: eomList ) { 00122 eom->applyRightBoundaryCondition( m ); 00123 } 00124 } 00125 00126 for( auto eom: eomList ) { 00127 eom->integStep_CalcEvolutionRHS( m ); 00128 } 00129 } 00130 00131 /** Executed at the end of each integration step. 00132 * Finalizes the step (e.g., post processing of EoM) and waits exchange of the 00133 * boundaries. 00134 */ 00135 void integStep_End( Int m1, Int m ) 00136 { 00137 for( auto eom: eomList ) { 00138 eom->integStep_Finalize( m1, m ); 00139 } 00140 gridDriver->waitExchange (); 00141 } 00142 00143 /** Executed at each checkpoint. Reports the integration time, writes 00144 * the output data and checks for eventual NaNs. 00145 */ 00146 void integStep_Checkpoint( Int m ) 00147 { 00148 reportIntegrationTime( m ); 00149 00150 output->write( m, cur_t, delta_t ); 00151 00152 // Check for NaNs in the central/lower region of the grid 00153 // 00154 Int nFrom = nGhost + nLen/10; 00155 Int nTo = nGhost + std::min( output->get_nOut(), nLen/2 ); 00156 for( auto eom: eomList ) { 00157 if ( eom->integStep_CheckForNaNs( m, nFrom, nTo ) ) { 00158 quit( m ); 00159 } 00160 } 00161 } 00162 00163 /** Evolves the state variables using the FT scheme. 00164 */ 00165 void integStep_Euler( Int mNext, Int m ); 00166 00167 /** Evolves the sate variables using the ICN averaging. 00168 */ 00169 void integStep_AvgICN( Int mNext, Int mMid, Int m ); 00170 00171 /** Time evolution using the iterated Crank-Nicolson (ICN). 00172 */ 00173 void integrate_AvgICN( int ICN_steps ); 00174 00175 /** Prepares the state variables for the intermediate MoL steps. 00176 */ 00177 void integStep_MoLInit( Int m1, Int m, Real beta_delta_t ); 00178 00179 /** The intermediate method of lines (MoL) steps. 00180 */ 00181 void integStep_MoL( Int m1, Int m, Real alpha ); 00182 00183 /** Time evolution using the method of lines (MoL). 00184 */ 00185 void integrate_MoL( const MoLDescriptor& MoL ); 00186 00187 public: 00188 00189 /** Constructs the MoL integrator as specified in the parameter file. 00190 */ 00191 MoLIntegrator( Parameters& params, UniformGrid& ug, GridOutputWriter& out ) 00192 : GridUser( ug ), output( &out ) 00193 { 00194 constantGF.reserve( 100 ); 00195 evolvedGF.reserve( 100 ); 00196 00197 constantGF.push_back( fld::r ); // Radial coordinate is kept constant by default 00198 00199 params.get( "integ.t_0", t_0, 0.0 ); 00200 params.get( "integ.t_1", t_1, 1.0 ); 00201 params.get( "integ.delta_t", delta_t, 0.1 ); 00202 params.get( "integ.dissip", dissip, 0.0 ); 00203 00204 std::string name = params.get( "integ.method", method, 0, knownMethods ); 00205 00206 // Fix the sign of delta_t 00207 // 00208 delta_t = t_1 < t_0 ? -abs(delta_t) : abs(delta_t); 00209 cur_t = t_0; // The initial `t` 00210 00211 // Calculate Kreiss-Oliger dissipation corrected delta r 00212 // 00213 dissip_delta_r = dissip / delta_r; 00214 00215 // Expected output file data size (in bytes) 00216 // 00217 const Real sz = ( ( t_1 - t_0 ) / delta_t / output->get_mSkip() + 1 ) 00218 * output->recordSizeInBytes(); 00219 00220 slog << "Integrator:" << std::endl << std::endl 00221 << " t_0 = " << t_0 << ", t_1 = " << t_1 << ", delta_t = " << delta_t 00222 << ", method = " << name << " (#" << method << ")" << std::endl 00223 << " Kreiss-Oliger dissipation = " << dissip << std::endl 00224 << " Expected output data size = " 00225 << round( ( sz < 1e3 ? sz : sz < 1e6 ? sz/1e3 : sz/1e6 ) * 10 ) / 10 00226 << ( sz < 1e3 ? "" : sz < 1e6 ? " kB" : " MB" ) 00227 << std::endl << std::endl; 00228 00229 // Register to receive the signals from the system 00230 // 00231 knownIntegrators.push_back( this ); 00232 if ( knownIntegrators.size () == 1 ) { 00233 signal( SIGINT, signalHandler ); 00234 signal( SIGTERM, signalHandler ); 00235 } 00236 } 00237 00238 /** Moves the integration final time to the current time effectively 00239 * stopping the integration. 00240 */ 00241 void quit( Int m = -1 ) 00242 { 00243 t_1 = m < 0 ? cur_t : GF( fld::t, m, nGhost ); 00244 slog << std::endl << "Quitting at t = " << t_1; 00245 } 00246 00247 /** Adds the EoM to the list of all evolved subsystems. 00248 */ 00249 void addToEvolution( IntegFace* eomInterface ) 00250 { 00251 eomList.push_back( eomInterface ); 00252 } 00253 00254 /** The given fields will be kept constant by the integrator. 00255 */ 00256 void keepConstant( const std::vector<Int>& gfs ) 00257 { 00258 constantGF.insert( constantGF.end(), gfs.begin(), gfs.end() ); 00259 } 00260 00261 /** Add the given functions to the evolution list. 00262 */ 00263 void keepEvolved( const std::vector<fld::EvolvedBy>& gfs ) 00264 { 00265 evolvedGF.insert( evolvedGF.end(), gfs.begin(), gfs.end() ); 00266 } 00267 00268 /** Returns the time step. 00269 */ 00270 Real dt () const { 00271 return delta_t; 00272 } 00273 00274 /** Integrates the evolution equations in the given grid realm. 00275 */ 00276 bool evolveEquations () 00277 { 00278 // Start timing the real-time 00279 rt_start = std::chrono::high_resolution_clock::now (); 00280 00281 switch( method ) 00282 { 00283 case 0: integrate_AvgICN ( 0 ); break; 00284 case 1: integrate_AvgICN ( 1 ); break; 00285 case 2: integrate_AvgICN ( 2 ); break; 00286 case 3: integrate_AvgICN ( 3 ); break; 00287 case 4: integrate_MoL ( ICN2 ); break; 00288 case 5: integrate_MoL ( ICN3 ); break; 00289 case 6: integrate_MoL ( RK1 ); break; 00290 case 7: integrate_MoL ( RK2 ); break; 00291 case 8: integrate_MoL ( RK3 ); break; 00292 case 9: integrate_MoL ( RK4 ); break; 00293 } 00294 00295 return true; 00296 } 00297 }; 00298 00299 ///////////////////////////////////////////////////////////////////////////////////////// 00300 // Definitions of MoLIntegrator static members 00301 ///////////////////////////////////////////////////////////////////////////////////////// 00302 00303 std::vector<MoLIntegrator*> MoLIntegrator::knownIntegrators; 00304 00305 std::map<std::string,int> MoLIntegrator::knownMethods = 00306 { 00307 { "Euler", 0 }, { "avgICN2", 1 }, 00308 { "avgICN3", 2 }, { "avgICN4", 3 }, 00309 { "ICN2", 4 }, { "ICN3", 5 }, 00310 { "RK1", 6 }, { "RK2", 7 }, 00311 { "RK3", 8 }, { "RK4", 9 } 00312 }; 00313 00314 ///////////////////////////////////////////////////////////////////////////////////////// 00315 // Definitions of MoLIntegrator methods 00316 ///////////////////////////////////////////////////////////////////////////////////////// 00317 00318 void MoLIntegrator::integStep_Euler( Int m1, Int m ) 00319 { 00320 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00321 { 00322 GF( fld::t, m1, n ) = GF( fld::t, m, n ) + delta_t; // Evolve time 00323 00324 for( auto f: constantGF ) 00325 { 00326 GF( f, m1, n ) = GF( f, m, n ); 00327 } 00328 00329 for( auto e: evolvedGF ) 00330 { 00331 GF( e.f, m1, n ) = GF( e.f, m, n ) + delta_t * GF( e.f_t, m, n ) 00332 - KreissOligerTermGF( e.f, delta_t ); 00333 } 00334 } 00335 } 00336 00337 void MoLIntegrator::integStep_AvgICN 00338 ( 00339 Int m2, // final time 00340 Int m1, // middle time 00341 Int m // initial time 00342 ) 00343 { 00344 const Real half_dt = delta_t / 2; 00345 00346 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00347 { 00348 GF( fld::t, m2, n ) = GF( fld::t, m, n ) + delta_t; // Evolve time 00349 00350 for( auto f: constantGF ) 00351 { 00352 GF( f, m2, n ) = GF( f, m, n ); 00353 } 00354 00355 for( auto e: evolvedGF ) 00356 { 00357 GF( e.f, m2, n ) = GF( e.f, m, n ) 00358 + half_dt * ( GF( e.f_t, m1, n) + GF( e.f_t, m, n ) ) 00359 - KreissOligerTermGF( e.f, delta_t ); 00360 } 00361 } 00362 } 00363 00364 /** Integrates the evolution equations using the iterated Crank-Nicolson (ICN) method. 00365 * If ICN_steps == 0, then the method effectively becomes FTCS. 00366 * See @cite Press:2007nrec, @cite Baumgarte:2010nr 00367 */ 00368 void MoLIntegrator::integrate_AvgICN( int ICN_steps ) 00369 { 00370 slog << "Employing Method of Lines, Iterative Crank-Nicolson (avg)" 00371 << ", N = " << ICN_steps << std::endl << std::endl; 00372 00373 cur_t = t_0; // The initial `t` 00374 00375 for( Int n = 0; n < nLen + 2 * nGhost; ++n ) { 00376 GF( fld::t, 0, n ) = t_0; 00377 } 00378 00379 Int mStep = 0; // The step counter (used to filter checkpoints) 00380 while( abs(cur_t) < abs(t_1) ) 00381 { 00382 for( Int m = 0; m < mLen && abs(cur_t) < abs(t_1); ++m ) 00383 { 00384 Int mNext = m + 1 >= mLen ? 0 : m + 1; 00385 00386 ///////////////////////////////////////////////////////////////////////////// 00387 // Predictor. Use FTCS for the first step. 00388 // Note that this will be the only step, if ICN_steps == 0. 00389 // 00390 Int m_1 = ICN_steps >= 1 ? mLen + 1 : mNext; 00391 00392 integStep_Begin ( m ); 00393 integStep_Euler ( m_1, m ); 00394 integStep_End ( m_1, m ); 00395 00396 ///////////////////////////////////////////////////////////////////////////// 00397 // Corrector iterations. There will be no iterations, if ICN_steps == 0. 00398 // 00399 for( int i = 0; i < ICN_steps; ++i ) 00400 { 00401 Int m_2 = ( i + 1 >= ICN_steps ? mNext : m_1 + 1 ); 00402 00403 integStep_Begin ( m_1 ); 00404 integStep_AvgICN ( m_2, m_1, m ); 00405 integStep_End ( m_2, m_1 ); 00406 00407 m_1 = m_2; 00408 } 00409 00410 if( ( mStep++ % output->get_mSkip () ) == 0 ) { 00411 integStep_Checkpoint( m ); 00412 } 00413 00414 cur_t = GF( fld::t, m, nGhost ); // The current time step 00415 } 00416 } 00417 } 00418 00419 void MoLIntegrator::integStep_MoL( Int m1, Int m, Real ac ) 00420 { 00421 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00422 { 00423 GF( fld::t, m1, n ) += ac * GF( fld::t, m, n ); // Evolve time 00424 00425 for( auto e: evolvedGF ) 00426 { 00427 GF( e.f, m1, n ) += ac * GF( e.f, m, n ); 00428 } 00429 } 00430 } 00431 00432 void MoLIntegrator::integStep_MoLInit( Int m1, Int m, Real bc_delta_t ) 00433 { 00434 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00435 { 00436 // There is no explicit dependence on t, and t_t(m,n) = 1. 00437 // Also, r is constant along the `line`, i.e., r_t(m,n) = 0. 00438 // 00439 GF( fld::t, m1, n ) = bc_delta_t; // Evolve time 00440 00441 for( auto f: constantGF ) 00442 { 00443 GF( f, m1, n ) = GF( f, m, n ); 00444 } 00445 00446 for( auto e: evolvedGF ) 00447 { 00448 GF( e.f, m1, n ) = bc_delta_t * GF( e.f_t, m, n ) 00449 - KreissOligerTermGF( e.f, bc_delta_t ); 00450 } 00451 } 00452 } 00453 00454 /** The method of lines (MoL) separates the time integration from the rest of 00455 * an evolution scheme. A `N`-step MoL method evolves the equation: 00456 * <pre> 00457 * partial_t Y = L[Y], </pre> 00458 * 00459 * using the following algorithm: 00460 * <pre> 00461 * Y^{(0)} = Y^{m}, 00462 * Y^{(i+1)} = Sum_{j=0}^{i} ( alpha_{ij} Y^{(j)} ) + <- integStep_MoL 00463 * Delta t beta_i L[Y^{(i)}], for i = 0, ..., N-1, <- integStep_MoLInit 00464 * Y^{m+1} = Y^{(N)}. </pre> 00465 * 00466 * The variables `Y^{(i)}`, `i = 0, ..., N`, are intermediate. 00467 * 00468 * The method is completely specified by `N`, and arrays alpha and beta 00469 * in the structure MoLDescriptor. 00470 */ 00471 void MoLIntegrator::integrate_MoL( const MoLDescriptor& MoL ) 00472 { 00473 slog << "Employing Method of Lines, " << MoL.name 00474 << ", N = " << MoL.N << std::endl << std::endl; 00475 00476 cur_t = t_0; // The initial `t` 00477 00478 for( Int n = 0; n < nLen + 2 * nGhost; ++n ) { 00479 GF( fld::t, 0, n ) = t_0; 00480 } 00481 00482 // The array mk[] contains redirections to grid rows (m's): 00483 // mk[0] points to Y^{m} == Y^{(0)}, 00484 // mk[1], ..., mk[Mol.N-1] point to the intermediate Y^{(1)}, ..., Y^{(N-1)}, 00485 // mk[MoL.N] points to Y^{m+1} == Y^{(N)}. 00486 // 00487 int mk[ 10 ]; // int mk[ MoL.N + 1 ]; // Note that N is at most 8 00488 for( int i = 0; i <= MoL.N; ++i ) { 00489 mk[i] = mLen + i; 00490 } 00491 00492 Int mStep = 0; // The step counter 00493 while( abs(cur_t) < abs(t_1) ) 00494 { 00495 for( Int m = 0; m < mLen && abs(cur_t) < abs(t_1); ++m ) 00496 { 00497 mk[0] = m; // Y^{(m)} 00498 mk[MoL.N] = m + 1 >= mLen ? 0 : m + 1; // Y^{(m+1)} 00499 00500 for( int i = 0; i < MoL.N; ++i ) 00501 { 00502 integStep_Begin( mk[i] ); 00503 00504 integStep_MoLInit( mk[i+1], mk[i], MoL.beta[i] * delta_t ); 00505 00506 for( int j = 0; j <= i; j++ ) { 00507 integStep_MoL( mk[i+1], mk[j], MoL.alpha[i][j] ); 00508 } 00509 00510 integStep_End( mk[i+1], mk[i] ); 00511 } 00512 00513 if( ( mStep++ % output->get_mSkip () ) == 0 ) { 00514 integStep_Checkpoint( m ); 00515 } 00516 00517 cur_t = GF( fld::t, m, nGhost ); // The current time step 00518 } 00519 } 00520 } 00521 /*@}*/ 00522 /////////////////////////////////////////////////////////////////////////////////////////