Nektar++
DriverParallelInTime.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File DriverParallelInTime.h
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 parallel-in-time solver
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_SOLVERUTILS_DRIVERPARALLELINTIME_H
36#define NEKTAR_SOLVERUTILS_DRIVERPARALLELINTIME_H
37
38#include <SolverUtils/Driver.h>
40
41namespace Nektar
42{
43namespace SolverUtils
44{
45
46/// Base class for the development of parallel-in-time solvers.
48{
49public:
50protected:
51 /// Constructor
55
56 /// Destructor
58
59 /// Virtual function for initialisation implementation.
61 std::ostream &out = std::cout) override;
62
63 /// Virtual function for solve implementation.
64 SOLVER_UTILS_EXPORT virtual void v_Execute(
65 std::ostream &out = std::cout) override;
66
68
70
72
74
76
78
80
82 const size_t iter, NekDouble fineSolveTime, NekDouble coarseSolveTime,
83 NekDouble restTime, NekDouble interTime, NekDouble commTime,
84 NekDouble predictorTime, NekDouble overheadTime);
85
86 void SetParallelInTimeEquationSystem(std::string AdvectiveType);
87
88 void GetParametersFromSession(void);
89
90 void InitialiseEqSystem(bool turnoff_output);
91
92 void AllocateMemory(void);
93
95
96 void PrintCoarseSolverInfo(std::ostream &out = std::cout);
97
98 void PrintFineSolverInfo(std::ostream &out = std::cout);
99
100 void PrintHeaderTitle1(const std::string &title);
101
102 void PrintHeaderTitle2(const std::string &title);
103
104 void PrintComputationalTime(const NekDouble time);
105
107 Array<OneD, Array<OneD, NekDouble>> &array, int &convergence);
108
111
113 int &convergence);
114
116
119
121
123
125 const Array<OneD, const Array<OneD, NekDouble>> &in);
126
128 const Array<OneD, const Array<OneD, NekDouble>> &in);
129
130 void UpdateSolution(const Array<OneD, const Array<OneD, NekDouble>> &in);
131
132 void EvaluateExactSolution(const NekDouble &time);
133
135
136 void SolutionConvergenceSummary(const NekDouble &time);
137
138 void UpdateErrorNorm(const bool normalized);
139
140 void PrintErrorNorm(const bool normalized);
141
142 void Interpolator(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
143 Array<OneD, Array<OneD, NekDouble>> &outarray);
144
146
147 void SpeedUpAnalysis();
148
149 void PrintSpeedUp(NekDouble fineSolveTime, NekDouble coarseSolveTime,
150 NekDouble restTime, NekDouble interTime,
151 NekDouble commTime, NekDouble predictorOverheadTime,
152 NekDouble overheadTime);
153
157
158 /// Timestep for fine solver.
160
161 /// Timestep for coarse solver.
163
164 /// Total time integration interval.
166
167 /// Time for chunks
169
170 /// Local time
172
173 /// Number of steps for info I/O
174 size_t m_infoSteps = 0;
175
176 /// Number of steps for checkpoint
177 size_t m_checkSteps = 0;
178
179 /// Number of steps for the fine solver
180 size_t m_fineSteps = 1;
181
182 /// Number of steps for the coarse solver
183 size_t m_coarseSteps = 1;
184
185 /// Number of time chunks
186 size_t m_numChunks = 1;
187
188 /// Rank in time
189 size_t m_chunkRank = 0;
190
191 /// Maximum number of parallel-in-time iteration
192 size_t m_iterMaxPIT = 0;
193
194 // Number of windows for parallel-in-time time iteration
195 size_t m_numWindowsPIT = 1;
196
197 /// Using exact solution to compute error norms
199
200 /// ParallelInTime tolerance
202
203 // Temporary storage for parallel-in-time iteration
208 size_t m_nVar;
209 std::shared_ptr<SolverUtils::UnsteadySystem> m_fineEqSys;
210 std::shared_ptr<SolverUtils::UnsteadySystem> m_coarseEqSys;
218};
219
220/// Interpolate from an expansion to an expansion
223
224} // namespace SolverUtils
225} // namespace Nektar
226
227#endif // NEKTAR_SOLVERUTILS_DRIVERPARALLELINTIME_H
#define SOLVER_UTILS_EXPORT
Base class for the development of solvers.
Definition: Driver.h:67
Base class for the development of parallel-in-time solvers.
NekDouble m_totalTime
Total time integration interval.
size_t m_coarseSteps
Number of steps for the coarse solver.
void PrintSpeedUp(NekDouble fineSolveTime, NekDouble coarseSolveTime, NekDouble restTime, NekDouble interTime, NekDouble commTime, NekDouble predictorOverheadTime, NekDouble overheadTime)
virtual SOLVER_UTILS_EXPORT ~DriverParallelInTime()=default
Destructor.
void SendSolutionToNextProc(Array< OneD, Array< OneD, NekDouble > > &array, int &convergence)
size_t m_iterMaxPIT
Maximum number of parallel-in-time iteration.
size_t m_fineSteps
Number of steps for the fine solver.
void SetParallelInTimeEquationSystem(std::string AdvectiveType)
void CopyToFinePhysField(const Array< OneD, const Array< OneD, NekDouble > > &in)
NekDouble m_coarseTimeStep
Timestep for coarse solver.
void RecvInitialConditionFromPreviousProc(Array< OneD, Array< OneD, NekDouble > > &array, int &convergence)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fineFields
Array< OneD, MultiRegions::ExpListSharedPtr > m_coarseFields
void CopyFromFinePhysField(Array< OneD, Array< OneD, NekDouble > > &out)
Array< OneD, Array< OneD, NekDouble > > m_tmpfine
void Interpolator(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
std::shared_ptr< SolverUtils::UnsteadySystem > m_coarseEqSys
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
void PrintHeaderTitle1(const std::string &title)
void PrintHeaderTitle2(const std::string &title)
NekDouble m_tolerPIT
ParallelInTime tolerance.
void PrintCoarseSolverInfo(std::ostream &out=std::cout)
NekDouble m_fineTimeStep
Timestep for fine solver.
void PrintFineSolverInfo(std::ostream &out=std::cout)
bool m_exactSolution
Using exact solution to compute error norms.
size_t m_checkSteps
Number of steps for checkpoint.
void SolutionConvergenceMonitoring(const NekDouble &time)
virtual SOLVER_UTILS_EXPORT void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
Array< OneD, Array< OneD, NekDouble > > m_tmpcoarse
std::shared_ptr< SolverUtils::UnsteadySystem > m_fineEqSys
Array< OneD, Array< OneD, NekDouble > > m_exactsoln
NekDouble EstimateCommunicationTime(Array< OneD, Array< OneD, NekDouble > > &buffer1, Array< OneD, Array< OneD, NekDouble > > &buffer2)
void SolutionConvergenceSummary(const NekDouble &time)
void CopyFromCoarsePhysField(Array< OneD, Array< OneD, NekDouble > > &out)
void UpdateSolution(const Array< OneD, const Array< OneD, NekDouble > > &in)
void CopyToCoarsePhysField(const Array< OneD, const Array< OneD, NekDouble > > &in)
void CopySolutionVector(const Array< OneD, const Array< OneD, NekDouble > > &in, Array< OneD, Array< OneD, NekDouble > > &out)
size_t m_infoSteps
Number of steps for info I/O.
virtual NekDouble v_ComputeSpeedUp(const size_t iter, NekDouble fineSolveTime, NekDouble coarseSolveTime, NekDouble restTime, NekDouble interTime, NekDouble commTime, NekDouble predictorTime, NekDouble overheadTime)
SOLVER_UTILS_EXPORT DriverParallelInTime(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
void InterpExp1ToExp2(const Array< OneD, MultiRegions::ExpListSharedPtr > exp1, Array< OneD, MultiRegions::ExpListSharedPtr > &exp2)
Interpolate from an expansion to an expansion.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble