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

Classes

struct  ConfigOption
 Represents a command-line configuration option. More...
 
struct  Field
 
class  FieldConvertComm
 
class  InputDat
 Input module for Xml files. More...
 
class  InputFld
 
class  InputModule
 Abstract base class for input modules. More...
 
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...
 
class  Iso
 
class  IsoVertex
 
class  Module
 
class  Octree
 
class  OutputFileBase
 Converter from fld to vtk. More...
 
class  OutputFld
 Output to fld format. More...
 
class  OutputInfo
 
class  OutputModule
 Abstract base class for output modules. More...
 
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  ProcessAverageFld
 This processing module averages the input fld files. 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
 
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  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  ProcessModule
 Abstract base class for processing modules. 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  ProcessPowerSpectrum
 This processing module outputs power spectrum at given regions for a 3DH1D fields. 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...
 
class  ProcessZeroHomogeneousPlane
 This processing module replaces all expansions by a single plane from 3DH1D fields, defined by the parameter planeid. More...
 
struct  TriFaceHash
 
struct  TriFaceIDs
 

Typedefs

typedef std::shared_ptr< FieldFieldSharedPtr
 
typedef std::shared_ptr< Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > InterpolatorSharedPtr
 
typedef std::pair< ModuleType, std::string > ModuleKey
 
typedef std::shared_ptr< InputModuleInputModuleSharedPtr
 
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 990 of file Field.hpp.

◆ InputModuleSharedPtr

Definition at line 294 of file Module.h.

◆ InterpolatorSharedPtr

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

◆ IsoSharedPtr

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

Definition at line 184 of file ProcessIsoContour.h.

◆ ModuleFactory

Definition at line 331 of file Module.h.

◆ ModuleKey

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

Definition at line 180 of file Module.h.

◆ ModuleSharedPtr

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

Definition at line 329 of file Module.h.

◆ TriFaceMap

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

Definition at line 93 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 74 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 64 of file Module.h.

◆ TecplotZoneType

Enumerator
eOrdered 
eFELineSeg 
eFETriangle 
eFEQuadrilateral 
eFETetrahedron 
eFEBrick 
eFEPolygon 
eFEPolyhedron 

Definition at line 44 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 47 of file Module.cpp.

48{
49 static ModuleFactory instance;
50 return instance;
51}
Provides a generic Factory class.
Definition: NekFactory.hpp:104

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 111 of file ProcessQualityMetric.cpp.

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

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 68 of file ProcessL2Criterion.cpp.

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

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

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

◆ operator!=()

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

Definition at line 963 of file ProcessIsoContour.cpp.

964{
965 return ((x.m_x - y.m_x) * (x.m_x - y.m_x) +
966 (x.m_y - y.m_y) * (x.m_y - y.m_y) +
967 (x.m_z - y.m_z) * (x.m_z - y.m_z) <
969 ? 0
970 : 1;
971}
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 56 of file Module.cpp.

57{
58 return os << ModuleTypeMap[rhs.first] << ": " << rhs.second;
59}
const std::string ModuleTypeMap[]
Definition: Module.h:72

References ModuleTypeMap.

◆ operator==() [1/2]

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

Definition at line 952 of file ProcessIsoContour.cpp.

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

◆ operator==() [2/2]

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

Definition at line 76 of file ProcessDisplacement.cpp.

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

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 116 of file Module.h.

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

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 98 of file Module.h.

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

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 305 of file ProcessIsoContour.cpp.

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

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 278 of file ProcessIsoContour.cpp.

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

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 125 of file OutputTecplot.cpp.

126{
127 if (data.size())
128 {
129 outfile.write(reinterpret_cast<char *>(&data[0]),
130 data.size() * sizeof(T));
131 }
132}

◆ 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 106 of file OutputTecplot.cpp.

107{
108 // Convert string to array of int32_t
109 for (std::string::size_type i = 0; i < data.size(); ++i)
110 {
111 char strChar = data[i];
112 NekInt32 strCharInt = strChar;
113 WriteStream(outfile, strCharInt);
114 }
115
116 // Now dump out zero character to terminate
117 WriteStream(outfile, 0);
118}
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 138 of file OutputTecplot.cpp.

139{
140 if (data.size())
141 {
142 outfile.write(reinterpret_cast<char *>(&data[0]),
143 data.size() * sizeof(T));
144 }
145}

◆ 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 93 of file OutputTecplot.cpp.

94{
95 T tmp = data;
96 outfile.write(reinterpret_cast<char *>(&tmp), sizeof(T));
97}

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 90 of file Module.h.

◆ ModuleTypeMap

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

Definition at line 72 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 48 of file OutputTecplot.cpp.

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