12#include <yaml-cpp/yaml.h>
41 const YAML::Node &node);
61 PetscErrorCode
init(
const MPI_Comm &world,
const YAML::Node &node);
70 PetscErrorCode
write();
244 const PetscFileMode &mode,
245 PetscViewer &viewer);
278 const std::string &filePath);
286 virtual PetscErrorCode
readTimeHDF5(
const std::string &filePath,
boundary::BoundaryBase, type::Boundary, and factory function.
Solve the incompressible Navier-Stokes equations with a projection method (Perot 1993).
Vec rhs2
Right-hand side vector of the Poisson system.
petibm::type::Solution solution
Data object holding the velocity and pressure fields.
petibm::type::LinSolver pSolver
Poisson linear solver.
PetscReal nu
Viscous diffusion coefficient.
PetscInt nstart
Starting time-step index.
PetscLogStage stageMonitor
Log stage for monitoring user-defined regions of the domain.
Vec bc1
Inhomogeneous boundary terms for the velocity system.
PetscLogStage stageSolvePoisson
Log stage for solving the Poisson system.
PetscErrorCode ioInitialData()
Read or write initial data.
YAML::Node config
YAML configuration settings.
std::vector< petibm::type::Probe > probes
Probes to monitor the solution.
petibm::type::Boundary bc
Information about the domain boundaries.
NavierStokesSolver()=default
Default constructor.
PetscInt nsave
Frequency at which the solution fields are written to files.
virtual PetscErrorCode updatePressure()
virtual PetscErrorCode assembleRHSVelocity()
Assemble the RHS vector of the velocity system.
PetscErrorCode advance()
Advance the solution by one time step.
virtual PetscErrorCode readRestartDataHDF5(const std::string &filePath)
Read data required to restart a simulation from a HDF5 file.
petibm::type::TimeIntegration convCoeffs
Time scheme for the convective terms.
virtual PetscErrorCode writeLinSolversInfo()
Write numbers of iterations and residuals of solvers to file.
PetscLogStage stageInitialize
Log stage for the initialization phase.
virtual PetscErrorCode readTimeHDF5(const std::string &filePath, PetscReal &t)
Read the time value from a HDF5 file.
bool finished()
Evaluate if the simulation is finished.
Vec rhs1
Right-hand side vector of the velocity system.
Vec dP
Pressure-correction vector.
virtual PetscErrorCode writeSolutionHDF5(const std::string &filePath)
Write the solution fields into a HDF5 file.
Mat A
Implicit operator for the velocity solver.
virtual PetscErrorCode createOperators()
Create operators.
petibm::type::LinSolver vSolver
Velocity linear solver.
PetscMPIInt commRank
Rank of the process in the MPI communicator.
PetscReal dt
Time-step size.
Mat N
Convective operator (matrix-free).
Mat LCorrection
Laplacian correction operator for boundary conditions.
std::vector< Vec > conv
Explicit convective terms.
virtual PetscErrorCode assembleRHSPoisson()
Assemble the RHS vector of the Poisson system.
MPI_Comm comm
MPI communicator.
virtual PetscErrorCode solvePoisson()
Solve the Poisson system.
virtual PetscErrorCode writeRestartDataHDF5(const std::string &filePath)
Write data required to restart a simulation into a HDF5 file.
PetscBool isRefP
True if we pin the pressure at a reference point.
~NavierStokesSolver()
Default destructor.
petibm::type::TimeIntegration diffCoeffs
Time scheme for the diffusion terms.
virtual PetscErrorCode applyDivergenceFreeVelocity()
PetscInt nrestart
Frequency at which data to restart are written to files.
Mat DCorrection
Divergence correction for boundary conditions.
PetscMPIInt commSize
Size of the MPI communicator.
PetscLogStage stageRHSVelocity
Log stage for assembling the RHS of the velocity system.
PetscErrorCode write()
Write solution and solver info to files.
PetscErrorCode init(const MPI_Comm &world, const YAML::Node &node)
Initialize the Navier-Stokes solver.
virtual PetscErrorCode setNullSpace()
Set Poisson nullspace or pin pressure at a reference point.
PetscLogStage stageWrite
Log stage when write the solution fields.
std::vector< Vec > diff
Explicit diffusion terms.
PetscLogStage stageUpdate
Log stage for updating field variables.
PetscInt ite
Time-step index.
petibm::type::Mesh mesh
Structured Cartesian mesh object.
PetscLogStage stageSolveVelocity
Log stage for solving the velocity system.
PetscLogStage stageRHSPoisson
Log stage for assembling the RHS of the Poisson system.
virtual PetscErrorCode solveVelocity()
Solve the velocity system.
virtual PetscErrorCode monitorProbes()
Monitor the solution at probes.
Mat BNG
Projection operator.
PetscInt nt
Number of time steps to compute.
virtual PetscErrorCode writeTimeHDF5(const PetscReal &t, const std::string &filePath)
Write the time value into a HDF5 file.
virtual PetscErrorCode createVectors()
Create vectors.
PetscErrorCode destroy()
Manually destroy data.
Mat D
Divergence operator.
virtual PetscErrorCode createPetscViewerASCII(const std::string &filePath, const PetscFileMode &mode, PetscViewer &viewer)
Create an ASCII PetscViewer.
Mat DBNG
Poisson operator.
PetscViewer solversViewer
ASCII PetscViewer object to output solvers info.
std::shared_ptr< boundary::BoundaryBase > Boundary
Type definition of petibm::type::Boundary.
std::shared_ptr< linsolver::LinSolverBase > LinSolver
Type definition of LinSolver.
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
std::shared_ptr< solution::SolutionBase > Solution
Type definition of solution object.
std::shared_ptr< timeintegration::TimeIntegrationBase > TimeIntegration
Definition of type::TimeIntegration.
Def. of LinSolverBase, LinSolver, and factory function.
Prototype of mesh::MeshBase, type::Mesh, and factory function.
Prototypes of factory functions for operators.
Prototype of probe classes, type::Probe, and factory function.
Definition of the class petibm::solution::SolutionBase, the type definition petibm::type::Solution,...
Definition of TimeIntegration related classes.