cdifs_obj Derived Type

type, public, extends(solver_obj) :: cdifs_obj

Constant Density Incompressible Flow Solver


Inherits

type~~cdifs_obj~~InheritsGraph type~cdifs_obj cdifs_obj type~bc_set bc_set type~cdifs_obj->type~bc_set bcs type~block_obj block_obj type~cdifs_obj->type~block_obj block type~bodyforce_obj bodyforce_obj type~cdifs_obj->type~bodyforce_obj bodyforce type~collision_obj collision_obj type~cdifs_obj->type~collision_obj collisions type~eulerian_obj_i eulerian_obj_i type~cdifs_obj->type~eulerian_obj_i maskV type~eulerian_obj_r eulerian_obj_r type~cdifs_obj->type~eulerian_obj_r V, P, dP, ibS, ibVF, ibF, ibN, Vold, resV, rhs, divu, Vm, srcV type~eulerian_set eulerian_set type~cdifs_obj->type~eulerian_set fields type~hdf5_obj hdf5_obj type~cdifs_obj->type~hdf5_obj hdf5 type~hypre_obj hypre_obj type~cdifs_obj->type~hypre_obj hypre, VFSolver type~marker_set marker_set type~cdifs_obj->type~marker_set IB type~monitor_set monitor_set type~cdifs_obj->type~monitor_set monitors, pmonitor type~op_obj op_obj type~cdifs_obj->type~op_obj op type~respart_set ResPart_set type~cdifs_obj->type~respart_set RP type~solver_obj solver_obj type~cdifs_obj->type~solver_obj type~bc_set->type~block_obj block type~bc_set->type~hdf5_obj hdf5 type~hashtbl_obj hashtbl_obj type~bc_set->type~hashtbl_obj tbl type~parallel_obj parallel_obj type~bc_set->type~parallel_obj parallel type~region_obj region_obj type~bc_set->type~region_obj region type~block_obj->type~hdf5_obj hdf5 MPI_Datatype MPI_Datatype type~block_obj->MPI_Datatype gc_slab_r, gc_slab_i type~axis_obj axis_obj type~block_obj->type~axis_obj axis, axis_partition type~block_obj->type~parallel_obj parallel type~collision_obj->type~block_obj cblock type~collision_obj->type~marker_set IB type~collision_obj->type~monitor_set monitors type~collision_obj->type~respart_set RP type~collision_obj->type~parallel_obj parallel type~parser_obj parser_obj type~collision_obj->type~parser_obj parser type~particle_set particle_set type~collision_obj->type~particle_set PP type~sllist_obj sllist_obj type~collision_obj->type~sllist_obj RPneighbors, PPneighbors, IBneighbors type~timer_obj timer_obj type~collision_obj->type~timer_obj timer type~eulerian_obj_base eulerian_obj_base type~eulerian_obj_i->type~eulerian_obj_base type~eulerian_obj_r->type~eulerian_obj_base type~eulerian_set->type~block_obj block type~eulerian_ptr eulerian_ptr type~eulerian_set->type~eulerian_ptr field type~eulerian_set->type~hashtbl_obj tbl type~eulerian_set->type~parallel_obj parallel type~hdf5_obj->type~hashtbl_obj tbl type~hdf5_obj->type~parallel_obj parallel type~hypre_obj->type~block_obj block type~hypre_obj->type~eulerian_obj_i irow MPI_Comm MPI_Comm type~hypre_obj->MPI_Comm comm c_ptr c_ptr type~hypre_obj->c_ptr p_rhs_values, p_sol_values, p_values, p_rows, p_cols, p_tmpi type~hypre_obj->type~parallel_obj parallel type~marker_set->type~bc_set bcs type~marker_set->type~monitor_set monitors type~marker_set->type~op_obj op type~lagrangian_set lagrangian_set type~marker_set->type~lagrangian_set type~marker_set->type~parser_obj parser type~marker_set->type~timer_obj timer type~monitor_set->type~hashtbl_obj tbl type~monitor_obj monitor_obj type~monitor_set->type~monitor_obj m type~monitor_set->type~parallel_obj parallel type~op_obj->type~block_obj block type~op_obj->type~parallel_obj parallel type~respart_set->type~bc_set bcs type~respart_set->type~marker_set ib type~respart_set->type~monitor_set monitors type~respart_set->type~op_obj op type~respart_set->type~lagrangian_set type~respart_set->type~parser_obj parser type~respart_set->type~timer_obj timer type~solver_obj->type~parallel_obj parallel type~solver_obj->type~parser_obj parser type~solver_obj->type~timer_obj timer type~eulerian_obj_base->type~block_obj block type~eulerian_obj_base->type~parallel_obj parallel type~eulerian_ptr->type~eulerian_obj_base p type~hashtbl_obj->type~sllist_obj vec type~lagrangian_set->type~block_obj block type~lagrangian_set->MPI_Datatype MPI_TYPE type~lagrangian_set->type~parallel_obj parallel type~lagrangian_obj lagrangian_obj type~lagrangian_set->type~lagrangian_obj p, sample type~column_obj column_obj type~monitor_obj->type~column_obj col type~parallel_obj->MPI_Datatype REAL_SP, REAL_DP, REAL_WP, COMPLEX_SP, COMPLEX_DP, COMPLEX_WP, INTEGER, INT8, LOGICAL type~communicators communicators type~parallel_obj->type~communicators comm type~patch patch type~parallel_obj->type~patch rank type~entry_obj entry_obj type~parser_obj->type~entry_obj entries type~particle_set->type~monitor_set monitors type~particle_set->type~op_obj op type~particle_set->type~lagrangian_set type~particle_set->type~parser_obj parser type~particle_set->type~timer_obj timer type~region_obj->type~block_obj region type~region_obj->type~hashtbl_obj tbl type~bc_obj bc_obj type~region_obj->type~bc_obj BC type~sllist_obj->type~sllist_obj child type~timer_obj->type~hashtbl_obj tbl type~timer_obj->type~parallel_obj parallel type~timer_obj->type~parser_obj parser type~communicators->MPI_Comm w, g

