![]() |
Bimetric 3+1 toolkit for spherical symmetry
|
00001 /** 00002 * @file integrator.h 00003 * @brief Time evolution using 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 using the Method of Lines (MoL). 00015 * 00016 * The implemented methods are declared in MoL::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 MoL : 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 bool running; //!< Indicates if the integration should proceed 00042 Real t_0; //!< The initial `t`. 00043 Real t_1; //!< Integrate up to this `t`. 00044 Real delta_t; //!< The integration step 00045 Real dissip; //!< Kreiss-Oliger dissipation coefficient. 00046 Real dissip_delta_r; //!< `dissip / delta_r` 00047 Real cur_t; //!< Current time 00048 00049 /** A list of the implemented integration methods. 00050 */ 00051 static std::map<std::string,int> knownMethods; 00052 00053 /** Keep the list of all created integrator objects (used by signalHandler). 00054 */ 00055 static std::vector<MoL*> knownIntegrators; 00056 00057 /** Forward the signal to all integrators (signaling them to quit). 00058 */ 00059 static void signalHandler( int signum ) 00060 { 00061 std::cout << std::endl; 00062 std::cerr << "*** Signal " << signum << " received" << std::endl; 00063 00064 if( knownIntegrators.size () == 0 ) { 00065 exit( -1 ); 00066 } 00067 for( auto i: knownIntegrators ) { 00068 i->quit (); 00069 } 00070 } 00071 00072 /** Report the integration time and the estimated real-time to finish. 00073 */ 00074 void reportIntegrationTime( Int m ) 00075 { 00076 if ( mpiRank() > 1 ) { 00077 return; // only the first in rank can access the standard output 00078 } 00079 00080 Real tnow = GF( fld::t, m, nGhost ); // current integration time 00081 00082 auto elapsed = 1e-3 * std::chrono::duration_cast<std::chrono::milliseconds> 00083 ( std::chrono::high_resolution_clock::now () - rt_start ).count (); 00084 00085 auto rt_factor = elapsed / ( tnow - t_0 ); // real time/integration time 00086 auto rt_total = round( ( t_1 - t_0 ) * rt_factor ); // projected total time 00087 auto rt_left = round( ( t_1 - tnow ) * rt_factor ); // estimated time to cover 00088 00089 std::cout << " t = " << std::setw(6) << std::setprecision(2) << std::fixed 00090 << tnow << std::setw(-1) << std::setprecision(-1) << std::defaultfloat 00091 << ", real = " << round( elapsed ) 00092 << " s (" << round( 100 * ( tnow - t_0 ) / ( t_1 - t_0 ) ) 00093 << "% of " << rt_total 00094 << " s), ETA " << rt_left << " s" 00095 << " \r"; 00096 std::flush( std::cout ); 00097 } 00098 00099 /** Executed at the beginning of each integration step. 00100 */ 00101 void integStep_Begin( Int m ) 00102 { 00103 /// - Prepare the intermediate variables. 00104 /// 00105 for( auto eom: eomList ) { 00106 eom->integStep_Prepare( m ); 00107 } 00108 00109 /// - Exchange boundaries with the neighboring MPI ranks and/or fix the 00110 /// left/right boundary conditions (depending on our MPI rank). 00111 /// Quit if MPI has failed or we received abort message. 00112 /// 00113 if( gridDriver->mpiSize() > 1 00114 && ! gridDriver->exchangeBoundaries( m ) ) 00115 { 00116 running = false; 00117 return; 00118 } 00119 if ( gridDriver->isFirstInRank () ) { 00120 for( auto eom: eomList ) { 00121 eom->applyLeftBoundaryCondition( m ); 00122 } 00123 } 00124 if ( gridDriver->isLastInRank () ) { 00125 for( auto eom: eomList ) { 00126 eom->applyRightBoundaryCondition( m ); 00127 } 00128 } 00129 00130 /// - Calculate the RHS of time evolution equations. 00131 /// 00132 for( auto eom: eomList ) { 00133 eom->integStep_CalcEvolutionRHS( m ); 00134 } 00135 } 00136 00137 /** Executed at the end of each integration step. 00138 * Returns `false` if we need to break the integration loop. 00139 */ 00140 void integStep_End( Int m1, Int m ) 00141 { 00142 /// - Finalize the integration step (e.g., doing diagnostic and post 00143 /// processing of EoM like calculating the constraint violations). 00144 /// 00145 for( auto eom: eomList ) { 00146 eom->integStep_Finalize( m1, m ); 00147 } 00148 00149 /// - Wait on finishing the exchange of the boundaries (when using MPI). 00150 /// 00151 if ( ! gridDriver->waitExchange () ) { 00152 quit( m ); // Quit if MPI has failed 00153 } 00154 } 00155 00156 /** Executed at each checkpoint. 00157 */ 00158 void integStep_Checkpoint( Int m ) 00159 { 00160 /// - Report the elapsed (both integration and real) time. 00161 /// 00162 reportIntegrationTime( m ); 00163 00164 /// - Check for eventual NaNs in the central/lower region of the grid. 00165 /// 00166 for( auto eom: eomList ) 00167 { 00168 if ( ! eom->integStep_Diagnostics( m, 00169 output->get_nFrom (), output->get_nTo () ) ) { 00170 if ( running ) { // quit only once 00171 quit( m ); 00172 } 00173 } 00174 } 00175 00176 /// - Send the grid slice to the output. 00177 /// 00178 output->write( m, cur_t, delta_t ); 00179 } 00180 00181 /** Evolves the state variables using the FT scheme. 00182 */ 00183 void integStep_Euler( Int mNext, Int m ); 00184 00185 /** Evolves the sate variables using the ICN averaging. 00186 */ 00187 void integStep_AvgICN( Int mNext, Int mMid, Int m ); 00188 00189 /** Time evolution using the iterated Crank-Nicolson (ICN). 00190 */ 00191 void integrate_AvgICN( int ICN_steps ); 00192 00193 /** Prepares the state variables for the intermediate MoL steps. 00194 */ 00195 void integStep_MoLInit( Int m1, Int m, Real beta_delta_t ); 00196 00197 /** The intermediate method of lines (MoL) steps. 00198 */ 00199 void integStep_MoL( Int m1, Int m, Real alpha ); 00200 00201 /** Time evolution using the method of lines (MoL). 00202 */ 00203 void integrate_MoL( const MoLDescriptor& MoL ); 00204 00205 public: 00206 00207 /** Constructs the MoL integrator as specified in the parameter file. 00208 */ 00209 MoL( Parameters& params, UniformGrid& ug, GridOutputWriter& out ) 00210 : GridUser( ug ), output( &out ) 00211 { 00212 constantGF.reserve( 64 ); 00213 evolvedGF.reserve( 128 ); 00214 00215 constantGF.push_back( fld::r ); // Radial coordinate is kept constant by default 00216 00217 params.get( "integ.t_0", t_0, 0.0 ); 00218 params.get( "integ.t_1", t_1, 1.0 ); 00219 params.get( "integ.delta_t", delta_t, 0.1 ); 00220 params.get( "integ.dissip", dissip, 0.0 ); 00221 00222 std::string name = params.get( "integ.method", method, 0, knownMethods ); 00223 00224 // Fix the sign of delta_t 00225 // 00226 delta_t = t_1 < t_0 ? -abs(delta_t) : abs(delta_t); 00227 00228 cur_t = t_0; // The initial `t` 00229 running = false; // Stopped integration by default (will be enabled later on) 00230 00231 // Calculate Kreiss-Oliger dissipation corrected delta r 00232 // 00233 dissip_delta_r = dissip / delta_r; 00234 00235 // Expected output file data size (in bytes) 00236 // 00237 const Real sz = ( ( t_1 - t_0 ) / delta_t / output->get_mSkip() + 1 ) 00238 * output->recordSizeInBytes(); 00239 00240 slog << "Integrator:" << std::endl << std::endl 00241 << " t_0 = " << t_0 << ", t_1 = " << t_1 << ", delta_t = " << delta_t 00242 << ", method = " << name << " (#" << method << ")" << std::endl 00243 << " Kreiss-Oliger dissipation = " << dissip << std::endl 00244 << " Expected output data size = " 00245 << round( ( sz < 1e3 ? sz : sz < 1e6 ? sz/1e3 : sz/1e6 ) * 10 ) / 10 00246 << ( sz < 1e3 ? "" : sz < 1e6 ? " kB" : " MB" ) 00247 << ", GF count = " << output->GFCount () 00248 << std::endl << std::endl; 00249 00250 // Register to receive the signals from the system 00251 // 00252 knownIntegrators.push_back( this ); 00253 if ( knownIntegrators.size () == 1 ) { 00254 signal( SIGINT, signalHandler ); 00255 signal( SIGTERM, signalHandler ); 00256 } 00257 } 00258 00259 /** The destructor (here only responsible for cleaning up standard output). 00260 */ 00261 ~MoL() 00262 { 00263 if( mpiRank() == 0 ) { 00264 std::cout << std::endl << std::endl; 00265 } 00266 } 00267 00268 /** Moves the integration final time to the current time and marks 00269 * the integration stopped. 00270 */ 00271 void quit( Int m = -1 ) 00272 { 00273 running = false; // Stop the integration 00274 00275 t_1 = m < 0 ? cur_t : GF( fld::t, m, nGhost ); 00276 00277 if( mpiRank() == 0 ) std::cout << std::endl; 00278 slog << "*** Quitting at t = " << t_1 << std::endl; 00279 00280 /// - In case of NaNs, signal the other MPI ranks that we are quitting. 00281 /// @todo Distribute the abort msg to all ranks instead of calling MPI_Abort. 00282 /// 00283 gridDriver->abortExchange (); 00284 } 00285 00286 /** Adds the EoM to the list of all evolved subsystems. 00287 */ 00288 void addToEvolution( IntegFace* eomInterface ) 00289 { 00290 eomList.push_back( eomInterface ); 00291 } 00292 00293 /** The given fields will be kept constant by the integrator. 00294 */ 00295 void keepConstant( const std::vector<Int>& gfs ) 00296 { 00297 constantGF.insert( constantGF.end(), gfs.begin(), gfs.end() ); 00298 } 00299 00300 /** Add the given functions to the evolution list. 00301 */ 00302 void keepEvolved( const std::vector<fld::EvolvedBy>& gfs ) 00303 { 00304 evolvedGF.insert( evolvedGF.end(), gfs.begin(), gfs.end() ); 00305 } 00306 00307 /** Returns the time step. 00308 */ 00309 Real dt () const { 00310 return delta_t; 00311 } 00312 00313 /** Integrates the evolution equations in the given grid realm. 00314 */ 00315 bool evolveEquations () 00316 { 00317 #if 0 00318 if( delta_r != r(0,nGhost+1) - r(0,nGhost) ) { 00319 slog << " Error: delta_r = " << delta_r << " != grid spacing = " 00320 << r(0,nGhost+1) - r(0,nGhost) 00321 << ", diff = " << delta_r - ( r(0,nGhost+1) - r(0,nGhost) ) 00322 << std::endl; 00323 } 00324 #endif 00325 00326 output->writeHeader (); 00327 00328 // Start timing the real-time 00329 rt_start = std::chrono::high_resolution_clock::now (); 00330 00331 switch( method ) 00332 { 00333 case 0: integrate_AvgICN ( 0 ); break; 00334 case 1: integrate_AvgICN ( 1 ); break; 00335 case 2: integrate_AvgICN ( 2 ); break; 00336 case 3: integrate_AvgICN ( 3 ); break; 00337 case 4: integrate_MoL ( ICN2 ); break; 00338 case 5: integrate_MoL ( ICN3 ); break; 00339 case 6: integrate_MoL ( RK1 ); break; 00340 case 7: integrate_MoL ( RK2 ); break; 00341 case 8: integrate_MoL ( RK3 ); break; 00342 case 9: integrate_MoL ( RK4 ); break; 00343 } 00344 00345 return true; 00346 } 00347 }; 00348 00349 ///////////////////////////////////////////////////////////////////////////////////////// 00350 // Definitions of MoL static members 00351 ///////////////////////////////////////////////////////////////////////////////////////// 00352 00353 std::vector<MoL*> MoL::knownIntegrators; 00354 00355 std::map<std::string,int> MoL::knownMethods = 00356 { 00357 { "Euler", 0 }, { "avgICN2", 1 }, 00358 { "avgICN3", 2 }, { "avgICN4", 3 }, 00359 { "ICN2", 4 }, { "ICN3", 5 }, 00360 { "RK1", 6 }, { "RK2", 7 }, 00361 { "RK3", 8 }, { "RK4", 9 } 00362 }; 00363 00364 ///////////////////////////////////////////////////////////////////////////////////////// 00365 // Definitions of MoL methods 00366 ///////////////////////////////////////////////////////////////////////////////////////// 00367 00368 void MoL::integStep_Euler( Int m1, Int m ) 00369 { 00370 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00371 { 00372 GF( fld::t, m1, n ) = GF( fld::t, m, n ) + delta_t; // Evolve time 00373 00374 for( auto f: constantGF ) 00375 { 00376 GF( f, m1, n ) = GF( f, m, n ); 00377 } 00378 00379 for( auto e: evolvedGF ) 00380 { 00381 GF( e.f, m1, n ) = GF( e.f, m, n ) + delta_t * GF( e.f_t, m, n ) 00382 - KreissOligerTerm( e.f, delta_t ); 00383 } 00384 } 00385 } 00386 00387 void MoL::integStep_AvgICN 00388 ( 00389 Int m2, // final time 00390 Int m1, // middle time 00391 Int m // initial time 00392 ) 00393 { 00394 const Real half_dt = delta_t / 2; 00395 00396 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00397 { 00398 GF( fld::t, m2, n ) = GF( fld::t, m, n ) + delta_t; // Evolve time 00399 00400 for( auto f: constantGF ) 00401 { 00402 GF( f, m2, n ) = GF( f, m, n ); 00403 } 00404 00405 for( auto e: evolvedGF ) 00406 { 00407 GF( e.f, m2, n ) = GF( e.f, m, n ) 00408 + half_dt * ( GF( e.f_t, m1, n) + GF( e.f_t, m, n ) ) 00409 - KreissOligerTerm( e.f, delta_t ); 00410 } 00411 } 00412 } 00413 00414 /** Integrates the evolution equations using the iterated Crank-Nicolson (ICN) method. 00415 * If ICN_steps == 0, then the method effectively becomes FTCS. 00416 * See @cite Press:2007nrec, @cite Baumgarte:2010nr 00417 */ 00418 void MoL::integrate_AvgICN( int ICN_steps ) 00419 { 00420 slog << "Employing Method of Lines, Iterative Crank-Nicolson (avg)" 00421 << ", N = " << ICN_steps << std::endl << std::endl; 00422 00423 cur_t = t_0; // The initial `t` 00424 running = true; // Enable integration 00425 00426 for( Int n = 0; n < nTotal; ++n ) { 00427 GF( fld::t, 0, n ) = t_0; 00428 } 00429 00430 Int mStep = 0; // The step counter (used to filter checkpoints) 00431 while( running && abs(cur_t) < abs(t_1) ) 00432 { 00433 for( Int m = 0; running && m < mLen && abs(cur_t) < abs(t_1); ++m ) 00434 { 00435 Int mNext = m + 1 >= mLen ? 0 : m + 1; 00436 00437 ///////////////////////////////////////////////////////////////////////////// 00438 // Predictor. Use FTCS for the first step. 00439 // Note that this will be the only step, if ICN_steps == 0. 00440 // 00441 Int m_1 = ICN_steps >= 1 ? mLen + 1 : mNext; 00442 00443 integStep_Begin ( m ); if( ! running ) break; 00444 integStep_Euler ( m_1, m ); 00445 integStep_End ( m_1, m ); 00446 00447 ///////////////////////////////////////////////////////////////////////////// 00448 // Corrector iterations. There will be no iterations, if ICN_steps == 0. 00449 // 00450 for( int i = 0; running && i < ICN_steps; ++i ) 00451 { 00452 Int m_2 = ( i + 1 >= ICN_steps ? mNext : m_1 + 1 ); 00453 00454 integStep_Begin ( m_1 ); if( ! running ) break; 00455 integStep_AvgICN ( m_2, m_1, m ); 00456 integStep_End ( m_2, m_1 ); 00457 00458 m_1 = m_2; 00459 } 00460 00461 if( ( mStep++ % output->get_mSkip () ) == 0 ) { 00462 integStep_Checkpoint( m ); 00463 } 00464 00465 cur_t = GF( fld::t, m, nGhost ); // The current time step 00466 } 00467 } 00468 } 00469 00470 void MoL::integStep_MoL( Int m1, Int m, Real ac ) 00471 { 00472 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00473 { 00474 GF( fld::t, m1, n ) += ac * GF( fld::t, m, n ); // Evolve time 00475 00476 for( auto e: evolvedGF ) 00477 { 00478 GF( e.f, m1, n ) += ac * GF( e.f, m, n ); 00479 } 00480 } 00481 } 00482 00483 void MoL::integStep_MoLInit( Int m1, Int m, Real bc_delta_t ) 00484 { 00485 OMP_parallel_for( Int n = nGhost; n < nGhost + nLen; ++n ) 00486 { 00487 // There is no explicit dependence on t, and t_t(m,n) = 1. 00488 // Also, r is constant along the `line`, i.e., r_t(m,n) = 0. 00489 // 00490 GF( fld::t, m1, n ) = bc_delta_t; // Evolve time 00491 00492 for( auto f: constantGF ) 00493 { 00494 GF( f, m1, n ) = GF( f, m, n ); 00495 } 00496 00497 for( auto e: evolvedGF ) 00498 { 00499 GF( e.f, m1, n ) = bc_delta_t * GF( e.f_t, m, n ) 00500 - KreissOligerTerm( e.f, bc_delta_t ); 00501 } 00502 } 00503 } 00504 00505 /** The method of lines (MoL) separates the time integration from the rest of 00506 * an evolution scheme. A `N`-step MoL method evolves the equation: 00507 * <pre> 00508 * partial_t Y = L[Y], </pre> 00509 * 00510 * using the following algorithm: 00511 * <pre> 00512 * Y^{(0)} = Y^{m}, 00513 * Y^{(i+1)} = Sum_{j=0}^{i} ( alpha_{ij} Y^{(j)} ) + <- integStep_MoL 00514 * Delta t beta_i L[Y^{(i)}], for i = 0, ..., N-1, <- integStep_MoLInit 00515 * Y^{m+1} = Y^{(N)}. </pre> 00516 * 00517 * The variables `Y^{(i)}`, `i = 0, ..., N`, are intermediate. 00518 * 00519 * The method is completely specified by `N`, and arrays alpha and beta 00520 * in the structure MoLDescriptor. 00521 */ 00522 void MoL::integrate_MoL( const MoLDescriptor& MoL ) 00523 { 00524 slog << "Employing Method of Lines, " << MoL.name 00525 << ", N = " << MoL.N << std::endl << std::endl; 00526 00527 cur_t = t_0; // The initial `t` 00528 running = true; // Enable integration 00529 00530 for( Int n = 0; n < nTotal; ++n ) { 00531 GF( fld::t, 0, n ) = t_0; 00532 } 00533 00534 // The array mk[] contains redirections to grid rows (m's): 00535 // mk[0] points to Y^{m} == Y^{(0)}, 00536 // mk[1], ..., mk[Mol.N-1] point to the intermediate Y^{(1)}, ..., Y^{(N-1)}, 00537 // mk[MoL.N] points to Y^{m+1} == Y^{(N)}. 00538 // 00539 int mk[ 10 ]; // int mk[ MoL.N + 1 ]; // Note that N is at most 8 00540 for( int i = 0; i <= MoL.N; ++i ) { 00541 mk[i] = mLen + i; 00542 } 00543 00544 Int mStep = 0; // The step counter 00545 while( running && abs(cur_t) < abs(t_1) ) 00546 { 00547 for( Int m = 0; running && m < mLen && abs(cur_t) < abs(t_1); ++m ) 00548 { 00549 mk[0] = m; // Y^{(m)} 00550 mk[MoL.N] = m + 1 >= mLen ? 0 : m + 1; // Y^{(m+1)} 00551 00552 for( int i = 0; running && i < MoL.N; ++i ) 00553 { 00554 integStep_Begin( mk[i] ); if( ! running ) break; 00555 00556 integStep_MoLInit( mk[i+1], mk[i], MoL.beta[i] * delta_t ); 00557 00558 for( int j = 0; j <= i; j++ ) { 00559 integStep_MoL( mk[i+1], mk[j], MoL.alpha[i][j] ); 00560 } 00561 00562 integStep_End( mk[i+1], mk[i] ); 00563 } 00564 00565 if( ( mStep++ % output->get_mSkip () ) == 0 ) { 00566 integStep_Checkpoint( m ); 00567 } 00568 00569 cur_t = GF( fld::t, m, nGhost ); // The current time step 00570 } 00571 } 00572 } 00573 /*@}*/ 00574 /////////////////////////////////////////////////////////////////////////////////////////