Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Static Protected Attributes | 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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::DriverAdaptive:
Collaboration graph
[legend]

Static Public Member Functions

static DriverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 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)
 Constructor. More...
 
virtual SOLVER_UTILS_EXPORT ~DriverAdaptive ()
 Destructor. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (std::ostream &out=std::cout)
 Second-stage initialisation. More...
 
virtual SOLVER_UTILS_EXPORT void v_Execute (std::ostream &out=std::cout)
 Virtual function for solve implementation. More...
 
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...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Driver
 Driver (const LibUtilities::SessionReaderSharedPtr pSession)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT
Array< OneD, NekDouble
v_GetRealEvl (void)
 
virtual SOLVER_UTILS_EXPORT
Array< OneD, NekDouble
v_GetImagEvl (void)
 

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
 

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,
EquationSystemSharedPtr
GetEqu ()
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
GetRealEvl (void)
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
GetImagEvl (void)
 
- 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
 I the Coupling between SFD and arnoldi. More...
 
Array< OneD,
EquationSystemSharedPtr
m_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 48 of file DriverAdaptive.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::DriverAdaptive::DriverAdaptive ( const LibUtilities::SessionReaderSharedPtr  pSession)
protected

Constructor.

Definition at line 58 of file DriverAdaptive.cpp.

60  : Driver(pSession)
61 {
62 }
Driver(const LibUtilities::SessionReaderSharedPtr pSession)
Initialises EquationSystem class members.
Definition: Driver.cpp:78
Nektar::SolverUtils::DriverAdaptive::~DriverAdaptive ( )
protectedvirtual

Destructor.

Definition at line 67 of file DriverAdaptive.cpp.

68 {
69 }

Member Function Documentation

static DriverSharedPtr Nektar::SolverUtils::DriverAdaptive::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 54 of file DriverAdaptive.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

56  {
57  DriverSharedPtr p =
59  p->InitObject();
60  return p;
61  }
boost::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:52
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void Nektar::SolverUtils::DriverAdaptive::ReplaceExpansion ( Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
std::map< int, int >  deltaP 
)
protected

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

Parameters
fieldsInput fields.
deltaPMap of polynomial order expansions

Definition at line 385 of file DriverAdaptive.cpp.

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

Referenced by v_Execute().

