Abstract Base Class of a linear interpolation object.
More...
#include <lininterp.h>
Abstract Base Class of a linear interpolation object.
- See also
- Miscellaneous functions, petibm::type::LinInterp, petibm::misc::createLinInterp
Definition at line 27 of file lininterp.h.
◆ LinInterpBase()
petibm::misc::LinInterpBase::LinInterpBase |
( |
| ) |
|
|
default |
◆ ~LinInterpBase()
petibm::misc::LinInterpBase::~LinInterpBase |
( |
| ) |
|
◆ destroy()
PetscErrorCode petibm::misc::LinInterpBase::destroy |
( |
| ) |
|
Manually destroy data in the object.
Definition at line 59 of file lininterp.cpp.
◆ getBLGridlineIndices()
PetscErrorCode petibm::misc::LinInterpBase::getBLGridlineIndices |
( |
const type::Mesh & |
mesh, |
|
|
const type::Field & |
field |
|
) |
| |
|
protected |
Get the gridline indices of the front-bottom-left neighbor.
- Parameters
-
mesh | [in] Cartesian mesh object |
field | [in] Index of the field to monitor |
- Returns
- PetscErrorCode
Definition at line 94 of file lininterp.cpp.
◆ getBoxCoords()
PetscErrorCode petibm::misc::LinInterpBase::getBoxCoords |
( |
const type::Mesh & |
mesh, |
|
|
const type::Field & |
field |
|
) |
| |
|
protected |
Get the coordinates of the neighbors.
- Parameters
-
mesh | [in] Cartesian mesh object |
field | [in] Index of the field to monitor |
- Returns
- PetscErrorCode
Definition at line 112 of file lininterp.cpp.
◆ init()
Initialize the linear interpolation object.
- Parameters
-
comm | [in] MPI communicator |
point | [in] Coordinates of the target point to interpolate |
mesh | [in] Cartesian mesh object |
field | [in] Index of the field to monitor |
- Returns
- PetscErrorCode
Definition at line 70 of file lininterp.cpp.
◆ interpolate()
virtual PetscErrorCode petibm::misc::LinInterpBase::interpolate |
( |
const DM & |
da, |
|
|
const Vec & |
vec, |
|
|
PetscReal & |
v |
|
) |
| |
|
pure virtual |
Interpolate the value at the target point.
- Parameters
-
da | [in] Parallel layout of the field to interpolate |
vec | [in] Vector with the field data |
val | [out] Interpolated value |
- Returns
- PetscErrorCode
Implemented in petibm::misc::TriLinInterp, and petibm::misc::BiLinInterp.
◆ bl
Coordinates of the front-bottom-left neighbor.
Definition at line 53 of file lininterp.h.
◆ comm
MPI_Comm petibm::misc::LinInterpBase::comm |
|
protected |
◆ commRank
PetscMPIInt petibm::misc::LinInterpBase::commRank |
|
protected |
Rank of the local process in the communicator.
Definition at line 68 of file lininterp.h.
◆ commSize
PetscMPIInt petibm::misc::LinInterpBase::commSize |
|
protected |
Size of the MPI communicator.
Definition at line 65 of file lininterp.h.
◆ idxDirs
Gridline indices of the front-bottom-left neighbor.
Definition at line 59 of file lininterp.h.
◆ target
Coordinates of the point to interpolate.
Definition at line 50 of file lininterp.h.
◆ tr
Coordinates of the back-top-right neighbor.
Definition at line 56 of file lininterp.h.
The documentation for this class was generated from the following files: