MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
VolumeData< nDim, SolverType > Class Template Reference

Class to handle sampling of volume data. More...

#include <samplingdata.h>

Inheritance diagram for VolumeData< nDim, SolverType >:
[legend]
Collaboration diagram for VolumeData< nDim, SolverType >:
[legend]

Public Types

using Base = SamplingData< nDim, SolverType >
 

Public Member Functions

 VolumeData (SolverType &solver_)
 
virtual ~VolumeData ()
 
MString getBaseName () const override
 Return base name of derived class, e.g. used for naming properties and output file names. More...
 
void indexXD (const MInt index, const MInt *const nPoints, MInt *indexXD)
 Convert 1D index into 2D/3D index. More...
 
MBool isMpiRoot () const
 
MPI_Comm mpiComm () const
 
SolverTypesolver ()
 Return reference to solver. More...
 
const SolverTypesolver () const
 
MInt solverId () const
 
- Public Member Functions inherited from SamplingData< nDim, SolverType >
 SamplingData (SolverType &solver)
 
void init ()
 
void save (MBool finalTimestep)
 
void setInputOutputProperties ()
 
MBool enabled () const
 
virtual MInt loadInputFile (const MString NotUsed(inputFileName), const MInt NotUsed(fileNo), std::vector< MFloat > &NotUsed(coordinates))
 
virtual MString getBaseName () const
 Return base name of derived class, e.g. used for naming properties and output file names. More...
 

Public Attributes

MBool m_enabled
 
std::vector< SamplingDataSeriesm_timeSeries
 Holds all properties and buffers for each input file. More...
 

Protected Member Functions

MBool generatePoints () const override
 Always generate points. More...
 
MBool hasInputDataFile (const MInt fileNo) const override
 Points are generated not read from file, just check for property 'volumeDataShape_*'. More...
 
MString getInputDataFile (const MInt NotUsed(fileNo)) const
 
void readAdditionalProperties (const MInt fileNo) override
 Read additional properties for volume sampling. More...
 
MInt noSamplingPoints (const MInt fileNo) override
 Return number of volume sampling points. More...
 
void getSamplingPoint (const MInt fileNo, const MInt pointId, MFloat *const coordinates) override
 Return sampling point coordinates. More...
 
- Protected Member Functions inherited from SamplingData< nDim, SolverType >
void init (SamplingDataSeries &timeSeries)
 
void save (SamplingDataSeries &timeSeries, MBool finalTimeStep)
 Saves sampling data Saves all states of conservative variables at specified points. Also writes sampling data files in specified intervals. More...
 
virtual MBool hasInputDataFile (const MInt fileNo) const
 Check if an input file exists, might be overwritten in derived class if points are generated. More...
 
virtual MBool generatePoints () const
 Default behaviour for point generation which can be overwritten in derived class. More...
 
virtual MInt noSamplingPoints (const MInt NotUsed(fileNo))
 Return the number of sampling points which are generated by the derived class (no input file) More...
 
virtual void getSamplingPoint (const MInt NotUsed(fileNo), const MInt NotUsed(pointId), MFloat *const NotUsed(coordinates))
 Return the coordinates of the sampling point with the given id when generating points. More...
 
void indexXD (const MInt index, const MInt *const nPoints, MInt *indexXD)
 Convert 1D index into 2D/3D index. More...
 
void saveSamplingPointCoordinates (SamplingDataSeries &timeSeries)
 Save point coordinates of time series. More...
 
virtual void readAdditionalProperties (const MInt NotUsed(fileNo))
 Read additional properties which are required by the derived class. More...
 
MInt solverId () const
 
MBool isMpiRoot () const
 
MPI_Comm mpiComm () const
 
MInt noVars () const
 Return total number of sampling variables. More...
 
MString getFileBaseName () const
 Return base name of files for used sampling feature. More...
 
SolverTypesolver ()
 Return reference to solver. More...
 
const SolverTypesolver () const
 

Private Attributes

std::vector< MStringm_volumeShape {}
 
