My Project
programmer's documentation
Functions
cs_lagr_prototypes.h File Reference
#include "cs_base.h"
#include "cs_mesh.h"
#include "cs_mesh_quantities.h"
#include "cs_mesh_bad_cells.h"
#include "cs_domain.h"
#include "cs_lagr.h"
#include "cs_lagr_tracking.h"
#include "cs_lagr_stat.h"
Include dependency graph for cs_lagr_prototypes.h:

Go to the source code of this file.

Functions

void cs_user_lagr_ef (cs_real_t dt_p, const cs_real_t taup[], const cs_real_3_t tlag[], const cs_real_3_t piil[], const cs_real_33_t bx[], const cs_real_t tsfext[], const cs_real_33_t vagaus[], const cs_real_3_t gradpr[], const cs_real_33_t gradvf[], cs_real_t rho_p[], cs_real_3_t fextla[])
 User definition of an external force field acting on the particles. More...
 
void cs_user_lagr_in (cs_lagr_particle_set_t *particles, const cs_lagr_injection_set_t *zis, const cs_lnum_t particle_range[2], const cs_lnum_t particle_face_id[], const cs_real_t visc_length[])
 User modification of newly injected particles. More...
 
void cs_user_lagr_model (void)
 
void cs_user_lagr_rt (cs_lnum_t id_p, cs_real_t re_p, cs_real_t uvwr, cs_real_t rho_f, cs_real_t rho_p, cs_real_t nu_f, cs_real_t taup[], const cs_real_t dt[])
 Modification of the calculation of the particle relaxation time with respect to the chosen formulation for the drag coefficient. More...
 
void cs_user_lagr_rt_t (cs_lnum_t id_p, cs_real_t re_p, cs_real_t uvwr, cs_real_t rho_f, cs_real_t rho_p, cs_real_t nu_f, cs_real_t cp_f, cs_real_t k_f, cs_real_t tauc[], const cs_real_t dt[])
 Modification of the calculation of the thermal relaxation time of the particles with respect to the chosen formulation of the Nusselt number. More...
 
void cs_user_lagr_imposed_motion (const cs_real_t coords[3], const cs_real_t dt, cs_real_t disp[3])
 Impose the motion of a particle flagged CS_LAGR_PART_IMPOSED_MOTION. More...
 
void cs_user_lagr_extra_operations (const cs_real_t dt[])
 User function (non-mandatory intervention) More...
 
void cs_user_lagr_sde (const cs_real_t dt[], cs_real_t taup[], cs_real_3_t tlag[], cs_real_t tempct[])
 User integration of the SDE for the user-defined variables. More...
 
void cs_user_lagr_boundary_conditions (const int itypfb[])
 Define particle boundary conditions. More...
 
void cs_lagr_user_boundary_interaction (cs_lagr_particle_set_t *particles, cs_lnum_t p_id, cs_lnum_t face_id, const cs_real_t face_norm[3], const cs_real_t c_intersect[3], cs_real_t t_intersect, int b_zone_id, int *event_flag, cs_lagr_tracking_state_t *tracking_state)
 Handling of a particle interaction with a boundary of type CS_LAGR_BC_USER. More...
 
void cs_user_lagr_volume_conditions (void)
 Define particle volume conditions. More...
 
void cs_lagr_user_internal_interaction (cs_lagr_particle_set_t *particles, cs_lnum_t p_id, cs_lnum_t face_id, const cs_real_t face_norm[3], const cs_real_t c_intersect[3], cs_real_t t_intersect, cs_lagr_tracking_state_t *tracking_state)
 Handling of a particle interaction with a interior face of type CS_LAGR_BC_USER. More...
 

Function Documentation

◆ cs_lagr_user_boundary_interaction()

void cs_lagr_user_boundary_interaction ( cs_lagr_particle_set_t particles,
cs_lnum_t  p_id,
cs_lnum_t  face_id,
const cs_real_t  face_norm[3],
const cs_real_t  c_intersect[3],
cs_real_t  t_intersect,
int  b_zone_id,
int *  event_flag,
cs_lagr_tracking_state_t tracking_state 
)

Handling of a particle interaction with a boundary of type CS_LAGR_BC_USER.

Parameters
[in,out]particlespointer to particle set
[in]p_idparticle id
[in]face_idboundary face id
[in]face_normunit face (or face subdivision) normal
[in]c_intersectcoordinates of intersection with the face
[in]t_intersectrelative distance (in [0, 1]) of the intersection point with the face relative to the initial trajectory segment
[in]b_zone_idboundary zone id of the matching face
[in,out]event_flagevent flag in case events are available
[in,out]tracking_stateparticle tracking state

