Nektar++
Classes | Typedefs | Enumerations | Functions | Variables
Nektar::FieldUtils Namespace Reference

Classes

struct  Field
 
class  FieldConvertComm
 
class  InputDat
 Input module for Xml files. More...
 
class  InputFld
 
class  InputNek5000
 
class  InputPts
 
class  InputSemtex
 
class  InputXml
 
class  Interpolator
 A class that contains algorithms for interpolation between pts fields, expansions and different meshes. More...
 
struct  ConfigOption
 Represents a command-line configuration option. More...
 
class  Module
 
class  InputModule
 Abstract base class for input modules. More...
 
class  ProcessModule
 Abstract base class for processing modules. More...
 
class  OutputModule
 Abstract base class for output modules. More...
 
class  Octree
 
class  OutputFileBase
 Converter from fld to vtk. More...
 
class  OutputFld
 Output to fld format. More...
 
class  OutputInfo
 
class  OutputPts
 Converter from fld to pts. More...
 
class  OutputStdOut
 
class  OutputTecplot
 Tecplot output class. More...
 
class  OutputTecplotBinary
 Tecplot output class, specifically for binary field output. More...
 
class  OutputVtk
 Converter from fld to vtk. More...
 
class  OutputVtkBase
 Converter from fld to vtk. More...
 
class  OutputXml
 Converter from fld to vtk. More...
 
class  ProcessAddCompositeID
 This processing module adds a fld with the composite ID. More...
 
class  ProcessAddFld
 This processing module scales the input fld file. More...
 
class  ProcessBodyFittedVelocity
 This processing module calculates the wall shear stress and adds it as an extra-field to the output file, and writes it to a surface output file. More...
 
class  ProcessBoundaryExtract
 This processing module sets up for the boundary field to be extracted. More...
 
class  ProcessC0Projection
 This processing module calculates the Q Criterion and adds it as an extra-field to the output file. More...
 
class  ProcessCFL
 This processing module calculates the CFL and adds it as an extra-field to the output file. More...
 
class  ProcessCombineAvg
 This processing module combines two fld files containing average fields. More...
 
class  ProcessConcatenateFld
 This processing module sets up for the boundary field to be extracted. More...
 
class  ProcessCreateExp
 This processing module scales the input fld file. More...
 
class  ProcessDeform
 
struct  TriFaceIDs
 
struct  TriFaceHash
 
class  ProcessDisplacement
 
class  ProcessDOF
 This processing module calculates the number of DOF. More...
 
class  ProcessEquiSpacedOutput
 This processing module interpolates one field to another. More...
 
class  ProcessFieldFromString
 This processing module adds a new field from a string definition. More...
 
class  ProcessGrad
 This processing module calculates the gradient and adds it as an extra-field to the output file. More...
 
class  ProcessHalfModeToFourier
 This processing modifies Fourier Half Mode to a Fourier expansion. More...
 
class  ProcessHomogeneousPlane
 This processing module replaces all expansions by a single plane from 3DH1D fields, defined by the parameter planeid. More...
 
class  ProcessHomogeneousStretch
 This processing module stretches the homogeneous direction of a 3DH1D expansion by an integer factor. More...
 
class  ProcessInnerProduct
 This processing module computes the inner product between two fields. More...
 
class  ProcessInterpField
 This processing module interpolates one field to another. More...
 
class  ProcessInterpPointDataToFld
 This processing module interpolates one field to another. More...
 
class  ProcessInterpPoints
 This processing module interpolates one field to another. More...
 
class  ProcessInterpPtsToPts
 This processing module interpolates one field to another. More...
 
class  Iso
 
class  IsoVertex
 
class  ProcessIsoContour
 This processing module extracts an isocontour. More...
 
class  ProcessJacobianEnergy
 This processing module scales the input fld file. More...
 
class  ProcessL2Criterion
 This processing module calculates the Lambda 2 Criterion and adds it as an extra-field to the output file. More...
 
class  ProcessMapping
 This processing module scales the input fld file. More...
 
class  ProcessMean
 This processing module computes the mean of each field. More...
 
class  ProcessMeanMode
 This processing module replaces all expansions by the mean mode from 3DH1D fields. More...
 
class  ProcessMultiShear
 This processing module calculates the shear stress metrics and writes it to a surface output file. More...
 
class  ProcessNumModes
 This processing module determine the number of modes and adds it as an extra-field to the output file. More...
 
class  ProcessPhiFromFile
 
class  ProcessPointDataToFld
 This processing module interpolates one field to another. More...
 
class  ProcessPrintFldNorms
 This processing module prints the L2 and LInf norms of the variables in the field. More...
 
class  ProcessQCriterion
 This processing module calculates the Q Criterion and adds it as an extra-field to the output file. More...
 
class  ProcessQualityMetric
 This processing module scales the input fld file. More...
 
class  ProcessRemoveField
 This processing module adds a new field from a string definition. More...
 
class  ProcessScaleInFld
 This processing module scales the input fld file. More...
 
class  ProcessScalGrad
 This processing module calculates the scalar gradient field and writes it to a surface output file. More...
 
