8.3
general documentation
cs_equation_iterative_solve.cpp File Reference

This file gathers functions that solve advection diffusion equations with source terms for one time step for a scalar, vector or tensor variable. More...

#include "cs_defs.h"
#include <stdarg.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <assert.h>
#include <float.h>
#include <math.h>
#include <mpi.h>
#include "bft_mem.h"
#include "bft_error.h"
#include "bft_printf.h"
#include "cs_array.h"
#include "cs_balance.h"
#include "cs_blas.h"
#include "cs_convection_diffusion.h"
#include "cs_dispatch.h"
#include "cs_field.h"
#include "cs_field_pointer.h"
#include "cs_halo.h"
#include "cs_log.h"
#include "cs_math.h"
#include "cs_mesh.h"
#include "cs_gradient.h"
#include "cs_mesh_quantities.h"
#include "cs_multigrid.h"
#include "cs_parameters.h"
#include "cs_porous_model.h"
#include "cs_prototypes.h"
#include "cs_timer.h"
#include "cs_parall.h"
#include "cs_matrix_building.h"
#include "cs_matrix_default.h"
#include "cs_sles.h"
#include "cs_sles_default.h"
#include "cs_equation_iterative_solve.h"
+ Include dependency graph for cs_equation_iterative_solve.cpp:

Functions

void cs_equation_iterative_solve_scalar (int idtvar, int iterns, int f_id, const char *name, int iescap, int imucpp, cs_real_t normp, cs_equation_param_t *eqp, const cs_real_t pvara[], const cs_real_t pvark[], const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_viscm[], const cs_real_t b_viscm[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], const cs_real_t rovsdt[], cs_real_t smbrp[], cs_real_t pvar[], cs_real_t dpvar[], const cs_real_t xcpp[], cs_real_t eswork[])
 Solve an advection diffusion equation with source terms for one time step for the variable $ a $. More...
 