Components

Type Visibility Attributes Name Initial
real(kind=wp), public, allocatable :: DIV(:,:,:,:,:)

Divergence operator for velocity

type(marker_set), public :: IB

Immersed boundaries

real(kind=wp), public, allocatable :: LAP(:,:,:,:,:,:)

Viscous Laplacian for velocity

type(eulerian_obj_r), public :: P

Fluid pressure

type(ResPart_set), public :: RP

Resolved particles

type(eulerian_obj_r), public :: V(3)

Fluid velocity fields

type(hypre_obj), public :: VFSolver

HYPRE Solver for the volume-fraction equation

type(eulerian_obj_r), public :: Vm(3)

Velocities at mid points

type(eulerian_obj_r), public :: Vold(3)

Old fluid velocity (at n)

type(bc_set), public :: bcs

Boundary conditions

type(block_obj), public :: block

Block information

class(bodyforce_obj), public, allocatable :: bodyforce

Bodyforce

type(collision_obj), public :: collisions

Utility that handles collisions

type(eulerian_obj_r), public :: dP

Pressure correction

type(eulerian_obj_r), public :: divu

Divergence field

type(eulerian_set), public :: fields

Eulerian data container

real(kind=wp), public :: gravity(3) = 0.0_wp

Gravity

type(hdf5_obj), public, pointer :: hdf5 => null()

Associated hdf5 object

type(hypre_obj), public :: hypre

HYPRE Solver for the pressure-Poisson equation

type(eulerian_obj_r), public :: ibF(3)

Immersed boundary force

type(eulerian_obj_r), public :: ibN(3)

IB normals field

type(eulerian_obj_r), public :: ibS

Surface density function

type(eulerian_obj_r), public :: ibVF

Solid volume fraction

type(eulerian_obj_i), public :: maskV(3)

Velocity field masks

integer, public :: maxitp = 1

Collision substeps

type(monitor_set), public :: monitors

Monitors to print to stdout and files at end of each time step

character(len=:), public, allocatable :: name

