mus_particle_config_module Module

mus_particle_config_module contains routines for loading the data for LBM-DEM simulations of particles in a flow from the lua script. All particle data is read from the "particles" table in the lua script. This table should look as follows: particles = { nParticles = 1, -- Total number of particles in the simulation -- particle kind, can be: -- * 'MEM' for fully-resolved particles using the Momentum-Exchange Method -- * 'DPS' for four-way coupled unresolved particles based on the -- Volume-Averaged Navier-Stokes equations -- * 'DPS_twoway' for two-way coupled unresolved particles, neglecting the -- effect of volume fraction -- * 'DPS_oneway' for one-way coupled unresolved particles. That is, particles -- do not affect the flow kind = 'DPS',

-- Describe how particles should interact with boundaries. They can either bounce -- off them ('wall') or have periodic boundaries ('periodic'). Currently only -- interactions with prismatic boundaries are supported. If not specified, -- particles will be removed from the domain once they hit a boundary of the -- fluid domain (also works for arbitrary geometries). boundaries = { domain_bnd = { 0.0, L, 0.0, W, 0.0, H }, bnd_kind = {'wall', 'wall', 'wall', 'wall', 'wall', 'wall'} },

-- Number of DEM subcycles per LBM time step to use nDEMsubcycles = 50,

-- How often to log the particle data (1 = log every LBM time step) particleLogInterval = 1,

-- Size of the buffers in the particle communication routines. Should be -- sufficient to hold all particles anticipated to be on a domain boundary -- at any time during the simulation. Safest option is to make this equal -- to or greater than nParticles. particleBufferSize = 100,

-- Collision time (in physical units) of particle-particle collisions particle_collision_time = 0.1*dt,

-- Collision tolerance, i.e. what should the gap between particles be -- before we call it a collision particle_collision_tol = 0.0,

-- Size of the gap before capping the calculation of lubrication forces critical_gap_Flub = 0.01,

-- Particle positions. Should contain nParticles tables of the form -- {x,y,z,rx,ry,rz} with x, y, z the translational positions and -- rx, ry, rz the rotation angles about the corresponding axes. position = { { x_p_phy, y_p_phy, z_p_phy, 0.0, 0.0, 0.0} },

-- Particle initial velocities. Should contain nParticles tables of the form -- {ux, uy, uz, urx, ury, urz} velocity = { {0.0, 0.0, 0.0, 0.0, 0.0, 0.0} },

-- External forces on particles. Should contain nParticles tables of the form -- {Fx, Fy, Fz, Frx, Fry, Frz} force = { {0.0, 0.0, -F_gravity + F_buoyancy, 0.0, 0.0, 0.0} },

-- Particle radii, should have nParticles entries radius = { 0.5*Dia_p_phy },

-- Particle masses, should have nParticles entries mass = { m_p_phy } }



Functions

public function check_particle_scheme_kind_compatibility(particle_kind, scheme_kind) result(is_compatible)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: particle_kind
character(len=*), intent(in) :: scheme_kind

Return Value logical


Subroutines

public subroutine mus_load_particle_creator_timing(conf, p_thandle, particle_creator, Nparticles, chunkSize, scheme, geometry, myRank)

Arguments

Type IntentOptional Attributes Name
type(flu_State) :: conf

configuration

integer, intent(in) :: p_thandle

Handle to particle table

type(mus_particle_creator_type), intent(inout) :: particle_creator

Particle creator object

integer, intent(in) :: Nparticles

Number of particles in lua table

integer, intent(in) :: chunkSize

Size of the chunks of particle data to read from the lua script

type(mus_scheme_type), intent(in) :: scheme

Scheme

type(mus_geom_type), intent(in) :: geometry

Geometry

integer, intent(in) :: myRank

This MPI process rank

public subroutine mus_load_particle_interpolation(conf, p_thandle, layout, interpolator)

Arguments

Type IntentOptional Attributes Name
type(flu_State) :: conf

configuration

integer, intent(in) :: p_thandle

Handle to parent table

character(len=labelLen) :: layout

Layout label, e.g. d3q19, d2q9

type(mus_particle_interpolator_type), intent(inout) :: interpolator

tracker to initialize values of

public subroutine mus_load_debugtracker(conf, p_thandle, tracker)

Arguments

Type IntentOptional Attributes Name
type(flu_State) :: conf

configuration

integer, intent(in) :: p_thandle

Handle to parent table

type(mus_particle_debugtracking_type), intent(inout) :: tracker

tracker to initialize values of

public subroutine mus_load_particle_boundaries(conf, p_thandle, bndData)

Arguments

Type IntentOptional Attributes Name
type(flu_State) :: conf

configuration

integer, intent(in) :: p_thandle

Handle to parent table

type(mus_particle_boundarydata_type), intent(inout) :: bndData

Particle boundary data

public subroutine mus_load_particlekind(particle_kind, conf)

Get the particle kind from the configuration

Read more…

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(inout) :: particle_kind

particleGroup to add particles loaded from the lua script to

type(flu_State) :: conf

configuration

public subroutine mus_load_particle_collisions(particleGroup, conf, p_thandle)

Arguments

Type IntentOptional Attributes Name
type(mus_particle_group_type), intent(inout) :: particleGroup

particleGroup to add particles loaded from the lua script to

type(flu_State) :: conf

configuration

integer, intent(in) :: p_thandle

Handle to particle table (which should be opened before calling this routine)

public subroutine mus_load_particles(particleGroup, particle_kind, conf, chunkSize, scheme, geometry, myRank)

Routine to load the particle data from the musubi.lua file meant to be called from within mus_load_config after the lua file has already been opened

Arguments

Type IntentOptional Attributes Name
type(mus_particle_group_type), intent(inout) :: particleGroup