388 {
389  int nExp, nDim;
390 
391  // Get field definitions
392  std::vector<LibUtilities::FieldDefinitionsSharedPtr> fielddefs =
393  fields[0]->GetFieldDefinitions();
394 
395  int expDim = 2;
396 
397  if (fielddefs[0]->m_numHomogeneousDir == 1)
398  {
399  nDim = 3;
400  }
401  else
402  {
403  nDim = 2;
404  }
405 
406  // Add variables to field definition
407  for (int i = 0; i < fielddefs.size(); ++i)
408  {
409  for (int j = 0; j < fields.num_elements(); ++j)
410  {
411  fielddefs[i]->m_fields.push_back(m_session->GetVariable(j));
412  }
413  }
414 
415  // Get tinyxml objects
416  TiXmlElement *exp_tag = m_session->GetElement("NEKTAR/EXPANSIONS");
417 
418  // Clear current expansions
419  exp_tag->Clear();
420 
421  // Write new expansion information
422  for (int f = 0; f < fielddefs.size(); ++f)
423  {
424  nExp = fielddefs[f]->m_elementIDs.size();
425 
426  //---------------------------------------------
427  // Write ELEMENTS
428  TiXmlElement *elemTag = new TiXmlElement("ELEMENTS");
429  exp_tag->LinkEndChild(elemTag);
430 
431  // Write FIELDS
432  std::string fieldsString;
433  {
434  std::stringstream fieldsStringStream;
435  bool first = true;
436  for (std::vector<int>::size_type i = 0;
437  i < fielddefs[f]->m_fields.size(); i++)
438  {
439  if (!first)
440  {
441  fieldsStringStream << ",";
442  }
443  fieldsStringStream << fielddefs[f]->m_fields[i];
444  first = false;
445  }
446  fieldsString = fieldsStringStream.str();
447  }
448  elemTag->SetAttribute("FIELDS", fieldsString);
449 
450  // Write SHAPE
451  std::string shapeString;
452  {
453  std::stringstream shapeStringStream;
454  shapeStringStream
455  << LibUtilities::ShapeTypeMap[fielddefs[f]->m_shapeType];
456  if (fielddefs[f]->m_numHomogeneousDir == 1)
457  {
458  shapeStringStream << "-HomogenousExp1D";
459  }
460  else if (fielddefs[f]->m_numHomogeneousDir == 2)
461  {
462  shapeStringStream << "-HomogenousExp2D";
463  }
464 
465  shapeString = shapeStringStream.str();
466  }
467  elemTag->SetAttribute("SHAPE", shapeString);
468 
469  // Write BASIS
470  std::string basisString;
471  {
472  std::stringstream basisStringStream;
473  bool first = true;
474  for (std::vector<LibUtilities::BasisType>::size_type i = 0;
475  i < fielddefs[f]->m_basis.size(); i++)
476  {
477  if (!first)
478  {
479  basisStringStream << ",";
480  }
481  basisStringStream
482  << LibUtilities::BasisTypeMap[fielddefs[f]->m_basis[i]];
483  first = false;
484  }
485  basisString = basisStringStream.str();
486  }
487  elemTag->SetAttribute("BASIS", basisString);
488 
489  // Write homogeneuous length details
490  if (fielddefs[f]->m_numHomogeneousDir)
491  {
492  std::string homoLenString;
493  {
494  std::stringstream homoLenStringStream;
495  bool first = true;
496  for (int i = 0; i < fielddefs[f]->m_numHomogeneousDir; ++i)
497  {
498  if (!first)
499  {
500  homoLenStringStream << ",";
501  }
502  homoLenStringStream
503  << fielddefs[f]->m_homogeneousLengths[i];
504  first = false;
505  }
506  homoLenString = homoLenStringStream.str();
507  }
508  elemTag->SetAttribute("HOMOGENEOUSLENGTHS", homoLenString);
509  }
510 
511  // Write homogeneuous planes/lines details
512  if (fielddefs[f]->m_numHomogeneousDir)
513  {
514  if (fielddefs[f]->m_homogeneousYIDs.size() > 0)
515  {
516  std::string homoYIDsString;
517  {
518  std::stringstream homoYIDsStringStream;
519  bool first = true;
520  for (int i = 0; i < fielddefs[f]->m_homogeneousYIDs.size();
521  i++)
522  {
523  if (!first)
524  {
525  homoYIDsStringStream << ",";
526  }
527  homoYIDsStringStream
528  << fielddefs[f]->m_homogeneousYIDs[i];
529  first = false;
530  }
531  homoYIDsString = homoYIDsStringStream.str();
532  }
533  elemTag->SetAttribute("HOMOGENEOUSYIDS", homoYIDsString);
534  }
535 
536  if (fielddefs[f]->m_homogeneousZIDs.size() > 0)
537  {
538  std::string homoZIDsString;
539  {
540  std::stringstream homoZIDsStringStream;
541  bool first = true;
542  for (int i = 0; i < fielddefs[f]->m_homogeneousZIDs.size();
543  i++)
544  {
545  if (!first)
546  {
547  homoZIDsStringStream << ",";
548  }
549  homoZIDsStringStream
550  << fielddefs[f]->m_homogeneousZIDs[i];
551  first = false;
552  }
553  homoZIDsString = homoZIDsStringStream.str();
554  }
555  elemTag->SetAttribute("HOMOGENEOUSZIDS", homoZIDsString);
556  }
557  }
558 
559  // Write NUMMODESPERDIR
560  std::string numModesString;
561  {
562  std::stringstream numModesStringStream;
563 
564  numModesStringStream << "MIXORDER:";
565  bool first = true;
566  int eId;
567  Array<OneD, int> order(nDim, 0);
568  for (int n = 0; n < nExp; n++)
569  {
570  eId = fielddefs[f]->m_elementIDs[n];
571 
572  for (int i = 0; i < expDim; i++)
573  {
574  order[i] = deltaP[eId];
575  }
576 
577  for (int i = 0; i < nDim; i++)
578  {
579  if (fielddefs[f]->m_uniOrder)
580  {
581  order[i] += fielddefs[f]->m_numModes[i];
582  }
583  else
584  {
585  order[i] += fielddefs[f]->m_numModes[n * nDim + i];
586  }
587 
588  if (!first)
589  {
590  numModesStringStream << ",";
591  }
592 
593  numModesStringStream << order[i];
594  first = false;
595  }
596  }
597  numModesString = numModesStringStream.str();
598  }
599  elemTag->SetAttribute("NUMMODESPERDIR", numModesString);
600 
601  // Write IDs. Should ideally look at ways of compressing this stream if
602  // just sequential.
603  elemTag->SetAttribute(
604  "ID", ParseUtils::GenerateSeqString(fielddefs[f]->m_elementIDs));
605  }
606 }
const char *const BasisTypeMap[]
Definition: Foundations.hpp:47
static std::string GenerateSeqString(const std::vector< unsigned int > &elmtids)
Definition: ParseUtils.hpp:159
const char *const ShapeTypeMap[]
Definition: ShapeType.hpp:66
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
void Nektar::SolverUtils::DriverAdaptive::v_Execute ( std::ostream &  out = std::cout)
protectedvirtual

