PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
misc.h
Go to the documentation of this file.
1
8#pragma once
9
10#include <algorithm>
11#include <cmath>
12#include <functional>
13#include <vector>
14
15#include <petscdmda.h>
16#include <petscsys.h>
17
18#include <petibm/mesh.h>
19#include <petibm/type.h>
20
35namespace petibm
36{
41namespace misc
42{
60PetscErrorCode checkPeriodicBC(const type::IntVec2D &bcTypes,
61 type::BoolVec2D &periodic);
62
73PetscErrorCode checkBoundaryProc(const DM &da, const type::IntVec1D &n,
74 const type::BCLoc &loc, PetscBool &onThisProc);
75
86PetscErrorCode getGhostPointList(const type::Mesh &mesh,
87 const type::Field &field,
88 const type::BCLoc &loc,
89 type::GhostPointsList &points);
90
99PetscErrorCode getPerpendAxes(const PetscInt &self, type::IntVec1D &pAxes);
100
132PetscErrorCode getGhostTargetStencil(const type::IntVec1D &n,
133 const type::BCLoc &loc,
134 const type::IntVec1D &pIdx,
135 MatStencil &ghost, MatStencil &target);
136
148inline PetscErrorCode stretchGrid(const PetscReal &bg, const PetscReal &ed,
149 const PetscInt &n, const PetscReal &r,
150 type::RealVec1D &dL)
151{
152 PetscFunctionBeginUser;
153
154 dL.resize(n);
155
156 // calculate the size of the first cell
157 dL[0] = (ed - bg) * (r - 1.0) / (std::pow(r, n) - 1.0);
158
159 // dL[i] = dL[i-1] * r
160 for (auto it = dL.begin() + 1; it < dL.end(); ++it) *it = *(it - 1) * r;
161
162 PetscFunctionReturn(0);
163} // stretchGrid
164
171{
172 const PetscInt &bg, &ed;
173};
174
184inline PetscErrorCode doubleLoops(
185 const LoopBound &bound1, const LoopBound &bound2,
186 const std::function<PetscErrorCode(const PetscInt &, const PetscInt &)> &f)
187{
188 PetscErrorCode ierr;
189
190 PetscFunctionBeginUser;
191
192 for (PetscInt i = bound1.bg; i < bound1.ed; i++)
193 for (PetscInt j = bound2.bg; j < bound2.ed; j++)
194 {
195 ierr = f(i, j); CHKERRQ(ierr);
196 }
197
198 PetscFunctionReturn(0);
199} // doubleLoops
200
211inline PetscErrorCode tripleLoops(
212 const LoopBound &bound1, const LoopBound &bound2, const LoopBound &bound3,
213 const std::function<PetscErrorCode(const PetscInt &, const PetscInt &,
214 const PetscInt &)> &f)
215{
216 using namespace std::placeholders;
217
218 PetscErrorCode ierr;
219
220 PetscFunctionBeginUser;
221
222 for (PetscInt i = bound1.bg; i < bound1.ed; ++i)
223 {
224 ierr = doubleLoops(bound2, bound3, std::bind(f, i, _1, _2));
225 CHKERRQ(ierr);
226 }
227
228 PetscFunctionReturn(0);
229} // tripleLoops
230
231} // end of namespace misc
232
233} // end of namespace petibm
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
PetscErrorCode getGhostTargetStencil(const type::IntVec1D &n, const type::BCLoc &loc, const type::IntVec1D &pIdx, MatStencil &ghost, MatStencil &target)
Get the stencils of a desired ghost point and its corresponding boundary point.
Definition: misc.cpp:226
PetscErrorCode checkBoundaryProc(const DM &da, const type::IntVec1D &n, const type::BCLoc &loc, PetscBool &onThisProc)
Check if a boundary is on this process.
Definition: misc.cpp:85
PetscErrorCode getPerpendAxes(const PetscInt &self, type::IntVec1D &pAxes)
An utility to get the perpendicular axes of a desired axis.
Definition: misc.cpp:199
PetscErrorCode doubleLoops(const LoopBound &bound1, const LoopBound &bound2, const std::function< PetscErrorCode(const PetscInt &, const PetscInt &)> &f)
A helper function to carry out a double loop on a given function.
Definition: misc.h:184
PetscErrorCode tripleLoops(const LoopBound &bound1, const LoopBound &bound2, const LoopBound &bound3, const std::function< PetscErrorCode(const PetscInt &, const PetscInt &, const PetscInt &)> &f)
A helper function to carry out a triple loop on a given function.
Definition: misc.h:211
PetscErrorCode getGhostPointList(const type::Mesh &mesh, const type::Field &field, const type::BCLoc &loc, type::GhostPointsList &points)
Get a list of ghost points on a desired boundary.
Definition: misc.cpp:129
PetscErrorCode stretchGrid(const PetscReal &bg, const PetscReal &ed, const PetscInt &n, const PetscReal &r, type::RealVec1D &dL)
Calculate and return cell sizes of stretched grid in one direction.
Definition: misc.h:148
PetscErrorCode checkPeriodicBC(const type::IntVec2D &bcTypes, type::BoolVec2D &periodic)
Check if there is any periodic boundary condition and if these periodic BCs make sense.
Definition: misc.cpp:19
std::vector< IntVec1D > IntVec2D
2D std::vector holding PetscInt.
Definition: type.h:135
Field
Legal fields.
Definition: type.h:80
std::vector< BoolVec1D > BoolVec2D
2D std::vector holding PetscBool.
Definition: type.h:149
std::map< MatStencil, GhostPointInfo > GhostPointsList
A map between MatStencil and GhostPointInfo.
Definition: type.h:177
std::vector< PetscInt > IntVec1D
1D std::vector holding PetscInt.
Definition: type.h:133
BCLoc
Location of a boundary.
Definition: type.h:108
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
A helper struct to make looping function easier.
Definition: misc.h:171
const PetscInt & bg
Definition: misc.h:172
const PetscInt & ed
Definition: misc.h:172
Definition of user-defined types for convenience.