◆ cs_lagr_user_internal_interaction()

void cs_lagr_user_internal_interaction ( cs_lagr_particle_set_t particles,
cs_lnum_t  p_id,
cs_lnum_t  face_id,
const cs_real_t  face_norm[3],
const cs_real_t  c_intersect[3],
cs_real_t  t_intersect,
cs_lagr_tracking_state_t tracking_state 
)

Handling of a particle interaction with a interior face of type CS_LAGR_BC_USER.

Parameters
[in,out]particlespointer to particle set
[in]p_idparticle id
[in]face_idinterior face id
[in]face_normunit face (or face subdivision) normal
[in]c_intersectcoordinates of intersection with the face
[in]t_intersectrelative distance (in [0, 1]) of the intersection point with the face relative to the initial trajectory segment
[in,out]tracking_stateparticle tracking state

◆ cs_user_lagr_boundary_conditions()

void cs_user_lagr_boundary_conditions ( const int  bc_type[])

Define particle boundary conditions.

This is used for the definition of inlet and other boundaries, based on predefined boundary zones (cs_zone_t).

Parameters
[in]bc_typetype of the boundary faces

◆ cs_user_lagr_ef()

void cs_user_lagr_ef ( cs_real_t  dt_p,
const cs_real_t  taup[],
const cs_real_3_t  tlag[],
const cs_real_3_t  piil[],
const cs_real_33_t  bx[],
const cs_real_t  tsfext[],
const cs_real_33_t  vagaus[],
const cs_real_3_t  gradpr[],
const cs_real_33_t  gradvf[],
cs_real_t  rho_p[],
cs_real_3_t  fextla[] 
)

User definition of an external force field acting on the particles.

It must be prescribed in every cell and be homogeneous to gravity (m/s^2) By default gravity and drag force are the only forces acting on the particles (the gravity components gx gy gz are assigned in the GUI or in usipsu)

Parameters
[in]dt_ptime step (for the cell)
[in]taupparticle relaxation time
[in]tlagrelaxation time for the flow
[in]piilterm in the integration of the sde
[in]bxcharacteristics of the turbulence
[in]tsfextinfos for the return coupling
[in]vagausGaussian random variables
[in]gradprpressure gradient
[in]gradvfgradient of the flow velocity
[in,out]rho_pparticle density
[out]fextlauser external force field (m/s^2)$

◆ cs_user_lagr_extra_operations()

void cs_user_lagr_extra_operations ( const cs_real_t  dt[])

User function (non-mandatory intervention)

User-defined modifications on the variables at the end of the Lagrangian time step and calculation of user-defined additional statistics on the particles.

User-defined modifications on the variables at the end of the Lagrangian time step and calculation of user-defined additional statistics on the particles.

Parameters
[in]dttime step (per cell)

◆ cs_user_lagr_imposed_motion()

void cs_user_lagr_imposed_motion ( const cs_real_t  coords[3],
cs_real_t  dt,
cs_real_t  disp[3] 
)

Impose the motion of a particle flagged CS_LAGR_PART_IMPOSED_MOTION.

User-defined modifications on the particle position and its velocity.

Parameters
[in]coordsold particle coordinates
[in]dttime step (per particle)
[out]dispparticle dispacement

◆ cs_user_lagr_in()

void cs_user_lagr_in ( cs_lagr_particle_set_t particles,
const cs_lagr_injection_set_t zis,
const cs_lnum_t  particle_range[2],
const cs_lnum_t  particle_face_id[],
const cs_real_t  visc_length[] 
)

User modification of newly injected particles.

This function is called after the initialization of the new particles in order to modify them according to new particle profiles (injection profiles, position of the injection point, statistical weights, correction of the diameter if the standard-deviation option is activated).

This function is called for each injection zone and class. Particles with ids between pset->n_particles and n_elts are initialized but may be modidied by this function.

Parameters
[in,out]particlesparticle set
[in]zisinjection data for this set
[in]particle_rangestart and past-the-end ids of new particles for this zone and class
[in]particle_face_idface ids of new particles if zone is a boundary, NULL otherwise
[in]visc_lengthviscous layer thickness (size: number of mesh boundary faces)

This function is called after the initialization of the new particles in order to modify them according to new particle profiles (injection profiles, position of the injection point, statistical weights, correction of the diameter if the standard-deviation option is activated).

This function is called for each injection zone and set. Particles with ids between pset->n_particles and n_elts are initialized but may be modidied by this function.