Solver's name.

type(op_obj), public :: op

Differential Operators

character(len=str8), public, allocatable :: output_var(:)

Names of variables to output

character(len=str8), public, allocatable :: outputs(:)

List of outputs (one output can contain multiple variables)

real(kind=wp), public, allocatable :: pGRAD(:,:,:,:,:)

Gradient operator for pressure

type(parallel_obj), public, pointer :: parallel => null()

Associated parallel structure

type(parser_obj), public, pointer :: parser => null()

Associated parser

type(monitor_set), public :: pmonitor

Monitor for the pressure-Poisson equation, prints at each substep

type(eulerian_obj_r), public :: resV(3)

Velocity residual

real(kind=wp), public :: rho = 1.0_wp

Fluid density

type(eulerian_obj_r), public :: rhs

Right-hand side

logical, public :: singularVF = .false.

Flag for singular volume fraction Poisson equation

type(eulerian_obj_r), public :: srcV(3)

Momentum source

integer, public :: subit = 2

Solver sub-iteration

type(timer_obj), public, pointer :: timer => null()

Associated timer

logical, public :: use_IB = .false.

Use Immersed Boundaries

logical, public :: use_RP = .false.

Use Resolved Particles

logical, public :: use_col = .false.

Use collisions

real(kind=wp), public :: visc = 1.0_wp

Fluid viscosity


Type-Bound Procedures

procedure, public :: AdvanceSolution => cdifs_obj_AdvanceSolution

  • interface

    private module subroutine cdifs_obj_AdvanceSolution(this)

    Advances solution from n to n+1.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

procedure, public :: ComputeSolidVF => cdifs_obj_ComputeSolidVF

  • private subroutine cdifs_obj_ComputeSolidVF(this)

    Computes the solid volume fraction (ibVF) by solving a Poisson equation. This subroutine computes the IB normals fields, applies approriate boundary conditions, and takes the divergence, which forms the right-hand side of the poisson equation. It then uses a HYPRE solver to solve this equation. In case where the Poisson operator is singular and using resolved particles, this subroutine will apply a correction that sets the solid volume fraction at particle centroids.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

procedure, public :: Finalize => cdifs_obj_Final

  • private subroutine cdifs_obj_Final(this)

    Finalizes solver and frees memory.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

procedure, public :: Initialize => cdifs_obj_Init

  • private subroutine cdifs_obj_Init(this, timer, parallel, parser)

    Initializes the CDIFS solver.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

    type(timer_obj), intent(in), target :: timer

    Timer utility

    type(parallel_obj), intent(in), target :: parallel

    Parallel machinery

    type(parser_obj), intent(in), target :: parser

    Parser for input file

procedure, public :: LinkHDF5Object => cdifs_obj_LinkHDF5Object

  • private subroutine cdifs_obj_LinkHDF5Object(this, hdf5)

    Links HDF5 Object.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

    type(hdf5_obj), intent(in), target :: hdf5

    HDF5 utility

procedure, public :: Monitor => cdifs_obj_Monitor

  • interface

    private module subroutine cdifs_obj_Monitor(this)

    Analyzes data and posts to monitors.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

procedure, public :: PrepareSolver => cdifs_obj_PrepareSolver

  • interface

    private module subroutine cdifs_obj_PrepareSolver(this)

    Prepares data for run.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

procedure, public :: WriteOutputData => cdifs_obj_WriteOutputData

  • interface

    private module subroutine cdifs_obj_WriteOutputData(this)

    Processes and outputs single-precision data for visualization.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

procedure, public :: WriteRestartData => cdifs_obj_WriteRestartData

  • interface

    private module subroutine cdifs_obj_WriteRestartData(this)

    Writes restart data to disk.

    Arguments

    Type IntentOptional Attributes Name
    class(cdifs_obj), intent(inout) :: this

    CDIFS solver