std::vector< std::vector< MFloat > > m_volumeParameterFloat {}
 
std::vector< std::vector< MInt > > m_volumeParameterInt {}
 

Additional Inherited Members

- Protected Attributes inherited from SamplingData< nDim, SolverType >
MBool m_enabled = false
 
MBool m_interpolatePointData
 
std::vector< SamplingDataSeriesm_timeSeries {}
 Holds all properties and buffers for each input file. More...
 
std::vector< MIntm_solverSamplingVarIds {}
 List of sampling variables. More...
 
std::vector< MIntm_noSolverSamplingVars {}
 Number of variables for each sampling variable. More...
 
std::vector< std::vector< MString > > m_solverSamplingVarNames {}
 List of variable names for each sampling variable. More...
 

Detailed Description

template<MInt nDim, class SolverType>
class VolumeData< nDim, SolverType >

Definition at line 220 of file samplingdata.h.

Member Typedef Documentation

◆ Base

template<MInt nDim, class SolverType >
using VolumeData< nDim, SolverType >::Base = SamplingData<nDim, SolverType>

Definition at line 222 of file samplingdata.h.

Constructor & Destructor Documentation

◆ VolumeData()

template<MInt nDim, class SolverType >
VolumeData< nDim, SolverType >::VolumeData ( SolverType solver_)
inline

Definition at line 233 of file samplingdata.h.

This base class is responsible for the sampling data feature that provides the required general metho...
Definition: samplingdata.h:35

◆ ~VolumeData()

template<MInt nDim, class SolverType >
virtual VolumeData< nDim, SolverType >::~VolumeData ( )
inlinevirtual

Definition at line 234 of file samplingdata.h.

234{};

Member Function Documentation

◆ generatePoints()

template<MInt nDim, class SolverType >
MBool VolumeData< nDim, SolverType >::generatePoints ( ) const
inlineoverrideprotectedvirtual

Reimplemented from SamplingData< nDim, SolverType >.

Definition at line 240 of file samplingdata.h.

240{ return true; }

◆ getBaseName()

template<MInt nDim, class SolverType >
MString VolumeData< nDim, SolverType >::getBaseName ( ) const
inlineoverridevirtual

Reimplemented from SamplingData< nDim, SolverType >.

Definition at line 236 of file samplingdata.h.

236{ return "volume"; }

◆ getInputDataFile()

template<MInt nDim, class SolverType >
MString VolumeData< nDim, SolverType >::getInputDataFile ( const MInt   NotUsedfileNo) const
inlineprotected

Definition at line 253 of file samplingdata.h.

253 {
254 return ""; // No input file, points are generated
255 }

◆ getSamplingPoint()

template<MInt nDim, class SolverType >
void VolumeData< nDim, SolverType >::getSamplingPoint ( const MInt  fileNo,
const MInt  pointId,
MFloat *const  coordinates 
)
overrideprotected

Definition at line 1390 of file samplingdata.h.

