programmer's documentation
Functions/Subroutines
coditv.f90 File Reference

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

Functions/Subroutines

subroutine coditv (idtvar, ivar, iconvp, idiffp, ndircp, imrgra, nswrsp, nswrgp, imligp, ircflp, ivisep, ischcp, isstpp, iescap, idftnp, iswdyp, iwarnp, blencp, epsilp, epsrsp, epsrgp, climgp, relaxp, thetap, pvara, pvark, coefav, coefbv, cofafv, cofbfv, flumas, flumab, viscfm, viscbm, viscfs, viscbs, secvif, secvib, icvflb, icvfli, fimp, smbrp, pvar, eswork)
 

Detailed Description

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!

Function/Subroutine Documentation

subroutine coditv ( integer  idtvar,
integer  ivar,
integer  iconvp,
integer  idiffp,
integer  ndircp,
integer  imrgra,
integer  nswrsp,
integer  nswrgp,
integer  imligp,
integer  ircflp,
integer  ivisep,
integer  ischcp,
integer  isstpp,
integer  iescap,
integer  idftnp,
integer  iswdyp,
integer  iwarnp,
double precision  blencp,
double precision  epsilp,
double precision  epsrsp,
double precision  epsrgp,
double precision  climgp,
double precision  relaxp,
double precision  thetap,
double precision, dimension(3,ncelet)  pvara,
double precision, dimension(3,ncelet)  pvark,
double precision, dimension(3,ndimfb)  coefav,
double precision, dimension(3,3,ndimfb)  coefbv,
double precision, dimension(3,ndimfb)  cofafv,
double precision, dimension(3,3,ndimfb)  cofbfv,
double precision, dimension(nfac)  flumas,
double precision, dimension(nfabor)  flumab,
double precision, dimension(*)  viscfm,
double precision, dimension(nfabor)  viscbm,
double precision, dimension(*)  viscfs,
double precision, dimension(nfabor)  viscbs,
double precision, dimension(nfac)  secvif,
double precision, dimension(nfabor)  secvib,
integer  icvflb,
integer, dimension(nfabor)  icvfli,
double precision, dimension(3,3,ncelet)  fimp,
double precision, dimension(3,ncelet)  smbrp,
double precision, dimension(3,ncelet)  pvar,
double precision, dimension(3,ncelet)  eswork 
)
Parameters
[in]idtvarindicator of the temporal scheme
[in]ivarindex of the current variable
[in]iconvpindicator
  • 1 convection,
  • 0 otherwise
[in]idiffpindicator
  • 1 diffusion,
  • 0 otherwise
[in]ndircpindicator (0 if the diagonal is stepped aside)
[in]imrgraindicateur
  • 0 iterative gradient
  • 1 least squares gradient
[in]nswrspnumber of reconstruction sweeps for the Right Hand Side
[in]nswrgpnumber of reconstruction sweeps for the gradients
[in]imligpclipping gradient method
  • < 0 no clipping
  • = 0 thank to neighboring gradients
  • = 1 thank to the mean gradient
[in]ircflpindicator
  • 1 flux reconstruction,
  • 0 otherwise
[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]ischcpindicator
  • 1 centered
  • 0 2nd order
[in]isstppindicator
  • 1 without slope test
  • 0 with slope test
[in]iescapcompute the predictor indicator if 1
[in]idftnpindicator
  • 1 the diffusivity is scalar
  • 6 the diffusivity is a symmetric tensor
[in]iswdypindicator
  • 0 no dynamic relaxation
  • 1 dynamic relaxation depending on $ \delta \vect{\varia}^k $
  • 2 dynamic relaxation depending on $ \delta \vect{\varia}^k $ and $ \delta \vect{\varia}^{k-1} $
[in]iwarnpverbosity
[in]blencpfraction of upwinding
[in]epsilpprecision pour resol iter
[in]epsrsprelative precision for the iterative process
[in]epsrgprelative precision for the gradient reconstruction
[in]climgpclipping coefficient for the computation of the gradient
[in]relaxpcoefficient of relaxation
[in]thetapweighting coefficient for the theta-schema,
  • thetap = 0: explicit scheme
  • thetap = 0.5: time-centered scheme (mix between Crank-Nicolson and Adams-Bashforth)
  • thetap = 1: implicit scheme
[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]flumasmass flux at interior faces
[in]flumabmass flux at boundary faces
[in]viscfm$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the matrix
[in]viscbm$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the matrix
[in]viscfs$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]viscbs$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the r.h.s.
[in]secvifsecondary viscosity at interior faces
[in]secvibsecondary viscosity at boundary faces
[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)