PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
NavierStokesSolver Class Reference

Solve the incompressible Navier-Stokes equations with a projection method (Perot 1993). More...

#include <navierstokes.h>

Inheritance diagram for NavierStokesSolver:
[legend]

Public Member Functions

 NavierStokesSolver ()=default
 Default constructor. More...
 
 NavierStokesSolver (const MPI_Comm &world, const YAML::Node &node)
 Constructor. Initialize the Navier-Stokes solver. More...
 
 ~NavierStokesSolver ()
 Default destructor. More...
 
PetscErrorCode destroy ()
 Manually destroy data. More...
 
PetscErrorCode init (const MPI_Comm &world, const YAML::Node &node)
 Initialize the Navier-Stokes solver. More...
 
PetscErrorCode ioInitialData ()
 Read or write initial data. More...
 
PetscErrorCode advance ()
 Advance the solution by one time step. More...
 
PetscErrorCode write ()
 Write solution and solver info to files. More...
 
bool finished ()
 Evaluate if the simulation is finished. More...
 

Protected Member Functions

virtual PetscErrorCode assembleRHSVelocity ()
 Assemble the RHS vector of the velocity system. More...
 
virtual PetscErrorCode solveVelocity ()
 Solve the velocity system. More...
 
virtual PetscErrorCode assembleRHSPoisson ()
 Assemble the RHS vector of the Poisson system. More...
 
virtual PetscErrorCode solvePoisson ()
 Solve the Poisson system. More...
 
virtual PetscErrorCode applyDivergenceFreeVelocity ()
 
virtual PetscErrorCode updatePressure ()
 
virtual PetscErrorCode createOperators ()
 Create operators. More...
 
virtual PetscErrorCode createVectors ()
 Create vectors. More...
 
virtual PetscErrorCode setNullSpace ()
 Set Poisson nullspace or pin pressure at a reference point. More...
 
virtual PetscErrorCode createPetscViewerASCII (const std::string &filePath, const PetscFileMode &mode, PetscViewer &viewer)
 Create an ASCII PetscViewer. More...
 
virtual PetscErrorCode writeSolutionHDF5 (const std::string &filePath)
 Write the solution fields into a HDF5 file. More...
 
virtual PetscErrorCode writeRestartDataHDF5 (const std::string &filePath)
 Write data required to restart a simulation into a HDF5 file. More...
 
virtual PetscErrorCode readRestartDataHDF5 (const std::string &filePath)
 Read data required to restart a simulation from a HDF5 file. More...
 
virtual PetscErrorCode writeLinSolversInfo ()
 Write numbers of iterations and residuals of solvers to file. More...
 
virtual PetscErrorCode writeTimeHDF5 (const PetscReal &t, const std::string &filePath)
 Write the time value into a HDF5 file. More...
 
virtual PetscErrorCode readTimeHDF5 (const std::string &filePath, PetscReal &t)
 Read the time value from a HDF5 file. More...
 
virtual PetscErrorCode monitorProbes ()
 Monitor the solution at probes. More...
 

Protected Attributes

MPI_Comm comm
 MPI communicator. More...
 
PetscMPIInt commSize
 Size of the MPI communicator. More...
 
PetscMPIInt commRank
 Rank of the process in the MPI communicator. More...
 
YAML::Node config
 YAML configuration settings. More...
 
petibm::type::Mesh mesh
 Structured Cartesian mesh object. More...
 
petibm::type::Boundary bc
 Information about the domain boundaries. More...
 
petibm::type::Solution solution
 Data object holding the velocity and pressure fields. More...
 
petibm::type::TimeIntegration convCoeffs
 Time scheme for the convective terms. More...
 
petibm::type::TimeIntegration diffCoeffs
 Time scheme for the diffusion terms. More...
 
petibm::type::LinSolver vSolver
 Velocity linear solver. More...
 
petibm::type::LinSolver pSolver
 Poisson linear solver. More...
 
std::vector< petibm::type::Probeprobes
 Probes to monitor the solution. More...
 
PetscReal dt
 Time-step size. More...
 
PetscInt ite
 Time-step index. More...
 
PetscReal t
 Time value. More...
 
PetscInt nstart
 Starting time-step index. More...
 
