PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
ibpm.h
Go to the documentation of this file.
1
10#pragma once
11
12#include <petibm/bodypack.h>
13
14#include "../navierstokes/navierstokes.h"
15
23{
24public:
26 IBPMSolver() = default;
27
34 IBPMSolver(const MPI_Comm &world, const YAML::Node &node);
35
38
40 PetscErrorCode destroy();
41
48 PetscErrorCode init(const MPI_Comm &world, const YAML::Node &node);
49
51
53
55 PetscErrorCode write();
56
58
59protected:
62
64 Vec P;
65
67 IS isDE[2];
68
70 PetscLogStage stageIntegrateForces;
71
73 PetscViewer forcesViewer;
74
76 virtual PetscErrorCode assembleRHSPoisson();
77
79 virtual PetscErrorCode createOperators();
80
82 virtual PetscErrorCode createVectors();
83
85 virtual PetscErrorCode setNullSpace();
86
92 virtual PetscErrorCode writeSolutionHDF5(const std::string &filePath);
93
99 virtual PetscErrorCode writeRestartDataHDF5(const std::string &filePath);
100
106 virtual PetscErrorCode readRestartDataHDF5(const std::string &filePath);
107
109 virtual PetscErrorCode writeForcesASCII();
110
111}; // IBPMSolver
body::BodyPackBase, type::BodyPack, and factory function.
Immersed-boundary method proposed by Taira and Colonius (2007).
Definition: ibpm.h:23
virtual PetscErrorCode writeRestartDataHDF5(const std::string &filePath)
Write data required to restart a simulation into a HDF5 file.
Definition: ibpm.cpp:338
virtual PetscErrorCode createOperators()
Create operators.
Definition: ibpm.cpp:100
PetscErrorCode init(const MPI_Comm &world, const YAML::Node &node)
Initialize the IBPM solver.
Definition: ibpm.cpp:54
~IBPMSolver()
Default destructor.
Definition: ibpm.cpp:23
virtual PetscErrorCode writeForcesASCII()
Write the forces acting on the bodies into an ASCII file.
Definition: ibpm.cpp:387
PetscLogStage stageIntegrateForces
Log stage for integrating the Lagrangian forces.
Definition: ibpm.h:70
virtual PetscErrorCode createVectors()
Create vectors.
Definition: ibpm.cpp:206
virtual PetscErrorCode readRestartDataHDF5(const std::string &filePath)
Read data required to restart a simulation from a HDF5 file.
Definition: ibpm.cpp:357
PetscErrorCode write()
Write solution, forces, and solvers info to files.
Definition: ibpm.cpp:85
PetscErrorCode destroy()
Manually destroy data.
Definition: ibpm.cpp:37
IBPMSolver()=default
Default constructor.
virtual PetscErrorCode assembleRHSPoisson()
Assemble the RHS vector of the Poisson system.
Definition: ibpm.cpp:286
virtual PetscErrorCode writeSolutionHDF5(const std::string &filePath)
Write the solution fields into a HDF5 file.
Definition: ibpm.cpp:316
virtual PetscErrorCode setNullSpace()
Set Poisson nullspace or pin pressure at a reference point.
Definition: ibpm.cpp:242
IS isDE[2]
Global index sets for pressure field and Lagrangian forces.
Definition: ibpm.h:67
petibm::type::BodyPack bodies
Pack of immersed bodies.
Definition: ibpm.h:61
Vec P
PETSc Vec object with pressure field and Lagrangian forces.
Definition: ibpm.h:64
PetscViewer forcesViewer
ASCII PetscViewer object to output the forces.
Definition: ibpm.h:73
Solve the incompressible Navier-Stokes equations with a projection method (Perot 1993).
Definition: navierstokes.h:30
PetscErrorCode ioInitialData()
Read or write initial data.
PetscErrorCode advance()
Advance the solution by one time step.
bool finished()
Evaluate if the simulation is finished.
std::shared_ptr< body::BodyPackBase > BodyPack
Definition of type::BodyPack.
Definition: bodypack.h:262