Nektar++
DriverParareal.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File DriverParareal.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Driver class for the parareal solver
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <iomanip>
36
39#include <boost/format.hpp>
40
41namespace Nektar::SolverUtils
42{
43std::string DriverParareal::className =
48
49/**
50 *
51 */
55 : DriverParallelInTime(pSession, pGraph)
56{
57}
58
59/**
60 *
61 */
62void DriverParareal::v_InitObject(std::ostream &out)
63{
69}
70
71/**
72 *
73 */
74void DriverParareal::v_Execute([[maybe_unused]] std::ostream &out)
75{
76 // Timing.
78 NekDouble totalTime = 0.0, predictorTime = 0.0, coarseSolveTime = 0.0,
79 fineSolveTime = 0.0, correctionTime = 0.0;
80
81 // Get and assert parameters from session file.
83
84 // Initialie time step parameters.
89
90 // Pre-solve for one time-step to initialize preconditioner.
92
93 // Start iteration windows.
94 m_comm->GetTimeComm()->Block();
96 for (size_t w = 0; w < m_numWindowsPIT; w++)
97 {
98 timer.Start();
99 // Initialize time for the current window.
101
102 // Print window number.
103 PrintHeader((boost::format("WINDOWS #%1%") % (w + 1)).str(), '*');
104
105 // Update coarse initial condition.
107
108 // Run predictor.
109 for (size_t i = 0; i < m_nVar; ++i)
110 {
111 RecvFromPreviousProc(m_EqSys[m_coarseLevel]->UpdatePhysField(i));
112 }
113 if (m_chunkRank > 0)
114 {
116 }
118 for (size_t i = 0; i < m_nVar; ++i)
119 {
120 SendToNextProc(m_EqSys[m_coarseLevel]->UpdatePhysField(i));
121 }
122
123 // Interpolate coarse solution.
125
126 // Compute exact solution, if necessary.
127 if (m_exactSolution)
128 {
130 }
131 timer.Stop();
132 predictorTime += timer.Elapsed().count();
133 totalTime += timer.Elapsed().count();
134
135 // Solution convergence monitoring.
136 timer.Start();
139 timer.Stop();
140 totalTime += timer.Elapsed().count();
141 if (m_chunkRank == m_numChunks - 1 &&
142 m_comm->GetSpaceComm()->GetRank() == 0)
143 {
144 std::cout << "Total Computation Time : " << totalTime << "s"
145 << std::endl
146 << std::flush;
147 }
148
149 // Start Parareal iteration.
150 size_t iter = 1;
151 int convergenceCurr = false;
152 int convergencePrev = (m_chunkRank == 0);
153 while (iter <= m_iterMaxPIT && !convergenceCurr)
154 {
155 // Use previous parareal solution as "exact solution", if necessary.
156 timer.Start();
157 if (!m_exactSolution)
158 {
160 }
161 timer.Stop();
162 totalTime += timer.Elapsed().count();
163
164 // Calculate fine solution (parallel-in-time).
165 timer.Start();
168 timer.Stop();
169 fineSolveTime += timer.Elapsed().count();
170 totalTime += timer.Elapsed().count();
171
172 // Compute F -> F - Gold
173 timer.Start();
175 timer.Stop();
176 correctionTime += timer.Elapsed().count();
177 totalTime += timer.Elapsed().count();
178
179 // Receive coarse solution from previous processor.
180 timer.Start();
182 timer.Stop();
183 totalTime += timer.Elapsed().count();
184
185 // Calculate coarse solution (serial-in-time).
186 timer.Start();
189 iter);
190 timer.Stop();
191 coarseSolveTime += timer.Elapsed().count();
192 totalTime += timer.Elapsed().count();
193
194 // Compute F -> F + Gnew
195 timer.Start();
197 timer.Stop();
198 correctionTime += timer.Elapsed().count();
199 totalTime += timer.Elapsed().count();
200
201 // Solution convergence monitoring.
204 if (m_chunkRank == m_numChunks - 1 &&
205 m_comm->GetSpaceComm()->GetRank() == 0)
206 {
207 std::cout << "Total Computation Time : " << totalTime << "s"
208 << std::endl
209 << std::flush;
210 std::cout << " - Predictor Time : " << predictorTime << "s"
211 << std::endl
212 << std::flush;
213 std::cout << " - Coarse Solve Time : " << coarseSolveTime << "s"
214 << std::endl
215 << std::flush;
216 std::cout << " - Fine Solve Time : " << fineSolveTime << "s"
217 << std::endl
218 << std::flush;
219 std::cout << " - Correction Time : " << correctionTime << "s"
220 << std::endl
221 << std::flush;
222 }
223
224 // Check convergence of L2 error for each time chunk.
225 convergenceCurr = (vL2ErrorMax() < m_tolerPIT && convergencePrev) ||
226 (m_chunkRank + 1 == iter);
227
228 // Send solution to next processor.
229 timer.Start();
230 SendToNextProc(m_fineSolution, convergenceCurr);
231 timer.Stop();
232 totalTime += timer.Elapsed().count();
233
234 // Increment iteration index.
235 iter++;
236 }
237
238 // Copy converged check points.
240
241 // Write time chunk solution to files.
243
244 // Apply windowing.
245 timer.Start();
247 timer.Stop();
248 totalTime += timer.Elapsed().count();
249 }
250
251 m_comm->GetTimeComm()->Block();
252 PrintHeader("SUMMARY", '*');
255 if (m_chunkRank == m_numChunks - 1 &&
256 m_comm->GetSpaceComm()->GetRank() == 0)
257 {
258 std::cout << "Total Computation Time : " << totalTime << "s"
259 << std::endl
260 << std::flush;
261 std::cout << " - Predictor Time : " << predictorTime << "s" << std::endl
262 << std::flush;
263 std::cout << " - Coarse Solve Time : " << coarseSolveTime << "s"
264 << std::endl
265 << std::flush;
266 std::cout << " - Fine Solve Time : " << fineSolveTime << "s"
267 << std::endl
268 << std::flush;
269 std::cout << " - Correction Time : " << correctionTime << "s"
270 << std::endl
271 << std::flush;
272 }
273}
274
275/**
276 *
277 */
279{
280 // Allocate storage for Parareal solver.
284 for (size_t i = 0; i < m_nVar; ++i)
285 {
288 m_fineSolution[i] = m_EqSys[m_fineLevel]->UpdatePhysField(i);
291 ? m_EqSys[m_coarseLevel]->UpdatePhysField(i)
293 }
294}
295
296/**
297 *
298 */
300{
301 // Assert time-stepping parameters
302 ASSERTL0(
304 "Total number of fine step should be divisible by number of chunks.");
305
306 ASSERTL0(
308 "Total number of coarse step should be divisible by number of chunks.");
309
311 "Total number of fine step should be divisible by number of "
312 "windows times number of chunks.");
313
315 "Total number of coarse step should be divisible by number of "
316 "windows times number of chunks.");
317
320 "Fine and coarse total computational times do not match");
321
323 ->GetTimeIntegrationScheme()
324 ->GetNumIntegrationPhases() == 1,
325 "Only single step time-integration schemes currently supported "
326 "for Parareal");
327
329 ->GetTimeIntegrationScheme()
330 ->GetNumIntegrationPhases() == 1,
331 "Only single step time-integration schemes currently supported "
332 "for Parareal");
333
334 // Assert I/O parameters
335 if (m_EqSys[m_fineLevel]->GetInfoSteps())
336 {
337 ASSERTL0(m_nsteps[m_fineLevel] % (m_EqSys[m_fineLevel]->GetInfoSteps() *
339 0,
340 "number of IO_InfoSteps should divide number of fine steps "
341 "per time chunk");
342 }
343
344 if (m_EqSys[m_coarseLevel]->GetInfoSteps())
345 {
347 (m_EqSys[m_coarseLevel]->GetInfoSteps() * m_numChunks *
349 0,
350 "number of IO_InfoSteps should divide number of coarse steps "
351 "per time chunk");
352 }
353
354 if (m_EqSys[m_fineLevel]->GetCheckpointSteps())
355 {
357 (m_EqSys[m_fineLevel]->GetCheckpointSteps() *
359 0,
360 "number of IO_CheckSteps should divide number of fine steps "
361 "per time chunk");
362 }
363
364 if (m_EqSys[m_coarseLevel]->GetCheckpointSteps())
365 {
367 (m_EqSys[m_coarseLevel]->GetCheckpointSteps() *
369 0,
370 "number of IO_CheckSteps should divide number of coarse steps "
371 "per time chunk");
372 }
373}
374
375/**
376 *
377 */
379{
380 // Interpolate solution to fine field.
381 Interpolate(m_EqSys[timeLevel]->UpdateFields(),
384}
385
386/**
387 *
388 */
390{
391 // Restrict fine field to coarse solution.
392 Interpolate(m_EqSys[m_fineLevel]->UpdateFields(),
393 m_EqSys[timeLevel]->UpdateFields(), m_initialCondition,
395}
396
397/**
398 *
399 */
400void DriverParareal::UpdateSolution(const size_t timeLevel,
401 const NekDouble time, const size_t nstep,
402 const size_t wd, const size_t iter)
403{
404 // Number of checkpoint by chunk.
405 size_t nChkPts =
406 m_EqSys[timeLevel]->GetCheckpointSteps()
407 ? m_nsteps[timeLevel] / m_EqSys[timeLevel]->GetCheckpointSteps()
408 : 1;
409
410 // Checkpoint index.
411 size_t iChkPts = (m_chunkRank + wd * m_numChunks) * nChkPts + 1;
412
413 // Reinitialize check point number for each parallel-in-time
414 // iteration.
415 m_EqSys[timeLevel]->SetCheckpointNumber(iChkPts);
416
417 // Update parallel-in-time iteration number.
418 m_EqSys[timeLevel]->SetIterationNumberPIT(iter);
419
420 // Update parallel-in-time window number.
421 m_EqSys[timeLevel]->SetWindowNumberPIT(wd);
422
423 m_EqSys[timeLevel]->SetTime(time);
424 m_EqSys[timeLevel]->SetSteps(nstep);
425 m_EqSys[timeLevel]->DoSolve();
426}
427
428/**
429 *
430 */
432{
433 // Correct solution F -> F - Gold.
434 for (size_t i = 0; i < m_nVar; ++i)
435 {
437 m_coarseSolution[i], 1, m_fineSolution[i], 1);
438 }
439}
440
441/**
442 *
443 */
445{
446 // Interpolate coarse solution.
448
449 // Correct solution F -> F + Gnew.
450 for (size_t i = 0; i < m_nVar; ++i)
451 {
453 m_coarseSolution[i], 1, m_fineSolution[i], 1);
454 }
455}
456
457/**
458 *
459 */
461{
463 {
464 // Interpolate coarse solution to fine field.
465 Interpolate(m_EqSys[m_coarseLevel]->UpdateFields(),
466 m_EqSys[m_fineLevel]->UpdateFields(),
468 }
469}
470
471/**
472 *
473 */
475{
476 if (w == m_numWindowsPIT - 1)
477 {
478 // No windowing required for the last window.
479 return;
480 }
481
482 // Use last chunk solution as initial condition for the next
483 // window.
484 if (m_chunkRank == m_numChunks - 1)
485 {
486 for (size_t i = 0; i < m_nVar; ++i)
487 {
489 m_EqSys[m_fineLevel]->UpdatePhysField(i), 1,
490 m_initialCondition[i], 1);
491 }
492 }
493
494 // Broadcast I.C. for windowing.
495 for (size_t i = 0; i < m_nVar; ++i)
496 {
497 m_comm->GetTimeComm()->Bcast(m_initialCondition[i], m_numChunks - 1);
498 }
499}
500
501/**
502 *
503 */
504void DriverParareal::CopyConvergedCheckPoints(const size_t w, const size_t k)
505{
506 // Determine max number of iteration.
507 size_t kmax = k;
508 m_comm->GetTimeComm()->AllReduce(kmax, Nektar::LibUtilities::ReduceMax);
509
510 if (m_comm->GetSpaceComm()->GetRank() == 0)
511 {
512 for (size_t j = k; j < kmax; j++)
513 {
514 // Copy converged solution files from directory corresponding to
515 // iteration j - 1 to the directory corresponding to iteration j.
516
517 auto sessionName = m_EqSys[m_fineLevel]->GetSessionName();
518
519 // Input directory name.
520 std::string indir =
521 sessionName + "_" + std::to_string(j - 1) + ".pit";
522
523 /// Output directory name.
524 std::string outdir = sessionName + "_" + std::to_string(j) + ".pit";
525
526 for (size_t timeLevel = 0; timeLevel < m_nTimeLevel; timeLevel++)
527 {
528 // Number of checkpoint by chunk.
529 size_t nChkPts =
530 m_EqSys[timeLevel]->GetCheckpointSteps()
531 ? m_nsteps[timeLevel] /
532 m_EqSys[timeLevel]->GetCheckpointSteps()
533 : 0;
534
535 // Checkpoint index.
536 size_t iChkPts = (m_chunkRank + w * m_numChunks) * nChkPts;
537
538 for (size_t i = 1; i <= nChkPts; i++)
539 {
540 // Filename corresponding to checkpoint iChkPts.
541 std::string filename = sessionName + "_timeLevel" +
542 std::to_string(timeLevel) + "_" +
543 std::to_string(iChkPts + i) + ".chk";
544
545 // Intput full file name.
546 std::string infullname = indir + "/" + filename;
547
548 // Output full file name.
549 std::string outfullname = outdir + "/" + filename;
550
551 // Remove output file if already existing.
552 fs::remove_all(outfullname);
553
554 // Copy converged solution files.
555 fs::copy(infullname, outfullname);
556 }
557 }
558 }
559 }
560}
561
562/**
563 *
564 */
566{
567 PrintHeader("PRINT SOLUTION FILES", '-');
568
569 // Update field coefficients.
571
572 // Output solution files.
573 m_EqSys[m_fineLevel]->Output();
574}
575
576} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:80
Base class for the development of parallel-in-time solvers.
NekDouble m_totalTime
Total time integration interval.
void UpdateFieldCoeffs(const size_t timeLevel, const Array< OneD, const Array< OneD, NekDouble > > &in=NullNekDoubleArrayOfArray)
NekDouble m_chunkTime
Time integration interval per chunk.
void SendToNextProc(Array< OneD, Array< OneD, NekDouble > > &array, int &convergence)
size_t m_iterMaxPIT
Maximum number of parallel-in-time iteration.
Array< OneD, size_t > m_nsteps
Number of time steps for each time level.
SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
void EvaluateExactSolution(const size_t timeLevel, const NekDouble &time)
void RecvFromPreviousProc(Array< OneD, Array< OneD, NekDouble > > &array, int &convergence)
Array< OneD, std::shared_ptr< UnsteadySystem > > m_EqSys
Equation system to solve.
void PrintHeader(const std::string &title, const char c)
NekDouble m_tolerPIT
ParallelInTime tolerance.
void PrintSolverInfo(std::ostream &out=std::cout)
bool m_exactSolution
Using exact solution to compute error norms.
void CopyFromPhysField(const size_t timeLevel, Array< OneD, Array< OneD, NekDouble > > &out)
void CopyToPhysField(const size_t timeLevel, const Array< OneD, const Array< OneD, NekDouble > > &in)
Array< OneD, size_t > m_npts
Number of dof for each time level.
Array< OneD, Array< OneD, NekDouble > > m_exactsoln
void SolutionConvergenceSummary(const size_t timeLevel)
Array< OneD, NekDouble > m_timestep
Time step for each time level.
void SolutionConvergenceMonitoring(const size_t timeLevel, const size_t iter)
void Interpolate(const Array< OneD, MultiRegions::ExpListSharedPtr > &infield, const Array< OneD, MultiRegions::ExpListSharedPtr > &outfield, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SOLVER_UTILS_EXPORT void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
static constexpr size_t m_fineLevel
void UpdateInitialConditionFromSolver(const size_t timeLevel)
static constexpr size_t m_coarseLevel
static std::string className
Name of the class.
void CopyConvergedCheckPoints(const size_t w, const size_t k)
SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
void UpdateSolution(const size_t timeLevel, const NekDouble time, const size_t nstep, const size_t wd, const size_t iter)
Array< OneD, Array< OneD, NekDouble > > m_coarseSolution
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
SOLVER_UTILS_EXPORT DriverParareal(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.
void UpdateSolverInitialCondition(const size_t timeLevel)
Array< OneD, Array< OneD, NekDouble > > m_fineSolution
Array< OneD, Array< OneD, NekDouble > > m_initialCondition
def copy(self)
Definition: pycml.py:2663
std::shared_ptr< SessionReader > SessionReaderSharedPtr
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:66
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::vector< double > w(NPUPPER)
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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.
Definition: Vmath.hpp:220