PetscInt nt
 Number of time steps to compute. More...
 
PetscInt nsave
 Frequency at which the solution fields are written to files. More...
 
PetscInt nrestart
 Frequency at which data to restart are written to files. More...
 
PetscReal nu
 Viscous diffusion coefficient. More...
 
Mat L
 Laplacian operator. More...
 
Mat LCorrection
 Laplacian correction operator for boundary conditions. More...
 
Mat G
 Gradient operator. More...
 
Mat D
 Divergence operator. More...
 
Mat DCorrection
 Divergence correction for boundary conditions. More...
 
Mat N
 Convective operator (matrix-free). More...
 
Mat A
 Implicit operator for the velocity solver. More...
 
Mat BNG
 Projection operator. More...
 
Mat DBNG
 Poisson operator. More...
 
Vec dP
 Pressure-correction vector. More...
 
Vec bc1
 Inhomogeneous boundary terms for the velocity system. More...
 
Vec rhs1
 Right-hand side vector of the velocity system. More...
 
Vec rhs2
 Right-hand side vector of the Poisson system. More...
 
std::vector< Vec > conv
 Explicit convective terms. More...
 
std::vector< Vec > diff
 Explicit diffusion terms. More...
 
PetscBool isRefP
 True if we pin the pressure at a reference point. More...
 
PetscLogStage stageInitialize
 Log stage for the initialization phase. More...
 
PetscLogStage stageRHSVelocity
 Log stage for assembling the RHS of the velocity system. More...
 
PetscLogStage stageSolveVelocity
 Log stage for solving the velocity system. More...
 
PetscLogStage stageRHSPoisson
 Log stage for assembling the RHS of the Poisson system. More...
 
PetscLogStage stageSolvePoisson
 Log stage for solving the Poisson system. More...
 
PetscLogStage stageUpdate
 Log stage for updating field variables. More...
 
PetscLogStage stageWrite
 Log stage when write the solution fields. More...
 
PetscLogStage stageMonitor
 Log stage for monitoring user-defined regions of the domain. More...
 
PetscViewer solversViewer
 ASCII PetscViewer object to output solvers info. More...
 

Detailed Description

Solve the incompressible Navier-Stokes equations with a projection method (Perot 1993).

See also
Navier-Stokes solver

Definition at line 29 of file navierstokes.h.

Constructor & Destructor Documentation

◆ NavierStokesSolver() [1/2]

NavierStokesSolver::NavierStokesSolver ( )
default

Default constructor.

◆ NavierStokesSolver() [2/2]

NavierStokesSolver::NavierStokesSolver ( const MPI_Comm &  world,
const YAML::Node &  node 
)

Constructor. Initialize the Navier-Stokes solver.

Parameters
world[in] MPI communicator
node[in] YAML configuration settings

Definition at line 18 of file navierstokes.cpp.

◆ ~NavierStokesSolver()

NavierStokesSolver::~NavierStokesSolver ( )

Default destructor.

Definition at line 24 of file navierstokes.cpp.

Member Function Documentation

◆ advance()

PetscErrorCode NavierStokesSolver::advance ( )

Advance the solution by one time step.

Definition at line 240 of file navierstokes.cpp.

◆ applyDivergenceFreeVelocity()

PetscErrorCode NavierStokesSolver::applyDivergenceFreeVelocity ( )
protectedvirtual

Definition at line 583 of file navierstokes.cpp.

◆ assembleRHSPoisson()

PetscErrorCode NavierStokesSolver::assembleRHSPoisson ( )
protectedvirtual

Assemble the RHS vector of the Poisson system.

Reimplemented in IBPMSolver.

Definition at line 540 of file navierstokes.cpp.

◆ assembleRHSVelocity()

PetscErrorCode NavierStokesSolver::assembleRHSVelocity ( )
protectedvirtual

Assemble the RHS vector of the velocity system.

Reimplemented in DecoupledIBPMSolver.

Definition at line 432 of file navierstokes.cpp.

◆ createOperators()

PetscErrorCode NavierStokesSolver::createOperators ( )
protectedvirtual

Create operators.

Reimplemented in IBPMSolver.