Source Code

  type, extends(solver_obj):: cdifs_obj
    !> Constant Density Incompressible Flow Solver
    type(block_obj)                   :: block                                 !! Block information
    type(op_obj)                      :: op                                    !! Differential Operators
    type(bc_set)                      :: bcs                                   !! Boundary conditions
    type(hypre_obj)                   :: hypre                                 !! HYPRE Solver for the pressure-Poisson equation
    type(hypre_obj)                   :: VFSolver                              !! HYPRE Solver for the volume-fraction equation
    real(wp),             allocatable :: pGRAD(:,:,:,:,:)                      !! Gradient operator for pressure
    real(wp),             allocatable :: DIV(:,:,:,:,:)                        !! Divergence operator for velocity
    real(wp),             allocatable :: LAP(:,:,:,:,:,:)                      !! Viscous Laplacian for velocity
    type(Eulerian_set)                :: fields                                !! Eulerian data container
    type(Eulerian_obj_r)              :: V(3)                                  !! Fluid velocity fields
    type(Eulerian_obj_r)              :: P                                     !! Fluid pressure
    type(Eulerian_obj_r)              :: dP                                    !! Pressure correction
    logical                           :: use_RP         = .false.              !! Use Resolved Particles
    logical                           :: use_IB         = .false.              !! Use Immersed Boundaries
    type(ResPart_set)                 :: RP                                    !! Resolved particles
    type(marker_set)                  :: IB                                    !! Immersed boundaries
    type(Eulerian_obj_r)              :: ibS                                   !! Surface density function
    type(Eulerian_obj_r)              :: ibVF                                  !! Solid volume fraction
    type(Eulerian_obj_r)              :: ibF(3)                                !! Immersed boundary force
    type(Eulerian_obj_r)              :: ibN(3)                                !! IB normals field
    logical                           :: singularVF    = .false.               !! Flag for singular volume fraction Poisson equation
    type(Eulerian_obj_r)              :: Vold(3)                               !! Old fluid velocity (at n)
    type(Eulerian_obj_r)              :: resV(3)                               !! Velocity residual
    type(Eulerian_obj_r)              :: rhs                                   !! Right-hand side
    type(Eulerian_obj_r)              :: divu                                  !! Divergence field
    type(Eulerian_obj_r)              :: Vm(3)                                 !! Velocities at mid points
    class(bodyforce_obj), allocatable :: bodyforce                             !! Bodyforce
    type(Eulerian_obj_r)              :: srcV(3)                               !! Momentum source
    type(Eulerian_obj_i)              :: maskV(3)                              !! Velocity field masks
    type(hdf5_obj),           pointer :: hdf5          => null()               !! Associated hdf5 object
    type(monitor_set)                 :: monitors                              !! Monitors to print to stdout and files at end of each time step
    type(monitor_set)                 :: pmonitor                              !! Monitor for the pressure-Poisson equation, prints at each substep
    character(len=str8),  allocatable :: outputs(:)                            !! List of outputs (one output can contain multiple variables)
    character(len=str8),  allocatable :: output_var(:)                         !! Names of variables to output
    real(wp)                          :: rho            = 1.0_wp               !! Fluid density
    real(wp)                          :: visc           = 1.0_wp               !! Fluid viscosity
    real(wp)                          :: gravity(3)     = 0.0_wp               !! Gravity
    integer                           :: subit          = 2                    !! Solver sub-iteration
    logical                           :: use_col        = .false.              !! Use collisions
    type(collision_obj)               :: collisions                            !! Utility that handles collisions
    integer                           :: maxitp         = 1                    !! Collision substeps
    contains
      procedure :: Initialize         => cdifs_obj_Init
      procedure :: Finalize           => cdifs_obj_Final
      procedure :: PrepareSolver      => cdifs_obj_PrepareSolver
      procedure :: AdvanceSolution    => cdifs_obj_AdvanceSolution
      procedure :: Monitor            => cdifs_obj_Monitor
      procedure :: WriteRestartData   => cdifs_obj_WriteRestartData
      procedure :: WriteOutputData    => cdifs_obj_WriteOutputData
      procedure :: LinkHDF5Object     => cdifs_obj_LinkHDF5Object
      procedure :: ComputeSolidVF     => cdifs_obj_ComputeSolidVF
  end type cdifs_obj