Parameters
[in,out]particlesparticle set
[in]ziszone injection set data
[in]particle_rangestart and past-the-end ids of new particles for this zone and class
[in]particle_face_idface ids of new particles if zone is a boundary, NULL otherwise
[in]visc_lengthviscous layer thickness (size: number of mesh boundary faces)

◆ cs_user_lagr_model()

void cs_user_lagr_model ( void  )

◆ cs_user_lagr_rt()

void cs_user_lagr_rt ( cs_lnum_t  id_p,
cs_real_t  re_p,
cs_real_t  uvwr,
cs_real_t  rho_f,
cs_real_t  rho_p,
cs_real_t  nu_f,
cs_real_t  taup[],
const cs_real_t  dt[] 
)

Modification of the calculation of the particle relaxation time with respect to the chosen formulation for the drag coefficient.

This function is called in a loop on the particles, so be careful to avoid too costly operations.

This function is called in a loop on the particles, so be careful to avoid too costly operations.

           m   Cp
            p    p
  Tau = ---------------
     c          2
           PI d    h
               p    e

 Tau  : Thermal relaxation time (value to be computed)
    c

 m    : Particle mass
  p

 Cp   : Particle specific heat
   p

 d    : Particle diameter
  p

 h    : Coefficient of thermal exchange
  e

he coefficient of thermal exchange is calculated from a Nusselt number, itself evaluated by a correlation (Ranz-Marshall by default)

       h  d
        e  p
Nu = --------  = 2 + 0.55 Re **(0.5) Prt**(0.33)
      Lambda                p

Lambda : Thermal conductivity of the carrier field

Re     : Particle Reynolds number
  p

Prt    : Prandtl number
Parameters
[in]id_pparticle id
[in]re_pparticle Reynolds number
[in]uvwrrelative velocity of the particle (flow-seen velocity - part. velocity)
[in]rho_ffluid density at particle position
[in]rho_pparticle density
[in]nu_fkinematic viscosity of the fluid at particle position
[out]taupthermal relaxation time
[in]dttime step (per cell)

◆ cs_user_lagr_rt_t()

void cs_user_lagr_rt_t ( cs_lnum_t  id_p,
cs_real_t  re_p,
cs_real_t  uvwr,
cs_real_t  rho_f,
cs_real_t  rho_p,
cs_real_t  nu_f,
cs_real_t  cp_f,
cs_real_t  k_f,
cs_real_t  tauc[],
const cs_real_t  dt[] 
)

Modification of the calculation of the thermal relaxation time of the particles with respect to the chosen formulation of the Nusselt number.

This function is called in a loop on the particles, so be careful to avoid too costly operations.

Modification of the calculation of the thermal relaxation time of the particles with respect to the chosen formulation of the Nusselt number.

This function is called in a loop on the particles, so be careful to avoid too costly operations.

Parameters
[in]id_pparticle id
[in]re_pparticle Reynolds number
[in]uvwrrelative velocity of the particle (flow-seen velocity - part. velocity)
[in]rho_ffluid density at particle position
[in]rho_pparticle density
[in]nu_fkinematic viscosity of the fluid at particle position
[in]cp_fspecific heat of the fluid at particle position
[in]k_fdiffusion coefficient of the fluid at particle position
[out]taucthermal relaxation time
[in]dttime step (per cell)

◆ cs_user_lagr_sde()

void cs_user_lagr_sde ( const cs_real_t  dt[],
cs_real_t  taup[],
cs_real_3_t  tlag[],
cs_real_t  tempct[] 
)

User integration of the SDE for the user-defined variables.

The variables are constant by default. The SDE must be of the form:

\[ \frac{dT}{dt}=\frac{T - PIP}{Tca} \]

T: particle attribute representing the variable Tca: characteristic time for the sde to be prescribed in the array auxl1 PIP: coefficient of the SDE (pseudo RHS) to be prescribed in the array auxl2. If the chosen scheme is first order (nordre=1) then, at the first and only call pip is expressed as a function of the quantities of the previous time step (contained in the particle data). If the chosen scheme is second order (nordre=2) then, at the first call (nor=1) pip is expressed as a function of the quantities of the previous time step, and at the second passage (nor=2) pip is expressed as a function of the quantities of the current time step.

Parameters
[in]dttime step (per cell)
[in]taupparticle relaxation time
[in]tlagrelaxation time for the flow
[in]tempctcharacteristic thermal time and implicit source term of return coupling

◆ cs_user_lagr_volume_conditions()

void cs_user_lagr_volume_conditions ( void  )

Define particle volume conditions.

This is used definition of for injection based on predefined volume zones (cs_zone_t).

This is used for the definition of volume injections, based on predefined volume zones (cs_zone_t).