1390 {
1391 TRACE();
1392
1393 if(pointId < 0 || pointId > noSamplingPoints(fileNo)) {
1394 TERMM(1, "Invalid point id");
1395 }
1396
1397 // Generate coordinates depending on volume shape
1398 std::vector<MInt> indexNd(3);
1399 if(m_volumeShape[fileNo] == "box") {
1400 // Points for box are created with [x,y,z] = [i,j,k] ordering, i.e., k varies the fastest
1401 indexXD(pointId, &m_volumeParameterInt[fileNo][0], &indexNd[0]);
1402
1403 for(MInt i = 0; i < nDim; i++) {
1404 const MFloat minCoord = m_volumeParameterFloat[fileNo][i];
1405 const MFloat maxCoord = m_volumeParameterFloat[fileNo][nDim + i];
1406 const MFloat delta = (maxCoord - minCoord) / static_cast<MFloat>(m_volumeParameterInt[fileNo][i] - 1);
1407 coordinates[i] = minCoord + indexNd[i] * delta;
1408 }
1409 } else if(m_volumeShape[fileNo] == "cylinder") {
1410 /* TERMM(1, "TODO: test volume data sampling with cylinder shape generation"); */
1411 // ordering: z, r, phi
1412 // Note: center point (r=0) only appears once
1413 const MInt nPointsPerCirc = m_volumeParameterInt[fileNo][1] * m_volumeParameterInt[fileNo][2] + 1;
1414 indexNd[0] = floor(pointId / nPointsPerCirc); // z
1415
1416 const MInt tmpIndex = pointId - indexNd[0] * nPointsPerCirc;
1417 indexNd[1] = ceil(static_cast<MFloat>(tmpIndex) / static_cast<MFloat>(m_volumeParameterInt[fileNo][2]));
1418 indexNd[2] = std::max((tmpIndex - 1) % m_volumeParameterInt[fileNo][2], 0);
1419
1420 const MInt dir = m_volumeParameterInt[fileNo][3];
1421 const MInt dir2 = (dir == 0) ? 1 : 0;
1422 const MInt dir3 = 3 - dir - dir2;
1423
1424 const MFloat* startPos = &m_volumeParameterFloat[fileNo][0];
1425 const MFloat nz = static_cast<MFloat>(m_volumeParameterInt[fileNo][0]);
1426 const MFloat dist_z = static_cast<MFloat>(indexNd[0]) * (startPos[3] - startPos[dir]) / (nz - 1);
1427
1428 const MFloat r = m_volumeParameterFloat[fileNo][4];
1429 const MFloat nr = m_volumeParameterInt[fileNo][1];
1430 const MFloat dist_r = static_cast<MFloat>(indexNd[1]) * r / nr;
1431
1432 const MFloat nphi = m_volumeParameterInt[fileNo][2];
1433 const MFloat phi = static_cast<MFloat>(indexNd[2]) * 2.0 * PI / nphi;
1434
1435 coordinates[dir] = startPos[dir] + dist_z;
1436 coordinates[dir2] = startPos[dir2] + dist_r * sin(phi);
1437 coordinates[dir3] = startPos[dir3] + dist_r * cos(phi);
1438 } else {
1439 TERMM(1, "Unknown volume data shape: " + m_volumeShape[fileNo]);
1440 }
1441}
std::vector< std::vector< MInt > > m_volumeParameterInt
Definition: samplingdata.h:267
std::vector< std::vector< MFloat > > m_volumeParameterFloat
Definition: samplingdata.h:266
MInt noSamplingPoints(const MInt fileNo) override
Return number of volume sampling points.
std::vector< MString > m_volumeShape
Definition: samplingdata.h:265
void indexXD(const MInt index, const MInt *const nPoints, MInt *indexXD)
Convert 1D index into 2D/3D index.
int32_t MInt
Definition: maiatypes.h:62
double MFloat
Definition: maiatypes.h:52
T cos(const T a, const T b, const T x)
Cosine slope filter.
Definition: filter.h:125

◆ hasInputDataFile()

template<MInt nDim, class SolverType >
MBool VolumeData< nDim, SolverType >::hasInputDataFile ( const MInt  fileNo) const
inlineoverrideprotectedvirtual

Reimplemented from SamplingData< nDim, SolverType >.

Definition at line 243 of file samplingdata.h.

243 {
244 const MString propName = getBaseName() + "DataShape_" + std::to_string(fileNo);
245 if(Context::propertyExists(propName, solverId())) {
246 const MString shape = Context::getSolverProperty<MString>(propName, solverId(), AT_);
247 return shape != "";
248 } else {
249 return false;
250 }
251 }
static MBool propertyExists(const MString &name, MInt solver=m_noSolvers)
This function checks if a property exists in general.
Definition: context.cpp:494
MString getBaseName() const override
Return base name of derived class, e.g. used for naming properties and output file names.
Definition: samplingdata.h:236
MInt solverId() const
Definition: samplingdata.h:118
std::basic_string< char > MString
Definition: maiatypes.h:55

