8.1
general documentation
cs_equation_iterative_solve.h
Go to the documentation of this file.
1 #ifndef __CS_EQUATION_ITERATIVE_SOLVE_H__
2 #define __CS_EQUATION_ITERATIVE_SOLVE_H__
3 
4 /*============================================================================
5  * Scalar convection diffusion equation solver.
6  *============================================================================*/
7 
8 /*
9  This file is part of code_saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2023 EDF S.A.
12 
13  This program is free software; you can redistribute it and/or modify it under
14  the terms of the GNU General Public License as published by the Free Software
15  Foundation; either version 2 of the License, or (at your option) any later
16  version.
17 
18  This program is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21  details.
22 
23  You should have received a copy of the GNU General Public License along with
24  this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25  Street, Fifth Floor, Boston, MA 02110-1301, USA.
26 */
27 
28 /*----------------------------------------------------------------------------*/
29 
30 /*----------------------------------------------------------------------------
31  * Local headers
32  *----------------------------------------------------------------------------*/
33 
34 #include "cs_defs.h"
35 #include "cs_parameters.h"
36 
38 
39 /*=============================================================================
40  * Public function prototypes
41  *============================================================================*/
42 
43 /*----------------------------------------------------------------------------*/
49 /*----------------------------------------------------------------------------*/
50 
51 /*----------------------------------------------------------------------------*/
52 /*
53  * \brief Solve an advection diffusion equation with source
54  * terms for one time step for the variable \f$ a \f$.
55  *
56  * The equation reads:
57  *
58  * \f[
59  * f_s^{imp}(a^{n+1}-a^n)
60  * + \divs \left( a^{n+1} \rho \vect{u} - \mu \grad a^{n+1} \right)
61  * = Rhs
62  * \f]
63  *
64  * This equation is rewritten as:
65  *
66  * \f[
67  * f_s^{imp} \delta a
68  * + \divs \left( \delta a \rho \vect{u} - \mu \grad \delta a \right)
69  * = Rhs^1
70  * \f]
71  *
72  * where \f$ \delta a = a^{n+1} - a^n\f$ and
73  * \f$ Rhs^1 = Rhs - \divs( a^n \rho \vect{u} - \mu \grad a^n)\f$
74  *
75  *
76  * It is in fact solved with the following iterative process:
77  *
78  * \f[
79  * f_s^{imp} \delta a^k
80  * + \divs \left(\delta a^k \rho \vect{u}-\mu\grad\delta a^k \right)
81  * = Rhs^k
82  * \f]
83  *
84  * where \f$Rhs^k=Rhs-f_s^{imp}(a^k-a^n)
85  * - \divs \left( a^k\rho\vect{u}-\mu\grad a^k \right)\f$
86  *
87  * Be careful, it is forbidden to modify \f$ f_s^{imp} \f$ here!
88  *
89  * Please refer to the
90  * <a href="../../theory.pdf#codits"><b>codits</b></a> section of the
91  * theory guide for more informations.
92  *
93  * \param[in] idtvar indicator of the temporal scheme
94  * \param[in] iterns external sub-iteration number
95  * \param[in] f_id field id (or -1)
96  * \param[in] name associated name if f_id < 0, or NULL
97  * \param[in] iescap compute the predictor indicator if 1
98  * \param[in] imucpp indicator
99  * - 0 do not multiply the convectiv term by Cp
100  * - 1 do multiply the convectiv term by Cp
101  * \param[in] normp Reference norm to solve the system (optional)
102  * if negative: recomputed here
103  * \param[in] var_cal_opt pointer to a cs_var_cal_opt_t structure which
104  * contains variable calculation options
105  * \param[in] pvara variable at the previous time step
106  * \f$ a^n \f$
107  * \param[in] pvark variable at the previous sub-iteration
108  * \f$ a^k \f$.
109  * If you sub-iter on Navier-Stokes, then
110  * it allows to initialize by something else than
111  * pvara (usually pvar=pvara)
112  * \param[in] coefap boundary condition array for the variable
113  * (explicit part)
114  * \param[in] coefbp boundary condition array for the variable
115  * (implicit part)
116  * \param[in] cofafp boundary condition array for the diffusion
117  * of the variable (explicit part)
118  * \param[in] cofbfp boundary condition array for the diffusion
119  * of the variable (implicit part)
120  * \param[in] i_massflux mass flux at interior faces
121  * \param[in] b_massflux mass flux at boundary faces
122  * \param[in] i_viscm \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
123  * at interior faces for the matrix
124  * \param[in] b_viscm \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
125  * at boundary faces for the matrix
126  * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
127  * at interior faces for the r.h.s.
128  * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
129  * at boundary faces for the r.h.s.
130  * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
131  * \param[in] weighf internal face weight between cells i j in case
132  * of tensor diffusion
133  * \param[in] weighb boundary face weight for cells i in case
134  * of tensor diffusion
135  * \param[in] icvflb global indicator of boundary convection flux
136  * - 0 upwind scheme at all boundary faces
137  * - 1 imposed flux at some boundary faces
138  * \param[in] icvfli boundary face indicator array of convection flux
139  * - 0 upwind scheme
140  * - 1 imposed flux
141  * \param[in] rovsdt \f$ f_s^{imp} \f$
142  * \param[in] smbrp Right hand side \f$ Rhs^k \f$
143  * \param[in,out] pvar current variable
144  * \param[out] dpvar last variable increment
145  * \param[in] xcpp array of specific heat (Cp)
146  * \param[out] eswork prediction-stage error estimator
147  * (if iescap > 0)
148  */
149 /*----------------------------------------------------------------------------*/
150 
151 void
153  int iterns,
154  int f_id,
155  const char *name,
156  int iescap,
157  int imucpp,
158  cs_real_t normp,
160  const cs_real_t pvara[],
161  const cs_real_t pvark[],
162  const cs_real_t coefap[],
163  const cs_real_t coefbp[],
164  const cs_real_t cofafp[],
165  const cs_real_t cofbfp[],
166  const cs_real_t i_massflux[],
167  const cs_real_t b_massflux[],
168  const cs_real_t i_viscm[],
169  const cs_real_t b_viscm[],
170  const cs_real_t i_visc[],
171  const cs_real_t b_visc[],
172  cs_real_6_t viscel[],
173  const cs_real_2_t weighf[],
174  const cs_real_t weighb[],
175  int icvflb,
176  const int icvfli[],
177  const cs_real_t rovsdt[],
178  cs_real_t smbrp[],
179  cs_real_t pvar[],
180  cs_real_t dpvar[],
181  const cs_real_t xcpp[],
182  cs_real_t eswork[]);
183 
184 /*----------------------------------------------------------------------------*/
185 /*
186  * \brief This function solves an advection diffusion equation with source terms
187  * for one time step for the vector variable \f$ \vect{a} \f$.
188  *
189  * The equation reads:
190  *
191  * \f[
192  * \tens{f_s}^{imp}(\vect{a}^{n+1}-\vect{a}^n)
193  * + \divv \left( \vect{a}^{n+1} \otimes \rho \vect {u}
194  * - \mu \gradt \vect{a}^{n+1}\right)
195  * = \vect{Rhs}
196  * \f]
197  *
198  * This equation is rewritten as:
199  *
200  * \f[
201  * \tens{f_s}^{imp} \delta \vect{a}
202  * + \divv \left( \delta \vect{a} \otimes \rho \vect{u}
203  * - \mu \gradt \delta \vect{a} \right)
204  * = \vect{Rhs}^1
205  * \f]
206  *
207  * where \f$ \delta \vect{a} = \vect{a}^{n+1} - \vect{a}^n\f$ and
208  * \f$ \vect{Rhs}^1 = \vect{Rhs}
209  * - \divv \left( \vect{a}^n \otimes \rho \vect{u}
210  * - \mu \gradt \vect{a}^n \right)\f$
211  *
212  *
213  * It is in fact solved with the following iterative process:
214  *
215  * \f[
216  * \tens{f_s}^{imp} \delta \vect{a}^k
217  * + \divv \left( \delta \vect{a}^k \otimes \rho \vect{u}
218  * - \mu \gradt \delta \vect{a}^k \right)
219  * = \vect{Rhs}^k
220  * \f]
221  *
222  * where \f$ \vect{Rhs}^k = \vect{Rhs}
223  * - \tens{f_s}^{imp} \left(\vect{a}^k-\vect{a}^n \right)
224  * - \divv \left( \vect{a}^k \otimes \rho \vect{u}
225  * - \mu \gradt \vect{a}^k \right)\f$
226  *
227  * Be careful, it is forbidden to modify \f$ \tens{f_s}^{imp} \f$ here!
228  *
229  * \param[in] idtvar indicator of the temporal scheme
230  * \param[in] iterns external sub-iteration number
231  * \param[in] f_id field id (or -1)
232  * \param[in] name associated name if f_id < 0, or NULL
233  * \param[in] ivisep indicator to take \f$ \divv
234  * \left(\mu \gradt \transpose{\vect{a}} \right)
235  * -2/3 \grad\left( \mu \dive \vect{a} \right)\f$
236  * - 1 take into account,
237  * - 0 otherwise
238  * \param[in] iescap compute the predictor indicator if >= 1
239  * \param[in] var_cal_opt pointer to a cs_var_cal_opt_t structure which
240  * contains variable calculation options
241  * \param[in] pvara variable at the previous time step
242  * \f$ \vect{a}^n \f$
243  * \param[in] pvark variable at the previous sub-iteration
244  * \f$ \vect{a}^k \f$.
245  * If you sub-iter on Navier-Stokes, then
246  * it allows to initialize by something else than
247  * \c pvara (usually \c pvar= \c pvara)
248  * \param[in] coefav boundary condition array for the variable
249  * (explicit part)
250  * \param[in] coefbv boundary condition array for the variable
251  * (implicit part)
252  * \param[in] cofafv boundary condition array for the diffusion
253  * of the variable (Explicit part)
254  * \param[in] cofbfv boundary condition array for the diffusion
255  * of the variable (Implicit part)
256  * \param[in] i_massflux mass flux at interior faces
257  * \param[in] b_massflux mass flux at boundary faces
258  * \param[in] i_viscm \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
259  * at interior faces for the matrix
260  * \param[in] b_viscm \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
261  * at boundary faces for the matrix
262  * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
263  * at interior faces for the r.h.s.
264  * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
265  * at boundary faces for the r.h.s.
266  * \param[in] i_secvis secondary viscosity at interior faces
267  * \param[in] b_secvis secondary viscosity at boundary faces
268  * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
269  * \param[in] weighf internal face weight between cells i j in case
270  * of tensor diffusion
271  * \param[in] weighb boundary face weight for cells i in case
272  * of tensor diffusion
273  * \param[in] icvflb global indicator of boundary convection flux
274  * - 0 upwind scheme at all boundary faces
275  * - 1 imposed flux at some boundary faces
276  * \param[in] icvfli boundary face indicator array of convection flux
277  * - 0 upwind scheme
278  * - 1 imposed flux
279  * \param[in, out] fimp \f$ \tens{f_s}^{imp} \f$
280  * \param[in] smbrp Right hand side \f$ \vect{Rhs}^k \f$
281  * \param[in, out] pvar current variable
282  * \param[out] eswork prediction-stage error estimator
283  * (if iescap >= 0)
284  */
285 /*----------------------------------------------------------------------------*/
286 
287 void
289  int iterns,
290  int f_id,
291  const char *name,
292  int ivisep,
293  int iescap,
295  const cs_real_3_t pvara[],
296  const cs_real_3_t pvark[],
297  const cs_real_3_t coefav[],
298  const cs_real_33_t coefbv[],
299  const cs_real_3_t cofafv[],
300  const cs_real_33_t cofbfv[],
301  const cs_real_t i_massflux[],
302  const cs_real_t b_massflux[],
303  cs_real_t i_viscm[],
304  const cs_real_t b_viscm[],
305  const cs_real_t i_visc[],
306  const cs_real_t b_visc[],
307  const cs_real_t i_secvis[],
308  const cs_real_t b_secvis[],
309  cs_real_6_t viscel[],
310  const cs_real_2_t weighf[],
311  const cs_real_t weighb[],
312  int icvflb,
313  const int icvfli[],
314  cs_real_33_t fimp[],
315  cs_real_3_t smbrp[],
316  cs_real_3_t pvar[],
317  cs_real_3_t eswork[]);
318 
319 /*----------------------------------------------------------------------------*/
320 /*
321  * \brief This function solves an advection diffusion equation with source
322  * terms for one time step for the symmetric tensor variable
323  * \f$ \tens{\variat} \f$.
324  *
325  * The equation reads:
326  *
327  * \f[
328  * \tens{f_s}^{imp}(\tens{\variat}^{n+1}-\tens{\variat}^n)
329  * + \divt \left( \tens{\variat}^{n+1} \otimes \rho \vect {u}
330  * - \mu \gradtt \tens{\variat}^{n+1}\right)
331  * = \tens{Rhs}
332  * \f]
333  *
334  * This equation is rewritten as:
335  *
336  * \f[
337  * \tens{f_s}^{imp} \delta \tens{\variat}
338  * + \divt \left( \delta \tens{\variat} \otimes \rho \vect{u}
339  * - \mu \gradtt \delta \tens{\variat} \right)
340  * = \tens{Rhs}^1
341  * \f]
342  *
343  * where \f$ \delta \tens{\variat} = \tens{\variat}^{n+1} - \tens{\variat}^n\f$
344  * and \f$ \tens{Rhs}^1 = \tens{Rhs}
345  * - \divt \left( \tens{\variat}^n \otimes \rho \vect{u}
346  * - \mu \gradtt \tens{\variat}^n \right)\f$
347  *
348  *
349  * It is in fact solved with the following iterative process:
350  *
351  * \f[
352  * \tens{f_s}^{imp} \delta \tens{\variat}^k
353  * + \divt \left( \delta \tens{\variat}^k \otimes \rho \vect{u}
354  * - \mu \gradtt \delta \tens{\variat}^k \right)
355  * = \tens{Rhs}^k
356  * \f]
357  *
358  * where \f$ \tens{Rhs}^k = \tens{Rhs}
359  * - \tens{f_s}^{imp} \left(\tens{\variat}^k-\tens{\variat}^n \right)
360  * - \divt \left( \tens{\variat}^k \otimes \rho \vect{u}
361  * - \mu \gradtt \tens{\variat}^k \right)\f$
362  *
363  * Be careful, it is forbidden to modify \f$ \tens{f_s}^{imp} \f$ here!
364  *
365  * \param[in] idtvar indicator of the temporal scheme
366  * \param[in] f_id field id (or -1)
367  * \param[in] name associated name if f_id < 0, or NULL
368  * \param[in] var_cal_opt pointer to a cs_var_cal_opt_t structure which
369  * contains variable calculation options
370  * \param[in] pvara variable at the previous time step
371  * \f$ \vect{a}^n \f$
372  * \param[in] pvark variable at the previous sub-iteration
373  * \f$ \vect{a}^k \f$.
374  * If you sub-iter on Navier-Stokes, then
375  * it allows to initialize by something else than
376  * pvara (usually pvar=pvara)
377  * \param[in] coefats boundary condition array for the variable
378  * (Explicit part)
379  * \param[in] coefbts boundary condition array for the variable
380  * (Impplicit part)
381  * \param[in] cofafts boundary condition array for the diffusion
382  * of the variable (Explicit part)
383  * \param[in] cofbfts boundary condition array for the diffusion
384  * of the variable (Implicit part)
385  * \param[in] i_massflux mass flux at interior faces
386  * \param[in] b_massflux mass flux at boundary faces
387  * \param[in] i_viscm \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
388  * at interior faces for the matrix
389  * \param[in] b_viscm \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
390  * at boundary faces for the matrix
391  * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
392  * at interior faces for the r.h.s.
393  * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
394  * at boundary faces for the r.h.s.
395  * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
396  * \param[in] weighf internal face weight between cells i j in case
397  * of tensor diffusion
398  * \param[in] weighb boundary face weight for cells i in case
399  * of tensor diffusion
400  * \param[in] icvflb global indicator of boundary convection flux
401  * - 0 upwind scheme at all boundary faces
402  * - 1 imposed flux at some boundary faces
403  * \param[in] icvfli boundary face indicator array of convection flux
404  * - 0 upwind scheme
405  * - 1 imposed flux
406  * \param[in] fimp \f$ \tens{f_s}^{imp} \f$
407  * \param[in] smbrp Right hand side \f$ \vect{Rhs}^k \f$
408  * \param[in,out] pvar current variable
409  */
410 /*----------------------------------------------------------------------------*/
411 
412 void
414  int f_id,
415  const char *name,
417  const cs_real_6_t pvara[],
418  const cs_real_6_t pvark[],
419  const cs_real_6_t coefats[],
420  const cs_real_66_t coefbts[],
421  const cs_real_6_t cofafts[],
422  const cs_real_66_t cofbfts[],
423  const cs_real_t i_massflux[],
424  const cs_real_t b_massflux[],
425  const cs_real_t i_viscm[],
426  const cs_real_t b_viscm[],
427  const cs_real_t i_visc[],
428  const cs_real_t b_visc[],
429  cs_real_6_t viscel[],
430  const cs_real_2_t weighf[],
431  const cs_real_t weighb[],
432  int icvflb,
433  const int icvfli[],
434  const cs_real_66_t fimp[],
435  cs_real_6_t smbrp[],
436  cs_real_6_t pvar[]);
437 
439 
440 #endif /* __CS_EQUATION_ITERATIVE_SOLVE_H__ */
#define BEGIN_C_DECLS
Definition: cs_defs.h:514
double cs_real_t
Floating-point value.
Definition: cs_defs.h:319
cs_real_t cs_real_66_t[6][6]
6x6 matrix of floating-point values
Definition: cs_defs.h:344
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:334
cs_real_t cs_real_2_t[2]
vector of 2 floating-point values
Definition: cs_defs.h:333
cs_real_t cs_real_6_t[6]
vector of 6 floating-point values
Definition: cs_defs.h:336
#define END_C_DECLS
Definition: cs_defs.h:515
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition: cs_defs.h:343
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 symm...
Definition: cs_equation_iterative_solve.c:2063
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[])
Solve an advection diffusion equation with source terms for one time step for the variable .
Definition: cs_equation_iterative_solve.c:205
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 vect...
Definition: cs_equation_iterative_solve.c:1128
integer, dimension(:), allocatable, target icvfli
boundary convection flux indicator of a Rusanov or an analytical flux (some boundary contributions of...
Definition: cfpoin.f90:52
cs_equation_param_t * var_cal_opt
Definition: keywords.h:135
integer(c_int), pointer, save idtvar
option for a variable time step
Definition: optcal.f90:260
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition: cs_equation_param.h:192