class  ProcessStreamFunction
 This processing module calculates the stream function of a 2D field and adds it as an extra-field to the output file. More...
 
class  ProcessSurfDistance
 This processing module calculates the height of an element connected to a surface and adds it as an extra-field to the output file. More...
 
class  ProcessVelocityDivergence
 This processing module calculates the divergence of the velocity field and adds it as an extra-field to the output file. More...
 
class  ProcessVorticity
 This processing module calculates the vorticity and adds it as an extra-field to the output file. More...
 
class  ProcessWallNormalData
 This processing module calculates the wall shear stress and adds it as an extra-field to the output file, and writes it to a surface output file. More...
 
class  ProcessWSS
 This processing module calculates the wall shear stress and adds it as an extra-field to the output file, and writes it to a surface output file. More...
 

Typedefs

typedef std::shared_ptr< FieldFieldSharedPtr
 
typedef std::shared_ptr< Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > InterpolatorSharedPtr
 
typedef std::shared_ptr< InputModuleInputModuleSharedPtr
 
typedef std::pair< ModuleType, std::string > ModuleKey
 
typedef std::shared_ptr< ModuleModuleSharedPtr
 
typedef LibUtilities::NekFactory< ModuleKey, Module, FieldSharedPtrModuleFactory
 
typedef std::unordered_map< TriFaceIDs, int, TriFaceHashTriFaceMap
 
typedef std::shared_ptr< IsoIsoSharedPtr
 

Enumerations

enum  ModuleType { eInputModule , eProcessModule , eOutputModule , SIZE_ModuleType }
 
enum  ModulePriority {
  eCreateGraph , eCreateFieldData , eModifyFieldData , eCreateExp ,
  eFillExp , eModifyExp , eBndExtraction , eCreatePts ,
  eConvertExpToPts , eModifyPts , eOutput , SIZE_ModulePriority
}
 
enum  TecplotZoneType {
  eOrdered = 0 , eFELineSeg , eFETriangle , eFEQuadrilateral ,
  eFETetrahedron , eFEBrick , eFEPolygon , eFEPolyhedron
}
 

Functions

ModuleFactoryGetModuleFactory ()
 
std::ostream & operator<< (std::ostream &os, const ModuleKey &rhs)
 
template<typename T >
void swap_endian (T &u)
 Swap endian ordering of the input variable. More...
 
template<typename T >
void swap_endian (std::vector< T > &u)
 
template<typename T >
void WriteStream (std::ostream &outfile, T data)
 Helper function to write binary data to stream. More...
 
template<>
void WriteStream (std::ostream &outfile, std::string data)
 Specialisation of WriteStream to support writing std::string. More...
 
template<typename T >
void WriteStream (std::ostream &outfile, Array< OneD, T > data)
 Specialisation of WriteStream to support writing Nektar::Array datatype. More...
 
template<typename T >
void WriteStream (std::ostream &outfile, std::vector< T > data)
 Specialisation of WriteStream to support writing std::vector datatype. More...
 
bool operator== (TriFaceIDs const &p1, TriFaceIDs const &p2)
 
void TwoPairs (Array< OneD, NekDouble > &cx, Array< OneD, NekDouble > &cy, Array< OneD, NekDouble > &cz, int &pr)
 
void ThreeSimilar (const int i, const int j, const int k, int &pr, int &ps)
 
bool operator== (const IsoVertex &x, const IsoVertex &y)
 
bool operator!= (const IsoVertex &x, const IsoVertex &y)
 
void MatSymEVals (NekDouble d1, NekDouble d2, NekDouble d3, NekDouble a, NekDouble b, NekDouble c, NekDouble &l1, NekDouble &l2, NekDouble &l3)
 Calculates eigenvalues of a 3x3 Symmetric matrix. More...
 
vector< DNekMatMappingIdealToRef (SpatialDomains::GeometrySharedPtr geom, StdRegions::StdExpansionSharedPtr chi)
 

Variables

const std::string ModuleTypeMap [] = {"Input", "Process", "Output"}
 
const char *const ModulePriorityMap []
 
std::string TecplotZoneTypeMap []
 

Typedef Documentation

◆ FieldSharedPtr

typedef std::shared_ptr<Field> Nektar::FieldUtils::FieldSharedPtr

Definition at line 991 of file Field.hpp.

◆ InputModuleSharedPtr

Definition at line 285 of file Module.h.

◆ InterpolatorSharedPtr

Definition at line 108 of file FieldUtils/Interpolator.h.

◆ IsoSharedPtr

typedef std::shared_ptr<Iso> Nektar::FieldUtils::IsoSharedPtr

Definition at line 186 of file ProcessIsoContour.h.

◆ ModuleFactory

Definition at line 323 of file Module.h.

◆ ModuleKey

typedef std::pair<ModuleType, std::string> Nektar::FieldUtils::ModuleKey

Definition at line 317 of file Module.h.

◆ ModuleSharedPtr

typedef std::shared_ptr<Module> Nektar::FieldUtils::ModuleSharedPtr

Definition at line 321 of file Module.h.

