14#include "../navierstokes/navierstokes.h"
46 PetscErrorCode
init(
const MPI_Comm &world,
const YAML::Node &node);
54 PetscErrorCode
write();
body::BodyPackBase, type::BodyPack, and factory function.
Immersed-boundary method proposed by Li et. al. (2016).
Mat E
Regularization operator.
PetscLogStage stageRHSForces
Log stage for assembling the RHS of the forces system.
Mat BNH
Projection operator for the forces.
petibm::type::BodyPack bodies
Pack of immersed bodies.
virtual PetscErrorCode writeLinSolversInfo()
Write numbers of iterations and residuals of solvers to file.
virtual PetscErrorCode applyNoSlip()
Update the velocity to satisfy the no-slip condition.
PetscViewer forcesViewer
ASCII PetscViewer object to output the forces.
virtual PetscErrorCode writeRestartDataHDF5(const std::string &filePath)
Write data required to restart a simulation into a HDF5 file.
Vec f
Vector to hold the forces at time step n.
virtual PetscErrorCode readRestartDataHDF5(const std::string &filePath)
Read data required to restart a simulation from a HDF5 file.
PetscErrorCode init(const MPI_Comm &world, const YAML::Node &node)
Initialize the decoupled IBPM solver.
DecoupledIBPMSolver()=default
Default constructor.
virtual PetscErrorCode createExtraOperators()
Assemble additional operators.
PetscLogStage stageSolveForces
Log stage for solving the forces system.
virtual PetscErrorCode writeForcesASCII()
Write the forces acting on the bodies into an ASCII file.
virtual PetscErrorCode solveForces()
Solve the system for the boundary forces.
PetscErrorCode destroy()
Manually destroy data.
PetscLogStage stageIntegrateForces
Log stage for integrating the Lagrangian forces.
Vec rhsf
Right-hand side of the forces system.
~DecoupledIBPMSolver()
Default destructor.
Vec df
Force-increment vector.
petibm::type::LinSolver fSolver
Linear solver for the Lagrangian forces.
virtual PetscErrorCode createExtraVectors()
Create additional vectors.
virtual PetscErrorCode assembleRHSForces()
Assemble the RHS vector of the system for the Lagrangian forces.
PetscErrorCode advance()
Advance the solution by one time step.
PetscErrorCode write()
Write solution and solver info to files.
virtual PetscErrorCode assembleRHSVelocity()
Assemble the RHS vector of the velocity system.
virtual PetscErrorCode updateForces()
Update the Lagrangian forces.
Mat EBNH
Left-hand side operator of the system for the forces.
Solve the incompressible Navier-Stokes equations with a projection method (Perot 1993).
PetscErrorCode ioInitialData()
Read or write initial data.
bool finished()
Evaluate if the simulation is finished.
std::shared_ptr< body::BodyPackBase > BodyPack
Definition of type::BodyPack.
std::shared_ptr< linsolver::LinSolverBase > LinSolver
Type definition of LinSolver.