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

Surface data sampling class. Records all sampling variables on all surface elements and outputs additional geometric information for further postprocessing the data, e.g., using the Ffowcs Williams-Hawkings method for acoustic far-field predictions. More...

#include <samplingdata.h>

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

Public Types

using Base = SamplingData< nDim, SolverType >
 

Public Member Functions

 SurfaceData (SolverType &solver_)
 
virtual ~SurfaceData ()
 
void saveGeomInfo (const MString inputFileName, const MInt fileNo, MFloat *coordinates)
 Store geometric information for all elements of a surface. More...
 
MInt loadInputFile (const MString inputFileName, const MInt fileNo, std::vector< MFloat > &coordinates) override
 Load an input file for the surface data sampling feature. More...
 
MString getBaseName () const override
 Return base name of derived class, e.g. used for naming properties and output file names. 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...
 

Private Attributes

GeometryXD< nDim >::type * m_geometry = nullptr
 

Additional Inherited Members

- 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
 
- 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 SurfaceData< nDim, SolverType >

Definition at line 190 of file samplingdata.h.

Member Typedef Documentation

◆ Base

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

Definition at line 192 of file samplingdata.h.

Constructor & Destructor Documentation

◆ SurfaceData()

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

Definition at line 202 of file samplingdata.h.

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

◆ ~SurfaceData()

template<MInt nDim, class SolverType >
SurfaceData< nDim, SolverType >::~SurfaceData
virtual

Definition at line 1091 of file samplingdata.h.

1091 {
1092 // Free geometry object if existing
1093 if(m_geometry != nullptr) {
1094 // TODO labels:PP,IO fix deallocation
1095 // delete m_geometry;
1096 }
1097}
GeometryXD< nDim >::type * m_geometry
Definition: samplingdata.h:214

Member Function Documentation

◆ getBaseName()

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

Reimplemented from SamplingData< nDim, SolverType >.

Definition at line 210 of file samplingdata.h.

210{ return "surface"; }

◆ 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

◆ loadInputFile()

template<MInt nDim, class SolverType >
MInt SurfaceData< nDim, SolverType >::loadInputFile ( const MString  inputFileName,
const MInt  fileNo,
std::vector< MFloat > &  coordinates 
)
override

The surface is created as an geometry object. For all surface elements that make up a surface in 2D or 3D the centroid is computed in which the sampling variables are determined. Additional geometric information on the surface elements is stored in a separated data file.

Definition at line 1106 of file samplingdata.h.

1108 {
1109 TRACE();
1110 if(!isMpiRoot()) {
1111 TERMM(1, "SurfaceData::loadInputFile() should only be called by mpi root process!");
1112 }
1113
1114 // Create geometry object
1115 auto geometry = new typename GeometryXD<nDim>::type(0, inputFileName, MPI_COMM_SELF);
1116 // Store pointer in class member variable
1117 m_geometry = geometry;
1118
1119 // Determine number of geometry/surface elements
1120 MInt noElements = geometry->GetNoElements();
1121 coordinates.resize(nDim * noElements);
1122
1123 std::array<MFloat, nDim * nDim> vertex;
1124
1125 // Compute element centroids and store in coordinates
1126 for(MInt i = 0; i < noElements; i++) {
1127 geometry->elements[i].getVertices(&vertex[0]);
1128 geometry->elements[i].calcCentroid(&vertex[0], &coordinates[i * nDim]);
1129 }
1130
1131 // Write specific data for all elements into an output file
1132 saveGeomInfo(inputFileName, fileNo, &coordinates[0]);
1133
1134 return noElements;
1135}
MBool isMpiRoot() const
Definition: samplingdata.h:119
void saveGeomInfo(const MString inputFileName, const MInt fileNo, MFloat *coordinates)
Store geometric information for all elements of a surface.
int32_t MInt
Definition: maiatypes.h:62

◆ 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(); }

◆ saveGeomInfo()

template<MInt nDim, class SysEqn >
void SurfaceData< nDim, SysEqn >::saveGeomInfo ( const MString  inputFileName,
const MInt  fileNo,
MFloat coordinates 
)

This function calculates and saves specific information of 2D and 3D elements that make up a surface: vertices, length/surface area, centroids and surface normal vectors.