◆ TriFaceMap

typedef std::unordered_map<TriFaceIDs, int, TriFaceHash> Nektar::FieldUtils::TriFaceMap

Definition at line 95 of file ProcessDisplacement.cpp.

Enumeration Type Documentation

◆ ModulePriority

Enumerator
eCreateGraph 
eCreateFieldData 
eModifyFieldData 
eCreateExp 
eFillExp 
eModifyExp 
eBndExtraction 
eCreatePts 
eConvertExpToPts 
eModifyPts 
eOutput 
SIZE_ModulePriority 

Definition at line 76 of file Module.h.

◆ ModuleType

Denotes different types of mesh converter modules: so far only input, output and process modules are defined.

Enumerator
eInputModule 
eProcessModule 
eOutputModule 
SIZE_ModuleType 

Definition at line 66 of file Module.h.

67 {
72 };

◆ TecplotZoneType

Enumerator
eOrdered 
eFELineSeg 
eFETriangle 
eFEQuadrilateral 
eFETetrahedron 
eFEBrick 
eFEPolygon 
eFEPolyhedron 

Definition at line 46 of file OutputTecplot.h.

Function Documentation

◆ GetModuleFactory()

FIELD_UTILS_EXPORT ModuleFactory & Nektar::FieldUtils::GetModuleFactory ( )

Returns an instance of the module factory, held as a singleton.

Definition at line 49 of file Module.cpp.

50 {
51  static ModuleFactory instance;
52  return instance;
53 }
LibUtilities::NekFactory< ModuleKey, Module, FieldSharedPtr > ModuleFactory
Definition: Module.h:323

Referenced by Nektar::SolverUtils::FilterFieldConvert::CreateModules(), main(), Module_Create(), and Module_Register().

◆ MappingIdealToRef()

vector<DNekMat> Nektar::FieldUtils::MappingIdealToRef ( SpatialDomains::GeometrySharedPtr  geom,
StdRegions::StdExpansionSharedPtr  chi 
)
inline

Definition at line 115 of file ProcessQualityMetric.cpp.

