PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
lininterp.h
Go to the documentation of this file.
1
9#pragma once
10
11#include <petscksp.h>
12#include <petscsys.h>
13
14#include <petibm/mesh.h>
15#include <petibm/type.h>
16
17namespace petibm
18{
19namespace misc
20{
28{
29public:
31 LinInterpBase() = default;
32
35
37 PetscErrorCode destroy();
38
46 virtual PetscErrorCode interpolate(const DM &da, const Vec &vec, PetscReal &v) = 0;
47
48protected:
51
54
57
60
62 MPI_Comm comm;
63
65 PetscMPIInt commSize;
66
68 PetscMPIInt commRank;
69
78 PetscErrorCode init(const MPI_Comm &comm,
79 const type::RealVec1D &point,
80 const type::Mesh &mesh,
81 const type::Field &field);
82
89 PetscErrorCode getBLGridlineIndices(const type::Mesh &mesh,
90 const type::Field &field);
91
98 PetscErrorCode getBoxCoords(const type::Mesh &mesh,
99 const type::Field &field);
100
101}; // LinInterpBase
102
109{
110public:
119 TriLinInterp(const MPI_Comm &comm,
120 const type::RealVec1D &point,
121 const type::Mesh &mesh,
122 const type::Field &field);
123
125 ~TriLinInterp() = default;
126
128 PetscErrorCode interpolate(const DM &da, const Vec &vec, PetscReal &v);
129
130}; // TriLinInterp
131
138{
139public:
148 BiLinInterp(const MPI_Comm &comm,
149 const type::RealVec1D &point,
150 const type::Mesh &mesh,
151 const type::Field &field);
152
154 ~BiLinInterp() = default;
155
157 PetscErrorCode interpolate(const DM &da, const Vec &vec, PetscReal &v);
158
159}; // BiLinInterp
160
161} // end of namespace misc
162
163namespace type
164{
172typedef std::shared_ptr<misc::LinInterpBase> LinInterp;
173
174} // end of namespace type
175
176namespace misc
177{
191PetscErrorCode createLinInterp(const MPI_Comm &comm,
192 const type::RealVec1D &point,
193 const type::Mesh &mesh,
194 const type::Field &field,
195 type::LinInterp &interp);
196
197} // end of namespace misc
198
199} // end of namespace petibm
Class to perform a bilinear interpolation at a point in a 3D domain.
Definition: lininterp.h:138
PetscErrorCode interpolate(const DM &da, const Vec &vec, PetscReal &v)
Definition: lininterp.cpp:187
BiLinInterp(const MPI_Comm &comm, const type::RealVec1D &point, const type::Mesh &mesh, const type::Field &field)
Constructor. Initialize the linear interpolation object.
Definition: lininterp.cpp:178
~BiLinInterp()=default
Default destructor.
Abstract Base Class of a linear interpolation object.
Definition: lininterp.h:28
type::RealVec1D bl
Coordinates of the front-bottom-left neighbor.
Definition: lininterp.h:53
PetscMPIInt commSize
Size of the MPI communicator.
Definition: lininterp.h:65
type::RealVec1D target
Coordinates of the point to interpolate.
Definition: lininterp.h:50
type::RealVec1D tr
Coordinates of the back-top-right neighbor.
Definition: lininterp.h:56
MPI_Comm comm
MPI communicator.
Definition: lininterp.h:62
PetscErrorCode getBoxCoords(const type::Mesh &mesh, const type::Field &field)
Get the coordinates of the neighbors.
Definition: lininterp.cpp:112
PetscMPIInt commRank
Rank of the local process in the communicator.
Definition: lininterp.h:68
PetscErrorCode init(const MPI_Comm &comm, const type::RealVec1D &point, const type::Mesh &mesh, const type::Field &field)
Initialize the linear interpolation object.
Definition: lininterp.cpp:70
virtual PetscErrorCode interpolate(const DM &da, const Vec &vec, PetscReal &v)=0
Interpolate the value at the target point.
LinInterpBase()=default
Default constructor.
PetscErrorCode getBLGridlineIndices(const type::Mesh &mesh, const type::Field &field)
Get the gridline indices of the front-bottom-left neighbor.
Definition: lininterp.cpp:94
type::IntVec1D idxDirs
Gridline indices of the front-bottom-left neighbor.
Definition: lininterp.h:59
PetscErrorCode destroy()
Manually destroy data in the object.
Definition: lininterp.cpp:59
Class to perform a trilinear interpolation at a point in a 3D domain.
Definition: lininterp.h:109
PetscErrorCode interpolate(const DM &da, const Vec &vec, PetscReal &v)
Definition: lininterp.cpp:144
~TriLinInterp()=default
Default destructor.
TriLinInterp(const MPI_Comm &comm, const type::RealVec1D &point, const type::Mesh &mesh, const type::Field &field)
Constructor. Initialize the linear interpolation object.
Definition: lininterp.cpp:135
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
std::shared_ptr< misc::LinInterpBase > LinInterp
Type definition of LinInterp.
Definition: lininterp.h:172
PetscErrorCode createLinInterp(const MPI_Comm &comm, const type::RealVec1D &point, const type::Mesh &mesh, const type::Field &field, type::LinInterp &interp)
Factory function to create a linear interpolation object.
Definition: lininterp.cpp:16
Field
Legal fields.
Definition: type.h:80
std::vector< PetscInt > IntVec1D
1D std::vector holding PetscInt.
Definition: type.h:133
std::vector< PetscReal > RealVec1D
1D std::vector holding PetscReal.
Definition: type.h:140
Prototype of mesh::MeshBase, type::Mesh, and factory function.
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of user-defined types for convenience.