Virtual function for solve implementation.

Implements Nektar::SolverUtils::Driver.

Definition at line 79 of file DriverAdaptive.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTriangle, Nektar::GlobalMapping::Mapping::Load(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::Driver::m_session, Nektar::GlobalMapping::MappingSharedPtr, Nektar::LibUtilities::ReduceMax, ReplaceExpansion(), Nektar::SolverUtils::Driver::v_InitObject(), Vmath::Vmul(), Vmath::Vsub(), and Vmath::Vsum().

80 {
81  time_t starttime, endtime;
82  NekDouble CPUtime;
83 
84  m_equ[0]->PrintSummary(out);
85 
86  // First run using original order
87  time(&starttime);
88  m_equ[0]->DoInitialise();
89 
90  // Obtain initial time in case a restart was used
91  NekDouble startTime = m_equ[0]->GetFinalTime();
92  m_equ[0]->DoSolve();
93 
94  // Load session parameters and solver info
95  bool isHomogeneous1D;
96  int nRuns, minP, maxP, sensorVar;
97  NekDouble lowerTol, upperTol;
98  m_session->LoadParameter ("NumRuns", nRuns, 1);
99  m_session->LoadParameter ("AdaptiveMinModes", minP, 4);
100  m_session->LoadParameter ("AdaptiveMaxModes", maxP, 12);
101  m_session->LoadParameter ("AdaptiveLowerTolerance", lowerTol, 1e-8);
102  m_session->LoadParameter ("AdaptiveUpperTolerance", upperTol, 1e-6);
103  m_session->LoadParameter ("AdaptiveSensorVariable", sensorVar, 0);
104  m_session->MatchSolverInfo("Homogeneous", "1D", isHomogeneous1D, false);
105 
106  // Get number of elements and planes
107  int nExp, nPlanes;
108  if (isHomogeneous1D)
109  {
110  nExp = m_equ[0]->UpdateFields()[0]->GetPlane(0)->GetExpSize();
111  nPlanes = m_equ[0]->UpdateFields()[0]->GetZIDs().num_elements();
112  }
113  else
114  {
115  nExp = m_equ[0]->UpdateFields()[0]->GetExpSize();
116  nPlanes = 1;
117  }
118 
119  int nFields = m_equ[0]->UpdateFields().num_elements();
120  int numSteps = m_session->GetParameter("NumSteps");
121  NekDouble period = m_session->GetParameter("TimeStep") * numSteps;
122 
123  Array<OneD, NekDouble> coeffs, phys, physReduced, tmpArray;
124 
125  // Get mapping
128  m_equ[0]->UpdateFields());
129 
130  // Adaptive loop
132  for (int i = 1; i < nRuns; i++)
133  {
134  // Get field expansions
135  Array<OneD, MultiRegions::ExpListSharedPtr> fields =
136  m_equ[0]->UpdateFields();
137 
138  // Determine the change to be applied in the order
139  map<int, int> deltaP;
140  int offset = 0;
141  for (int n = 0; n < nExp; n++)
142  {
143  offset = fields[sensorVar]->GetPhys_Offset(n);
144  Exp = fields[sensorVar]->GetExp(n);
145 
146  int P = Exp->GetBasis(0)->GetNumModes();
147  int Q = Exp->GetBasis(1)->GetNumModes();
148  int qa = Exp->GetBasis(0)->GetNumPoints();
149  int qb = Exp->GetBasis(1)->GetNumPoints();
150 
151  // Declare orthogonal basis.
152  LibUtilities::PointsKey pa(qa, Exp->GetBasis(0)->GetPointsType());
153  LibUtilities::PointsKey pb(qb, Exp->GetBasis(1)->GetPointsType());
154 
156  switch (Exp->GetGeom()->GetShapeType())
157  {
159  {
160  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P - 1,
161  pa);
162  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, Q - 1,
163  pb);
164  OrthoExp = MemoryManager<
165  StdRegions::StdQuadExp>::AllocateSharedPtr(Ba, Bb);
166  break;
167  }
169  {
170  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, P - 1,
171  pa);
172  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B, Q - 1,
173  pb);
174  OrthoExp =
176  Ba, Bb);
177  break;
178  }
179  default:
180  ASSERTL0(false, "Shape not supported.");
181  break;
182  }
183 
184  int nq = OrthoExp->GetTotPoints();
185 
186  NekDouble error = 0;
187  NekDouble fac = 0;
188  NekDouble tmp = 0;
189 
190  coeffs = Array<OneD, NekDouble>(OrthoExp->GetNcoeffs());
191  physReduced = Array<OneD, NekDouble>(OrthoExp->GetTotPoints());
192  tmpArray = Array<OneD, NekDouble>(OrthoExp->GetTotPoints(), 0.0);
193 
194  // Refinement based only on one variable
195  for (int plane = 0; plane < nPlanes; plane++)
196  {
197  if (isHomogeneous1D)
198  {
199  phys =
200  fields[sensorVar]->GetPlane(plane)->GetPhys() + offset;
201  }
202  else
203  {
204  phys = fields[sensorVar]->GetPhys() + offset;
205  }
206 
207  // Project solution to lower order
208  OrthoExp->FwdTrans(phys, coeffs);
209  OrthoExp->BwdTrans(coeffs, physReduced);
210 
211  // Calculate error =||phys-physReduced||^2 / ||phys||^2
212  Vmath::Vsub(nq, phys, 1, physReduced, 1, tmpArray, 1);
213  Vmath::Vmul(nq, tmpArray, 1, tmpArray, 1, tmpArray, 1);
214  tmp = Exp->Integral(tmpArray);
215 
216  Vmath::Vmul(nq, phys, 1, phys, 1, tmpArray, 1);
217  fac = Exp->Integral(tmpArray);
218 
219  tmp = abs(tmp / fac);
220 
221  if (tmp != tmp)
222  {
223  ASSERTL0(false,
224  "Adaptive procedure encountered NaN value.");
225  }
226 
227  // Get maximum value along planes
228  error = (tmp > error) ? tmp : error;
229  }
230 
231  // Combine planes from different processes
232  m_session->GetComm()->GetColumnComm()->AllReduce(
233  error, LibUtilities::ReduceMax);
234 
235  // Override tolerances if function is defined
236  if (m_session->DefinesFunction("AdaptiveLowerTolerance"))
237  {
238  int nq = Exp->GetTotPoints();
239 
240  // Obtain points from the the element
241  Array<OneD, NekDouble> xc0, xc1, xc2;
242  xc0 = Array<OneD, NekDouble>(nq, 0.0);
243  xc1 = Array<OneD, NekDouble>(nq, 0.0);
244  xc2 = Array<OneD, NekDouble>(nq, 0.0);
245  Exp->GetCoords(xc0, xc1, xc2);
246 
247  // Evaluate function from session file
248  Array<OneD, NekDouble> tolerance(nq, 0.0);
250  m_session->GetFunction("AdaptiveLowerTolerance", 0);
251  ffunc->Evaluate(xc0, xc1, xc2, tolerance);
252  lowerTol = Vmath::Vsum(nq, tolerance, 1) / nq;
253  }
254 
255  if (m_session->DefinesFunction("AdaptiveUpperTolerance"))
256  {
257  int nq = Exp->GetTotPoints();
258 
259  // Obtain points from the the element
260  Array<OneD, NekDouble> xc0, xc1, xc2;
261  xc0 = Array<OneD, NekDouble>(nq, 0.0);
262  xc1 = Array<OneD, NekDouble>(nq, 0.0);
263  xc2 = Array<OneD, NekDouble>(nq, 0.0);
264  Exp->GetCoords(xc0, xc1, xc2);
265 
266  // Evaluate function from session file
267  Array<OneD, NekDouble> tolerance(nq, 0.0);
269  m_session->GetFunction("AdaptiveUpperTolerance", 0);
270  ffunc->Evaluate(xc0, xc1, xc2, tolerance);
271  upperTol = Vmath::Vsum(nq, tolerance, 1) / nq;
272  }
273 
274  // Determine what to do with the polynomial order
275  if ((error > upperTol) && (P < maxP))
276  {
277  deltaP[Exp->GetGeom()->GetGlobalID()] = 1;
278  }
279  else if ((error < lowerTol) && P > minP)
280  {
281  deltaP[Exp->GetGeom()->GetGlobalID()] = -1;
282  }
283  else
284  {
285  deltaP[Exp->GetGeom()->GetGlobalID()] = 0;
286  }
287  }
288 
289  // Write new expansion section to the session reader
290  ReplaceExpansion(fields, deltaP);
291 
292  // Reset GlobalLinSys Manager to avoid using too much memory
293  //
294  // @todo This could be made better by replacing individual matrices
295  // within the linear system.
296  if (LibUtilities::NekManager<MultiRegions::GlobalLinSysKey,
297  MultiRegions::GlobalLinSys>::
298  PoolCreated(std::string("GlobalLinSys")))
299  {
300  LibUtilities::NekManager<MultiRegions::GlobalLinSysKey,
301  MultiRegions::GlobalLinSys>::
302  ClearManager(std::string("GlobalLinSys"));
303  }
304 
305  int chkNumber = m_equ[0]->GetCheckpointNumber();
306  int chkSteps = m_equ[0]->GetCheckpointSteps();
307 
308  // Initialise driver again
310 
311  // Update mapping (must be before m_equ[0]->DoInitialise();)
312  mapping->ReplaceField(m_equ[0]->UpdateFields());
313 
314  // Set chkSteps to zero to avoid writing initial condition
315  m_equ[0]->SetCheckpointSteps(0);
316 
317  // Initialise equation
318  m_equ[0]->DoInitialise();
319  m_equ[0]->SetInitialStep(i * numSteps);
320  m_equ[0]->SetSteps(i * numSteps + numSteps);
321  m_equ[0]->SetTime(startTime + i * period);
322  m_equ[0]->SetBoundaryConditions(startTime + i * period);
323  m_equ[0]->SetCheckpointNumber(chkNumber);
324  m_equ[0]->SetCheckpointSteps(chkSteps);
325 
326  // Project solution to new expansion
327  for (int n = 0; n < nFields; n++)
328  {
329  m_equ[0]->UpdateFields()[n]->ExtractCoeffsToCoeffs(
330  fields[n], fields[n]->GetCoeffs(),
331  m_equ[0]->UpdateFields()[n]->UpdateCoeffs());
332  m_equ[0]->UpdateFields()[n]->BwdTrans_IterPerExp(
333  m_equ[0]->UpdateFields()[n]->GetCoeffs(),
334  m_equ[0]->UpdateFields()[n]->UpdatePhys());
335  }
336 
337  // Solve equation
338  m_equ[0]->DoSolve();
339  }
340 
341  time(&endtime);
342 
343  m_equ[0]->Output();
344 
345  if (m_comm->GetRank() == 0)
346  {
347  CPUtime = difftime(endtime, starttime);
348  cout << "-------------------------------------------" << endl;
349  cout << "Total Computation Time = " << CPUtime << "s" << endl;
350  cout << "-------------------------------------------" << endl;
351  }
352 
353  // Evaluate and output computation time and solution accuracy.
354  // The specific format of the error output is essential for the
355  // regression tests to work.
356 
357  // Evaluate L2 Error
358  for (int i = 0; i < m_equ[0]->GetNvariables(); ++i)
359  {
360  Array<OneD, NekDouble> exactsoln(m_equ[0]->GetTotPoints(), 0.0);
361 
362  // Evaluate "ExactSolution" function, or zero array
363  m_equ[0]->EvaluateExactSolution(i, exactsoln, m_equ[0]->GetFinalTime());
364 
365  NekDouble vL2Error = m_equ[0]->L2Error (i, exactsoln);
366  NekDouble vLinfError = m_equ[0]->LinfError(i, exactsoln);
367 
368  if (m_comm->GetRank() == 0)
369  {
370  out << "L 2 error (variable " << m_equ[0]->GetVariable(i)
371  << ") : " << vL2Error << endl;
372  out << "L inf error (variable " << m_equ[0]->GetVariable(i)
373  << ") : " << vLinfError << endl;
374  }
375  }
376 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
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.
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Principle Orthogonal Functions .
Definition: BasisType.h:47
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout)
Definition: Driver.cpp:92
Principle Orthogonal Functions .
Definition: BasisType.h:46
double NekDouble
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
GLOBAL_MAPPING_EXPORT typedef boost::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:51
boost::shared_ptr< Equation > EquationSharedPtr
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:329
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:723
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
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
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
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:169
void Nektar::SolverUtils::DriverAdaptive::v_InitObject ( std::ostream &  out = std::cout)
protectedvirtual

Second-stage initialisation.

Reimplemented from Nektar::SolverUtils::Driver.

Definition at line 74 of file DriverAdaptive.cpp.

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

75 {
77 }
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout)
Definition: Driver.cpp:92

Friends And Related Function Documentation

friend class MemoryManager< DriverAdaptive >
friend

Definition at line 51 of file DriverAdaptive.h.

Member Data Documentation

string Nektar::SolverUtils::DriverAdaptive::className
static
Initial value:

Name of the class.

Definition at line 64 of file DriverAdaptive.h.

string Nektar::SolverUtils::DriverAdaptive::driverLookupId
staticprotected
Initial value:

Definition at line 84 of file DriverAdaptive.h.