grans_obj Derived Type

type, public, extends(solver_obj) :: grans_obj

Granular flow solver with resolved particles


Inherits

type~~grans_obj~~InheritsGraph type~grans_obj grans_obj type~bc_set bc_set type~grans_obj->type~bc_set bcs type~block_obj block_obj type~grans_obj->type~block_obj block type~collision_obj collision_obj type~grans_obj->type~collision_obj collisions type~eulerian_obj_r eulerian_obj_r type~grans_obj->type~eulerian_obj_r ibVF, PVF, ibS, Fp, ibF, ibN, rhs type~eulerian_set eulerian_set type~grans_obj->type~eulerian_set fields type~hdf5_obj hdf5_obj type~grans_obj->type~hdf5_obj hdf5 type~hypre_obj hypre_obj type~grans_obj->type~hypre_obj VFSolver type~marker_set marker_set type~grans_obj->type~marker_set IB type~monitor_set monitor_set type~grans_obj->type~monitor_set monitors type~op_obj op_obj type~grans_obj->type~op_obj op type~particle_set particle_set type~grans_obj->type~particle_set PP type~respart_set ResPart_set type~grans_obj->type~respart_set RP type~solver_obj solver_obj type~grans_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~particle_set PP 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~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_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 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~eulerian_obj_i eulerian_obj_i type~hypre_obj->type~eulerian_obj_i irow 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~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~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_obj_i->type~eulerian_obj_base 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~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
type(eulerian_obj_r), public :: Fp(3)

Particle Eulerian forcing

type(marker_set), public :: IB

Immersed solids (walls)

type(particle_set), public :: PP

Point particles

type(eulerian_obj_r), public :: PVF

Particle volume fraction

type(ResPart_set), public :: RP

Resolved particles

type(hypre_obj), public :: VFSolver

HYPRE Solver for the volume-fraction equation

integer, public :: VF_MaxIT

VF iterations

real(kind=wp), public :: VF_MaxTol

VF maximum relative tolerance

integer, public :: VF_it

Sub-iterations for VF solves

real(kind=wp), public :: VF_rel

Relative error at end of VF solves

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

Name of VF solver to use

type(bc_set), public :: bcs

Boundary conditions

type(block_obj), public :: block

Block information

type(collision_obj), public :: collisions

Collision utility

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(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 Area Density

type(eulerian_obj_r), public :: ibVF

Solid volume fraction (walls)

real(kind=wp), public :: intRHS

Magnitude of RHS in VF equation

type(monitor_set), public :: monitors

Monitors to print to stdout and files

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)

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

Associated parallel structure

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

Associated parser

real(kind=wp), public :: rel_error_max = 1.0e-6_wp

Maximum relative error

type(eulerian_obj_r), public :: rhs

Right-hand side

logical, public :: singularVF = .false.

Flag for singular volume fraction Poisson equation

integer, public :: subit = 1

Solver sub-iteration

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

Associated timer

logical, public :: use_IB = .false.

Use Immersed Boundaries

logical, public :: use_PP = .false.

Use Point Particles

logical, public :: use_RP = .false.

Use Resolved Particles

logical, public :: use_col = .false.

Collision flag


Type-Bound Procedures

procedure, public :: AdvanceSolution => grans_obj_AdvanceSolution

  • interface

    private module subroutine grans_obj_AdvanceSolution(this)

    Advances solution from n to n+1.

    Arguments

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

    GRANS solver

procedure, public :: ComputeSolidVF => grans_obj_ComputeSolidVF

  • private subroutine grans_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(grans_obj), intent(inout) :: this

    GRANS solver

procedure, public :: Finalize => grans_obj_Final

  • private subroutine grans_obj_Final(this)

    Finalizes solver and frees memory.

    Arguments

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

    GRANS solver

procedure, public :: Initialize => grans_obj_Init

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

    Initializes the GRANS solver.

    Arguments

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

    GRANS 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 => grans_obj_LinkHDF5Object

  • private subroutine grans_obj_LinkHDF5Object(this, hdf5)

    Links HDF5 Object.

    Arguments

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

    GRANS solver

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

    HDF5 utility

