Peano
|
Namespaces | |
namespace | internal |
namespace | musclhancock |
namespace | riemann |
namespace | rusanov |
namespace | tests |
Functions | |
void | applyBoundaryConditions (std::function< void(const double *__restrict__ Qinside, double *__restrict__ Qoutside, const tarch::la::Vector< Dimensions, double > &faceCentre, const tarch::la::Vector< Dimensions, double > &volumeH, double t, double dt, int normal) > boundaryCondition, const tarch::la::Vector< Dimensions, double > &faceCentre, const tarch::la::Vector< Dimensions, double > &patchSize, double t, double dt, int numberOfVolumesPerAxisInPatch, int overlap, int unknowns, int faceNumber, double *__restrict__ Q) |
Apply boundary conditions. | |
static tarch::la::Vector< 2, double > | getVolumeSize (const tarch::la::Vector< 2, double > &h, int numberOfVolumesPerAxisInPatch) |
We need this routine within vectorised and GPUised code. | |
static tarch::la::Vector< 3, double > | getVolumeSize (const tarch::la::Vector< 3, double > &h, int numberOfVolumesPerAxisInPatch) |
static tarch::la::Vector< 2, double > | getFaceSize (const tarch::la::Vector< 2, double > &h, int numberOfVolumesPerAxisInPatch) |
static tarch::la::Vector< 3, double > | getFaceSize (const tarch::la::Vector< 3, double > &h, int numberOfVolumesPerAxisInPatch) |
static tarch::la::Vector< 2, double > | getVolumeCentre (const tarch::la::Vector< 2, double > &x, const tarch::la::Vector< 2, double > &h, int numberOfVolumesPerAxisInPatch, const tarch::la::Vector< 2, int > &index) |
In ExaHyPE's Finite Volume setup, a cell hosts a patch of Finite Volumes. | |
static tarch::la::Vector< 3, double > | getVolumeCentre (const tarch::la::Vector< 3, double > &x, const tarch::la::Vector< 3, double > &h, int numberOfVolumesPerAxisInPatch, const tarch::la::Vector< 3, int > &index) |
static tarch::la::Vector< 2, double > | getFaceCentre (const tarch::la::Vector< 2, double > &x, const tarch::la::Vector< 2, double > &h, int numberOfVolumesPerAxisInPatch, int overlap, int normal, const tarch::la::Vector< 2, int > &index) |
static tarch::la::Vector< 3, double > | getFaceCentre (const tarch::la::Vector< 3, double > &x, const tarch::la::Vector< 3, double > &h, int numberOfVolumesPerAxisInPatch, int overlap, int normal, const tarch::la::Vector< 3, int > &index) |
double | getVolumeLength (const tarch::la::Vector< 2, double > &h, int numberOfVolumesPerAxisInPatch) |
With GCC 10, it was impossible to return/copy the vector class. | |
double | getVolumeLength (const tarch::la::Vector< 3, double > &h, int numberOfVolumesPerAxisInPatch) |
std::string | plotVolume (const double *__restrict__ Q, int unknowns) |
Helper routine that I need in the log statements. | |
void | validatePatch (const double *__restrict__ Q, int unknowns, int auxiliaryVariables, int numberOfVolumesPerAxisInPatch, int haloSize, const std::string &location="", bool triggerNonCriticalAssertion=true, double *minValues=nullptr, double *maxValues=nullptr) |
Just runs over the patch and ensures that no entry is non or infinite. | |
std::string | plotPatch (const double *__restrict__ Q, int unknowns, int auxiliaryVariables, int numberOfVolumesPerAxisInPatch, int haloSize, bool prettyPrint=false) |
Plot patch. | |
std::string | plotPatchOverlap (const double *__restrict__ Q, int unknowns, int auxiliaryVariables, int numberOfGridCellsPerPatchPerAxis, int haloSize, int normal, bool prettyPrint=false) |
void | copyHalfOfHalo (int unknownsPlusAuxiliaryVariables, int numberOfGridCellsPerPatchPerAxis, int haloSize, int normal, bool isRightLayer, const double *__restrict__ srcQ, double *__restrict__ destQ) |
A face always holds a left and a right overlap. | |
void | mapInnerNeighbourVoxelAlongBoundayOntoAuxiliaryVariable (double *__restrict__ Q, int unknowns, int auxiliaryVariables, int numberOfVolumesPerAxisInPatch, int fromIndex, int toIndex) |
Simple extrapolation within patch. | |
template<int SourceIndex, int DestIndex, typename Particle > | |
void | projectValueOntoParticle_piecewiseConstant (const peano4::datamanagement::CellMarker &marker, int voxelsPerAxis, int unknownsPerVoxel, const double *__restrict__ Q, Particle &particle) |
Project one quantity from the patch data onto the particle. | |
template<typename Particle > | |
void | projectAllValuesOntoParticle_piecewiseConstant (const peano4::datamanagement::CellMarker &marker, int voxelsPerAxis, int unknownsPerVoxel, const double *__restrict__ Q, Particle &particle) |
Project all values from solution onto particle. | |
template<int SourceIndex, int DestIndex, typename Particle > | |
void | projectValueOntoParticle_piecewiseLinear (const peano4::datamanagement::CellMarker &marker, int voxelsPerAxis, int unknownsPerVoxel, const double *__restrict__ Q, Particle &particle) |
template<typename Particle > | |
void | projectAllValuesOntoParticle_piecewiseLinear (const peano4::datamanagement::CellMarker &marker, int voxelsPerAxis, int unknownsPerVoxel, const double *__restrict__ Q, Particle &particle) |
Project all values from solution onto particle and use a piecewise linear interpolation. | |
template<typename Particle > | |
void | projectAllValuesOntoParticle_piecewiseLinear_explicit_Euler (const peano4::datamanagement::CellMarker &marker, int voxelsPerAxis, int unknownsPerVoxel, const double *__restrict__ Q, double timeStepSize, tarch::la::Vector< Dimensions, int > indices, tarch::la::Vector< Dimensions, int > factors, Particle &particle) |
Project all values from solution onto particle and use a piecewise linear interpolation and advance particle along an explicit Euler. | |
template<int SourceIndex, int... SourceIndices, typename Particle > | |
void | projectValuesOntoParticle_piecewiseConstant (const peano4::datamanagement::CellMarker &marker, int voxelsPerAxis, int unknownsPerVoxel, const double *__restrict__ Q, Particle &particle) |
Map multiple variables from input field onto particle. | |
template<int SourceIndex, int... SourceIndices, typename Particle > | |
void | projectValuesOntoParticle_piecewiseLinear (const peano4::datamanagement::CellMarker &marker, int voxelsPerAxis, int unknownsPerVoxel, const double *__restrict__ Q, Particle &particle) |
void exahype2::fv::applyBoundaryConditions | ( | std::function< void(const double *__restrict__ Qinside, double *__restrict__ Qoutside, const tarch::la::Vector< Dimensions, double > &faceCentre, const tarch::la::Vector< Dimensions, double > &volumeH, double t, double dt, int normal) > | boundaryCondition, |
const tarch::la::Vector< Dimensions, double > & | faceCentre, | ||
const tarch::la::Vector< Dimensions, double > & | patchSize, | ||
double | t, | ||
double | dt, | ||
int | numberOfVolumesPerAxisInPatch, | ||
int | overlap, | ||
int | unknowns, | ||
int | faceNumber, | ||
double *__restrict__ | Q ) |
Apply boundary conditions.
Works only for a halo size of 1.
faceNumber | Is usually taken from marker.getSelectedFaceNumber() and is thus a number between 0 and 2d-1. |
Definition at line 10 of file BoundaryConditions.cpp.
References _log, dfore, getVolumeSize(), logDebug, logTraceInWith4Arguments, logTraceOut, and tarch::la::multiplyComponents().
void exahype2::fv::copyHalfOfHalo | ( | int | unknownsPlusAuxiliaryVariables, |
int | numberOfGridCellsPerPatchPerAxis, | ||
int | haloSize, | ||
int | normal, | ||
bool | isRightLayer, | ||
const double *__restrict__ | srcQ, | ||
double *__restrict__ | destQ ) |
A face always holds a left and a right overlap.
This routine copies over half of the overlap
isRightLayer |
Definition at line 15 of file PatchUtils.cpp.
References dfore, and toolbox::blockstructured::serialiseVoxelIndexInOverlap().
Referenced by exahype2::dg::copyOneSideOfFaceProjection().
|
static |
Definition at line 128 of file PatchUtils.h.
References getFaceSize().
|
static |
Definition at line 155 of file PatchUtils.h.
References getVolumeSize().
|
static |
Definition at line 59 of file PatchUtils.h.
References getVolumeSize().
Referenced by getFaceCentre().
|
static |
Definition at line 69 of file PatchUtils.h.
References getVolumeSize().
|
static |
In ExaHyPE's Finite Volume setup, a cell hosts a patch of Finite Volumes.
When we iterate over these volumes, we typically have to know the centre of the volume.
I use this routine in a lot of the loops that have to be vectorised. Therefore, it is absolutely essential that the compiler can inline them. Inlining works properly with most compilers if and only if the definition is available in the header. This is the reason why this implementation ended up here and not in the cpp file.
x | Centre of the cell |
h | Size of the cell |
index | Index of Finite Volume (in lexicographic ordering) |
Definition at line 94 of file PatchUtils.h.
References getVolumeSize().
Referenced by exahype2::fv::musclhancock::internal::computeTimeDerivative_LoopBody(), and exahype2::fv::musclhancock::internal::updateSolutionwithNCPandSource_LoopBody().
|
static |
Definition at line 111 of file PatchUtils.h.
References getVolumeSize().
double exahype2::fv::getVolumeLength | ( | const tarch::la::Vector< 2, double > & | h, |
int | numberOfVolumesPerAxisInPatch ) |
With GCC 10, it was impossible to return/copy the vector class.
We almost never need it however, as we work with cubes. This specialisation thus does the job.
Definition at line 38 of file PatchUtils.cpp.
Referenced by exahype2::fv::internal::projectValueOntoParticle_piecewiseLinear().
double exahype2::fv::getVolumeLength | ( | const tarch::la::Vector< 3, double > & | h, |
int | numberOfVolumesPerAxisInPatch ) |
Definition at line 46 of file PatchUtils.cpp.
References assertion2.
|
static |
We need this routine within vectorised and GPUised code.
Therefore, it has to be in the header (as link-time optimisation is quite tricky for many compilers). To avoid the implied multiple symbol errors, I declare the function as static. Now the compilers plot warning, but everything seems to work.
Definition at line 24 of file PatchUtils.h.
References assertion2.
Referenced by applyBoundaryConditions(), exahype2::fv::musclhancock::internal::computeQonFace_LoopBody(), exahype2::fv::musclhancock::internal::computeTimeDerivative_LoopBody(), getFaceCentre(), getFaceSize(), getFaceSize(), getVolumeCentre(), getVolumeCentre(), and exahype2::fv::musclhancock::internal::updateSolutionwithNCPandSource_LoopBody().
|
static |
Definition at line 43 of file PatchUtils.h.
void exahype2::fv::mapInnerNeighbourVoxelAlongBoundayOntoAuxiliaryVariable | ( | double *__restrict__ | Q, |
int | unknowns, | ||
int | auxiliaryVariables, | ||
int | numberOfVolumesPerAxisInPatch, | ||
int | fromIndex, | ||
int | toIndex ) |
Simple extrapolation within patch.
The routine accepts a NxNxN patch. It runs over the layer
A typical usage first adds the header to the solver, so the present postprocessing routines are known to the compiler:
thesolver.add_user_action_set_includes( """ #include "exahype2/fv/PostprocessingKernels.h" """)
After that, it adds the actual postprocessing call. In the example below, we take the first unknown (index 0) and map it onto the first auxiliary variable (index {{NUMBER_OF_UNKNOWNS}}).
thesolver.postprocess_updated_patch = """ ::exahype2::fv::mapInnerNeighbourVoxelAlongBoundayOntoAuxiliaryVariable( targetPatch, {{NUMBER_OF_VOLUMES_PER_AXIS}}, {{NUMBER_OF_UNKNOWNS}}, {{NUMBER_OF_AUXILIARY_VARIABLES}}, 0, {{NUMBER_OF_UNKNOWNS}} ); """
Definition at line 8 of file PostprocessingKernels.cpp.
References dfore.
std::string exahype2::fv::plotPatch | ( | const double *__restrict__ | Q, |
int | unknowns, | ||
int | auxiliaryVariables, | ||
int | numberOfVolumesPerAxisInPatch, | ||
int | haloSize, | ||
bool | prettyPrint = false ) |
Plot patch.
Usually used for debugging.
Assumes all data are held as AoS.
Definition at line 108 of file PatchUtils.cpp.
References tarch::la::count(), dfor, and peano4::utils::dLinearised().
Referenced by exahype2::fd::plotPatch().
std::string exahype2::fv::plotPatchOverlap | ( | const double *__restrict__ | Q, |
int | unknowns, | ||
int | auxiliaryVariables, | ||
int | numberOfGridCellsPerPatchPerAxis, | ||
int | haloSize, | ||
int | normal, | ||
bool | prettyPrint = false ) |
Definition at line 179 of file PatchUtils.cpp.
References assertion.
Referenced by exahype2::fd::plotPatchOverlap().
Helper routine that I need in the log statements.
Applies to AoS data only.
Definition at line 56 of file PatchUtils.cpp.
Referenced by exahype2::fd::plotGridCell().
void exahype2::fv::projectAllValuesOntoParticle_piecewiseConstant | ( | const peano4::datamanagement::CellMarker & | marker, |
int | voxelsPerAxis, | ||
int | unknownsPerVoxel, | ||
const double *__restrict__ | Q, | ||
Particle & | particle ) |
Project all values from solution onto particle.
This routine assumes that the particle hosts exactly the same number of unknown as you store within the solution, i.e. within each entry of Q. The count refers to the unknowns plus auxiliary variables within ExaHyPE. If you host a PDE with 10 unknowns and 3 material parameters, then this routine works if and only if the particle carries 13 unknowns.
The solution assumes that the Q represents a piecewise constant solution and consequently maps the data piecewise constant. This holds for Finite Volumes, e.g., where each unkonwn within Q stores data in the AoS format and represents all the values within that finite volume. As the FV scheme is 0th order, the data within the finite volume can be considered to be piecewise constant, and the piecewise constant interpolation hence makes sense. You might however decide that you prefer a higher order interpolation scheme. In this case, consult projectAllValuesOntoParticle_piecewiseLinear(), e.g.
Some codes do not track all the unknowns via the tracers. In this case, you have to use specialised routines to map the solution onto the particle such as projectValueOntoParticle_piecewiseConstant() or projectValuesOntoParticle_piecewiseConstant().
The usage of this routine within the exaype2.tracers.FiniteVolumesTracing is trivial: Pass in
as argument, and all is done. The routine here is a template, but the template type will be induced by C++ automatically through the function parameters.
Definition at line 108 of file Tracer.cpph.
References exahype2::fv::internal::projectValueOntoParticle_piecewiseConstant().
void exahype2::fv::projectAllValuesOntoParticle_piecewiseLinear | ( | const peano4::datamanagement::CellMarker & | marker, |
int | voxelsPerAxis, | ||
int | unknownsPerVoxel, | ||
const double *__restrict__ | Q, | ||
Particle & | particle ) |
Project all values from solution onto particle and use a piecewise linear interpolation.
Definition at line 139 of file Tracer.cpph.
References assertion, and exahype2::fv::internal::projectValueOntoParticle_piecewiseLinear().
void exahype2::fv::projectAllValuesOntoParticle_piecewiseLinear_explicit_Euler | ( | const peano4::datamanagement::CellMarker & | marker, |
int | voxelsPerAxis, | ||
int | unknownsPerVoxel, | ||
const double *__restrict__ | Q, | ||
double | timeStepSize, | ||
tarch::la::Vector< Dimensions, int > | indices, | ||
tarch::la::Vector< Dimensions, int > | factors, | ||
Particle & | particle ) |
Project all values from solution onto particle and use a piecewise linear interpolation and advance particle along an explicit Euler.
This is a natural extension of projectAllValuesOntoParticle_piecewiseLinear() which accepts two additional arguments. The first one is the time step size, and the second one is a set of Dimensions integers. The routine takes the input data and maps it onto the particle attributes. After that, it assumes that the indices indicies[0], indicies[1] and indicies[2] actually denote the particle's velocities and it updates the velocity according to an explicit Euler. We also provide a set of scale factor for velocity in case the user would like to tune the velocity accordingly. e.g. if the velocity is the inverse of the indiced quantities, enter {-1,-1,-1}.
To use this extended mapping, you have to use the Python API similar to the following snippet:
Compared to other variants, we have replaced the function call with projectAllValuesOntoParticle_piecewiseLinear_explicit_Euler. This kernel now accepts a few additional parameters besides
Definition at line 157 of file Tracer.cpph.
References assertion, and exahype2::fv::internal::projectValueOntoParticle_piecewiseLinear().
void exahype2::fv::projectValueOntoParticle_piecewiseConstant | ( | const peano4::datamanagement::CellMarker & | marker, |
int | voxelsPerAxis, | ||
int | unknownsPerVoxel, | ||
const double *__restrict__ | Q, | ||
Particle & | particle ) |
Project one quantity from the patch data onto the particle.
The quantities of interest are given as template parameters, i.e. through SourceIndex and DestIndex. This sounds weird (and likely is a strange choice), but some of ExaHyPE's projection routines want to have a plain function call to project the solution onto a particle, and adding the source and destination index as template argument allows us to have exactly this behaviour.
The usage of this routine within the exaype2.tracers.FiniteVolumesTracing is trivial: Pass in
if you want to map the fourth entry of the solution onto the third one of the tracer, e.g.
marker | Cell marker describing the cell's geometric/spatial properties. |
voxelsPerAxis | Describe patch. We assume that the pathch has not halo. |
unknownsPerVoxel | |
Q | Voxel field, i.e. actual patch data. Has the dimensions \( voxelsPerAxis^d \cdot unknownsPerVoxel \). |
particleX | Position of particle. |
unknown | Which unknown from data field to pick. |
Definition at line 92 of file Tracer.cpph.
References exahype2::fv::internal::projectValueOntoParticle_piecewiseConstant().
void exahype2::fv::projectValueOntoParticle_piecewiseLinear | ( | const peano4::datamanagement::CellMarker & | marker, |
int | voxelsPerAxis, | ||
int | unknownsPerVoxel, | ||
const double *__restrict__ | Q, | ||
Particle & | particle ) |
Definition at line 126 of file Tracer.cpph.
References exahype2::fv::internal::projectValueOntoParticle_piecewiseLinear().
void exahype2::fv::projectValuesOntoParticle_piecewiseConstant | ( | const peano4::datamanagement::CellMarker & | marker, |
int | voxelsPerAxis, | ||
int | unknownsPerVoxel, | ||
const double *__restrict__ | Q, | ||
Particle & | particle ) |
Map multiple variables from input field onto particle.
The usage of this routine within the exaype2.tracers.FiniteVolumesTracing is simplistic due to the variadic template:
maps the arguments 4, 5 and 9 from the solution field onto the tracer. Hence, the function implicitly assumes that the particle has (at least) three unkowns.
Definition at line 6 of file Tracer.cpph.
References exahype2::fv::internal::projectValuesOntoParticle_piecewiseConstant().
void exahype2::fv::projectValuesOntoParticle_piecewiseLinear | ( | const peano4::datamanagement::CellMarker & | marker, |
int | voxelsPerAxis, | ||
int | unknownsPerVoxel, | ||
const double *__restrict__ | Q, | ||
Particle & | particle ) |
Definition at line 19 of file Tracer.cpph.
References exahype2::fv::internal::projectValuesOntoParticle_piecewiseLinear().
void exahype2::fv::validatePatch | ( | const double *__restrict__ | Q, |
int | unknowns, | ||
int | auxiliaryVariables, | ||
int | numberOfVolumesPerAxisInPatch, | ||
int | haloSize, | ||
const std::string & | location = "", | ||
bool | triggerNonCriticalAssertion = true, | ||
double * | minValues = nullptr, | ||
double * | maxValues = nullptr ) |
Just runs over the patch and ensures that no entry is non or infinite.
In ExaHyPE, we primarily work with split approaches. That is, diagonal halo entries are never initialised properly: We can copy over the face-connected data, but we lack information through the diagonals. This routine takes this into account when it validates the entries.
Assumes all data are held as AoS.
location | String that tells system from where this routine got called |
minValues | Is either a nullptr or it points to a double array with exactly unknowns+auxiliaryVariables entries |
Definition at line 263 of file PatchUtils.cpp.
References _log, assertion, assertion12, dfor, peano4::utils::dLinearised(), logTraceInWith5Arguments, logTraceOut, and nonCriticalAssertion11.
Referenced by exahype2::fd::validatePatch().