Nektar++
DriverParallelInTime.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File DriverParallelInTime.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 parallel-in-time solver
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <iomanip>
36
38#include <LocalRegions/HexExp.h>
40#include <LocalRegions/PyrExp.h>
42#include <LocalRegions/SegExp.h>
43#include <LocalRegions/TetExp.h>
44#include <LocalRegions/TriExp.h>
46#include <boost/format.hpp>
47
48namespace Nektar::SolverUtils
49{
50
51/**
52 *
53 */
57 : Driver(pSession, pGraph)
58{
59}
60
61/**
62 *
63 */
65{
66 try
67 {
68 // Retrieve the type of evolution operator to use
70 m_session->GetSolverInfoAsEnum<EvolutionOperatorType>(
71 "EvolutionOperator");
72
73 m_nTimeLevel = 2; // Only two time levels currently implemented.
74
76
77 // Set the AdvectiveType tag and create EquationSystem objects.
78 switch (m_EvolutionOperator)
79 {
80 case eNonlinear:
82 break;
83 case eDirect:
85 break;
86 case eAdjoint:
88 break;
89 case eSkewSymmetric:
90 SetParallelInTimeEquationSystem("SkewSymmetric");
91 break;
92 default:
93 ASSERTL0(false, "Unrecognised evolution operator.");
94 }
95
96 // Set pointers.
98 for (size_t timeLevel = 0; timeLevel < m_nTimeLevel; timeLevel++)
99 {
100 m_EqSys[timeLevel] =
101 std::dynamic_pointer_cast<SolverUtils::UnsteadySystem>(
102 m_equ[timeLevel]);
103 }
104
105 // Set time communication parameters.
106 m_numChunks = m_comm->GetTimeComm()->GetSize();
107 m_chunkRank = m_comm->GetTimeComm()->GetRank();
108 }
109 catch (int e)
110 {
111 ASSERTL0(e == -1, "No such class class defined.");
112 out << "An error occurred during driver initialisation." << std::endl;
113 }
114}
115
116/**
117 *
118 */
119void DriverParallelInTime::v_Execute([[maybe_unused]] std::ostream &out)
120{
121}
122
123/**
124 * Set the ParallelInTime (coarse solver) session file
125 */
127 std::string AdvectiveType)
128{
129 // Retrieve the equation system to solve.
130 ASSERTL0(m_session->DefinesSolverInfo("EqType"),
131 "EqType SolverInfo tag must be defined.");
132 std::string vEquation = m_session->DefinesSolverInfo("SolverType")
133 ? m_session->GetSolverInfo("SolverType")
134 : m_session->GetSolverInfo("EqType");
135
136 // Check such a module exists for this equation.
137 ASSERTL0(GetEquationSystemFactory().ModuleExists(vEquation),
138 "EquationSystem '" + vEquation +
139 "' is not defined.\n"
140 "Ensure equation name is correct and module is compiled.\n");
141
142 // Set fine parallel-in-time solver.
143 m_session->SetTag("AdvectiveType", AdvectiveType);
144 m_session->SetTag("ParallelInTimeSolver", "TimeLevel0");
146 m_graph);
147
148 // Define argument for the coarse parallel-in-time solver.
149 int npx = m_session->DefinesCmdLineArgument("npx")
150 ? m_session->GetCmdLineArgument<int>("npx")
151 : 1;
152 int npy = m_session->DefinesCmdLineArgument("npy")
153 ? m_session->GetCmdLineArgument<int>("npy")
154 : 1;
155 int npz = m_session->DefinesCmdLineArgument("npz")
156 ? m_session->GetCmdLineArgument<int>("npz")
157 : 1;
158 int nsz = m_session->DefinesCmdLineArgument("nsz")
159 ? m_session->GetCmdLineArgument<int>("nsz")
160 : 1;
161 int npt = m_session->DefinesCmdLineArgument("npt")
162 ? m_session->GetCmdLineArgument<int>("npt")
163 : 1;
164
165 // Convert into string.
166 std::string npx_string = std::to_string(npx);
167 std::string npy_string = std::to_string(npy);
168 std::string npz_string = std::to_string(npz);
169 std::string nsz_string = std::to_string(nsz);
170 std::string npt_string = std::to_string(npt);
171 std::string driver_string = "Driver=" + m_session->GetSolverInfo("Driver");
172
173 // use-opt-file
174 bool useOptFile = m_session->DefinesCmdLineArgument("use-opt-file");
175 std::string optfilename = useOptFile ? m_session->GetFilenames()[0] : "";
176
177 char *argv[] = {const_cast<char *>("Solver"), // this is just a place holder
178 const_cast<char *>("--solverinfo"),
179 const_cast<char *>(driver_string.c_str()),
180 const_cast<char *>("--npx"),
181 const_cast<char *>(npx_string.c_str()),
182 const_cast<char *>("--npy"),
183 const_cast<char *>(npy_string.c_str()),
184 const_cast<char *>("--npz"),
185 const_cast<char *>(npz_string.c_str()),
186 const_cast<char *>("--nsz"),
187 const_cast<char *>(nsz_string.c_str()),
188 const_cast<char *>("--npt"),
189 const_cast<char *>(npt_string.c_str()),
190 const_cast<char *>("-f"),
191 const_cast<char *>("--use-opt-file"),
192 const_cast<char *>(optfilename.c_str()),
193 nullptr};
194
195 size_t argc = useOptFile ? 16 : 14;
196
197 // Get list of session file names.
198 std::vector<std::string> sessionFileNames;
199 for (auto &filename : m_session->GetFilenames())
200 {
201 // Remove optfile name, if necessary.
202 if (filename.substr(filename.find_last_of(".") + 1) != "opt")
203 {
204 sessionFileNames.push_back(filename);
205 }
206 }
207
208 // Set session for coarse solver.
209 for (size_t timeLevel = 1; timeLevel < m_nTimeLevel; timeLevel++)
210 {
212 argc, argv, sessionFileNames, m_session->GetComm(), timeLevel);
213
214 // Set graph for coarse solver.
217
218 // Set BndRegionOrdering (necessary for DG with periodic BC) FIXME
219 graph->SetBndRegionOrdering(m_graph->GetBndRegionOrdering());
220
221 // Set CompositeOrdering (necessary for DG with periodic BC) FIXME
222 graph->SetCompositeOrdering(m_graph->GetCompositeOrdering());
223
224 // Retrieve the equation system to solve.
225 ASSERTL0(session->DefinesSolverInfo("EqType"),
226 "EqType SolverInfo tag must be defined.");
227 auto vEquation = session->DefinesSolverInfo("SolverType")
228 ? session->GetSolverInfo("SolverType")
229 : session->GetSolverInfo("EqType");
230
231 // Check such a module exists for this equation.
232 ASSERTL0(
233 GetEquationSystemFactory().ModuleExists(vEquation),
234 "EquationSystem '" + vEquation +
235 "' is not defined.\n"
236 "Ensure equation name is correct and module is compiled.\n");
237
238 // Set coarse parallel-in-time solver.
239 session->SetTag("AdvectiveType", AdvectiveType);
240 session->SetTag("ParallelInTimeSolver",
241 "TimeLevel" + std::to_string(timeLevel));
243 vEquation, session, graph);
244 }
245}
246
247/**
248 *
249 */
251{
252 // Parallel-in-Time iteration parameters.
253 m_tolerPIT = m_session->DefinesParameter("PITToler")
254 ? m_session->GetParameter("PITToler")
255 : 1.0E-16;
256 m_iterMaxPIT = m_session->DefinesParameter("PITIterMax")
257 ? m_session->GetParameter("PITIterMax")
258 : m_numChunks;
259 m_numWindowsPIT = m_session->DefinesParameter("NumWindows")
260 ? m_session->GetParameter("NumWindows")
261 : 1;
262
263 // Other parameters.
264 m_exactSolution = m_session->DefinesParameter("ExactSolution")
265 ? m_session->GetParameter("ExactSolution")
266 : 0;
267}
268
269/**
270 *
271 */
273{
274 m_nVar = m_EqSys[0]->GetNvariables();
275
276 // Initialize fine solver.
277 if (turnoff_output)
278 {
279 for (size_t timeLevel = 0; timeLevel < m_nTimeLevel; timeLevel++)
280 {
281 m_EqSys[timeLevel]->SetInfoSteps(0);
282 m_EqSys[timeLevel]->SetCheckpointSteps(0);
283 }
284 }
285 m_EqSys[0]->DoInitialise(true);
286
287 // Initialize coarse solver(s).
288 for (size_t timeLevel = 1; timeLevel < m_nTimeLevel; timeLevel++)
289 {
290 m_EqSys[timeLevel]->DoInitialise(false);
291 }
292
293 // Initialize time stepping parameters.
297 for (size_t timeLevel = 0; timeLevel < m_nTimeLevel; timeLevel++)
298 {
299 m_timestep[timeLevel] =
300 m_EqSys[timeLevel]->EquationSystem::GetTimeStep();
301 m_nsteps[timeLevel] = m_EqSys[timeLevel]->EquationSystem::GetSteps();
302 m_npts[timeLevel] = m_EqSys[timeLevel]->EquationSystem::GetNpoints();
303 }
304
305 // Initialize errors.
307 for (size_t i = 0; i < m_nVar; ++i)
308 {
310 }
313}
314
315/**
316 *
317 */
319{
320 if (m_chunkRank == 0 && m_comm->GetSpaceComm()->GetRank() == 0)
321 {
322 for (size_t timeLevel = 0; timeLevel < m_nTimeLevel; timeLevel++)
323 {
324 std::cout << std::string(71, '=') << std::endl << std::flush;
325 std::cout << "=========================== TIME LEVEL " +
326 std::to_string(timeLevel) +
327 " INFO "
328 "========================="
329 << std::endl
330 << std::flush;
331
332 m_EqSys[timeLevel]->PrintSummary(out);
333
334 std::cout << std::endl << std::flush;
335 }
336 }
337}
338
339/**
340 *
341 */
342void DriverParallelInTime::PrintHeader(const std::string &title, const char c)
343{
344 if (m_chunkRank == m_numChunks - 1 &&
345 m_comm->GetSpaceComm()->GetRank() == 0)
346 {
347 std::cout << std::endl;
348 std::cout << std::string(43, c) << std::endl << std::flush;
349 std::cout << title << std::endl << std::flush;
350 std::cout << std::string(43, c) << std::endl << std::flush;
351 }
352}
353
354/**
355 *
356 */
358 Array<OneD, Array<OneD, NekDouble>> &array, int &convergence)
359{
360 if (m_chunkRank > 0)
361 {
362 if (!convergence)
363 {
364 m_comm->GetTimeComm()->Recv(m_chunkRank - 1, convergence);
365 for (size_t i = 0; i < m_nVar; ++i)
366 {
367 m_comm->GetTimeComm()->Recv(m_chunkRank - 1, array[i]);
368 }
369 }
370 }
371}
372
373/**
374 *
375 */
377{
378 if (m_chunkRank > 0)
379 {
380 m_comm->GetTimeComm()->Recv(m_chunkRank - 1, array);
381 }
382}
383
384/**
385 *
386 */
388 Array<OneD, Array<OneD, NekDouble>> &array, int &convergence)
389{
390 if (m_chunkRank < m_numChunks - 1)
391 {
392 m_comm->GetTimeComm()->Send(m_chunkRank + 1, convergence);
393 for (size_t i = 0; i < m_nVar; ++i)
394 {
395 m_comm->GetTimeComm()->Send(m_chunkRank + 1, array[i]);
396 }
397 }
398}
399
400/**
401 *
402 */
404{
405 if (m_chunkRank < m_numChunks - 1)
406 {
407 m_comm->GetTimeComm()->Send(m_chunkRank + 1, array);
408 }
409}
410
411/**
412 *
413 */
415 const Array<OneD, const Array<OneD, NekDouble>> &in,
417{
418 for (size_t i = 0; i < m_nVar; ++i)
419 {
420 Vmath::Vcopy(in[i].size(), in[i], 1, out[i], 1);
421 }
422}
423
424/**
425 *
426 */
428 const size_t timeLevel, Array<OneD, Array<OneD, NekDouble>> &out)
429{
430 for (size_t i = 0; i < m_nVar; ++i)
431 {
432 m_EqSys[timeLevel]->CopyFromPhysField(i, out[i]);
433 }
434}
435
436/**
437 *
438 */
440 const size_t timeLevel, const Array<OneD, const Array<OneD, NekDouble>> &in)
441{
442 for (size_t i = 0; i < m_nVar; ++i)
443 {
444 m_EqSys[timeLevel]->CopyToPhysField(i, in[i]);
445 m_EqSys[timeLevel]->UpdateFields()[i]->SetPhysState(true);
446 }
447}
448
449/**
450 *
451 */
453 const size_t timeLevel, const Array<OneD, const Array<OneD, NekDouble>> &in)
454{
455 for (size_t i = 0; i < m_nVar; ++i)
456 {
458 {
459 m_EqSys[timeLevel]->CopyToPhysField(i, in[i]);
460 }
461 m_EqSys[timeLevel]->UpdateFields()[i]->FwdTrans(
462 m_EqSys[timeLevel]->UpdateFields()[i]->GetPhys(),
463 m_EqSys[timeLevel]->UpdateFields()[i]->UpdateCoeffs());
464 }
465}
466
467/**
468 *
469 */
471 const NekDouble &time)
472{
473 for (size_t i = 0; i < m_nVar; ++i)
474 {
475 m_EqSys[timeLevel]->EvaluateExactSolution(i, m_exactsoln[i], time);
476 }
477}
478
479/**
480 *
481 */
483 const size_t iter)
484{
485 PrintHeader((boost::format("ITERATION %1%") % iter).str(), '-');
486 UpdateErrorNorm(timeLevel, true);
487 PrintErrorNorm(timeLevel, true);
488}
489
490/**
491 *
492 */
494{
495 UpdateErrorNorm(timeLevel, false);
496 PrintErrorNorm(timeLevel, false);
497}
498
499/**
500 *
501 */
502void DriverParallelInTime::UpdateErrorNorm(const size_t timeLevel,
503 const bool normalized)
504{
505 for (size_t i = 0; i < m_nVar; ++i)
506 {
507 m_vL2Errors[i] =
508 m_EqSys[timeLevel]->L2Error(i, m_exactsoln[i], normalized);
509 m_vLinfErrors[i] = m_EqSys[timeLevel]->LinfError(i, m_exactsoln[i]);
510 }
511}
512
513/**
514 *
515 */
516void DriverParallelInTime::PrintErrorNorm(const size_t timeLevel,
517 const bool normalized)
518{
519 if (m_chunkRank == m_numChunks - 1 &&
520 m_comm->GetSpaceComm()->GetRank() == 0)
521 {
522 for (size_t i = 0; i < m_nVar; ++i)
523 {
524 if (normalized)
525 {
526 std::cout << "L2 error (variable "
527 << m_EqSys[timeLevel]->GetVariable(i)
528 << ") : " << m_vL2Errors[i] << std::endl
529 << std::flush;
530 std::cout << "Linf error (variable "
531 << m_EqSys[timeLevel]->GetVariable(i)
532 << ") : " << m_vLinfErrors[i] << std::endl
533 << std::flush;
534 }
535 else
536 {
537 std::cout << "L 2 error (variable "
538 << m_EqSys[timeLevel]->GetVariable(i)
539 << ") : " << m_vL2Errors[i] << std::endl
540 << std::flush;
541 std::cout << "L inf error (variable "
542 << m_EqSys[timeLevel]->GetVariable(i)
543 << ") : " << m_vLinfErrors[i] << std::endl
544 << std::flush;
545 }
546 }
547 }
548}
549
550/**
551 *
552 */
554{
555 NekDouble L2Error = 0.0;
556 for (size_t i = 0; i < m_nVar; ++i)
557 {
558 L2Error = std::max(L2Error, m_vL2Errors[i]);
559 }
560 return L2Error;
561}
562
563/**
564 *
565 */
569{
570 if (m_numChunks == 1)
571 {
572 return 0.0;
573 }
574 else
575 {
576 // Average communication time over niter iteration.
577 size_t niter = 20;
579 for (size_t n = 0; n <= niter; n++)
580 {
581 if (n == 1)
582 {
583 timer.Start(); // Ignore the first iteration
584 }
585
586 if (m_chunkRank == 0)
587 {
588 for (size_t i = 0; i < buffer1.size(); ++i)
589 {
590 m_comm->GetTimeComm()->Send(m_numChunks - 1, buffer1[i]);
591 }
592 }
593
594 if (m_chunkRank == m_numChunks - 1)
595 {
596 for (size_t i = 0; i < buffer2.size(); ++i)
597 {
598 m_comm->GetTimeComm()->Recv(0, buffer2[i]);
599 }
600 }
601 }
602 timer.Stop();
603 return timer.Elapsed().count() / niter;
604 }
605}
606
607/**
608 *
609 */
613 const Array<OneD, Array<OneD, NekDouble>> &inarray,
615{
616 if (infield.size() != outfield.size())
617 {
618 NEKERROR(ErrorUtil::efatal, "not the same number of variables")
619 }
620
621 for (size_t n = 0; n < infield.size(); ++n)
622 {
623 // Interpolation from infield -> outfield assuming that infield and
624 // outfield are the same explists, but at potentially different
625 // polynomial orders.
626 if (infield[n]->GetExpSize() != outfield[n]->GetExpSize())
627 {
628 NEKERROR(ErrorUtil::efatal, "not the same mesh")
629 }
630
631 // Assign input/output array.
632 auto inphys = (inarray == NullNekDoubleArrayOfArray)
633 ? infield[n]->UpdatePhys()
634 : inarray[n];
635 auto outphys = (outarray == NullNekDoubleArrayOfArray)
636 ? outfield[n]->UpdatePhys()
637 : outarray[n];
638
639 // If same polynomial orders, simply copy solution.
640 if (infield[n]->GetTotPoints() == outfield[n]->GetTotPoints())
641 {
642 Vmath::Vcopy(infield[n]->GetTotPoints(), inphys, 1, outphys, 1);
643 }
644 // If different polynomial orders, interpolate solution.
645 else
646 {
647 for (size_t i = 0; i < infield[n]->GetExpSize(); ++i)
648 {
649 // Get the elements.
650 auto inElmt = infield[n]->GetExp(i);
651 auto outElmt = outfield[n]->GetExp(i);
652
653 // Get the offset of elements in the storage arrays.
654 size_t inoffset = infield[n]->GetPhys_Offset(i);
655 size_t outoffset = outfield[n]->GetPhys_Offset(i);
656
657 // Transform solution from physical to coefficient space.
658 Array<OneD, NekDouble> incoeff(inElmt->GetNcoeffs());
659 inElmt->FwdTrans(inphys + inoffset, incoeff);
660
661 // Interpolate elements.
663 if (inElmt->DetShapeType() == LibUtilities::Seg)
664 {
665 expPtr = std::make_shared<StdRegions::StdSegExp>(
667 inElmt->GetBasis(0)->GetBasisType(),
668 inElmt->GetBasis(0)->GetNumModes(),
669 outElmt->GetBasis(0)->GetPointsKey()));
670 }
671 else if (inElmt->DetShapeType() == LibUtilities::Quad)
672 {
673 expPtr = std::make_shared<StdRegions::StdQuadExp>(
675 inElmt->GetBasis(0)->GetBasisType(),
676 inElmt->GetBasis(0)->GetNumModes(),
677 outElmt->GetBasis(0)->GetPointsKey()),
679 inElmt->GetBasis(1)->GetBasisType(),
680 inElmt->GetBasis(1)->GetNumModes(),
681 outElmt->GetBasis(1)->GetPointsKey()));
682 }
683 else if (inElmt->DetShapeType() == LibUtilities::Tri)
684 {
685 expPtr = std::make_shared<StdRegions::StdTriExp>(
687 inElmt->GetBasis(0)->GetBasisType(),
688 inElmt->GetBasis(0)->GetNumModes(),
689 outElmt->GetBasis(0)->GetPointsKey()),
691 inElmt->GetBasis(1)->GetBasisType(),
692 inElmt->GetBasis(1)->GetNumModes(),
693 outElmt->GetBasis(1)->GetPointsKey()));
694 }
695 else if (inElmt->DetShapeType() == LibUtilities::Hex)
696 {
697 expPtr = std::make_shared<StdRegions::StdHexExp>(
699 inElmt->GetBasis(0)->GetBasisType(),
700 inElmt->GetBasis(0)->GetNumModes(),
701 outElmt->GetBasis(0)->GetPointsKey()),
703 inElmt->GetBasis(1)->GetBasisType(),
704 inElmt->GetBasis(1)->GetNumModes(),
705 outElmt->GetBasis(1)->GetPointsKey()),
707 inElmt->GetBasis(2)->GetBasisType(),
708 inElmt->GetBasis(2)->GetNumModes(),
709 outElmt->GetBasis(2)->GetPointsKey()));
710 }
711 else if (inElmt->DetShapeType() == LibUtilities::Prism)
712 {
713 expPtr = std::make_shared<StdRegions::StdPrismExp>(
715 inElmt->GetBasis(0)->GetBasisType(),
716 inElmt->GetBasis(0)->GetNumModes(),
717 outElmt->GetBasis(0)->GetPointsKey()),
719 inElmt->GetBasis(1)->GetBasisType(),
720 inElmt->GetBasis(1)->GetNumModes(),
721 outElmt->GetBasis(1)->GetPointsKey()),
723 inElmt->GetBasis(2)->GetBasisType(),
724 inElmt->GetBasis(2)->GetNumModes(),
725 outElmt->GetBasis(2)->GetPointsKey()));
726 }
727 else if (inElmt->DetShapeType() == LibUtilities::Pyr)
728 {
729 expPtr = std::make_shared<StdRegions::StdPyrExp>(
731 inElmt->GetBasis(0)->GetBasisType(),
732 inElmt->GetBasis(0)->GetNumModes(),
733 outElmt->GetBasis(0)->GetPointsKey()),
735 inElmt->GetBasis(1)->GetBasisType(),
736 inElmt->GetBasis(1)->GetNumModes(),
737 outElmt->GetBasis(1)->GetPointsKey()),
739 inElmt->GetBasis(2)->GetBasisType(),
740 inElmt->GetBasis(2)->GetNumModes(),
741 outElmt->GetBasis(2)->GetPointsKey()));
742 }
743 else if (inElmt->DetShapeType() == LibUtilities::Tet)
744 {
745 expPtr = std::make_shared<StdRegions::StdTetExp>(
747 inElmt->GetBasis(0)->GetBasisType(),
748 inElmt->GetBasis(0)->GetNumModes(),
749 outElmt->GetBasis(0)->GetPointsKey()),
751 inElmt->GetBasis(1)->GetBasisType(),
752 inElmt->GetBasis(1)->GetNumModes(),
753 outElmt->GetBasis(1)->GetPointsKey()),
755 inElmt->GetBasis(2)->GetBasisType(),
756 inElmt->GetBasis(2)->GetNumModes(),
757 outElmt->GetBasis(2)->GetPointsKey()));
758 }
759
760 // Transform solution from coefficient to physical space.
761 Array<OneD, NekDouble> tmp = outphys + outoffset;
762 expPtr->BwdTrans(incoeff, tmp);
763 }
764 }
765 }
766}
767
768} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
Describes the specification for a Basis.
Definition: Basis.h:45
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
Base class for the development of solvers.
Definition: Driver.h:65
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:83
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:80
SpatialDomains::MeshGraphSharedPtr m_graph
MeshGraph object.
Definition: Driver.h:89
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:98
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:92
Array< OneD, NekDouble > m_vL2Errors
Storage for parallel-in-time iteration.
void UpdateFieldCoeffs(const size_t timeLevel, const Array< OneD, const Array< OneD, NekDouble > > &in=NullNekDoubleArrayOfArray)
void SendToNextProc(Array< OneD, Array< OneD, NekDouble > > &array, int &convergence)
size_t m_iterMaxPIT
Maximum number of parallel-in-time iteration.
void PrintErrorNorm(const size_t timeLevel, const bool normalized)
void SetParallelInTimeEquationSystem(std::string AdvectiveType)
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)
SOLVER_UTILS_EXPORT void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
Array< OneD, size_t > m_npts
Number of dof for each time level.
Array< OneD, Array< OneD, NekDouble > > m_exactsoln
NekDouble EstimateCommunicationTime(Array< OneD, Array< OneD, NekDouble > > &buffer1, Array< OneD, Array< OneD, NekDouble > > &buffer2)
void SolutionConvergenceSummary(const size_t timeLevel)
Array< OneD, NekDouble > m_timestep
Time step for each time level.
void CopySolutionVector(const Array< OneD, const Array< OneD, NekDouble > > &in, Array< OneD, Array< OneD, NekDouble > > &out)
void UpdateErrorNorm(const size_t timeLevel, const bool normalized)
void SolutionConvergenceMonitoring(const size_t timeLevel, const size_t iter)
SOLVER_UTILS_EXPORT DriverParallelInTime(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.
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)
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
Definition: MeshGraph.cpp:115
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static DomainRangeShPtr NullDomainRangeShPtr
Definition: DomainRange.h:65
EquationSystemFactory & GetEquationSystemFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825