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...
 
SOLVER_UTILS_EXPORT ~DriverAdaptive () override
 Destructor. More...
 
SOLVER_UTILS_EXPORT void v_InitObject (std::ostream &out=std::cout) override
 Virtual function for initialisation implementation. More...
 
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 45 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 61 of file DriverAdaptive.cpp.

64 : Driver(pSession, pGraph)
65{
66}
Driver(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Initialises EquationSystem class members.
Definition: Driver.cpp:74

◆ ~DriverAdaptive()

Nektar::SolverUtils::DriverAdaptive::~DriverAdaptive ( )
overrideprotected

Destructor.

Definition at line 71 of file DriverAdaptive.cpp.

72{
73}

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 51 of file DriverAdaptive.h.

54 {
57 p->InitObject();
58 return p;
59 }
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:52

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 447 of file DriverAdaptive.cpp.

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

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 86 of file DriverAdaptive.cpp.

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

79{
81}

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:197
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:65

Name of the class.

Definition at line 62 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 80 of file DriverAdaptive.h.