void cs_equation_iterative_solve_vector (int idtvar, int iterns, int f_id, const char *name, int ivisep, int iescap, cs_equation_param_t *eqp, const cs_real_t pvara[][3], const cs_real_t pvark[][3], const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_viscm[], const cs_real_t b_viscm[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t i_secvis[], const cs_real_t b_secvis[], cs_real_t viscel[][6], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_t fimp[][3][3], cs_real_t smbrp[][3], cs_real_t pvar[][3], cs_real_t eswork[][3])
 This function solves an advection diffusion equation with source terms for one time step for the vector variable $ \vect{a} $. More...
 
void cs_equation_iterative_solve_tensor (int idtvar, int f_id, const char *name, cs_equation_param_t *eqp, const cs_real_t pvara[][6], const cs_real_t pvark[][6], const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_viscm[], const cs_real_t b_viscm[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t viscel[][6], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_t fimp[][6][6], cs_real_t smbrp[][6], cs_real_t pvar[][6])
 This function solves an advection diffusion equation with source terms for one time step for the symmetric tensor variable $ \tens{\variat} $. More...
 

Detailed Description

This file gathers functions that solve advection diffusion equations with source terms for one time step for a scalar, vector or tensor variable.

Function Documentation

◆ cs_equation_iterative_solve_scalar()

void cs_equation_iterative_solve_scalar ( int  idtvar,
int  iterns,
int  f_id,
const char *  name,
int  iescap,
int  imucpp,
cs_real_t  normp,
cs_equation_param_t eqp,
const cs_real_t  pvara[],
const cs_real_t  pvark[],
const cs_field_bc_coeffs_t bc_coeffs,
const cs_real_t  i_massflux[],
const cs_real_t  b_massflux[],
const cs_real_t  i_viscm[],
const cs_real_t  b_viscm[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
cs_real_6_t  viscel[],
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
int  icvflb,
const int  icvfli[],
const cs_real_t  rovsdt[],
cs_real_t  smbrp[],
cs_real_t  pvar[],
cs_real_t  dpvar[],
const cs_real_t  xcpp[],
cs_real_t  eswork[] 
)

Solve an advection diffusion equation with source terms for one time step for the variable $ a $.

The equation reads:

\[
f_s^{imp}(a^{n+1}-a^n)
+ \divs \left( a^{n+1} \rho \vect{u} - \mu \grad a^{n+1} \right)
= Rhs
\]

This equation is rewritten as:

\[
f_s^{imp} \delta a
+ \divs \left( \delta a \rho \vect{u} - \mu \grad \delta a \right)
= Rhs^1
\]

where $ \delta a = a^{n+1} - a^n$ and $ Rhs^1 = Rhs - \divs( a^n \rho \vect{u} - \mu \grad a^n)$

It is in fact solved with the following iterative process:

\[
f_s^{imp} \delta a^k
+ \divs \left(\delta a^k \rho \vect{u}-\mu\grad\delta a^k \right)
= Rhs^k
\]

where $Rhs^k=Rhs-f_s^{imp}(a^k-a^n)
- \divs \left( a^k\rho\vect{u}-\mu\grad a^k \right)$

Be careful, it is forbidden to modify $ f_s^{imp} $ here!

Please refer to the codits section of the theory guide for more informations.

Parameters
[in]idtvarindicator of the temporal scheme
[in]iternsexternal sub-iteration number
[in]f_idfield id (or -1)
[in]nameassociated name if f_id < 0, or nullptr
[in]iescapcompute the predictor indicator if 1
[in]imucppindicator
  • 0 do not multiply the convectiv term by Cp
  • 1 do multiply the convectiv term by Cp
[in]normpReference norm to solve the system (optional) if negative: recomputed here
[in]eqppointer to a cs_equation_param_t structure which contains variable calculation options
[in]pvaravariable at the previous time step $ a^n $
[in]pvarkvariable at the previous sub-iteration $ a^k $. If you sub-iter on Navier-Stokes, then it allows to initialize by something else than pvara (usually pvar=pvara)
[in]bc_coeffsboundary condition structure for the variable
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]i_viscm$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the matrix
[in]b_viscm$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the matrix
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in]rovsdt$ f_s^{imp} $
[in]smbrpRight hand side $ Rhs^k $
[in,out]pvarcurrent variable
[out]dpvarlast variable increment
[in]xcpparray of specific heat (Cp)
[out]esworkprediction-stage error estimator (if iescap > 0)

◆ cs_equation_iterative_solve_tensor()

void cs_equation_iterative_solve_tensor ( int  idtvar,
int  f_id,
const char *  name,
cs_equation_param_t eqp,
const cs_real_t  pvara[][6],
const cs_real_t  pvark[][6],
const cs_field_bc_coeffs_t bc_coeffs_ts,
const cs_real_t  i_massflux[],
const cs_real_t  b_massflux[],
const cs_real_t  i_viscm[],
const cs_real_t  b_viscm[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
cs_real_t  viscel[][6],
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
int  icvflb,
const int  icvfli[],
cs_real_t  fimp[][6][6],
cs_real_t  smbrp[][6],
cs_real_t  pvar[][6] 
)

This function solves an advection diffusion equation with source terms for one time step for the symmetric tensor variable $ \tens{\variat} $.

The equation reads:

\[
\tens{f_s}^{imp}(\tens{\variat}^{n+1}-\tens{\variat}^n)
+ \divt \left( \tens{\variat}^{n+1} \otimes \rho \vect {u}
             - \mu \gradtt \tens{\variat}^{n+1}\right)
= \tens{Rhs}
\]

This equation is rewritten as:

\[
\tens{f_s}^{imp} \delta \tens{\variat}
+ \divt \left( \delta \tens{\variat} \otimes \rho \vect{u}
             - \mu \gradtt \delta \tens{\variat} \right)
= \tens{Rhs}^1
\]

where $ \delta \tens{\variat} = \tens{\variat}^{n+1} - \tens{\variat}^n$ and $ \tens{Rhs}^1 = \tens{Rhs}
- \divt \left( \tens{\variat}^n \otimes \rho \vect{u}
             - \mu \gradtt \tens{\variat}^n \right)$

It is in fact solved with the following iterative process:

\[
\tens{f_s}^{imp} \delta \tens{\variat}^k
+ \divt \left( \delta \tens{\variat}^k \otimes \rho \vect{u}
             - \mu \gradtt \delta \tens{\variat}^k \right)
= \tens{Rhs}^k
\]

where $ \tens{Rhs}^k = \tens{Rhs}
- \tens{f_s}^{imp} \left(\tens{\variat}^k-\tens{\variat}^n \right)
- \divt \left( \tens{\variat}^k \otimes \rho \vect{u}
             - \mu \gradtt \tens{\variat}^k \right)$

Be careful, it is forbidden to modify $ \tens{f_s}^{imp} $ here!

Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idfield id (or -1)
[in]nameassociated name if f_id < 0, or nullptr
[in]eqppointer to a cs_equation_param_t structure which contains variable calculation options
[in]pvaravariable at the previous time step $ \vect{a}^n $
[in]pvarkvariable at the previous sub-iteration $ \vect{a}^k $. If you sub-iter on Navier-Stokes, then it allows to initialize by something else than pvara (usually pvar=pvara)
[in]bc_coeffs_tsboundary condition structure for the variable
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]i_viscm$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the matrix
[in]b_viscm$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the matrix
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in,out]fimp$ \tens{f_s}^{imp} $
[in,out]smbrpRight hand side $ \vect{Rhs}^k $
[in,out]pvarcurrent variable

◆ cs_equation_iterative_solve_vector()

void cs_equation_iterative_solve_vector ( int  idtvar,
int  iterns,
int  f_id,
const char *  name,
int  ivisep,
int  iescap,
cs_equation_param_t eqp,
const cs_real_t  pvara[][3],
const cs_real_t  pvark[][3],
const cs_field_bc_coeffs_t bc_coeffs_v,
const cs_real_t  i_massflux[],
const cs_real_t  b_massflux[],
const cs_real_t  i_viscm[],
const cs_real_t  b_viscm[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  i_secvis[],
const cs_real_t  b_secvis[],
cs_real_t  viscel[][6],
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
int  icvflb,
const int  icvfli[],
cs_real_t  fimp[][3][3],
cs_real_t  smbrp[][3],
cs_real_t  pvar[][3],
cs_real_t  eswork[][3] 
)

This function solves an advection diffusion equation with source terms for one time step for the vector variable $ \vect{a} $.

The equation reads:

\[
\tens{f_s}^{imp}(\vect{a}^{n+1}-\vect{a}^n)
+ \divv \left( \vect{a}^{n+1} \otimes \rho \vect {u}
             - \mu \gradt \vect{a}^{n+1}\right)
= \vect{Rhs}
\]

This equation is rewritten as:

\[
\tens{f_s}^{imp} \delta \vect{a}
+ \divv \left( \delta \vect{a} \otimes \rho \vect{u}
             - \mu \gradt \delta \vect{a} \right)
= \vect{Rhs}^1
\]

where $ \delta \vect{a} = \vect{a}^{n+1} - \vect{a}^n$ and $ \vect{Rhs}^1 = \vect{Rhs}
- \divv \left( \vect{a}^n \otimes \rho \vect{u}
             - \mu \gradt \vect{a}^n \right)$

It is in fact solved with the following iterative process:

\[
\tens{f_s}^{imp} \delta \vect{a}^k
+ \divv \left( \delta \vect{a}^k \otimes \rho \vect{u}
             - \mu \gradt \delta \vect{a}^k \right)
= \vect{Rhs}^k
\]

where $ \vect{Rhs}^k = \vect{Rhs}
- \tens{f_s}^{imp} \left(\vect{a}^k-\vect{a}^n \right)
- \divv \left( \vect{a}^k \otimes \rho \vect{u}
             - \mu \gradt \vect{a}^k \right)$

Be careful, it is forbidden to modify $ \tens{f_s}^{imp} $ here!

Parameters
[in]idtvarindicator of the temporal scheme
[in]iternsexternal sub-iteration number
[in]f_idfield id (or -1)
[in]nameassociated name if f_id < 0, or nullptr
[in]ivisepindicator to take $ \divv
                              \left(\mu \gradt \transpose{\vect{a}} \right)
                              -2/3 \grad\left( \mu \dive \vect{a} \right)$
  • 1 take into account,
  • 0 otherwise
[in]iescapcompute the predictor indicator if >= 1
[in]eqppointer to a cs_equation_param_t structure which contains variable calculation options
[in]pvaravariable at the previous time step $ \vect{a}^n $
[in]pvarkvariable at the previous sub-iteration $ \vect{a}^k $. If you sub-iter on Navier-Stokes, then it allows to initialize by something else than pvara (usually pvar= pvara)
[in]bc_coeffs_vboundary condition structure for the variable
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]i_viscm$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the matrix
[in]b_viscm$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the matrix
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the r.h.s.
[in]i_secvissecondary viscosity at interior faces
[in]b_secvissecondary viscosity at boundary faces
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in,out]fimp$ \tens{f_s}^{imp} $
[in,out]smbrpRight hand side $ \vect{Rhs}^k $
[in,out]pvarcurrent variable
[out]esworkprediction-stage error estimator (if iescap >= 0)