My Project
programmer's documentation
|
User subroutines for additional right-hand side source terms. More...
Functions/Subroutines | |
subroutine | ustsnv (nvar, nscal, ncepdp, ncesmp, ivar, icepdc, icetsm, itypsm, dt, ckupdc, smacel, crvexp, crvimp) |
Additional right-hand side source terms for velocity components equation (Navier-Stokes) More... | |
subroutine | ustssc (nvar, nscal, ncepdp, ncesmp, iscal, icepdc, icetsm, itypsm, dt, ckupdc, smacel, crvexp, crvimp) |
User subroutine. More... | |
subroutine | ustsvv (nvar, nscal, ncepdp, ncesmp, iscal, icepdc, icetsm, itypsm, dt, ckupdc, smacel, crvexp, crvimp) |
Additional right-hand side source terms for vectorial equations (user vectors and specific physics vectors). More... | |
subroutine | cs_user_turbulence_source_terms (nvar, nscal, ncepdp, ncesmp, f_id, icepdc, icetsm, itypsm, ckupdc, smacel, crvexp, crvimp) |
Additional right-hand side source terms for turbulence models. More... | |
subroutine | cs_user_turbulence_source_terms2 (nvar, nscal, ncepdp, ncesmp, f_id, icepdc, icetsm, itypsm, ckupdc, smacel, crvexp, crvimp) |
Additional right-hand side source terms for turbulence models and irijco =1. More... | |
User subroutines for additional right-hand side source terms.
See cs_user_source_terms and Examples of data settings for source terms with scalar in a channel for examples.
subroutine cs_user_turbulence_source_terms | ( | integer | nvar, |
integer | nscal, | ||
integer | ncepdp, | ||
integer | ncesmp, | ||
integer | f_id, | ||
integer, dimension(ncepdp) | icepdc, | ||
integer, dimension(ncesmp) | icetsm, | ||
integer, dimension(ncesmp,nvar) | itypsm, | ||
double precision, dimension(6,ncepdp) | ckupdc, | ||
double precision, dimension(ncesmp,nvar) | smacel, | ||
double precision, dimension(ncelet) | crvexp, | ||
double precision, dimension(ncelet) | crvimp | ||
) |
Additional right-hand side source terms for turbulence models.
The additional source term is decomposed into an explicit part (crvexp) and an implicit part (crvimp) that must be provided here. The resulting equations solved by the code are:
where is the turbulence field of index f_id
Note that crvexp, crvimp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:
The crvexp, crvimp arrays are already initialized to 0 before entering the routine. It is not needed to do it in the routine (waste of CPU time).
For stability reasons, Code_Saturne will not add -crvimp directly to the diagonal of the matrix, but Max(-crvimp,0). This way, the crvimp term is treated implicitely only if it strengthens the diagonal of the matrix. However, when using the second-order in time scheme, this limitation cannot be done anymore and -crvimp is added directly. The user should therefore test the negativity of crvimp by himself.
When using the second-order in time scheme, one should supply:
The selection of cells where to apply the source terms is based on a getcel command. For more info on the syntax of the getcel command, refer to the user manual or to the comments on the similar command getfbr in the routine cs_user_boundary_conditions.
[in] | nvar | total number of variables |
[in] | nscal | total number of scalars |
[in] | ncepdp | number of cells with head loss terms |
[in] | ncesmp | number of cells with mass source terms |
[in] | f_id | field index of the current turbulent variable |
[in] | icepdc | index number of cells with head loss terms |
[in] | icetsm | index number of cells with mass source terms |
[in] | itypsm | type of mass source term for each variable (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[in] | ckupdc | head loss coefficient |
[in] | smacel | value associated to each variable in the mass source terms or mass rate (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[out] | crvexp | explicit part of the source term |
[out] | crvimp | implicit part of the source term |
subroutine cs_user_turbulence_source_terms2 | ( | integer | nvar, |
integer | nscal, | ||
integer | ncepdp, | ||
integer | ncesmp, | ||
integer | f_id, | ||
integer, dimension(ncepdp) | icepdc, | ||
integer, dimension(ncesmp) | icetsm, | ||
integer, dimension(ncesmp,nvar) | itypsm, | ||
double precision, dimension(6,ncepdp) | ckupdc, | ||
double precision, dimension(ncesmp,nvar) | smacel, | ||
double precision, dimension(6,ncelet) | crvexp, | ||
double precision, dimension(6,6,ncelet) | crvimp | ||
) |
Additional right-hand side source terms for turbulence models and irijco =1.
The additional source term is decomposed into an explicit part (crvexp) and an implicit part (crvimp) that must be provided here. The resulting equations solved by the code are:
where is the turbulence field of index f_id
Note that crvexp, crvimp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:
The crvexp, crvimp arrays are already initialized to 0 before entering the routine. It is not needed to do it in the routine (waste of CPU time).
For stability reasons, Code_Saturne will not add -crvimp directly to the diagonal of the matrix, but Max(-crvimp,0). This way, the crvimp term is treated implicitely only if it strengthens the diagonal of the matrix. However, when using the second-order in time scheme, this limitation cannot be done anymore and -crvimp is added directly. The user should therefore test the negativity of crvimp by himself.
When using the second-order in time scheme, one should supply:
The selection of cells where to apply the source terms is based on a getcel command. For more info on the syntax of the getcel command, refer to the user manual or to the comments on the similar command getfbr in the routine cs_user_boundary_conditions.
[in] | nvar | total number of variables |
[in] | nscal | total number of scalars |
[in] | ncepdp | number of cells with head loss terms |
[in] | ncesmp | number of cells with mass source terms |
[in] | f_id | field index of the current turbulent variable |
[in] | icepdc | index number of cells with head loss terms |
[in] | icetsm | index number of cells with mass source terms |
[in] | itypsm | type of mass source term for each variable (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[in] | ckupdc | head loss coefficient |
[in] | smacel | value associated to each variable in the mass source terms or mass rate (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[out] | crvexp | explicit part of the source term |
[out] | crvimp | implicit part of the source term |
subroutine ustsnv | ( | integer | nvar, |
integer | nscal, | ||
integer | ncepdp, | ||
integer | ncesmp, | ||
integer | ivar, | ||
integer, dimension(ncepdp) | icepdc, | ||
integer, dimension(ncesmp) | icetsm, | ||
integer, dimension(ncesmp,nvar) | itypsm, | ||
double precision, dimension(ncelet) | dt, | ||
double precision, dimension(6,ncepdp) | ckupdc, | ||
double precision, dimension(ncesmp,nvar) | smacel, | ||
double precision, dimension(3,ncelet) | crvexp, | ||
double precision, dimension(3,3,ncelet) | crvimp | ||
) |
Additional right-hand side source terms for velocity components equation (Navier-Stokes)
The additional source term is decomposed into an explicit part (crvexp
) and an implicit part (crvimp
) that must be provided here. The resulting equation solved by the code for a velocity is:
Note that crvexp
and crvimp
are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:
The crvexp
and crvimp
arrays are already initialized to 0 before entering the the routine. It is not needed to do it in the routine (waste of CPU time).
crvexp(i, iel)
+ vel(j, iel)* crvimp(j, i, iel).For stability reasons, Code_Saturne will not add -crvimp directly to the diagonal of the matrix, but Max(-crvimp,0). This way, the crvimp term is treated implicitely only if it strengthens the diagonal of the matrix. However, when using the second-order in time scheme, this limitation cannot be done anymore and -crvimp is added directly. The user should therefore test the negativity of crvimp by himself.
When using the second-order in time scheme, one should supply:
The selection of cells where to apply the source terms is based on a getcel command. For more info on the syntax of the getcel command, refer to the user manual or to the comments on the similar command getfbr in the routine cs_user_boundary_conditions.
[in] | nvar | total number of variables |
[in] | nscal | total number of scalars |
[in] | ncepdp | number of cells with head loss terms |
[in] | ncesmp | number of cells with mass source terms |
[in] | ivar | index number of the current variable |
[in] | icepdc | index number of cells with head loss terms |
[in] | icetsm | index number of cells with mass source terms |
[in] | itypsm | type of mass source term for each variable (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[in] | dt | time step (per cell) |
[in] | ckupdc | head loss coefficient |
[in] | smacel | value associated to each variable in the mass source terms or mass rate (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[out] | crvexp | explicit part of the source term |
[out] | crvimp | implicit part of the source term |
subroutine ustssc | ( | integer | nvar, |
integer | nscal, | ||
integer | ncepdp, | ||
integer | ncesmp, | ||
integer | iscal, | ||
integer, dimension(ncepdp) | icepdc, | ||
integer, dimension(ncesmp) | icetsm, | ||
integer, dimension(ncesmp,nvar) | itypsm, | ||
double precision, dimension(ncelet) | dt, | ||
double precision, dimension(6,ncepdp) | ckupdc, | ||
double precision, dimension(ncesmp,nvar) | smacel, | ||
double precision, dimension(ncelet) | crvexp, | ||
double precision, dimension(ncelet) | crvimp | ||
) |
User subroutine.
Additional right-hand side source terms for scalar equations (user scalars and specific physics scalars).
The routine is called for each scalar, user or specific physisc. It is therefore necessary to test the value of the scalar number iscal to separate the treatments of the different scalars (if (iscal.eq.p) then ....).
The additional source term is decomposed into an explicit part (crvexp) and an implicit part (crvimp) that must be provided here. The resulting equation solved by the code for a scalar f is:
Note that crvexp and crvimp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:
The crvexp and crvimp arrays are already initialized to 0 before entering the the routine. It is not needed to do it in the routine (waste of CPU time).
For stability reasons, Code_Saturne will not add -crvimp directly to the diagonal of the matrix, but Max(-crvimp,0). This way, the crvimp term is treated implicitely only if it strengthens the diagonal of the matrix. However, when using the second-order in time scheme, this limitation cannot be done anymore and -crvimp is added directly. The user should therefore test the negativity of crvimp by himself.
When using the second-order in time scheme, one should supply:
The selection of cells where to apply the source terms is based on a getcel command. For more info on the syntax of the getcel command, refer to the user manual or to the comments on the similar command getfbr in the routine cs_user_boundary_conditions. WARNING: If scalar is the temperature, the resulting equation solved by the code is:
rho*Cp*volume*dT/dt + .... = crvimp*T + crvexp
Note that crvexp and crvimp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:
In case of a complex, non-linear source term, say F(f), for scalar f, the easiest method is to implement the source term explicitely.
df/dt = .... + F(f(n)) where f(n) is the value of f at time tn, the beginning of the time step.
This yields : crvexp = volume*F(f(n)) crvimp = 0
However, if the source term is potentially steep, this fully explicit method will probably generate instabilities. It is therefore wiser to partially implicit the term by writing:
df/dt = .... + dF/df*f(n+1) - dF/df*f(n) + F(f(n))
This yields: crvexp = volume*( F(f(n)) - dF/df*f(n) ) crvimp = volume*dF/df
[in] | nvar | total number of variables |
[in] | nscal | total number of scalars |
[in] | ncepdp | number of cells with head loss terms |
[in] | ncesmp | number of cells with mass source terms |
[in] | iscal | index number of the current scalar |
[in] | icepdc | index number of cells with head loss terms |
[in] | icetsm | index number of cells with mass source terms |
[in] | itypsm | type of mass source term for each variable \paramsee in |
[in] | dt | time step (per cell) |
[in] | ckupdc | head loss coefficient |
[in] | smacel | value associated to each variable in the mass |
[in] | source | terms or mass rate (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[out] | crvexp | explicit part of the source term |
[out] | crvimp | implicit part of the source term |
subroutine ustsvv | ( | integer | nvar, |
integer | nscal, | ||
integer | ncepdp, | ||
integer | ncesmp, | ||
integer | iscal, | ||
integer, dimension(ncepdp) | icepdc, | ||
integer, dimension(ncesmp) | icetsm, | ||
integer, dimension(ncesmp,nvar) | itypsm, | ||
double precision, dimension(ncelet) | dt, | ||
double precision, dimension(6,ncepdp) | ckupdc, | ||
double precision, dimension(ncesmp,nvar) | smacel, | ||
double precision, dimension(3,ncelet) | crvexp, | ||
double precision, dimension(3,3,ncelet) | crvimp | ||
) |
Additional right-hand side source terms for vectorial equations (user vectors and specific physics vectors).
The routine is called for each vector, user or specific physisc. It is therefore necessary to test the value of the vector number iscal to separate the treatments of the different vectors (if (iscal.eq.p) then ....).
The additional source term is decomposed into an explicit part (crvexp) and an implicit part (crvimp) that must be provided here. The resulting equation solved by the code for a vector f is:
Note that crvexp and crvimp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:
The crvexp and crvimp arrays are already initialized to 0 before entering the the routine. It is not needed to do it in the routine (waste of CPU time).
For stability reasons, Code_Saturne will not add -crvimp directly to the diagonal of the matrix, but Max(-crvimp,0). This way, the crvimp term is treated implicitely only if it strengthens the diagonal of the matrix. However, when using the second-order in time scheme, this limitation cannot be done anymore and -crvimp is added directly. The user should therefore test the negativity of crvimp by himself.
When using the second-order in time scheme, one should supply:
The selection of cells where to apply the source terms is based on a getcel command. For more info on the syntax of the getcel command, refer to the user manual or to the comments on the similar command getfbr in the routine cs_user_boundary_conditions. WARNING: If scalar is the temperature, the resulting equation solved by the code is:
rho*Cp*volume*dT/dt + .... = crvimp*T + crvexp
Note that crvexp and crvimp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:
In case of a complex, non-linear source term, say F(f), for scalar f, the easiest method is to implement the source term explicitely.
df/dt = .... + F(f(n)) where f(n) is the value of f at time tn, the beginning of the time step.
This yields : crvexp = volume*F(f(n)) crvimp = 0
However, if the source term is potentially steep, this fully explicit method will probably generate instabilities. It is therefore wiser to partially implicit the term by writing:
df/dt = .... + dF/df*f(n+1) - dF/df*f(n) + F(f(n))
This yields: crvexp = volume*( F(f(n)) - dF/df*f(n) ) crvimp = volume*dF/df
[in] | nvar | total number of variables |
[in] | nscal | total number of scalars |
[in] | ncepdp | number of cells with head loss terms |
[in] | ncesmp | number of cells with mass source terms |
[in] | iscal | index number of the current scalar |
[in] | icepdc | index number of cells with head loss terms |
[in] | icetsm | index number of cells with mass source terms |
[in] | itypsm | type of mass source term for each variable \paramsee in |
[in] | dt | time step (per cell) |
[in] | ckupdc | head loss coefficient |
[in] | smacel | value associated to each variable in the mass |
[in] | source | terms or mass rate (see Examples of data settings for mass source terms (cs_user_mass_source_terms.f90)) |
[out] | crvexp | explicit part of the source term |
[out] | crvimp | implicit part of the source term |