procedure, public :: Monitor => grans_obj_Monitor

  • interface

    private module subroutine grans_obj_Monitor(this)

    Analyzes data and posts to monitors.

    Arguments

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

    GRANS solver

procedure, public :: PrepareSolver => grans_obj_PrepareSolver

  • interface

    private module subroutine grans_obj_PrepareSolver(this)

    Prepares data for run.

    Arguments

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

    GRANS solver

procedure, public :: WriteOutputData => grans_obj_WriteOutputData

  • interface

    private module subroutine grans_obj_WriteOutputData(this)

    Processes and outputs single-precision data for visualization.

    Arguments

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

    GRANS solver

procedure, public :: WriteRestartData => grans_obj_WriteRestartData

  • interface

    private module subroutine grans_obj_WriteRestartData(this)

    Writes restart data to disk.

    Arguments

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

    GRANS solver

Source Code

  type, extends(solver_obj):: grans_obj
    !> Granular flow solver with resolved particles
    type(block_obj)                  :: block                                  !! Block information
    type(op_obj)                     :: op                                     !! Differential Operators
    type(bc_set)                     :: bcs                                    !! Boundary conditions
    type(hypre_obj)                  :: VFSolver                               !! HYPRE Solver for the volume-fraction equation
    logical                          :: singularVF    = .false.                !! Flag for singular volume fraction Poisson equation
    type(particle_set)               :: PP                                     !! Point particles
    type(ResPart_set)                :: RP                                     !! Resolved particles
    type(marker_set)                 :: IB                                     !! Immersed solids (walls)
    type(Eulerian_set)               :: fields                                 !! Eulerian data container
    type(Eulerian_obj_r)             :: ibVF                                   !! Solid volume fraction (walls)
    type(Eulerian_obj_r)             :: PVF                                    !! Particle volume fraction
    type(Eulerian_obj_r)             :: ibS                                    !! Surface Area Density
    type(Eulerian_obj_r)             :: Fp(3)                                  !! Particle Eulerian forcing
    type(Eulerian_obj_r)             :: ibF(3)                                 !! Immersed boundary force
    type(Eulerian_obj_r)             :: ibN(3)                                 !! IB normals field
    type(Eulerian_obj_r)             :: rhs                                    !! Right-hand side
    real(wp)                         :: gravity(3)    = 0.0_wp                 !! Gravity
    integer                          :: subit         = 1                      !! Solver sub-iteration
    real(wp)                         :: rel_error_max = 1.0e-6_wp              !! Maximum relative error
    logical                          :: use_RP        = .false.                !! Use Resolved Particles
    logical                          :: use_PP        = .false.                !! Use Point Particles
    logical                          :: use_IB        = .false.                !! Use Immersed Boundaries
    logical                          :: use_col       = .false.                !! Collision flag
    type(collision_obj)              :: collisions                             !! Collision utility
    type(hdf5_obj),          pointer :: hdf5         => null()                 !! Associated hdf5 object
    type(monitor_set)                :: monitors                               !! Monitors to print to stdout and files
    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
    character(len=:),    allocatable :: VF_solver                              !! Name of VF solver to use
    real(wp)                         :: VF_MaxTol                              !! VF maximum relative tolerance
    integer                          :: VF_MaxIT                               !! VF iterations
    integer                          :: VF_it                                  !! Sub-iterations for VF solves
    real(wp)                         :: VF_rel                                 !! Relative error at end of VF solves
    real(wp)                         :: intRHS                                 !! Magnitude of RHS in VF equation
    contains
      procedure :: Initialize           => grans_obj_Init
      procedure :: Finalize             => grans_obj_Final
      procedure :: PrepareSolver        => grans_obj_PrepareSolver
      procedure :: AdvanceSolution      => grans_obj_AdvanceSolution
      procedure :: Monitor              => grans_obj_Monitor
      procedure :: WriteRestartData     => grans_obj_WriteRestartData
      procedure :: WriteOutputData      => grans_obj_WriteOutputData
      procedure :: LinkHDF5Object       => grans_obj_LinkHDF5Object
      procedure :: ComputeSolidVF       => grans_obj_ComputeSolidVF
  end type grans_obj