Nektar++
DriverAdaptive.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: DriverAdaptive.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 adaptive solver
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iomanip>
36 
37 #include <GlobalMapping/Mapping.h>
41 #include <StdRegions/StdHexExp.h>
42 #include <StdRegions/StdPrismExp.h>
43 #include <StdRegions/StdPyrExp.h>
44 #include <StdRegions/StdQuadExp.h>
45 #include <StdRegions/StdTetExp.h>
46 #include <StdRegions/StdTriExp.h>
47 
48 using namespace std;
49 
50 namespace Nektar
51 {
52 namespace SolverUtils
53 {
54 
55 string DriverAdaptive::className = GetDriverFactory().RegisterCreatorFunction(
56  "Adaptive", DriverAdaptive::create);
57 string DriverAdaptive::driverLookupId =
58  LibUtilities::SessionReader::RegisterEnumValue("Driver", "Adaptive", 0);
59 
60 /**
61  *
62  */
63 DriverAdaptive::DriverAdaptive(
66  : Driver(pSession, pGraph)
67 {
68 }
69 
70 /**
71  *
72  */
74 {
75 }
76 
77 /**
78  *
79  */
80 void DriverAdaptive::v_InitObject(ostream &out)
81 {
83 }
84 
85 void DriverAdaptive::v_Execute(ostream &out)
86 {
88  NekDouble CPUtime;
89 
90  m_equ[0]->PrintSummary(out);
91 
92  // First run using original order
93  timer.Start();
94  m_equ[0]->DoInitialise();
95 
96  // Obtain initial time in case a restart was used
97  NekDouble startTime = m_equ[0]->GetFinalTime();
98  m_equ[0]->DoSolve();
99 
100  // Load session parameters and solver info
101  bool isHomogeneous1D;
102  int nRuns, minP, maxP, sensorVar;
103  NekDouble lowerTol, upperTol;
104  m_session->LoadParameter("NumRuns", nRuns, 1);
105  m_session->LoadParameter("AdaptiveMinModes", minP, 4);
106  m_session->LoadParameter("AdaptiveMaxModes", maxP, 12);
107  m_session->LoadParameter("AdaptiveLowerTolerance", lowerTol, 1e-8);
108  m_session->LoadParameter("AdaptiveUpperTolerance", upperTol, 1e-6);
109  m_session->LoadParameter("AdaptiveSensorVariable", sensorVar, 0);
110  m_session->MatchSolverInfo("Homogeneous", "1D", isHomogeneous1D, false);
111 
112  // Get number of elements and planes
113  int nExp, nPlanes;
114  if (isHomogeneous1D)
115  {
116  nExp = m_equ[0]->UpdateFields()[0]->GetPlane(0)->GetExpSize();
117  nPlanes = m_equ[0]->UpdateFields()[0]->GetZIDs().size();
118  }
119  else
120  {
121  nExp = m_equ[0]->UpdateFields()[0]->GetExpSize();
122  nPlanes = 1;
123  }
124  int expdim = m_equ[0]->UpdateFields()[0]->GetGraph()->GetMeshDimension();
125 
126  int nFields = m_equ[0]->UpdateFields().size();
127  int numSteps = m_session->GetParameter("NumSteps");
128  NekDouble period = m_session->GetParameter("TimeStep") * numSteps;
129 
130  Array<OneD, NekDouble> coeffs, phys, physReduced, tmpArray;
131 
132  // Get mapping
134  mapping = GlobalMapping::Mapping::Load(m_session, m_equ[0]->UpdateFields());
135 
136  // Adaptive loop
137  Array<OneD, int> P(expdim);
138  Array<OneD, int> numPoints(expdim);
141  for (int i = 1; i < nRuns; i++)
142  {
143  // Get field expansions
145  m_equ[0]->UpdateFields();
146 
147  // Determine the change to be applied in the order
148  map<int, int> deltaP;
149  int offset = 0;
150  for (int n = 0; n < nExp; n++)
151  {
152  offset = fields[sensorVar]->GetPhys_Offset(n);
153  Exp = fields[sensorVar]->GetExp(n);
154 
155  for (int k = 0; k < expdim; ++k)
156  {
157  P[k] = Exp->GetBasis(k)->GetNumModes();
158  numPoints[k] = Exp->GetBasis(k)->GetNumPoints();
159  ptsKey[k] = LibUtilities::PointsKey(
160  numPoints[k], Exp->GetBasis(k)->GetPointsType());
161  }
162 
163  // Declare orthogonal basis.
165  switch (Exp->GetGeom()->GetShapeType())
166  {
168  {
170  ptsKey[0]);
172  ptsKey[1]);
173  OrthoExp = MemoryManager<
174  StdRegions::StdQuadExp>::AllocateSharedPtr(Ba, Bb);
175  break;
176  }
178  {
180  ptsKey[0]);
182  ptsKey[1]);
183  OrthoExp =
185  Ba, Bb);
186  break;
187  }
189  {
191  ptsKey[0]);
193  ptsKey[1]);
195  ptsKey[2]);
196  OrthoExp =
198  Ba, Bb, Bc);
199  break;
200  }
202  {
204  ptsKey[0]);
206  ptsKey[1]);
208  ptsKey[2]);
209  OrthoExp =
211  Ba, Bb, Bc);
212  break;
213  }
215  {
217  ptsKey[0]);
219  ptsKey[1]);
221  ptsKey[2]);
222  OrthoExp = MemoryManager<
223  StdRegions::StdPrismExp>::AllocateSharedPtr(Ba, Bb, Bc);
224  break;
225  }
227  {
229  ptsKey[0]);
231  ptsKey[1]);
233  ptsKey[2]);
234  OrthoExp =
236  Ba, Bb, Bc);
237  break;
238  }
239  default:
240  ASSERTL0(false, "Shape not supported.");
241  break;
242  }
243 
244  int nq = OrthoExp->GetTotPoints();
245 
246  NekDouble error = 0;
247  NekDouble fac = 0;
248  NekDouble tmp = 0;
249 
250  coeffs = Array<OneD, NekDouble>(OrthoExp->GetNcoeffs());
251  physReduced = Array<OneD, NekDouble>(OrthoExp->GetTotPoints());
252  tmpArray = Array<OneD, NekDouble>(OrthoExp->GetTotPoints(), 0.0);
253 
254  // Refinement based only on one variable
255  for (int plane = 0; plane < nPlanes; plane++)
256  {
257  if (isHomogeneous1D)
258  {
259  phys =
260  fields[sensorVar]->GetPlane(plane)->GetPhys() + offset;
261  }
262  else
263  {
264  phys = fields[sensorVar]->GetPhys() + offset;
265  }
266 
267  // Project solution to lower order
268  OrthoExp->FwdTrans(phys, coeffs);
269  OrthoExp->BwdTrans(coeffs, physReduced);
270 
271  // Calculate error =||phys-physReduced||^2 / ||phys||^2
272  Vmath::Vsub(nq, phys, 1, physReduced, 1, tmpArray, 1);
273  Vmath::Vmul(nq, tmpArray, 1, tmpArray, 1, tmpArray, 1);
274  tmp = Exp->Integral(tmpArray);
275 
276  Vmath::Vmul(nq, phys, 1, phys, 1, tmpArray, 1);
277  fac = Exp->Integral(tmpArray);
278 
279  tmp = abs(tmp / fac);
280 
281  if (tmp != tmp)
282  {
283  ASSERTL0(false,
284  "Adaptive procedure encountered NaN value.");
285  }
286 
287  // Get maximum value along planes
288  error = (tmp > error) ? tmp : error;
289  }
290 
291  // Combine planes from different processes
292  m_session->GetComm()->GetColumnComm()->AllReduce(
293  error, LibUtilities::ReduceMax);
294 
295  // Override tolerances if function is defined
296  if (m_session->DefinesFunction("AdaptiveLowerTolerance"))
297  {
298  int nq = Exp->GetTotPoints();
299 
300  // Obtain points from the the element
301  Array<OneD, NekDouble> xc0, xc1, xc2;
302  xc0 = Array<OneD, NekDouble>(nq, 0.0);
303  xc1 = Array<OneD, NekDouble>(nq, 0.0);
304  xc2 = Array<OneD, NekDouble>(nq, 0.0);
305  Exp->GetCoords(xc0, xc1, xc2);
306 
307  // Evaluate function from session file
308  Array<OneD, NekDouble> tolerance(nq, 0.0);
310  m_session->GetFunction("AdaptiveLowerTolerance", 0);
311  ffunc->Evaluate(xc0, xc1, xc2, tolerance);
312  lowerTol = Vmath::Vsum(nq, tolerance, 1) / nq;
313  }
314 
315  if (m_session->DefinesFunction("AdaptiveUpperTolerance"))
316  {
317  int nq = Exp->GetTotPoints();
318 
319  // Obtain points from the the element
320  Array<OneD, NekDouble> xc0, xc1, xc2;
321  xc0 = Array<OneD, NekDouble>(nq, 0.0);
322  xc1 = Array<OneD, NekDouble>(nq, 0.0);
323  xc2 = Array<OneD, NekDouble>(nq, 0.0);
324  Exp->GetCoords(xc0, xc1, xc2);
325 
326  // Evaluate function from session file
327  Array<OneD, NekDouble> tolerance(nq, 0.0);
329  m_session->GetFunction("AdaptiveUpperTolerance", 0);
330  ffunc->Evaluate(xc0, xc1, xc2, tolerance);
331  upperTol = Vmath::Vsum(nq, tolerance, 1) / nq;
332  }
333 
334  // Determine what to do with the polynomial order
335  if ((error > upperTol) && (P[0] < maxP))
336  {
337  deltaP[Exp->GetGeom()->GetGlobalID()] = 1;
338  }
339  else if ((error < lowerTol) && P[0] > minP)
340  {
341  deltaP[Exp->GetGeom()->GetGlobalID()] = -1;
342  }
343  else
344  {
345  deltaP[Exp->GetGeom()->GetGlobalID()] = 0;
346  }
347  }
348 
349  // Write new expansion section to the session reader and re-read graph.
350  ReplaceExpansion(fields, deltaP);
351  m_graph->ReadExpansionInfo();
352 
353  // Reset GlobalLinSys Manager to avoid using too much memory
354  //
355  // @todo This could be made better by replacing individual matrices
356  // within the linear system.
359  PoolCreated(std::string("GlobalLinSys")))
360  {
363  ClearManager(std::string("GlobalLinSys"));
364  }
365 
366  int chkNumber = m_equ[0]->GetCheckpointNumber();
367  int chkSteps = m_equ[0]->GetCheckpointSteps();
368 
369  // Initialise driver again
371 
372  // Update mapping (must be before m_equ[0]->DoInitialise();)
373  mapping->ReplaceField(m_equ[0]->UpdateFields());
374 
375  // Set chkSteps to zero to avoid writing initial condition
376  m_equ[0]->SetCheckpointSteps(0);
377 
378  // Initialise equation
379  m_equ[0]->DoInitialise();
380  m_equ[0]->SetInitialStep(i * numSteps);
381  m_equ[0]->SetSteps(i * numSteps + numSteps);
382  m_equ[0]->SetTime(startTime + i * period);
383  m_equ[0]->SetBoundaryConditions(startTime + i * period);
384  m_equ[0]->SetCheckpointNumber(chkNumber);
385  m_equ[0]->SetCheckpointSteps(chkSteps);
386 
387  // Project solution to new expansion
388  for (int n = 0; n < nFields; n++)
389  {
390  m_equ[0]->UpdateFields()[n]->ExtractCoeffsToCoeffs(
391  fields[n], fields[n]->GetCoeffs(),
392  m_equ[0]->UpdateFields()[n]->UpdateCoeffs());
393  m_equ[0]->UpdateFields()[n]->BwdTrans(
394  m_equ[0]->UpdateFields()[n]->GetCoeffs(),
395  m_equ[0]->UpdateFields()[n]->UpdatePhys());
396  }
397 
398  // Solve equation
399  m_equ[0]->DoSolve();
400  }
401 
402  timer.Stop();
403 
404  m_equ[0]->Output();
405 
406  if (m_comm->GetRank() == 0)
407  {
408  CPUtime = timer.Elapsed().count();
409  cout << "-------------------------------------------" << endl;
410  cout << "Total Computation Time = " << CPUtime << "s" << endl;
411  cout << "-------------------------------------------" << endl;
412  }
413 
414  // Evaluate and output computation time and solution accuracy.
415  // The specific format of the error output is essential for the
416  // regression tests to work.
417 
418  // Evaluate L2 Error
419  for (int i = 0; i < m_equ[0]->GetNvariables(); ++i)
420  {
421  Array<OneD, NekDouble> exactsoln(m_equ[0]->GetTotPoints(), 0.0);
422 
423  // Evaluate "ExactSolution" function, or zero array
424  m_equ[0]->EvaluateExactSolution(i, exactsoln, m_equ[0]->GetFinalTime());
425 
426  NekDouble vL2Error = m_equ[0]->L2Error(i, exactsoln);
427  NekDouble vLinfError = m_equ[0]->LinfError(i, exactsoln);
428 
429  if (m_comm->GetRank() == 0)
430  {
431  out << "L 2 error (variable " << m_equ[0]->GetVariable(i)
432  << ") : " << vL2Error << endl;
433  out << "L inf error (variable " << m_equ[0]->GetVariable(i)
434  << ") : " << vLinfError << endl;
435  }
436  }
437 }
438 
439 /**
440  * @brief Update EXPANSIONS tag inside XML schema to reflect new polynomial
441  * order distribution.
442  *
443  * @param fields Input fields.
444  * @param deltaP Map of polynomial order expansions
445  */
447  Array<OneD, MultiRegions::ExpListSharedPtr> &fields, map<int, int> deltaP)
448 {
449  int nExp, nDim;
450  int expdim = m_equ[0]->UpdateFields()[0]->GetGraph()->GetMeshDimension();
451 
452  // Get field definitions
453  std::vector<LibUtilities::FieldDefinitionsSharedPtr> fielddefs =
454  fields[0]->GetFieldDefinitions();
455 
456  if (fielddefs[0]->m_numHomogeneousDir == 1)
457  {
458  nDim = expdim + 1;
459  }
460  else
461  {
462  nDim = expdim;
463  }
464 
465  // Add variables to field definition
466  for (int i = 0; i < fielddefs.size(); ++i)
467  {
468  for (int j = 0; j < fields.size(); ++j)
469  {
470  fielddefs[i]->m_fields.push_back(m_session->GetVariable(j));
471  }
472  }
473 
474  // Get tinyxml objects
475  TiXmlElement *exp_tag = m_session->GetElement("NEKTAR/EXPANSIONS");
476 
477  // Clear current expansions
478  exp_tag->Clear();
479 
480  // Write new expansion information
481  for (int f = 0; f < fielddefs.size(); ++f)
482  {
483  nExp = fielddefs[f]->m_elementIDs.size();
484 
485  //---------------------------------------------
486  // Write ELEMENTS
487  TiXmlElement *elemTag = new TiXmlElement("ELEMENTS");
488  exp_tag->LinkEndChild(elemTag);
489 
490  // Write FIELDS
491  std::string fieldsString;
492  {
493  std::stringstream fieldsStringStream;
494  bool first = true;
495  for (std::vector<int>::size_type i = 0;
496  i < fielddefs[f]->m_fields.size(); i++)
497  {
498  if (!first)
499  {
500  fieldsStringStream << ",";
501  }
502  fieldsStringStream << fielddefs[f]->m_fields[i];
503  first = false;
504  }
505  fieldsString = fieldsStringStream.str();
506  }
507  elemTag->SetAttribute("FIELDS", fieldsString);
508 
509  // Write SHAPE
510  std::string shapeString;
511  {
512  std::stringstream shapeStringStream;
513  shapeStringStream
514  << LibUtilities::ShapeTypeMap[fielddefs[f]->m_shapeType];
515  if (fielddefs[f]->m_numHomogeneousDir == 1)
516  {
517  shapeStringStream << "-HomogenousExp1D";
518  }
519  else if (fielddefs[f]->m_numHomogeneousDir == 2)
520  {
521  shapeStringStream << "-HomogenousExp2D";
522  }
523 
524  shapeString = shapeStringStream.str();
525  }
526  elemTag->SetAttribute("SHAPE", shapeString);
527 
528  // Write BASIS
529  std::string basisString;
530  {
531  std::stringstream basisStringStream;
532  bool first = true;
533  for (std::vector<LibUtilities::BasisType>::size_type i = 0;
534  i < fielddefs[f]->m_basis.size(); i++)
535  {
536  if (!first)
537  {
538  basisStringStream << ",";
539  }
540  basisStringStream
541  << LibUtilities::BasisTypeMap[fielddefs[f]->m_basis[i]];
542  first = false;
543  }
544  basisString = basisStringStream.str();
545  }
546  elemTag->SetAttribute("BASIS", basisString);
547 
548  // Write homogeneuous length details
549  if (fielddefs[f]->m_numHomogeneousDir)
550  {
551  std::string homoLenString;
552  {
553  std::stringstream homoLenStringStream;
554  bool first = true;
555  for (int i = 0; i < fielddefs[f]->m_numHomogeneousDir; ++i)
556  {
557  if (!first)
558  {
559  homoLenStringStream << ",";
560  }
561  homoLenStringStream
562  << fielddefs[f]->m_homogeneousLengths[i];
563  first = false;
564  }
565  homoLenString = homoLenStringStream.str();
566  }
567  elemTag->SetAttribute("HOMOGENEOUSLENGTHS", homoLenString);
568  }
569 
570  // Write homogeneuous planes/lines details
571  if (fielddefs[f]->m_numHomogeneousDir)
572  {
573  if (fielddefs[f]->m_homogeneousYIDs.size() > 0)
574  {
575  std::string homoYIDsString;
576  {
577  std::stringstream homoYIDsStringStream;
578  bool first = true;
579  for (int i = 0; i < fielddefs[f]->m_homogeneousYIDs.size();
580  i++)
581  {
582  if (!first)
583  {
584  homoYIDsStringStream << ",";
585  }
586  homoYIDsStringStream
587  << fielddefs[f]->m_homogeneousYIDs[i];
588  first = false;
589  }
590  homoYIDsString = homoYIDsStringStream.str();
591  }
592  elemTag->SetAttribute("HOMOGENEOUSYIDS", homoYIDsString);
593  }
594 
595  if (fielddefs[f]->m_homogeneousZIDs.size() > 0)
596  {
597  std::string homoZIDsString;
598  {
599  std::stringstream homoZIDsStringStream;
600  bool first = true;
601  for (int i = 0; i < fielddefs[f]->m_homogeneousZIDs.size();
602  i++)
603  {
604  if (!first)
605  {
606  homoZIDsStringStream << ",";
607  }
608  homoZIDsStringStream
609  << fielddefs[f]->m_homogeneousZIDs[i];
610  first = false;
611  }
612  homoZIDsString = homoZIDsStringStream.str();
613  }
614  elemTag->SetAttribute("HOMOGENEOUSZIDS", homoZIDsString);
615  }
616  }
617 
618  // Write NUMMODESPERDIR
619  std::string numModesString;
620  {
621  std::stringstream numModesStringStream;
622 
623  numModesStringStream << "MIXORDER:";
624  bool first = true;
625  int eId;
626  Array<OneD, int> order(nDim, 0);
627  for (int n = 0; n < nExp; n++)
628  {
629  eId = fielddefs[f]->m_elementIDs[n];
630 
631  for (int i = 0; i < expdim; i++)
632  {
633  order[i] = deltaP[eId];
634  }
635 
636  for (int i = 0; i < nDim; i++)
637  {
638  if (fielddefs[f]->m_uniOrder)
639  {
640  order[i] += fielddefs[f]->m_numModes[i];
641  }
642  else
643  {
644  order[i] += fielddefs[f]->m_numModes[n * nDim + i];
645  }
646 
647  if (!first)
648  {
649  numModesStringStream << ",";
650  }
651 
652  numModesStringStream << order[i];
653  first = false;
654  }
655  }
656  numModesString = numModesStringStream.str();
657  }
658  elemTag->SetAttribute("NUMMODESPERDIR", numModesString);
659 
660  // Write IDs. Should ideally look at ways of compressing this stream if
661  // just sequential.
662  elemTag->SetAttribute(
663  "ID", ParseUtils::GenerateSeqString(fielddefs[f]->m_elementIDs));
664  }
665 }
666 
667 } // namespace SolverUtils
668 } // namespace Nektar
#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:270
Describes the specification for a Basis.
Definition: Basis.h:50
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
Defines a specification for a set of points.
Definition: Points.h:59
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
A global linear system.
Definition: GlobalLinSys.h:72
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
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) override
Second-stage initialisation.
virtual SOLVER_UTILS_EXPORT void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
virtual SOLVER_UTILS_EXPORT ~DriverAdaptive()
Destructor.
Base class for the development of solvers.
Definition: Driver.h:67
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout)
Definition: Driver.cpp:87
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
SpatialDomains::MeshGraphSharedPtr m_graph
MeshGraph object.
Definition: Driver.h:94
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:97
Class representing a prismatic element in reference space.
Definition: StdPrismExp.h:49
GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:50
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:79
const char *const BasisTypeMap[]
Definition: Foundations.hpp:46
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
@ 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
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:64
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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:209
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:895
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:419
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298