programmer's documentation
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
Data setting for drift scalars

Introduction

This page provides an example of code blocks that may be used to perform a calculation with drift scalars.

Physical properties

Energy balance

Local variables to be added

The following local variables need to be defined for the examples in this section:

integer ivart, iel, ifac
integer ipcvis
integer ipcvsl, iscal, iflid, iscdri
integer f_id, keydri, nfld, keysca
double precision rho, viscl
double precision diamp, rhop, cuning
double precision xrtp, xk, xeps, beta1
character*80 fname
double precision, dimension(:), pointer :: taup
double precision, dimension(:), pointer :: taufpt
double precision, dimension(:), pointer :: crom

Initialization and finalization

The following initialization block needs to be added for the following examples:

ipcvis = ipproc(iviscl)
! Key id for drift scalar
call field_get_key_id("drift_scalar_model", keydri)
! Key id for scalar id
call field_get_key_id("scalar_id", keysca)
! Number of fields

In theory Fortran 95 deallocates locally-allocated arrays automatically, but deallocating arrays in a symmetric manner to their allocation is good practice, and it avoids using a different logic for C and Fortran.

Body

This example set the scalar laminar diffusivity (for Brownian motion) to take thermophoresis into account.

Here is the corresponding code:

! Loop over fields which are scalar with a drift
do iflid = 0, nfld-1
call field_get_key_int(iflid, keydri, iscdri)
! We only handle here scalar with a drift
if (btest(iscdri, drift_scalar_add_drift_flux)) then
! Position of variables, coefficients
! -----------------------------------
! Index of the scalar
call field_get_key_int(iflid, keysca, iscal)
! --- Number of the thermal variable
if (iscalt.gt.0) then
ivart = isca(iscalt)
else
write(nfecra,9010) iscalt
call csexit(1)
endif
! --- Rank of scalar's diffusivity (Brownian motion)
! in 'propce', physical properties at element centers: 'ipcvsl'
if (ivisls(iscal).gt.0) then
ipcvsl = ipproc(ivisls(iscal))
else
ipcvsl = 0
endif
! --- Coefficients of drift scalar CHOSEN BY THE USER
! Values given here are fictitious
! diamp: is the diameter of the particle class
! cuning: is the Cuningham correction factor
! rhop: particle density
diamp = 1.d-4
cuning = 1.d0
rhop = 1.d4
! Name of the scalar with a drift
call field_get_name(iflid, fname)
! Index of the corresponding relaxation time (taup)
call field_get_id('drift_tau_'//trim(fname), f_id)
call field_get_val_s(f_id, taup)
! Index of the corresponding interaction time particle--eddies (taufpt)
if (btest(iscdri, drift_scalar_turbophoresis)) then
call field_get_id('drift_turb_tau_'//trim(fname), f_id)
call field_get_val_s(f_id, taufpt)
endif
! Computation of the relaxation time of the particles
!----------------------------------------------------
if (diamp.le.1.d-6) then
! Cuningham's correction for submicronic particules
do iel = 1, ncel
taup(iel) = cuning*diamp**2*rhop/(18.d0*propce(iel,ipcvis))
enddo
else
do iel = 1, ncel
taup(iel) = diamp**2*rhop/(18.d0*propce(iel,ipcvis))
enddo
endif
! Compute the interaction time particle--eddies (tau_fpt)
!--------------------------------------------------------
if (btest(iscdri, drift_scalar_turbophoresis)) then
! k-epsilon or v2-f models
if (itytur.eq.2 .or. itytur.eq.5) then
do iel = 1, ncel
xk = rtp(iel,ik)
xeps = rtp(iel,iep)
taufpt(iel) = (3.d0/2.d0)*(cmu/sigmas(iscal))*xk/xeps
enddo
! Rij-epsilon models
else if (itytur.eq.3) then
beta1 = 0.5d0+3.d0/(4.d0*xkappa)
do iel = 1, ncel
xk = 0.5d0*( rtp(iel,ir11) &
+rtp(iel,ir22) &
+rtp(iel,ir33) )
xeps = rtp(iel,iep)
taufpt(iel) = xk/xeps/beta1
enddo
! k-omega models
else if (iturb.eq.60) then
do iel = 1, ncel
xk = rtp(iel,ik)
xeps = cmu*xk*rtp(iel,iomg)
taufpt(iel) = (3.d0/2.d0)*(cmu/sigmas(iscal))*xk/xeps
enddo
endif
endif
! Brownian diffusion at cell centers
!-----------------------------------
! --- Stop if the diffusivity is not variable
if (ipcvsl.le.0) then
write(nfecra,1010) iscal, iscal, ivisls(iscal)
call csexit(1)
endif
! Homogeneous to a dynamic viscosity
do iel = 1, ncel
xrtp = rtp(iel,ivart)
rho = crom(iel)
viscl = propce(iel, ipcvis)
propce(iel,ipcvsl) = rho*kboltz*xrtp*cuning/(3.d0*pi*diamp*viscl)
enddo
endif ! --- Tests on drift scalar
enddo