Definition at line 317 of file navierstokes.cpp.

◆ createPetscViewerASCII()

PetscErrorCode NavierStokesSolver::createPetscViewerASCII ( const std::string &  filePath,
const PetscFileMode &  mode,
PetscViewer &  viewer 
)
protectedvirtual

Create an ASCII PetscViewer.

Parameters
filePath[in] Path of the file to write in
mode[in] File mode
viewer[out] PetscViewer object
Returns
PetscErrorCode

Definition at line 749 of file navierstokes.cpp.

◆ createVectors()

PetscErrorCode NavierStokesSolver::createVectors ( )
protectedvirtual

Create vectors.

Reimplemented in IBPMSolver.

Definition at line 368 of file navierstokes.cpp.

◆ destroy()

PetscErrorCode NavierStokesSolver::destroy ( )

Manually destroy data.

Definition at line 38 of file navierstokes.cpp.

◆ finished()

bool NavierStokesSolver::finished ( )

Evaluate if the simulation is finished.

Definition at line 311 of file navierstokes.cpp.

◆ init()

PetscErrorCode NavierStokesSolver::init ( const MPI_Comm &  world,
const YAML::Node &  node 
)

Initialize the Navier-Stokes solver.

Create a structured Cartesian mesh and write it into a HDF5 file. Create the initial field solution. Create the linear solvers (along with their vectors and operators). Create probes to monitor the solution in some user-defined regions. Record some simulation parameters.

Parameters
world[in] MPI communicator
node[in] YAML configuration settings
Returns
PetscErrorCode

Definition at line 93 of file navierstokes.cpp.

◆ ioInitialData()

PetscErrorCode NavierStokesSolver::ioInitialData ( )

Read or write initial data.

Definition at line 207 of file navierstokes.cpp.

◆ monitorProbes()

PetscErrorCode NavierStokesSolver::monitorProbes ( )
protectedvirtual

Monitor the solution at probes.

Definition at line 840 of file navierstokes.cpp.

◆ readRestartDataHDF5()

PetscErrorCode NavierStokesSolver::readRestartDataHDF5 ( const std::string &  filePath)
protectedvirtual

Read data required to restart a simulation from a HDF5 file.

Parameters
filePath[in] Path of the file to read from
Returns
PetscErrorCode

Reimplemented in DecoupledIBPMSolver, and IBPMSolver.

Definition at line 691 of file navierstokes.cpp.

◆ readTimeHDF5()

PetscErrorCode NavierStokesSolver::readTimeHDF5 ( const std::string &  filePath,
PetscReal &  t 
)
protectedvirtual

Read the time value from a HDF5 file.

Parameters
filePath[in] Path of the file to read from
t[out] Time
Returns
PetscErrorCode

Definition at line 818 of file navierstokes.cpp.

◆ setNullSpace()

PetscErrorCode NavierStokesSolver::setNullSpace ( )
protectedvirtual

Set Poisson nullspace or pin pressure at a reference point.

Reimplemented in IBPMSolver.

Definition at line 395 of file navierstokes.cpp.

◆ solvePoisson()

PetscErrorCode NavierStokesSolver::solvePoisson ( )
protectedvirtual

Solve the Poisson system.

Definition at line 566 of file navierstokes.cpp.

◆ solveVelocity()

PetscErrorCode NavierStokesSolver::solveVelocity ( )
protectedvirtual

Solve the velocity system.

Definition at line 524 of file navierstokes.cpp.

◆ updatePressure()

PetscErrorCode NavierStokesSolver::updatePressure ( )
protectedvirtual

Definition at line 601 of file navierstokes.cpp.

◆ write()

PetscErrorCode NavierStokesSolver::write ( )

Write solution and solver info to files.

Definition at line 269 of file navierstokes.cpp.

◆ writeLinSolversInfo()

PetscErrorCode NavierStokesSolver::writeLinSolversInfo ( )
protectedvirtual

Write numbers of iterations and residuals of solvers to file.

Reimplemented in DecoupledIBPMSolver.

Definition at line 766 of file navierstokes.cpp.

◆ writeRestartDataHDF5()

PetscErrorCode NavierStokesSolver::writeRestartDataHDF5 ( const std::string &  filePath)
protectedvirtual