117 {
118  vector<DNekMat> ret;
119 
120  if (geom->GetShapeType() == LibUtilities::eQuadrilateral)
121  {
122  vector<Array<OneD, NekDouble>> xy;
123  for (int i = 0; i < geom->GetNumVerts(); i++)
124  {
125  Array<OneD, NekDouble> loc(2);
126  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
127  p->GetCoords(loc);
128  xy.push_back(loc);
129  }
130 
131  Array<OneD, const LibUtilities::BasisSharedPtr> b = chi->GetBase();
132  Array<OneD, NekDouble> u = b[0]->GetZ();
133  Array<OneD, NekDouble> v = b[1]->GetZ();
134 
135  for (int j = 0; j < b[1]->GetNumPoints(); j++)
136  {
137  for (int i = 0; i < b[0]->GetNumPoints(); i++)
138  {
139  NekDouble a1 = 0.5 * (1.0 - u[i]), a2 = 0.5 * (1.0 + u[i]);
140  NekDouble b1 = 0.5 * (1.0 - v[j]), b2 = 0.5 * (1.0 + v[j]);
141  DNekMat dxdz(2, 2, 1.0, eFULL);
142 
143  dxdz(0, 0) = 0.5 * (-b1 * xy[0][0] + b1 * xy[1][0] +
144  b2 * xy[2][0] - b2 * xy[3][0]);
145  dxdz(1, 0) = 0.5 * (-b1 * xy[0][1] + b1 * xy[1][1] +
146  b2 * xy[2][1] - b2 * xy[3][1]);
147 
148  dxdz(0, 1) = 0.5 * (-a1 * xy[0][0] - a2 * xy[1][0] +
149  a2 * xy[2][0] + a1 * xy[3][0]);
150  dxdz(1, 1) = 0.5 * (-a1 * xy[0][1] - a2 * xy[1][1] +
151  a2 * xy[2][1] + a1 * xy[3][1]);
152 
153  dxdz.Invert();
154  ret.push_back(dxdz);
155  }
156  }
157  }
158  else if (geom->GetShapeType() == LibUtilities::eTriangle)
159  {
160  vector<Array<OneD, NekDouble>> xy;
161  for (int i = 0; i < geom->GetNumVerts(); i++)
162  {
163  Array<OneD, NekDouble> loc(2);
164  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
165  p->GetCoords(loc);
166  xy.push_back(loc);
167  }
168 
169  Array<OneD, const LibUtilities::BasisSharedPtr> b = chi->GetBase();
170  Array<OneD, NekDouble> u = b[0]->GetZ();
171  Array<OneD, NekDouble> v = b[1]->GetZ();
172 
173  for (int i = 0; i < b[0]->GetNumPoints(); i++)
174  {
175  for (int j = 0; j < b[1]->GetNumPoints(); j++)
176  {
177  DNekMat dxdz(2, 2, 1.0, eFULL);
178  dxdz(0, 0) = -xy[0][0] / 2.0 + xy[1][0] / 2.0;
179 
180  dxdz(0, 1) = -xy[0][0] / 2.0 + xy[2][0] / 2.0;
181 
182  dxdz(1, 0) = -xy[0][1] / 2.0 + xy[1][1] / 2.0;
183 
184  dxdz(1, 1) = -xy[0][1] / 2.0 + xy[2][1] / 2.0;
185 
186  dxdz.Invert();
187  ret.push_back(dxdz);
188  }
189  }
190  }
191  else if (geom->GetShapeType() == LibUtilities::eTetrahedron)
192  {
193  vector<Array<OneD, NekDouble>> xyz;
194  for (int i = 0; i < geom->GetNumVerts(); i++)
195  {
196  Array<OneD, NekDouble> loc(3);
197  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
198  p->GetCoords(loc);
199  xyz.push_back(loc);
200  }
201 
202  Array<OneD, const LibUtilities::BasisSharedPtr> b = chi->GetBase();
203  Array<OneD, NekDouble> u = b[0]->GetZ();
204  Array<OneD, NekDouble> v = b[1]->GetZ();
205  Array<OneD, NekDouble> z = b[2]->GetZ();
206 
207  for (int i = 0; i < b[0]->GetNumPoints(); i++)
208  {
209  for (int j = 0; j < b[1]->GetNumPoints(); j++)
210  {
211  for (int k = 0; k < b[2]->GetNumPoints(); k++)
212  {
213  DNekMat dxdz(3, 3, 1.0, eFULL);
214  dxdz(0, 0) = -xyz[0][0] / 2.0 + xyz[1][0] / 2.0;
215 
216  dxdz(0, 1) = -xyz[0][0] / 2.0 + xyz[2][0] / 2.0;
217 
218  dxdz(0, 2) = -xyz[0][0] / 2.0 + xyz[3][0] / 2.0;
219 
220  dxdz(1, 0) = -xyz[0][1] / 2.0 + xyz[1][1] / 2.0;
221 
222  dxdz(1, 1) = -xyz[0][1] / 2.0 + xyz[2][1] / 2.0;
223 
224  dxdz(1, 2) = -xyz[0][1] / 2.0 + xyz[3][1] / 2.0;
225 
226  dxdz(2, 0) = -xyz[0][2] / 2.0 + xyz[1][2] / 2.0;
227 
228  dxdz(2, 1) = -xyz[0][2] / 2.0 + xyz[2][2] / 2.0;
229 
230  dxdz(2, 2) = -xyz[0][2] / 2.0 + xyz[3][2] / 2.0;
231 
232  dxdz.Invert();
233  ret.push_back(dxdz);
234  }
235  }
236  }
237  }
238  else if (geom->GetShapeType() == LibUtilities::ePrism)
239  {
240  vector<Array<OneD, NekDouble>> xyz;
241  for (int i = 0; i < geom->GetNumVerts(); i++)
242  {
243  Array<OneD, NekDouble> loc(3);
244  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
245  p->GetCoords(loc);
246  xyz.push_back(loc);
247  }
248 
249  Array<OneD, const LibUtilities::BasisSharedPtr> b = chi->GetBase();
250  Array<OneD, NekDouble> eta1 = b[0]->GetZ();
251  Array<OneD, NekDouble> eta2 = b[1]->GetZ();
252  Array<OneD, NekDouble> eta3 = b[2]->GetZ();
253 
254  for (int k = 0; k < b[2]->GetNumPoints(); k++)
255  {
256  for (int j = 0; j < b[1]->GetNumPoints(); j++)
257  {
258  for (int i = 0; i < b[0]->GetNumPoints(); i++)
259  {
260  NekDouble xi1 = 0.5 * (1 + eta1[i]) * (1 - eta3[k]) - 1.0;
261  NekDouble a2 = 0.5 * (1 + xi1);
262  NekDouble b1 = 0.5 * (1 - eta2[j]),
263  b2 = 0.5 * (1 + eta2[j]);
264  NekDouble c1 = 0.5 * (1 - eta3[k]),
265  c2 = 0.5 * (1 + eta3[k]);
266 
267  DNekMat dxdz(3, 3, 1.0, eFULL);
268 
269  dxdz(0, 0) = 0.5 * (-b1 * xyz[0][0] + b1 * xyz[1][0] +
270  b2 * xyz[2][0] - b2 * xyz[3][0]);
271  dxdz(1, 0) = 0.5 * (-b1 * xyz[0][1] + b1 * xyz[1][1] +
272  b2 * xyz[2][1] - b2 * xyz[3][1]);
273  dxdz(2, 0) = 0.5 * (-b1 * xyz[0][2] + b1 * xyz[1][2] +
274  b2 * xyz[2][2] - b2 * xyz[3][2]);
275 
276  dxdz(0, 1) = 0.5 * ((a2 - c1) * xyz[0][0] - a2 * xyz[1][0] +
277  a2 * xyz[2][0] + (c1 - a2) * xyz[3][0] -
278  c2 * xyz[4][0] + c2 * xyz[5][0]);
279  dxdz(1, 1) = 0.5 * ((a2 - c1) * xyz[0][1] - a2 * xyz[1][1] +
280  a2 * xyz[2][1] + (c1 - a2) * xyz[3][1] -
281  c2 * xyz[4][1] + c2 * xyz[5][1]);
282  dxdz(2, 1) = 0.5 * ((a2 - c1) * xyz[0][2] - a2 * xyz[1][2] +
283  a2 * xyz[2][2] + (c1 - a2) * xyz[3][2] -
284  c2 * xyz[4][2] + c2 * xyz[5][2]);
285 
286  dxdz(0, 2) = 0.5 * (-b1 * xyz[0][0] - b2 * xyz[3][0] +
287  b1 * xyz[4][0] + b2 * xyz[5][0]);
288  dxdz(1, 2) = 0.5 * (-b1 * xyz[0][1] - b2 * xyz[3][1] +
289  b1 * xyz[4][1] + b2 * xyz[5][1]);
290  dxdz(2, 2) = 0.5 * (-b1 * xyz[0][2] - b2 * xyz[3][2] +
291  b1 * xyz[4][2] + b2 * xyz[5][2]);
292 
293  dxdz.Invert();
294  ret.push_back(dxdz);
295  }
296  }
297  }
298  }
299  else if (geom->GetShapeType() == LibUtilities::eHexahedron)
300  {
301  vector<Array<OneD, NekDouble>> xyz;
302  for (int i = 0; i < geom->GetNumVerts(); i++)
303  {
304  Array<OneD, NekDouble> loc(3);
305  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
306  p->GetCoords(loc);
307  xyz.push_back(loc);
308  }
309 
310  Array<OneD, const LibUtilities::BasisSharedPtr> b = chi->GetBase();
311  Array<OneD, NekDouble> eta1 = b[0]->GetZ();
312  Array<OneD, NekDouble> eta2 = b[1]->GetZ();
313  Array<OneD, NekDouble> eta3 = b[2]->GetZ();
314 
315  for (int k = 0; k < b[2]->GetNumPoints(); k++)
316  {
317  for (int j = 0; j < b[1]->GetNumPoints(); j++)
318  {
319  for (int i = 0; i < b[0]->GetNumPoints(); i++)
320  {
321  NekDouble a1 = 0.5 * (1 - eta1[i]);
322  NekDouble a2 = 0.5 * (1 + eta1[i]);
323  NekDouble b1 = 0.5 * (1 - eta2[j]),
324  b2 = 0.5 * (1 + eta2[j]);
325  NekDouble c1 = 0.5 * (1 - eta3[k]),
326  c2 = 0.5 * (1 + eta3[k]);
327 
328  DNekMat dxdz(3, 3, 1.0, eFULL);
329 
330  dxdz(0, 0) =
331  -0.5 * b1 * c1 * xyz[0][0] + 0.5 * b1 * c1 * xyz[1][0] +
332  0.5 * b2 * c1 * xyz[2][0] - 0.5 * b2 * c1 * xyz[3][0] -
333  0.5 * b1 * c2 * xyz[5][0] + 0.5 * b1 * c2 * xyz[5][0] +
334  0.5 * b2 * c2 * xyz[6][0] - 0.5 * b2 * c2 * xyz[7][0];
335  dxdz(1, 0) =
336  -0.5 * b1 * c1 * xyz[0][1] + 0.5 * b1 * c1 * xyz[1][1] +
337  0.5 * b2 * c1 * xyz[2][1] - 0.5 * b2 * c1 * xyz[3][1] -
338  0.5 * b1 * c2 * xyz[5][1] + 0.5 * b1 * c2 * xyz[5][1] +
339  0.5 * b2 * c2 * xyz[6][1] - 0.5 * b2 * c2 * xyz[7][1];
340  dxdz(2, 0) =
341  -0.5 * b1 * c1 * xyz[0][2] + 0.5 * b1 * c1 * xyz[1][2] +
342  0.5 * b2 * c1 * xyz[2][2] - 0.5 * b2 * c1 * xyz[3][2] -
343  0.5 * b1 * c2 * xyz[5][2] + 0.5 * b1 * c2 * xyz[5][2] +
344  0.5 * b2 * c2 * xyz[6][2] - 0.5 * b2 * c2 * xyz[7][2];
345 
346  dxdz(0, 1) =
347  -0.5 * a1 * c1 * xyz[0][0] - 0.5 * a2 * c1 * xyz[1][0] +
348  0.5 * a2 * c1 * xyz[2][0] + 0.5 * a1 * c1 * xyz[3][0] -
349  0.5 * a1 * c2 * xyz[5][0] - 0.5 * a2 * c2 * xyz[5][0] +
350  0.5 * a2 * c2 * xyz[6][0] + 0.5 * a1 * c2 * xyz[7][0];
351  dxdz(1, 1) =
352  -0.5 * a1 * c1 * xyz[0][1] - 0.5 * a2 * c1 * xyz[1][1] +
353  0.5 * a2 * c1 * xyz[2][1] + 0.5 * a1 * c1 * xyz[3][1] -
354  0.5 * a1 * c2 * xyz[5][1] - 0.5 * a2 * c2 * xyz[5][1] +
355  0.5 * a2 * c2 * xyz[6][1] + 0.5 * a1 * c2 * xyz[7][1];
356  dxdz(2, 1) =
357  -0.5 * a1 * c1 * xyz[0][2] - 0.5 * a2 * c1 * xyz[1][2] +
358  0.5 * a2 * c1 * xyz[2][2] + 0.5 * a1 * c1 * xyz[3][2] -
359  0.5 * a1 * c2 * xyz[5][2] - 0.5 * a2 * c2 * xyz[5][2] +
360  0.5 * a2 * c2 * xyz[6][2] + 0.5 * a1 * c2 * xyz[7][2];
361 
362  dxdz(0, 0) =
363  -0.5 * b1 * a1 * xyz[0][0] - 0.5 * b1 * a2 * xyz[1][0] -
364  0.5 * b2 * a2 * xyz[2][0] - 0.5 * b2 * a1 * xyz[3][0] +
365  0.5 * b1 * a1 * xyz[5][0] + 0.5 * b1 * a2 * xyz[5][0] +
366  0.5 * b2 * a2 * xyz[6][0] + 0.5 * b2 * a1 * xyz[7][0];
367  dxdz(1, 0) =
368  -0.5 * b1 * a1 * xyz[0][1] - 0.5 * b1 * a2 * xyz[1][1] -
369  0.5 * b2 * a2 * xyz[2][1] - 0.5 * b2 * a1 * xyz[3][1] +
370  0.5 * b1 * a1 * xyz[5][1] + 0.5 * b1 * a2 * xyz[5][1] +
371  0.5 * b2 * a2 * xyz[6][1] + 0.5 * b2 * a1 * xyz[7][1];
372  dxdz(2, 0) =
373  -0.5 * b1 * a1 * xyz[0][2] - 0.5 * b1 * a2 * xyz[1][2] -
374  0.5 * b2 * a2 * xyz[2][2] - 0.5 * b2 * a1 * xyz[3][2] +
375  0.5 * b1 * a1 * xyz[5][2] + 0.5 * b1 * a2 * xyz[5][2] +
376  0.5 * b2 * a2 * xyz[6][2] + 0.5 * b2 * a1 * xyz[7][2];
377 
378  dxdz.Invert();
379  ret.push_back(dxdz);
380  }
381  }
382  }
383  }
384  else
385  {
386  ASSERTL0(false, "not coded");
387  }
388 
389  return ret;
390 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:59
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:50
double NekDouble

References ASSERTL0, Nektar::eFULL, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, CG_Iterations::loc, and CellMLToNektar.cellml_metadata::p.

Referenced by Nektar::FieldUtils::ProcessQualityMetric::GetQ().

◆ MatSymEVals()

void Nektar::FieldUtils::MatSymEVals ( NekDouble  d1,
NekDouble  d2,
NekDouble  d3,
NekDouble  a,
NekDouble  b,
NekDouble  c,
NekDouble l1,
NekDouble l2,
NekDouble l3 
)

Calculates eigenvalues of a 3x3 Symmetric matrix.

Parameters
d1,d2,d3- matrix diagonal entries at [0,0], [1,1] and [2,2]
a- matrix value at [0,1] and [1,0]
b- matrix value at [0,2] and [2,0]
c- matrix value at [1,2] and [2,1]
l1,l2,l3the computed eigenvalues, ordered l3 >= l2 >= l1

Definition at line 72 of file ProcessL2Criterion.cpp.

75 {
76  NekDouble p = a * a + b * b + c * c;
77  if (p == 0)
78  {
79  l1 = d1;
80  l2 = d2;
81  l3 = d3;
82  if (l1 > l3)
83  {
84  swap(l1, l3);
85  }
86  if (l1 > l2)
87  {
88  swap(l1, l2);
89  }
90  if (l2 > l3)
91  {
92  swap(l2, l3);
93  }
94  }
95  else
96  {
97  NekDouble q = (d1 + d2 + d3) / 3.0;
98  p = (d1 - q) * (d1 - q) + (d2 - q) * (d2 - q) + (d3 - q) * (d3 - q) +
99  2.0 * p;
100  p = sqrt(p / 6.0);
101  NekDouble r =
102  -0.5 *
103  (a * a * d3 - a * a * q - 2.0 * a * b * c + b * b * d2 - b * b * q +
104  c * c * d1 - c * c * q - d1 * d2 * d3 + d1 * d2 * q + d1 * d3 * q -
105  d1 * q * q + d2 * d3 * q - d2 * q * q - d3 * q * q + q * q * q) /
106  (p * p * p);
107 
108  NekDouble phi = 0;
109  if (r <= -1)
110  {
111  phi = M_PI / 3.0;
112  }
113  else if (r >= 1)
114  {
115  phi = 0.0;
116  }
117  else
118  {
119  phi = acos(r) / 3.0;
120  }
121 
122  // the eigenvalues satisfy eig3 >= eig2 >= eig1
123  l3 = q + 2.0 * p * cos(phi);
124  l1 = q + 2.0 * p * cos(phi + (2.0 * M_PI / 3.0));
125  // since trace(A) = eig1 + eig2 + eig3
126  l2 = 3.0 * q - l1 - l3;
127  }
128 }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References CellMLToNektar.cellml_metadata::p, and tinysimd::sqrt().

Referenced by Nektar::FieldUtils::ProcessL2Criterion::v_Process().

◆ operator!=()

bool Nektar::FieldUtils::operator!= ( const IsoVertex x,
const IsoVertex y 
)

Definition at line 961 of file ProcessIsoContour.cpp.

962 {
963  return ((x.m_x - y.m_x) * (x.m_x - y.m_x) +
964  (x.m_y - y.m_y) * (x.m_y - y.m_y) +
965  (x.m_z - y.m_z) * (x.m_z - y.m_z) <
967  ? 0
968  : 1;
969 }
static const NekDouble kNekZeroTol

◆ operator<<()

FIELD_UTILS_EXPORT std::ostream & Nektar::FieldUtils::operator<< ( std::ostream &  os,
const ModuleKey rhs 
)

Prints a given module key to a stream.

Definition at line 58 of file Module.cpp.

59 {
60  return os << ModuleTypeMap[rhs.first] << ": " << rhs.second;
61 }
const std::string ModuleTypeMap[]
Definition: Module.h:74

References ModuleTypeMap.

◆ operator==() [1/2]

bool Nektar::FieldUtils::operator== ( const IsoVertex x,
const IsoVertex y 
)

Definition at line 950 of file ProcessIsoContour.cpp.

951 {
952  return ((x.m_x - y.m_x) * (x.m_x - y.m_x) +
953  (x.m_y - y.m_y) * (x.m_y - y.m_y) +
954  (x.m_z - y.m_z) * (x.m_z - y.m_z) <
956  ? true
957  : false;
958 }

◆ operator==() [2/2]

bool Nektar::FieldUtils::operator== ( TriFaceIDs const &  p1,
TriFaceIDs const &  p2 
)

Definition at line 78 of file ProcessDisplacement.cpp.

79 {
80  std::vector<int> ids1(3), ids2(3);
81 
82  ids1[0] = p1.a;
83  ids1[1] = p1.b;
84  ids1[2] = p1.c;
85  ids2[0] = p2.a;
86  ids2[1] = p2.b;
87  ids2[2] = p2.c;
88 
89  std::sort(ids1.begin(), ids1.end());
90  std::sort(ids2.begin(), ids2.end());
91 
92  return ids1[0] == ids2[0] && ids1[1] == ids2[1] && ids1[2] == ids2[2];
93 }

References Nektar::FieldUtils::TriFaceIDs::a, Nektar::FieldUtils::TriFaceIDs::b, and Nektar::FieldUtils::TriFaceIDs::c.

◆ swap_endian() [1/2]

template<typename T >
void Nektar::FieldUtils::swap_endian ( std::vector< T > &  u)

Definition at line 118 of file Module.h.

119 {
120  size_t vecSize = u.size();
121  for (size_t i = 0; i < vecSize; ++i)
122  {
123  swap_endian(u[i]);
124  }
125 }
void swap_endian(std::vector< T > &u)
Definition: Module.h:118

References swap_endian().

◆ swap_endian() [2/2]

template<typename T >
void Nektar::FieldUtils::swap_endian ( T &  u)

Swap endian ordering of the input variable.

Definition at line 100 of file Module.h.

101 {
102  union
103  {
104  T u;
105  unsigned char u8[sizeof(T)];
106  } source, dest;
107 
108  source.u = u;
109 
110  for (size_t k = 0; k < sizeof(T); k++)
111  {
112  dest.u8[k] = source.u8[sizeof(T) - k - 1];
113  }
114 
115  u = dest.u;
116 }

Referenced by swap_endian(), Nektar::FieldUtils::InputNek5000::v_Process(), and Nektar::FieldUtils::InputSemtex::v_Process().

◆ ThreeSimilar()

void Nektar::FieldUtils::ThreeSimilar ( const int  i,
const int  j,
const int  k,
int &  pr,
int &  ps 
)

Definition at line 309 of file ProcessIsoContour.cpp.

312 {
313  switch (i + j + k)
314  {
315  case (3):
316  pr = 3;
317  ps = 4;
318  break;
319  case (4):
320  pr = 2;
321  ps = 4;
322  break;
323  case (5):
324  if (j == 1)
325  {
326  pr = 2;
327  ps = 3;
328  }
329  else
330  {
331  pr = 1;
332  ps = 4;
333  }
334  break;
335  case (6):
336  if (i == 0)
337  {
338  pr = 1;
339  ps = 3;
340  }
341  else
342  {
343  pr = 0;
344  ps = 4;
345  }
346  break;
347  case (7):
348  if (i == 0)
349  {
350  pr = 1;
351  ps = 2;
352  }
353  else
354  {
355  pr = 0;
356  ps = 3;
357  }
358  break;
359  case (8):
360  pr = 0;
361  ps = 2;
362  break;
363  case (9):
364  pr = 0;
365  ps = 1;
366  break;
367  default:
368  ASSERTL0(false, "Error in 5-point triangulation in ThreeSimilar");
369  break;
370  }
371 }

References ASSERTL0.

Referenced by Nektar::FieldUtils::ProcessIsoContour::ExtractContour().

◆ TwoPairs()

void Nektar::FieldUtils::TwoPairs ( Array< OneD, NekDouble > &  cx,
Array< OneD, NekDouble > &  cy,
Array< OneD, NekDouble > &  cz,
int &  pr 
)

Definition at line 282 of file ProcessIsoContour.cpp.

286 {
287  if (((cx[0] - cx[1]) == 0.0) && ((cy[0] - cy[1]) == 0.0) &&
288  ((cz[0] - cz[1]) == 0.0))
289  {
290  if (((cx[2] - cx[3]) == 0.0) && ((cy[2] - cy[3]) == 0.0) &&
291  ((cz[2] - cz[3]) == 0.0))
292  {
293  pr = 4;
294  }
295  else
296  {
297  pr = 3;
298  }
299  }
300  else
301  {
302  pr = 1;
303  }
304 }

Referenced by Nektar::FieldUtils::ProcessIsoContour::ExtractContour().

◆ WriteStream() [1/4]

template<typename T >
void Nektar::FieldUtils::WriteStream ( std::ostream &  outfile,
Array< OneD, T >  data 
)

Specialisation of WriteStream to support writing Nektar::Array datatype.

Definition at line 129 of file OutputTecplot.cpp.

130 {
131  if (data.size())
132  {
133  outfile.write(reinterpret_cast<char *>(&data[0]),
134  data.size() * sizeof(T));
135  }
136 }

◆ WriteStream() [2/4]

template<>
void Nektar::FieldUtils::WriteStream ( std::ostream &  outfile,
std::string  data 
)

Specialisation of WriteStream to support writing std::string.

Tecplot binary formats represent all strings by writing out their characters as 32-bit integers, followed by a 32-bit integer null (0) character to denote the end of the string.

Definition at line 110 of file OutputTecplot.cpp.

111 {
112  // Convert string to array of int32_t
113  for (std::string::size_type i = 0; i < data.size(); ++i)
114  {
115  char strChar = data[i];
116  NekInt32 strCharInt = strChar;
117  WriteStream(outfile, strCharInt);
118  }
119 
120  // Now dump out zero character to terminate
121  WriteStream(outfile, 0);
122 }
void WriteStream(std::ostream &outfile, std::vector< T > data)
Specialisation of WriteStream to support writing std::vector datatype.
std::int32_t NekInt32

References WriteStream().

◆ WriteStream() [3/4]

template<typename T >
void Nektar::FieldUtils::WriteStream ( std::ostream &  outfile,
std::vector< T >  data 
)

Specialisation of WriteStream to support writing std::vector datatype.

Definition at line 142 of file OutputTecplot.cpp.

143 {
144  if (data.size())
145  {
146  outfile.write(reinterpret_cast<char *>(&data[0]),
147  data.size() * sizeof(T));
148  }
149 }

◆ WriteStream() [4/4]

template<typename T >
void Nektar::FieldUtils::WriteStream ( std::ostream &  outfile,
data 
)

Helper function to write binary data to stream.

Definition at line 97 of file OutputTecplot.cpp.

98 {
99  T tmp = data;
100  outfile.write(reinterpret_cast<char *>(&tmp), sizeof(T));
101 }

Referenced by Nektar::FieldUtils::OutputTecplotBinary::v_WriteTecplotConnectivity(), Nektar::FieldUtils::OutputTecplotBinary::v_WriteTecplotHeader(), Nektar::FieldUtils::OutputTecplotBinary::v_WriteTecplotZone(), Nektar::FieldUtils::OutputTecplotBinary::WriteDoubleOrFloat(), and WriteStream().

Variable Documentation

◆ ModulePriorityMap

const char* const Nektar::FieldUtils::ModulePriorityMap[]
Initial value:
= {
"CreateGraph", "CreateFieldData", "ModifyFieldData", "CreateExp",
"FillExp", "ModifyExp", "BndExtraction", "CreatePts",
"ConvertExpToPts", "ModifyPts", "Output"}

Definition at line 92 of file Module.h.

◆ ModuleTypeMap

const std::string Nektar::FieldUtils::ModuleTypeMap[] = {"Input", "Process", "Output"}

Definition at line 74 of file Module.h.

Referenced by export_Module(), Module_Create(), Module_Register(), and operator<<().

◆ TecplotZoneTypeMap

std::string Nektar::FieldUtils::TecplotZoneTypeMap[]
Initial value:
= {"ORDERED", "LINESEG", "TRIANGLE",
"QUADRILATERAL", "TETRAHEDRON", "BRICK",
"POLYGON", "POLYHEDRON"}

Definition at line 52 of file OutputTecplot.cpp.

Referenced by Nektar::FieldUtils::OutputTecplot::v_WriteTecplotZone().