particleGroup to add particles loaded from the lua script to

character(len=*), intent(in) :: particle_kind

Kind of particle modelling to use

type(flu_State) :: conf

configuration

integer, intent(in) :: chunkSize

Size of the number of particles to be read in one chunk from the lua table

type(mus_scheme_type), intent(in) :: scheme

Scheme to determine if particles belong on this process

type(mus_geom_type), intent(in) :: geometry

Geometry to determine if particles belong on this process

integer, intent(in) :: myRank

This process rank

public subroutine load_particle_mem_data(conf, parent, particles, nparticles, chunksize)

Arguments

Type IntentOptional Attributes Name
type(flu_State), intent(in) :: conf
integer, intent(in), optional :: parent

handle to parent table if position, velocity tables etc are inside another table

type(dyn_particle_mem_array_type), intent(inout) :: particles

dynamic particle array to append particles read from file to

integer, intent(in) :: nparticles

total number of particles to read

integer, intent(in) :: chunksize

public subroutine load_particle_mem_creator_data(conf, parent, particle_creator, nparticles, chunksize, scheme, geometry, myrank)

Arguments

Type IntentOptional Attributes Name
type(flu_State), intent(in) :: conf
integer, intent(in), optional :: parent

handle to parent table if position, velocity tables etc are inside another table

type(mus_particle_creator_type), intent(inout) :: particle_creator

dynamic particle array to append particles read from file to

integer, intent(in) :: nparticles

total number of particles to read

integer, intent(in) :: chunksize
type(mus_scheme_type), intent(in) :: scheme

scheme

type(mus_geom_type), intent(in) :: geometry

geometry

integer, intent(in) :: myrank

this mpi process rank

public subroutine mus_load_particle_mem_data_chunk(conf, parent, particlebuffer, nchunkvals, key, istart)

Arguments

Type IntentOptional Attributes Name
type(flu_State), intent(in) :: conf
integer, intent(in), optional :: parent

handle to parent table if position, velocity tables etc are inside another table

type(mus_particle_MEM_type), intent(inout), allocatable :: particlebuffer(:)
integer, intent(out) :: nchunkvals
character(len=*), intent(in) :: key
integer, intent(in) :: istart

public subroutine load_particle_dps_data(conf, parent, particles, nparticles, chunksize)

Arguments

Type IntentOptional Attributes Name
type(flu_State), intent(in) :: conf
integer, intent(in), optional :: parent

handle to parent table if position, velocity tables etc are inside another table

type(dyn_particle_dps_array_type), intent(inout) :: particles

dynamic particle array to append particles read from file to

integer, intent(in) :: nparticles

total number of particles to read

integer, intent(in) :: chunksize

public subroutine load_particle_dps_creator_data(conf, parent, particle_creator, nparticles, chunksize, scheme, geometry, myrank)

Arguments

Type IntentOptional Attributes Name
type(flu_State), intent(in) :: conf
integer, intent(in), optional :: parent

handle to parent table if position, velocity tables etc are inside another table

type(mus_particle_creator_type), intent(inout) :: particle_creator

dynamic particle array to append particles read from file to

integer, intent(in) :: nparticles

total number of particles to read

integer, intent(in) :: chunksize
type(mus_scheme_type), intent(in) :: scheme

scheme

type(mus_geom_type), intent(in) :: geometry

geometry

integer, intent(in) :: myrank

this mpi process rank

public subroutine mus_load_particle_dps_data_chunk(conf, parent, particlebuffer, nchunkvals, key, istart)

Arguments

Type IntentOptional Attributes Name
type(flu_State), intent(in) :: conf
integer, intent(in), optional :: parent

handle to parent table if position, velocity tables etc are inside another table

type(mus_particle_DPS_type), intent(inout), allocatable :: particlebuffer(:)
integer, intent(out) :: nchunkvals
character(len=*), intent(in) :: key
integer, intent(in) :: istart

public subroutine mus_finalize_particleGroup(particleGroup)

Arguments

Type IntentOptional Attributes Name
type(mus_particle_group_type), intent(inout) :: particleGroup

public subroutine mus_particles_print_config(particleGroup, logUnit)

Arguments

Type IntentOptional Attributes Name
type(mus_particle_group_type), intent(in) :: particleGroup
integer, intent(in) :: logUnit

public subroutine mus_load_predefined_particleblob(conf, parent, blob_type, flag)

Arguments

Type IntentOptional Attributes Name
type(flu_State) :: conf

configuration

integer, intent(in) :: parent

Handle to parent table in which the "shape" table exists

character(len=LabelLen), intent(inout) :: blob_type

Output string containing the blob type

logical, intent(out) :: flag

Output flag set to TRUE when reading was successful, FALSE if not

public subroutine mus_load_particle_distribution(distribution, conf, parent)

Arguments

Type IntentOptional Attributes Name
type(mus_particle_blob_prob_type), intent(inout) :: distribution
type(flu_State) :: conf

configuration

integer, intent(in) :: parent

Handle to parent table in which the "shape" table exists

public subroutine mus_load_particleblob_cylinder(particleblob, conf, parent)

Arguments

Type IntentOptional Attributes Name
type(mus_particle_blob_cylinder_type), intent(inout) :: particleblob
type(flu_State) :: conf

configuration

integer, intent(in) :: parent

Handle to parent table in which the "shape" table exists

public subroutine mus_load_particleblob_prism(particleblob, conf, parent)

Arguments

Type IntentOptional Attributes Name
type(mus_particle_blob_prism_type), intent(inout) :: particleblob
type(flu_State) :: conf

configuration

integer, intent(in) :: parent

Handle to parent table in which the "shape" table exists