The center of the whole geometry is calculated for the determination of the outer surface normal vectors. There is also the option to manually set a geometric center in the property files, e.g., for the case of a line/plane to set a point inside the whole geometry given by multiple surfaces since the geometric center is on the line/plane itself and the direction of the normal vector cannot be determined. Note: this will only work for convex geometries, to make this work for general geometries the user could specify a set of points inside the geometry from which the nearest one is used for a surface element to determine the outer normal vector.

Definition at line 1152 of file samplingdata.h.

1152 {
1153 TRACE();
1154 using namespace maia::parallel_io;
1155
1156 const MInt noElements = m_geometry->GetNoElements();
1157 // Define output file
1158 std::ostringstream fileName;
1159 fileName << solver().outputDir() << "surface_geometryInfo_" << fileNo << ParallelIo::fileExt();
1160
1161 ParallelIo file(fileName.str(), PIO_REPLACE, MPI_COMM_SELF);
1162 file.setAttribute(inputFileName, "inputFileName");
1163
1164 // Get vertices of an element and calculate length/surface of element
1165 std::vector<MFloat> vertices(noElements * nDim * nDim);
1166 std::vector<MFloat> geomMeasures(noElements);
1167 std::vector<MFloat> normalVecs(noElements * nDim);
1168
1169 MFloat geomCenter[3]{};
1170
1171 // Check if geometry center is defined by user
1172 if(Context::propertyExists("surfaceDataGeometryCenter_" + std::to_string(fileNo), solverId())) {
1173 for(MInt i = 0; i < nDim; i++) {
1174 geomCenter[i] =
1175 Context::getSolverProperty<MFloat>("surfaceDataGeometryCenter_" + std::to_string(fileNo), solverId(), AT_, i);
1176 }
1177 } else {
1178 // Determine geometry center as average of all element centroids
1179 for(MInt i = 0; i < noElements; i++) {
1180 for(MInt dim = 0; dim < nDim; dim++) {
1181 geomCenter[dim] += coordinates[i * nDim + dim];
1182 }
1183 }
1184 for(MInt dim = 0; dim < nDim; dim++) {
1185 geomCenter[dim] /= static_cast<MFloat>(noElements);
1186 }
1187 }
1188
1189 m_log << "Surface geometry info: geomCenter =";
1190 for(MInt dim = 0; dim < nDim; dim++) {
1191 m_log << " " << geomCenter[dim];
1192 }
1193 m_log << std::endl;
1194
1195 MInt counter = 0;
1196 for(MInt i = 0; i < noElements; i++) {
1197 for(MInt j = 0; j < nDim; j++) {
1198 for(MInt k = 0; k < nDim; k++) {
1199 vertices[counter] = m_geometry->elements[i].m_vertices[j][k];
1200 counter++;
1201 }
1202 }
1203
1204 auto& element = m_geometry->elements[i];
1205 // Calc surface area/length of every element
1206 IF_CONSTEXPR(nDim == 3) {
1207 // Define two vectors for triangular surface in 3D
1208 const MFloat v12_x = element.m_vertices[0][0] - element.m_vertices[1][0];
1209 const MFloat v12_y = element.m_vertices[0][1] - element.m_vertices[1][1];
1210 const MFloat v12_z = element.m_vertices[0][2] - element.m_vertices[1][2];
1211 const MFloat v13_x = element.m_vertices[0][0] - element.m_vertices[2][0];
1212 const MFloat v13_y = element.m_vertices[0][1] - element.m_vertices[2][1];
1213 const MFloat v13_z = element.m_vertices[0][2] - element.m_vertices[2][2];
1214
1215 // Cross product
1216 const MFloat v_x = v12_y * v13_z - v12_z * v13_y;
1217 const MFloat v_y = v12_z * v13_x - v12_x * v13_z;
1218 const MFloat v_z = v12_x * v13_y - v12_y * v13_x;
1219 // Surface area is 1/2 the length of the cross product
1220 geomMeasures[i] = 0.5 * sqrt(pow(v_x, 2) + pow(v_y, 2) + pow(v_z, 2));
1221 }
1222 else {
1223 const MFloat p1_x = element.m_vertices[0][0];
1224 const MFloat p2_x = element.m_vertices[1][0];
1225 const MFloat p1_y = element.m_vertices[0][1];
1226 const MFloat p2_y = element.m_vertices[1][1];
1227 // Element length
1228 geomMeasures[i] = sqrt(pow(p1_x - p2_x, 2) + pow(p1_y - p2_y, 2));
1229 }
1230
1231 std::array<MFloat, nDim * nDim> vertex;
1232 m_geometry->elements[i].getVertices(&vertex[0]);
1233
1234 // Calculate surface normal vector
1235 m_geometry->elements[i].calcNormal(&vertex[0], &normalVecs[i * nDim]);
1236
1237 // Compute dot product of surface normal vector and the vector from the geometric center to the
1238 // surface centroid
1239 MFloat dotProduct = 0.0;
1240 for(MInt k = 0; k < nDim; k++) {
1241 const MFloat vecCenter = coordinates[i * nDim + k] - geomCenter[k];
1242 dotProduct += normalVecs[i * nDim + k] * vecCenter;
1243 }
1244 // If the dot product is positive the two vectors are pointing in the same 'direction', i.e.,
1245 // the angle between the vectors is less than 90deg; thus if not >0 change the sign of the
1246 // normal vector to point outside the geometry
1247 if(!(dotProduct > 0)) {
1248 for(MInt k = 0; k < nDim; k++) {
1249 normalVecs[i * nDim + k] = -normalVecs[i * nDim + k];
1250 }
1251 }
1252 }
1253
1254 const MFloat totalGeomMeasure = std::accumulate(geomMeasures.begin(), geomMeasures.end(), 0.0);
1255 const MString measure = (nDim == 3) ? "area" : "length";
1256 m_log << "Total segment " << measure << ": " << totalGeomMeasure << std::endl;
1257
1258 // TODO labels:PP for DEBUG
1259 std::stringstream dumpFileName;
1260 dumpFileName << solver().outputDir() << "surface_geometryInfo_" << fileNo << ".dump";
1261 // DEBUG output for plotting
1262 FILE* datei;
1263 datei = fopen(dumpFileName.str().c_str(), "w");
1264 /* datei = fopen("geometryInfo.dump", "a+"); */
1265 for(MInt i = 0; i < noElements; i++) {
1266 fprintf(datei, "%d", i);
1267
1268 for(MInt j = 0; j < nDim; j++) {
1269 fprintf(datei, " %.8f", coordinates[i * nDim + j]);
1270 }
1271
1272 for(MInt j = 0; j < nDim; j++) {
1273 fprintf(datei, " %.8f", normalVecs[i * nDim + j]);
1274 }
1275 fprintf(datei, "\n");
1276 }
1277 fclose(datei);
1278
1279 const MInt sizeArray1 = noElements * nDim;
1280 const MInt sizeArray2 = noElements * nDim * nDim;
1281
1282 // Write out geomInfo in output file
1283 file.defineArray(PIO_FLOAT, "vertices", sizeArray2);
1284 file.defineArray(PIO_FLOAT, "centroid", sizeArray1);
1285 file.defineArray(PIO_FLOAT, "normalVector", sizeArray1);
1286
1287 file.defineArray(PIO_FLOAT, "geomMeasure", noElements);
1288 IF_CONSTEXPR(nDim == 3) { file.setAttribute("area", "description", "geomMeasure"); }
1289 else {
1290 file.setAttribute("length", "description", "geomMeasure");
1291 }
1292
1293 file.setOffset(sizeArray2, 0);
1294 file.writeArray(&vertices[0], "vertices");
1295
1296 file.setOffset(sizeArray1, 0);
1297 file.writeArray(&coordinates[0], "centroid");
1298
1299 file.setOffset(sizeArray1, 0);
1300 file.writeArray(&normalVecs[0], "normalVector");
1301
1302 file.setOffset(noElements, 0);
1303 file.writeArray(&geomMeasures[0], "geomMeasure");
1304}
static MBool propertyExists(const MString &name, MInt solver=m_noSolvers)
This function checks if a property exists in general.
Definition: context.cpp:494
MFloat ** m_vertices
MInt solverId() const
Definition: samplingdata.h:118
SolverType & solver()
Return reference to solver.
Definition: samplingdata.h:129
InfoOutFile m_log
std::basic_string< char > MString
Definition: maiatypes.h:55
double MFloat
Definition: maiatypes.h:52
PARALLELIO_DEFAULT_BACKEND ParallelIo
Definition: parallelio.h:292

◆ 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_geometry

template<MInt nDim, class SolverType >
GeometryXD<nDim>::type* SurfaceData< nDim, SolverType >::m_geometry = nullptr
private

Definition at line 214 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.


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