Write data required to restart a simulation into a HDF5 file.

Parameters
filePath[in] Path of the file to write in
Returns
PetscErrorCode

Reimplemented in DecoupledIBPMSolver, and IBPMSolver.

Definition at line 637 of file navierstokes.cpp.

◆ writeSolutionHDF5()

PetscErrorCode NavierStokesSolver::writeSolutionHDF5 ( const std::string &  filePath)
protectedvirtual

Write the solution fields into a HDF5 file.

Parameters
filePath[in] Path of the file to write in
Returns
PetscErrorCode

Reimplemented in IBPMSolver.

Definition at line 618 of file navierstokes.cpp.

◆ writeTimeHDF5()

PetscErrorCode NavierStokesSolver::writeTimeHDF5 ( const PetscReal &  t,
const std::string &  filePath 
)
protectedvirtual

Write the time value into a HDF5 file.

Parameters
t[in] Time
filePath[in] Path of the file to write in
Returns
PetscErrorCode

Definition at line 797 of file navierstokes.cpp.

Member Data Documentation

◆ A

Mat NavierStokesSolver::A
protected

Implicit operator for the velocity solver.

Definition at line 155 of file navierstokes.h.

◆ bc

petibm::type::Boundary NavierStokesSolver::bc
protected

Information about the domain boundaries.

Definition at line 92 of file navierstokes.h.

◆ bc1

Vec NavierStokesSolver::bc1
protected

Inhomogeneous boundary terms for the velocity system.

Definition at line 167 of file navierstokes.h.

◆ BNG

Mat NavierStokesSolver::BNG
protected

Projection operator.

Definition at line 158 of file navierstokes.h.

◆ comm

MPI_Comm NavierStokesSolver::comm
protected

MPI communicator.

Definition at line 77 of file navierstokes.h.

◆ commRank

PetscMPIInt NavierStokesSolver::commRank
protected

Rank of the process in the MPI communicator.

Definition at line 83 of file navierstokes.h.

◆ commSize

PetscMPIInt NavierStokesSolver::commSize
protected

Size of the MPI communicator.

Definition at line 80 of file navierstokes.h.

◆ config

YAML::Node NavierStokesSolver::config
protected

YAML configuration settings.

Definition at line 86 of file navierstokes.h.

◆ conv

std::vector<Vec> NavierStokesSolver::conv
protected

Explicit convective terms.

Definition at line 176 of file navierstokes.h.

◆ convCoeffs

petibm::type::TimeIntegration NavierStokesSolver::convCoeffs
protected

Time scheme for the convective terms.

Definition at line 98 of file navierstokes.h.

◆ D

Mat NavierStokesSolver::D
protected

Divergence operator.

Definition at line 146 of file navierstokes.h.

◆ DBNG

Mat NavierStokesSolver::DBNG
protected

Poisson operator.

Definition at line 161 of file navierstokes.h.

◆ DCorrection

Mat NavierStokesSolver::DCorrection
protected

Divergence correction for boundary conditions.

Definition at line 149 of file navierstokes.h.

◆ diff

std::vector<Vec> NavierStokesSolver::diff
protected

Explicit diffusion terms.

Definition at line 179 of file navierstokes.h.

◆ diffCoeffs

petibm::type::TimeIntegration NavierStokesSolver::diffCoeffs
protected

Time scheme for the diffusion terms.

Definition at line 101 of file navierstokes.h.

◆ dP

Vec NavierStokesSolver::dP
protected

Pressure-correction vector.

Definition at line 164 of file navierstokes.h.

◆ dt

PetscReal NavierStokesSolver::dt
protected

Time-step size.

Definition at line 113 of file navierstokes.h.

◆ G

Mat NavierStokesSolver::G
protected

Gradient operator.

Definition at line 143 of file navierstokes.h.

◆ isRefP

PetscBool NavierStokesSolver::isRefP
protected

True if we pin the pressure at a reference point.

Definition at line 182 of file navierstokes.h.

◆ ite

PetscInt NavierStokesSolver::ite
protected

Time-step index.

Definition at line 116 of file navierstokes.h.

◆ L

Mat NavierStokesSolver::L
protected