◆ indexXD()

template<MInt nDim, class SolverType >
void SamplingData< nDim, SolverType >::indexXD ( const MInt  index,
const MInt *const  nPoints,
MInt indexXD 
)

Convert a one dimensional index into a 2D/3D index with the given number of points per direction (index of last dimension varies the fastest)

Definition at line 109 of file samplingdata.h.

1023 {
1024 const MInt noPoints2D = nPoints[0] * nPoints[1];
1025 const MInt maxNoPoints = (nDim == 2) ? noPoints2D : noPoints2D * nPoints[2];
1026 if(index < 0 || index >= maxNoPoints) {
1027 TERMM(1, "invalid index: " + std::to_string(index));
1028 }
1029
1030 // Last index (varies fastest)
1031 indexXD[nDim - 1] = index % nPoints[nDim - 1];
1032 // First/second index in 2D/3D
1033 indexXD[nDim - 2] = ((index - indexXD[nDim - 1]) / nPoints[nDim - 1]) % nPoints[nDim - 2];
1034 IF_CONSTEXPR(nDim == 3) {
1035 // First index in 3D
1036 indexXD[0] = (index - indexXD[1] * nPoints[1] - indexXD[2]) / (nPoints[2] * nPoints[1]);
1037 }
1038}
IdType index(const FloatType *const x, const IdType level)
Return Hilbert index for given location and level in 2D or 3D.
Definition: hilbert.h:165

◆ isMpiRoot()

template<MInt nDim, class SolverType >
MBool SamplingData< nDim, SolverType >::isMpiRoot ( ) const
inline

Definition at line 119 of file samplingdata.h.

119{ return m_solver.domainId() == 0; }
SolverType & m_solver
Definition: samplingdata.h:150

◆ mpiComm()

template<MInt nDim, class SolverType >
MPI_Comm SamplingData< nDim, SolverType >::mpiComm ( ) const
inline

Definition at line 120 of file samplingdata.h.

120{ return m_solver.mpiComm(); }

◆ noSamplingPoints()

template<MInt nDim, class SolverType >
MInt VolumeData< nDim, SolverType >::noSamplingPoints ( const MInt  fileNo)
overrideprotected

Definition at line 1367 of file samplingdata.h.

1367 {
1368 TRACE();
1369
1370 MInt noPoints = 0;
1371 if(m_volumeShape[fileNo] == "box") {
1372 noPoints = 1;
1373 for(MInt i = 0; i < nDim; i++) {
1374 noPoints *= m_volumeParameterInt[fileNo][i];
1375 }
1376 } else if(m_volumeShape[fileNo] == "cylinder") {
1377 // nZ * (nPhi * nR + 1) // +1 for center point
1378 noPoints =
1379 m_volumeParameterInt[fileNo][0] * (m_volumeParameterInt[fileNo][1] * m_volumeParameterInt[fileNo][2] + 1);
1380 } else {
1381 TERMM(1, "Unknown volume data shape: " + m_volumeShape[fileNo]);
1382 }
1383
1384 return noPoints;
1385}

◆ readAdditionalProperties()

template<MInt nDim, class SolverType >
void VolumeData< nDim, SolverType >::readAdditionalProperties ( const MInt  fileNo)
overrideprotected

Definition at line 1309 of file samplingdata.h.

1309 {
1310 TRACE();
1311
1312 const MString volumeShape =
1313 Context::getSolverProperty<MString>(getBaseName() + "DataShape_" + std::to_string(fileNo), solverId(), AT_);
1314
1315 MInt noInputValuesFloat = -1;
1316 MInt noInputValuesInt = -1;
1317
1318 if(volumeShape == "box") {
1319 // Box defined by two points and number of points in each spatial direction
1320 noInputValuesFloat = 2 * nDim;
1321 noInputValuesInt = nDim;
1322 } else if(volumeShape == "cylinder") {
1323 IF_CONSTEXPR(nDim == 2) { TERMM(1, "Cylinder shape not defined in 2D."); }
1324 // Cylinder defined by first center point, second center axial position and radius
1325 noInputValuesFloat = nDim + 2;
1326 // TODO labels:PP allow to generate non-equidistant points in radial direction to avoid that lots of
1327 // points are clustered in the center of the cylinder
1328 // Number of points in each direction: axial, radial, azimuthal; Cylinder axial direction
1329 // Note: given number of points in radial direction does not include center point
1330 noInputValuesInt = nDim + 1;
1331 } else {
1332 TERMM(1, "Unknown volume data shape: " + volumeShape);
1333 }
1334 m_volumeShape.push_back(volumeShape);
1335
1336 if(noInputValuesFloat < 1 || noInputValuesInt < 1) {
1337 TERMM(1, "Invalid number of input values for volume data shape " + volumeShape + ": "
1338 + std::to_string(noInputValuesFloat) + " " + std::to_string(noInputValuesInt));
1339 }
1340
1341 m_log << "Volume data #" << fileNo << ": " << volumeShape << " with parameters";
1342
1343 // Read float parameters
1344 std::vector<MFloat> paramFloat(noInputValuesFloat);
1345 for(MInt i = 0; i < noInputValuesFloat; i++) {
1346 paramFloat[i] = Context::getSolverProperty<MFloat>(getBaseName() + "DataParameterFloat_" + std::to_string(fileNo),
1347 solverId(), AT_, i);
1348 m_log << " " << paramFloat[i];
1349 }
1350 m_volumeParameterFloat.push_back(paramFloat);
1351
1352 // Read int parameters
1353 std::vector<MInt> paramInt(noInputValuesInt);
1354 for(MInt i = 0; i < noInputValuesInt; i++) {
1355 paramInt[i] = Context::getSolverProperty<MInt>(getBaseName() + "DataParameterInt_" + std::to_string(fileNo),
1356 solverId(), AT_, i);
1357 m_log << " " << paramInt[i];
1358 }
1359 m_volumeParameterInt.push_back(paramInt);
1360
1361 m_log << std::endl;
1362}
InfoOutFile m_log

◆ solver() [1/2]

template<MInt nDim, class SolverType >
SolverType & SamplingData< nDim, SolverType >::solver ( )
inline

Definition at line 129 of file samplingdata.h.

129{ return m_solver; }

◆ solver() [2/2]

template<MInt nDim, class SolverType >
const SolverType & SamplingData< nDim, SolverType >::solver ( ) const
inline

Definition at line 130 of file samplingdata.h.

130{ return m_solver; }

◆ solverId()

template<MInt nDim, class SolverType >
MInt SamplingData< nDim, SolverType >::solverId ( ) const
inline

Definition at line 118 of file samplingdata.h.

118{ return m_solver.solverId(); }

Member Data Documentation

◆ m_enabled

template<MInt nDim, class SolverType >
MBool SamplingData< nDim, SolverType >::m_enabled

Definition at line 134 of file samplingdata.h.

◆ m_timeSeries

template<MInt nDim, class SolverType >
std::vector<SamplingDataSeries> SamplingData< nDim, SolverType >::m_timeSeries

Definition at line 138 of file samplingdata.h.

◆ m_volumeParameterFloat

template<MInt nDim, class SolverType >
std::vector<std::vector<MFloat> > VolumeData< nDim, SolverType >::m_volumeParameterFloat {}
private

Definition at line 266 of file samplingdata.h.

◆ m_volumeParameterInt

template<MInt nDim, class SolverType >
std::vector<std::vector<MInt> > VolumeData< nDim, SolverType >::m_volumeParameterInt {}
private

Definition at line 267 of file samplingdata.h.

◆ m_volumeShape

template<MInt nDim, class SolverType >
std::vector<MString> VolumeData< nDim, SolverType >::m_volumeShape {}
private

Storage for volume parameters (shape; floating point values, e.g. coordinates; integer parameters, i.e. number of points)

Definition at line 265 of file samplingdata.h.


The documentation for this class was generated from the following file: