Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Static Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::SolverUtils::DriverAdaptive Class Reference

Base class for the adaptive polynomial order driver. More...

#include <DriverAdaptive.h>

Inheritance diagram for Nektar::SolverUtils::DriverAdaptive:
[legend]

Static Public Member Functions

static DriverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

SOLVER_UTILS_EXPORT DriverAdaptive (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Constructor. More...
 
virtual SOLVER_UTILS_EXPORT ~DriverAdaptive ()
 Destructor. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (std::ostream &out=std::cout) override
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_Execute (std::ostream &out=std::cout) override
 Virtual function for solve implementation. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Driver
 Driver (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (std::ostream &out=std::cout)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_Execute (std::ostream &out=std::cout)=0
 Virtual function for solve implementation. More...
 

Static Protected Attributes

static std::string driverLookupId
 
- Static Protected Attributes inherited from Nektar::SolverUtils::Driver
static std::string evolutionOperatorLookupIds []
 
static std::string evolutionOperatorDef
 
static std::string driverDefault
 

Private Member Functions

SOLVER_UTILS_EXPORT void ReplaceExpansion (Array< OneD, MultiRegions::ExpListSharedPtr > &fields, std::map< int, int > deltaP)
 Update EXPANSIONS tag inside XML schema to reflect new polynomial order distribution. More...
 

Friends

class MemoryManager< DriverAdaptive >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Driver
virtual ~Driver ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (std::ostream &out=std::cout)
 Initialise Object. More...
 
SOLVER_UTILS_EXPORT void Execute (std::ostream &out=std::cout)
 Execute driver. More...
 
SOLVER_UTILS_EXPORT Array< OneD, EquationSystemSharedPtrGetEqu ()
 
- Protected Attributes inherited from Nektar::SolverUtils::Driver
LibUtilities::CommSharedPtr m_comm
 Communication object. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 Session reader object. More...
 
LibUtilities::SessionReaderSharedPtr session_LinNS
 Coupling between SFD and arnoldi. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 MeshGraph object. More...
 
Array< OneD, EquationSystemSharedPtrm_equ
 Equation system to solve. More...
 
int m_nequ
 number of equations More...
 
enum EvolutionOperatorType m_EvolutionOperator
 Evolution Operator. More...
 

Detailed Description

Base class for the adaptive polynomial order driver.

Definition at line 47 of file DriverAdaptive.h.

Constructor & Destructor Documentation

◆ DriverAdaptive()

Nektar::SolverUtils::DriverAdaptive::DriverAdaptive ( const LibUtilities::SessionReaderSharedPtr  pSession,
const SpatialDomains::MeshGraphSharedPtr  pGraph 
)
protected

Constructor.

Definition at line 63 of file DriverAdaptive.cpp.

66 : Driver(pSession, pGraph)
67{
68}
Driver(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Initialises EquationSystem class members.
Definition: Driver.cpp:76

◆ ~DriverAdaptive()

Nektar::SolverUtils::DriverAdaptive::~DriverAdaptive ( )
protectedvirtual

Destructor.

Definition at line 73 of file DriverAdaptive.cpp.

74{
75}

Member Function Documentation

◆ create()

static DriverSharedPtr Nektar::SolverUtils::DriverAdaptive::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file DriverAdaptive.h.

56 {
59 p->InitObject();
60 return p;
61 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:54

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ ReplaceExpansion()

void Nektar::SolverUtils::DriverAdaptive::ReplaceExpansion ( Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
std::map< int, int >  deltaP 
)
private

Update EXPANSIONS tag inside XML schema to reflect new polynomial order distribution.

Parameters
fieldsInput fields.
deltaPMap of polynomial order expansions

Definition at line 449 of file DriverAdaptive.cpp.

451{
452 int nExp, nDim;
453 int expdim = m_equ[0]->UpdateFields()[0]->GetGraph()->GetMeshDimension();
454
455 // Get field definitions
456 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fielddefs =
457 fields[0]->GetFieldDefinitions();
458
459 if (fielddefs[0]->m_numHomogeneousDir == 1)
460 {
461 nDim = expdim + 1;
462 }
463 else
464 {
465 nDim = expdim;
466 }
467
468 // Add variables to field definition
469 for (int i = 0; i < fielddefs.size(); ++i)
470 {
471 for (int j = 0; j < fields.size(); ++j)
472 {
473 fielddefs[i]->m_fields.push_back(m_session->GetVariable(j));
474 }
475 }
476
477 // Get tinyxml objects
478 TiXmlElement *exp_tag = m_session->GetElement("NEKTAR/EXPANSIONS");
479
480 // Clear current expansions
481 exp_tag->Clear();
482
483 // Write new expansion information
484 for (int f = 0; f < fielddefs.size(); ++f)
485 {
486 nExp = fielddefs[f]->m_elementIDs.size();
487
488 //---------------------------------------------
489 // Write ELEMENTS
490 TiXmlElement *elemTag = new TiXmlElement("ELEMENTS");
491 exp_tag->LinkEndChild(elemTag);
492
493 // Write FIELDS
494 std::string fieldsString;
495 {
496 std::stringstream fieldsStringStream;
497 bool first = true;
498 for (std::vector<int>::size_type i = 0;
499 i < fielddefs[f]->m_fields.size(); i++)
500 {
501 if (!first)
502 {
503 fieldsStringStream << ",";
504 }
505 fieldsStringStream << fielddefs[f]->m_fields[i];
506 first = false;
507 }
508 fieldsString = fieldsStringStream.str();
509 }
510 elemTag->SetAttribute("FIELDS", fieldsString);
511
512 // Write SHAPE
513 std::string shapeString;
514 {
515 std::stringstream shapeStringStream;
516 shapeStringStream
517 << LibUtilities::ShapeTypeMap[fielddefs[f]->m_shapeType];
518 if (fielddefs[f]->m_numHomogeneousDir == 1)
519 {
520 shapeStringStream << "-HomogenousExp1D";
521 }
522 else if (fielddefs[f]->m_numHomogeneousDir == 2)
523 {
524 shapeStringStream << "-HomogenousExp2D";
525 }
526
527 shapeString = shapeStringStream.str();
528 }
529 elemTag->SetAttribute("SHAPE", shapeString);
530
531 // Write BASIS
532 std::string basisString;
533 {
534 std::stringstream basisStringStream;
535 bool first = true;
536 for (std::vector<LibUtilities::BasisType>::size_type i = 0;
537 i < fielddefs[f]->m_basis.size(); i++)
538 {
539 if (!first)
540 {
541 basisStringStream << ",";
542 }
543 basisStringStream
544 << LibUtilities::BasisTypeMap[fielddefs[f]->m_basis[i]];
545 first = false;
546 }
547 basisString = basisStringStream.str();
548 }
549 elemTag->SetAttribute("BASIS", basisString);
550
551 // Write homogeneuous length details
552 if (fielddefs[f]->m_numHomogeneousDir)
553 {
554 std::string homoLenString;
555 {
556 std::stringstream homoLenStringStream;
557 bool first = true;
558 for (int i = 0; i < fielddefs[f]->m_numHomogeneousDir; ++i)
559 {
560 if (!first)
561 {
562 homoLenStringStream << ",";
563 }
564 homoLenStringStream
565 << fielddefs[f]->m_homogeneousLengths[i];
566 first = false;
567 }
568 homoLenString = homoLenStringStream.str();
569 }
570 elemTag->SetAttribute("HOMOGENEOUSLENGTHS", homoLenString);
571 }
572
573 // Write homogeneuous planes/lines details
574 if (fielddefs[f]->m_numHomogeneousDir)
575 {
576 if (fielddefs[f]->m_homogeneousYIDs.size() > 0)
577 {
578 std::string homoYIDsString;
579 {
580 std::stringstream homoYIDsStringStream;
581 bool first = true;
582 for (int i = 0; i < fielddefs[f]->m_homogeneousYIDs.size();
583 i++)
584 {
585 if (!first)
586 {
587 homoYIDsStringStream << ",";
588 }
589 homoYIDsStringStream
590 << fielddefs[f]->m_homogeneousYIDs[i];
591 first = false;
592 }
593 homoYIDsString = homoYIDsStringStream.str();
594 }
595 elemTag->SetAttribute("HOMOGENEOUSYIDS", homoYIDsString);
596 }
597
598 if (fielddefs[f]->m_homogeneousZIDs.size() > 0)
599 {
600 std::string homoZIDsString;
601 {
602 std::stringstream homoZIDsStringStream;
603 bool first = true;
604 for (int i = 0; i < fielddefs[f]->m_homogeneousZIDs.size();
605 i++)
606 {
607 if (!first)
608 {
609 homoZIDsStringStream << ",";
610 }
611 homoZIDsStringStream
612 << fielddefs[f]->m_homogeneousZIDs[i];
613 first = false;
614 }
615 homoZIDsString = homoZIDsStringStream.str();
616 }
617 elemTag->SetAttribute("HOMOGENEOUSZIDS", homoZIDsString);
618 }
619 }
620
621 // Write NUMMODESPERDIR
622 std::string numModesString;
623 {
624 std::stringstream numModesStringStream;
625
626 numModesStringStream << "MIXORDER:";
627 bool first = true;
628 int eId;
629 Array<OneD, int> order(nDim, 0);
630 for (int n = 0; n < nExp; n++)
631 {
632 eId = fielddefs[f]->m_elementIDs[n];
633
634 for (int i = 0; i < expdim; i++)
635 {
636 order[i] = deltaP[eId];
637 }
638
639 for (int i = 0; i < nDim; i++)
640 {
641 if (fielddefs[f]->m_uniOrder)
642 {
643 order[i] += fielddefs[f]->m_numModes[i];
644 }
645 else
646 {
647 order[i] += fielddefs[f]->m_numModes[n * nDim + i];
648 }
649
650 if (!first)
651 {
652 numModesStringStream << ",";
653 }
654
655 numModesStringStream << order[i];
656 first = false;
657 }
658 }
659 numModesString = numModesStringStream.str();
660 }
661 elemTag->SetAttribute("NUMMODESPERDIR", numModesString);
662
663 // Write IDs. Should ideally look at ways of compressing this stream if
664 // just sequential.
665 elemTag->SetAttribute(
666 "ID", ParseUtils::GenerateSeqString(fielddefs[f]->m_elementIDs));
667 }
668}
static std::string GenerateSeqString(const std::vector< T > &v)
Generate a compressed comma-separated string representation of a vector of unsigned integers.
Definition: ParseUtils.h:72
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:85
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:79
const char *const BasisTypeMap[]
Definition: Foundations.hpp:46

References Nektar::LibUtilities::BasisTypeMap, Nektar::ParseUtils::GenerateSeqString(), Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::Driver::m_session, and Nektar::LibUtilities::ShapeTypeMap.

Referenced by v_Execute().

◆ v_Execute()

void Nektar::SolverUtils::DriverAdaptive::v_Execute ( std::ostream &  out = std::cout)
overrideprotectedvirtual

Virtual function for solve implementation.

Implements Nektar::SolverUtils::Driver.

Definition at line 88 of file DriverAdaptive.cpp.

89{
91 NekDouble CPUtime;
92
93 m_equ[0]->PrintSummary(out);
94
95 // First run using original order
96 timer.Start();
97 m_equ[0]->DoInitialise();
98
99 // Obtain initial time in case a restart was used
100 NekDouble startTime = m_equ[0]->GetFinalTime();
101 m_equ[0]->DoSolve();
102
103 // Load session parameters and solver info
104 bool isHomogeneous1D;
105 int nRuns, minP, maxP, sensorVar;
106 NekDouble lowerTol, upperTol;
107 m_session->LoadParameter("NumRuns", nRuns, 1);
108 m_session->LoadParameter("AdaptiveMinModes", minP, 4);
109 m_session->LoadParameter("AdaptiveMaxModes", maxP, 12);
110 m_session->LoadParameter("AdaptiveLowerTolerance", lowerTol, 1e-8);
111 m_session->LoadParameter("AdaptiveUpperTolerance", upperTol, 1e-6);
112 m_session->LoadParameter("AdaptiveSensorVariable", sensorVar, 0);
113 m_session->MatchSolverInfo("Homogeneous", "1D", isHomogeneous1D, false);
114
115 // Get number of elements and planes
116 int nExp, nPlanes;
117 if (isHomogeneous1D)
118 {
119 nExp = m_equ[0]->UpdateFields()[0]->GetPlane(0)->GetExpSize();
120 nPlanes = m_equ[0]->UpdateFields()[0]->GetZIDs().size();
121 }
122 else
123 {
124 nExp = m_equ[0]->UpdateFields()[0]->GetExpSize();
125 nPlanes = 1;
126 }
127 int expdim = m_equ[0]->UpdateFields()[0]->GetGraph()->GetMeshDimension();
128
129 int nFields = m_equ[0]->UpdateFields().size();
130 int numSteps = m_session->GetParameter("NumSteps");
131 NekDouble period = m_session->GetParameter("TimeStep") * numSteps;
132
133 Array<OneD, NekDouble> coeffs, phys, physReduced, tmpArray;
134
135 // Get mapping
137 mapping = GlobalMapping::Mapping::Load(m_session, m_equ[0]->UpdateFields());
138
139 // Adaptive loop
140 Array<OneD, int> P(expdim);
141 Array<OneD, int> numPoints(expdim);
142 Array<OneD, LibUtilities::PointsKey> ptsKey(expdim);
144 for (int i = 1; i < nRuns; i++)
145 {
146 // Get field expansions
147 Array<OneD, MultiRegions::ExpListSharedPtr> fields =
148 m_equ[0]->UpdateFields();
149
150 // Determine the change to be applied in the order
151 map<int, int> deltaP;
152 int offset = 0;
153 for (int n = 0; n < nExp; n++)
154 {
155 offset = fields[sensorVar]->GetPhys_Offset(n);
156 Exp = fields[sensorVar]->GetExp(n);
157
158 for (int k = 0; k < expdim; ++k)
159 {
160 P[k] = Exp->GetBasis(k)->GetNumModes();
161 numPoints[k] = Exp->GetBasis(k)->GetNumPoints();
162 ptsKey[k] = LibUtilities::PointsKey(
163 numPoints[k], Exp->GetBasis(k)->GetPointsType());
164 }
165
166 // Declare orthogonal basis.
168 switch (Exp->GetGeom()->GetShapeType())
169 {
171 {
172 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P[0] - 1,
173 ptsKey[0]);
174 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, P[1] - 1,
175 ptsKey[1]);
176 OrthoExp = MemoryManager<
177 StdRegions::StdQuadExp>::AllocateSharedPtr(Ba, Bb);
178 break;
179 }
181 {
182 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P[0] - 1,
183 ptsKey[0]);
184 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B, P[1] - 1,
185 ptsKey[1]);
186 OrthoExp =
188 Ba, Bb);
189 break;
190 }
192 {
193 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P[0] - 1,
194 ptsKey[0]);
195 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B, P[1] - 1,
196 ptsKey[1]);
197 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_C, P[2] - 1,
198 ptsKey[2]);
199 OrthoExp =
201 Ba, Bb, Bc);
202 break;
203 }
205 {
206 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P[0] - 1,
207 ptsKey[0]);
208 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, P[1] - 1,
209 ptsKey[1]);
210 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_C, P[2] - 1,
211 ptsKey[2]);
212 OrthoExp =
214 Ba, Bb, Bc);
215 break;
216 }
218 {
219 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P[0] - 1,
220 ptsKey[0]);
221 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, P[1] - 1,
222 ptsKey[1]);
223 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_B, P[2] - 1,
224 ptsKey[2]);
225 OrthoExp = MemoryManager<
226 StdRegions::StdPrismExp>::AllocateSharedPtr(Ba, Bb, Bc);
227 break;
228 }
230 {
231 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P[0] - 1,
232 ptsKey[0]);
233 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, P[1] - 1,
234 ptsKey[1]);
235 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_A, P[2] - 1,
236 ptsKey[2]);
237 OrthoExp =
239 Ba, Bb, Bc);
240 break;
241 }
242 default:
243 ASSERTL0(false, "Shape not supported.");
244 break;
245 }
246
247 int nq = OrthoExp->GetTotPoints();
248
249 NekDouble error = 0;
250 NekDouble fac = 0;
251 NekDouble tmp = 0;
252
253 coeffs = Array<OneD, NekDouble>(OrthoExp->GetNcoeffs());
254 physReduced = Array<OneD, NekDouble>(OrthoExp->GetTotPoints());
255 tmpArray = Array<OneD, NekDouble>(OrthoExp->GetTotPoints(), 0.0);
256
257 // Refinement based only on one variable
258 for (int plane = 0; plane < nPlanes; plane++)
259 {
260 if (isHomogeneous1D)
261 {
262 phys =
263 fields[sensorVar]->GetPlane(plane)->GetPhys() + offset;
264 }
265 else
266 {
267 phys = fields[sensorVar]->GetPhys() + offset;
268 }
269
270 // Project solution to lower order
271 OrthoExp->FwdTrans(phys, coeffs);
272 OrthoExp->BwdTrans(coeffs, physReduced);
273
274 // Calculate error =||phys-physReduced||^2 / ||phys||^2
275 Vmath::Vsub(nq, phys, 1, physReduced, 1, tmpArray, 1);
276 Vmath::Vmul(nq, tmpArray, 1, tmpArray, 1, tmpArray, 1);
277 tmp = Exp->Integral(tmpArray);
278
279 Vmath::Vmul(nq, phys, 1, phys, 1, tmpArray, 1);
280 fac = Exp->Integral(tmpArray);
281
282 tmp = abs(tmp / fac);
283
284 if (tmp != tmp)
285 {
286 ASSERTL0(false,
287 "Adaptive procedure encountered NaN value.");
288 }
289
290 // Get maximum value along planes
291 error = (tmp > error) ? tmp : error;
292 }
293
294 // Combine planes from different processes
295 m_session->GetComm()->GetColumnComm()->AllReduce(
297
298 // Override tolerances if function is defined
299 if (m_session->DefinesFunction("AdaptiveLowerTolerance"))
300 {
301 int nq = Exp->GetTotPoints();
302
303 // Obtain points from the the element
304 Array<OneD, NekDouble> xc0, xc1, xc2;
305 xc0 = Array<OneD, NekDouble>(nq, 0.0);
306 xc1 = Array<OneD, NekDouble>(nq, 0.0);
307 xc2 = Array<OneD, NekDouble>(nq, 0.0);
308 Exp->GetCoords(xc0, xc1, xc2);
309
310 // Evaluate function from session file
311 Array<OneD, NekDouble> tolerance(nq, 0.0);
313 m_session->GetFunction("AdaptiveLowerTolerance", 0);
314 ffunc->Evaluate(xc0, xc1, xc2, tolerance);
315 lowerTol = Vmath::Vsum(nq, tolerance, 1) / nq;
316 }
317
318 if (m_session->DefinesFunction("AdaptiveUpperTolerance"))
319 {
320 int nq = Exp->GetTotPoints();
321
322 // Obtain points from the the element
323 Array<OneD, NekDouble> xc0, xc1, xc2;
324 xc0 = Array<OneD, NekDouble>(nq, 0.0);
325 xc1 = Array<OneD, NekDouble>(nq, 0.0);
326 xc2 = Array<OneD, NekDouble>(nq, 0.0);
327 Exp->GetCoords(xc0, xc1, xc2);
328
329 // Evaluate function from session file
330 Array<OneD, NekDouble> tolerance(nq, 0.0);
332 m_session->GetFunction("AdaptiveUpperTolerance", 0);
333 ffunc->Evaluate(xc0, xc1, xc2, tolerance);
334 upperTol = Vmath::Vsum(nq, tolerance, 1) / nq;
335 }
336
337 // Determine what to do with the polynomial order
338 if ((error > upperTol) && (P[0] < maxP))
339 {
340 deltaP[Exp->GetGeom()->GetGlobalID()] = 1;
341 }
342 else if ((error < lowerTol) && P[0] > minP)
343 {
344 deltaP[Exp->GetGeom()->GetGlobalID()] = -1;
345 }
346 else
347 {
348 deltaP[Exp->GetGeom()->GetGlobalID()] = 0;
349 }
350 }
351
352 // Write new expansion section to the session reader and re-read graph.
353 ReplaceExpansion(fields, deltaP);
354 m_graph->ReadExpansionInfo();
355
356 // Reset GlobalLinSys Manager to avoid using too much memory
357 //
358 // @todo This could be made better by replacing individual matrices
359 // within the linear system.
360 if (LibUtilities::NekManager<MultiRegions::GlobalLinSysKey,
361 MultiRegions::GlobalLinSys>::
362 PoolCreated(std::string("GlobalLinSys")))
363 {
364 LibUtilities::NekManager<MultiRegions::GlobalLinSysKey,
365 MultiRegions::GlobalLinSys>::
366 ClearManager(std::string("GlobalLinSys"));
367 }
368
369 int chkNumber = m_equ[0]->GetCheckpointNumber();
370 int chkSteps = m_equ[0]->GetCheckpointSteps();
371
372 // Initialise driver again
374
375 // Update mapping (must be before m_equ[0]->DoInitialise();)
376 mapping->ReplaceField(m_equ[0]->UpdateFields());
377
378 // Set chkSteps to zero to avoid writing initial condition
379 m_equ[0]->SetCheckpointSteps(0);
380
381 // Initialise equation
382 m_equ[0]->DoInitialise();
383 m_equ[0]->SetInitialStep(i * numSteps);
384 m_equ[0]->SetSteps(i * numSteps + numSteps);
385 m_equ[0]->SetTime(startTime + i * period);
386 m_equ[0]->SetBoundaryConditions(startTime + i * period);
387 m_equ[0]->SetCheckpointNumber(chkNumber);
388 m_equ[0]->SetCheckpointSteps(chkSteps);
389
390 // Project solution to new expansion
391 for (int n = 0; n < nFields; n++)
392 {
393 m_equ[0]->UpdateFields()[n]->ExtractCoeffsToCoeffs(
394 fields[n], fields[n]->GetCoeffs(),
395 m_equ[0]->UpdateFields()[n]->UpdateCoeffs());
396 m_equ[0]->UpdateFields()[n]->BwdTrans(
397 m_equ[0]->UpdateFields()[n]->GetCoeffs(),
398 m_equ[0]->UpdateFields()[n]->UpdatePhys());
399 }
400
401 // Solve equation
402 m_equ[0]->DoSolve();
403 }
404
405 timer.Stop();
406
407 m_equ[0]->Output();
408
409 if (m_comm->GetRank() == 0)
410 {
411 CPUtime = timer.Elapsed().count();
412 cout << "-------------------------------------------" << endl;
413 cout << "Total Computation Time = " << CPUtime << "s" << endl;
414 cout << "-------------------------------------------" << endl;
415 }
416
417 // Evaluate and output computation time and solution accuracy.
418 // The specific format of the error output is essential for the
419 // regression tests to work.
420
421 // Evaluate L2 Error
422 for (int i = 0; i < m_equ[0]->GetNvariables(); ++i)
423 {
424 Array<OneD, NekDouble> exactsoln(m_equ[0]->GetTotPoints(), 0.0);
425
426 // Evaluate "ExactSolution" function, or zero array
427 m_equ[0]->EvaluateExactSolution(i, exactsoln, m_equ[0]->GetFinalTime());
428
429 NekDouble vL2Error = m_equ[0]->L2Error(i, exactsoln);
430 NekDouble vLinfError = m_equ[0]->LinfError(i, exactsoln);
431
432 if (m_comm->GetRank() == 0)
433 {
434 out << "L 2 error (variable " << m_equ[0]->GetVariable(i)
435 << ") : " << vL2Error << endl;
436 out << "L inf error (variable " << m_equ[0]->GetVariable(i)
437 << ") : " << vLinfError << endl;
438 }
439 }
440}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
static GLOBAL_MAPPING_EXPORT MappingSharedPtr Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Return a pointer to the mapping, creating it on first call.
Definition: Mapping.cpp:272
SOLVER_UTILS_EXPORT void ReplaceExpansion(Array< OneD, MultiRegions::ExpListSharedPtr > &fields, std::map< int, int > deltaP)
Update EXPANSIONS tag inside XML schema to reflect new polynomial order distribution.
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout)
Virtual function for initialisation implementation.
Definition: Driver.cpp:93
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:82
SpatialDomains::MeshGraphSharedPtr m_graph
MeshGraph object.
Definition: Driver.h:91
GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:53
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
@ P
Monomial polynomials .
Definition: BasisType.h:64
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:48
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:890
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.cpp:414
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298

References tinysimd::abs(), Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::Timer::Elapsed(), Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, Nektar::GlobalMapping::Mapping::Load(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::Driver::m_graph, Nektar::SolverUtils::Driver::m_session, Nektar::GlobalMapping::MappingSharedPtr, Nektar::LibUtilities::P, Nektar::LibUtilities::ReduceMax, ReplaceExpansion(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Nektar::SolverUtils::Driver::v_InitObject(), Vmath::Vmul(), Vmath::Vsub(), and Vmath::Vsum().

◆ v_InitObject()

void Nektar::SolverUtils::DriverAdaptive::v_InitObject ( std::ostream &  out = std::cout)
overrideprotectedvirtual

Virtual function for initialisation implementation.

Reimplemented from Nektar::SolverUtils::Driver.

Definition at line 80 of file DriverAdaptive.cpp.

81{
83}

References Nektar::SolverUtils::Driver::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< DriverAdaptive >

friend class MemoryManager< DriverAdaptive >
friend

Definition at line 1 of file DriverAdaptive.h.

Member Data Documentation

◆ className

string Nektar::SolverUtils::DriverAdaptive::className
static
Initial value:
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:67

Name of the class.

Definition at line 64 of file DriverAdaptive.h.

◆ driverLookupId

string Nektar::SolverUtils::DriverAdaptive::driverLookupId
staticprotected
Initial value:
=
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 83 of file DriverAdaptive.h.