Laplacian operator.

Definition at line 137 of file navierstokes.h.

◆ LCorrection

Mat NavierStokesSolver::LCorrection
protected

Laplacian correction operator for boundary conditions.

Definition at line 140 of file navierstokes.h.

◆ mesh

petibm::type::Mesh NavierStokesSolver::mesh
protected

Structured Cartesian mesh object.

Definition at line 89 of file navierstokes.h.

◆ N

Mat NavierStokesSolver::N
protected

Convective operator (matrix-free).

Definition at line 152 of file navierstokes.h.

◆ nrestart

PetscInt NavierStokesSolver::nrestart
protected

Frequency at which data to restart are written to files.

Definition at line 131 of file navierstokes.h.

◆ nsave

PetscInt NavierStokesSolver::nsave
protected

Frequency at which the solution fields are written to files.

Definition at line 128 of file navierstokes.h.

◆ nstart

PetscInt NavierStokesSolver::nstart
protected

Starting time-step index.

Definition at line 122 of file navierstokes.h.

◆ nt

PetscInt NavierStokesSolver::nt
protected

Number of time steps to compute.

Definition at line 125 of file navierstokes.h.

◆ nu

PetscReal NavierStokesSolver::nu
protected

Viscous diffusion coefficient.

Definition at line 134 of file navierstokes.h.

◆ probes

std::vector<petibm::type::Probe> NavierStokesSolver::probes
protected

Probes to monitor the solution.

Definition at line 110 of file navierstokes.h.

◆ pSolver

petibm::type::LinSolver NavierStokesSolver::pSolver
protected

Poisson linear solver.

Definition at line 107 of file navierstokes.h.

◆ rhs1

Vec NavierStokesSolver::rhs1
protected

Right-hand side vector of the velocity system.

Definition at line 170 of file navierstokes.h.

◆ rhs2

Vec NavierStokesSolver::rhs2
protected

Right-hand side vector of the Poisson system.

Definition at line 173 of file navierstokes.h.

◆ solution

petibm::type::Solution NavierStokesSolver::solution
protected

Data object holding the velocity and pressure fields.

Definition at line 95 of file navierstokes.h.

◆ solversViewer

PetscViewer NavierStokesSolver::solversViewer
protected

ASCII PetscViewer object to output solvers info.

Definition at line 209 of file navierstokes.h.

◆ stageInitialize

PetscLogStage NavierStokesSolver::stageInitialize
protected

Log stage for the initialization phase.

Definition at line 185 of file navierstokes.h.

◆ stageMonitor

PetscLogStage NavierStokesSolver::stageMonitor
protected

Log stage for monitoring user-defined regions of the domain.

Definition at line 206 of file navierstokes.h.

◆ stageRHSPoisson

PetscLogStage NavierStokesSolver::stageRHSPoisson
protected

Log stage for assembling the RHS of the Poisson system.

Definition at line 194 of file navierstokes.h.

◆ stageRHSVelocity

PetscLogStage NavierStokesSolver::stageRHSVelocity
protected

Log stage for assembling the RHS of the velocity system.

Definition at line 188 of file navierstokes.h.

◆ stageSolvePoisson

PetscLogStage NavierStokesSolver::stageSolvePoisson
protected

Log stage for solving the Poisson system.

Definition at line 197 of file navierstokes.h.

◆ stageSolveVelocity

PetscLogStage NavierStokesSolver::stageSolveVelocity
protected

Log stage for solving the velocity system.

Definition at line 191 of file navierstokes.h.

◆ stageUpdate

PetscLogStage NavierStokesSolver::stageUpdate
protected

Log stage for updating field variables.

Definition at line 200 of file navierstokes.h.

◆ stageWrite

PetscLogStage NavierStokesSolver::stageWrite
protected

Log stage when write the solution fields.

Definition at line 203 of file navierstokes.h.

◆ t

PetscReal NavierStokesSolver::t
protected

Time value.

Definition at line 119 of file navierstokes.h.

◆ vSolver

petibm::type::LinSolver NavierStokesSolver::vSolver
protected

Velocity linear solver.

Definition at line 104 of file navierstokes.h.


The documentation for this class was generated from the following files: