7.0
general documentation
cs_equation_iterative_solve.h File Reference
#include "cs_defs.h"
#include "cs_parameters.h"
+ Include dependency graph for cs_equation_iterative_solve.h:

Go to the source code of this file.

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_var_cal_opt_t *var_cal_opt, const cs_real_t pvara[], const cs_real_t pvark[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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[])
 This function solves 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_var_cal_opt_t *var_cal_opt, const cs_real_3_t pvara[], const cs_real_3_t pvark[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], 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 secvif[], const cs_real_t secvib[], cs_real_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_33_t fimp[], cs_real_3_t smbrp[], cs_real_3_t pvar[], cs_real_3_t eswork[])
 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_var_cal_opt_t *var_cal_opt, const cs_real_6_t pvara[], const cs_real_6_t pvark[], const cs_real_6_t coefats[], const cs_real_66_t coefbts[], const cs_real_6_t cofafts[], const cs_real_66_t cofbfts[], 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_66_t fimp[], cs_real_6_t smbrp[], cs_real_6_t pvar[])
 This function solves an advection diffusion equation with source terms for one time step for the symmetric tensor variable $ \tens{\variat} $. More...
 

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_var_cal_opt_t var_cal_opt,
const cs_real_t  pvara[],
const cs_real_t  pvark[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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[] 
)

This function solves 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!

Parameters
[in]idtvarindicator of the temporal scheme
[in]iternsexternal sub-iteration number
[in]f_idfield id (or -1)
[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]var_cal_optpointer to a cs_var_cal_opt_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]coefapboundary condition array for the variable (explicit part)
[in]coefbpboundary condition array for the variable (implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (implicit part)
[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
[in,out]dpvarlast variable increment
[in]xcpparray of specific heat (Cp)
[out]esworkprediction-stage error estimator (if iescap > 0)

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 NULL
[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]var_cal_optpointer to a cs_var_cal_opt_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]coefapboundary condition array for the variable (explicit part)
[in]coefbpboundary condition array for the variable (implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (implicit part)
[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
[in,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_var_cal_opt_t var_cal_opt,
const cs_real_6_t  pvara[],
const cs_real_6_t  pvark[],
const cs_real_6_t  coefats[],
const cs_real_66_t  coefbts[],
const cs_real_6_t  cofafts[],
const cs_real_66_t  cofbfts[],
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_66_t  fimp[],
cs_real_6_t  smbrp[],
cs_real_6_t  pvar[] 
)

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]var_cal_optpointer to a cs_var_cal_opt_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]coefatsboundary condition array for the variable (Explicit part)
[in]coefbtsboundary condition array for the variable (Impplicit part)
[in]cofaftsboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbftsboundary condition array for the diffusion of the variable (Implicit part)
[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]fimp$ \tens{f_s}^{imp} $
[in]smbrpRight hand side $ \vect{Rhs}^k $
[in,out]pvarcurrent variable

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 NULL
[in]var_cal_optpointer to a cs_var_cal_opt_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]coefatsboundary condition array for the variable (Explicit part)
[in]coefbtsboundary condition array for the variable (Impplicit part)
[in]cofaftsboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbftsboundary condition array for the diffusion of the variable (Implicit part)
[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]fimp$ \tens{f_s}^{imp} $
[in]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_var_cal_opt_t var_cal_opt,
const cs_real_3_t  pvara[],
const cs_real_3_t  pvark[],
const cs_real_3_t  coefav[],
const cs_real_33_t  coefbv[],
const cs_real_3_t  cofafv[],
const cs_real_33_t  cofbfv[],
const cs_real_t  i_massflux[],
const cs_real_t  b_massflux[],
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_6_t  viscel[],
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
int  icvflb,
const int  icvfli[],
cs_real_33_t  fimp[],
cs_real_3_t  smbrp[],
cs_real_3_t  pvar[],
cs_real_3_t  eswork[] 
)

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 NULL
[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]var_cal_optpointer to a cs_var_cal_opt_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]coefavboundary condition array for the variable (explicit part)
[in]coefbvboundary condition array for the variable (implicit part)
[in]cofafvboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfvboundary condition array for the diffusion of the variable (Implicit part)
[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]fimp$ \tens{f_s}^{imp} $
[in]smbrpRight hand side $ \vect{Rhs}^k $
[in,out]pvarcurrent variable
[out]esworkprediction-stage error estimator (if iescap > 0)

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 NULL
[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]var_cal_optpointer to a cs_var_cal_opt_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]coefavboundary condition array for the variable (explicit part)
[in]coefbvboundary condition array for the variable (implicit part)
[in]cofafvboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfvboundary condition array for the diffusion of the variable (Implicit part)
[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]fimp$ \tens{f_s}^{imp} $
[in]smbrpRight hand side $ \vect{Rhs}^k $
[in,out]pvarcurrent variable
[out]esworkprediction-stage error estimator (if iescap > 0)