9.1
general documentation
cs_balance.h
Go to the documentation of this file.
1#ifndef __CS_BALANCE_H__
2#define __CS_BALANCE_H__
3
4/*============================================================================
5 * Building of the right hand side for a transport of a field.
6 *============================================================================*/
7
8/*
9 This file is part of code_saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2025 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#include "base/cs_defs.h"
31
32/*----------------------------------------------------------------------------
33 * Standard C library headers
34 *----------------------------------------------------------------------------*/
35
36/*----------------------------------------------------------------------------
37 * Local headers
38 *----------------------------------------------------------------------------*/
39
40#include "base/cs_parameters.h"
41
42/*----------------------------------------------------------------------------*/
44
45/*=============================================================================
46 * Public function prototypes
47 *============================================================================*/
48
49/*----------------------------------------------------------------------------*/
53/*----------------------------------------------------------------------------*/
54
55void
57
58/*----------------------------------------------------------------------------*/
59/*
60 * \brief Wrapper to the function which adds the explicit part of the
61 * convection/diffusion terms of a transport equation of
62 * a scalar field \f$ \varia \f$.
63 *
64 * More precisely, the right hand side \f$ Rhs \f$ is updated as
65 * follows:
66 * \f[
67 * Rhs = Rhs - \sum_{\fij \in \Facei{\celli}} \left(
68 * \dot{m}_\ij \left( \varia_\fij - \varia_\celli \right)
69 * - \mu_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right)
70 * \f]
71 *
72 * Warning:
73 * - \f$ Rhs \f$ has already been initialized
74 * before calling cs_balance_scalar!
75 * - mind the minus sign
76 *
77 * Options for the convective scheme:
78 * - blencp = 0: upwind scheme for the advection
79 * - blencp = 1: no upwind scheme except in the slope test
80 * - ischcp = 0: SOLU
81 * - ischcp = 1: centered
82 * - ischcp = 2: SOLU with upwind gradient reconstruction
83 * - ischcp = 3: blending SOLU centered
84 * - ischcp = 4: NVD-TVD
85 * - imucpp = 0: do not multiply the convective part by \f$ C_p \f$
86 * - imucpp = 1: multiply the convective part by \f$ C_p \f$
87 *
88 * \param[in] idtvar indicator of the temporal scheme
89 * \param[in] f_id field id (or -1)
90 * \param[in] imucpp indicator
91 * - 0 do not multiply the convectiv term by Cp
92 * - 1 do multiply the convectiv term by Cp
93 * \param[in] imasac take mass accumulation into account?
94 * \param[in] inc indicator
95 * - 0 when solving an increment
96 * - 1 otherwise
97 * \param[in] eqp pointer to a cs_equation_param_t structure which
98 * contains variable calculation options
99 * \param[in] pvar solved variable (current time step)
100 * may be NULL if pvara != NULL
101 * \param[in] pvara solved variable (previous time step)
102 * may be NULL if pvar != NULL
103 * \param[in] bc_coeffs boundary condition structure for the variable
104 * \param[in] val_f boundary face value for gradient
105 * \param[in] flux boundary flux
106 * \param[in] i_massflux mass flux at interior faces
107 * \param[in] b_massflux mass flux at boundary faces
108 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
109 * at interior faces for the r.h.s.
110 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
111 * at boundary faces for the r.h.s.
112 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
113 * \param[in] xcpp array of specific heat (Cp)
114 * \param[in] weighf internal face weight between cells i j in case
115 * of tensor diffusion
116 * \param[in] weighb boundary face weight for cells i in case
117 * of tensor diffusion
118 * \param[in] icvflb global indicator of boundary convection flux
119 * - 0 upwind scheme at all boundary faces
120 * - 1 imposed flux at some boundary faces
121 * \param[in] icvfli boundary face indicator array of convection flux
122 * - 0 upwind scheme
123 * - 1 imposed flux
124 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
125 */
126/*----------------------------------------------------------------------------*/
127
128void
130 int f_id,
131 int imucpp,
132 int imasac,
133 int inc,
135 cs_real_t pvar[],
136 const cs_real_t pvara[],
137 const cs_field_bc_coeffs_t *bc_coeffs,
138 const cs_real_t i_massflux[],
139 const cs_real_t b_massflux[],
140 const cs_real_t i_visc[],
141 const cs_real_t b_visc[],
142 cs_real_6_t viscel[],
143 const cs_real_t xcpp[],
144 const cs_real_2_t weighf[],
145 const cs_real_t weighb[],
146 int icvflb,
147 const int icvfli[],
148 cs_real_t rhs[]);
149
151
152#if defined(__cplusplus)
153
154/*----------------------------------------------------------------------------*/
155/*
156 * \brief Wrapper to the function which adds the explicit part of the
157 * convection/diffusion terms of a transport equation of
158 * a scalar field \f$ \varia \f$.
159 *
160 * More precisely, the right hand side \f$ Rhs \f$ is updated as
161 * follows:
162 * \f[
163 * Rhs = Rhs - \sum_{\fij \in \Facei{\celli}} \left(
164 * \dot{m}_\ij \left( \varia_\fij - \varia_\celli \right)
165 * - \mu_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right)
166 * \f]
167 *
168 * Warning:
169 * - \f$ Rhs \f$ has already been initialized
170 * before calling cs_balance_scalar!
171 * - mind the minus sign
172 *
173 * Options for the convective scheme:
174 * - blencp = 0: upwind scheme for the advection
175 * - blencp = 1: no upwind scheme except in the slope test
176 * - ischcp = 0: SOLU
177 * - ischcp = 1: centered
178 * - ischcp = 2: SOLU with upwind gradient reconstruction
179 * - ischcp = 3: blending SOLU centered
180 * - ischcp = 4: NVD-TVD
181 * - imucpp = 0: do not multiply the convective part by \f$ C_p \f$
182 * - imucpp = 1: multiply the convective part by \f$ C_p \f$
183 *
184 * \param[in] idtvar indicator of the temporal scheme
185 * \param[in] f_id field id (or -1)
186 * \param[in] imucpp indicator
187 * - 0 do not multiply the convectiv term by Cp
188 * - 1 do multiply the convectiv term by Cp
189 * \param[in] imasac take mass accumulation into account?
190 * \param[in] inc indicator
191 * - 0 when solving an increment
192 * - 1 otherwise
193 * \param[in] eqp pointer to a cs_equation_param_t structure which
194 * contains variable calculation options
195 * \param[in] pvar solved variable (current time step)
196 * may be NULL if pvara != NULL
197 * \param[in] pvara solved variable (previous time step)
198 * may be NULL if pvar != NULL
199 * \param[in] bc_coeffs boundary condition structure for the variable
200 * \param[in] val_f boundary face value for gradient
201 * \param[in] flux boundary flux
202 * \param[in] i_massflux mass flux at interior faces
203 * \param[in] b_massflux mass flux at boundary faces
204 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
205 * at interior faces for the r.h.s.
206 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
207 * at boundary faces for the r.h.s.
208 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
209 * \param[in] xcpp array of specific heat (Cp)
210 * \param[in] weighf internal face weight between cells i j in case
211 * of tensor diffusion
212 * \param[in] weighb boundary face weight for cells i in case
213 * of tensor diffusion
214 * \param[in] icvflb global indicator of boundary convection flux
215 * - 0 upwind scheme at all boundary faces
216 * - 1 imposed flux at some boundary faces
217 * \param[in] icvfli boundary face indicator array of convection flux
218 * - 0 upwind scheme
219 * - 1 imposed flux
220 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
221 */
222/*----------------------------------------------------------------------------*/
223
224void
226 int f_id,
227 int imucpp,
228 int imasac,
229 int inc,
231 cs_real_t pvar[],
232 const cs_real_t pvara[],
233 const cs_field_bc_coeffs_t *bc_coeffs,
234 const cs_real_t val_f[],
235 const cs_real_t flux[],
236 const cs_real_t i_massflux[],
237 const cs_real_t b_massflux[],
238 const cs_real_t i_visc[],
239 const cs_real_t b_visc[],
240 cs_real_6_t viscel[],
241 const cs_real_t xcpp[],
242 const cs_real_2_t weighf[],
243 const cs_real_t weighb[],
244 int icvflb,
245 const int icvfli[],
246 cs_real_t rhs[]);
247
248/*----------------------------------------------------------------------------*/
249/*
250 * \brief Wrapper to the function which adds the explicit part of the
251 * convection/diffusion terms of a transport equation of
252 * a scalar field \f$ \varia \f$.
253 *
254 * More precisely, the right hand side \f$ Rhs \f$ is updated as
255 * follows:
256 * \f[
257 * Rhs = Rhs - \sum_{\fij \in \Facei{\celli}} \left(
258 * \dot{m}_\ij \left( \varia_\fij - \varia_\celli \right)
259 * - \mu_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right)
260 * \f]
261 *
262 * Warning:
263 * - \f$ Rhs \f$ has already been initialized
264 * before calling cs_balance_scalar!
265 * - mind the minus sign
266 *
267 * Options for the convective scheme:
268 * - blencp = 0: upwind scheme for the advection
269 * - blencp = 1: no upwind scheme except in the slope test
270 * - ischcp = 0: SOLU
271 * - ischcp = 1: centered
272 * - ischcp = 2: SOLU with upwind gradient reconstruction
273 * - ischcp = 3: blending SOLU centered
274 * - ischcp = 4: NVD-TVD
275 * - imucpp = 0: do not multiply the convective part by \f$ C_p \f$
276 * - imucpp = 1: multiply the convective part by \f$ C_p \f$
277 *
278 * \param[in] idtvar indicator of the temporal scheme
279 * \param[in] f_id field id (or -1)
280 * \param[in] imucpp indicator
281 * - 0 do not multiply the convectiv term by Cp
282 * - 1 do multiply the convectiv term by Cp
283 * \param[in] imasac take mass accumulation into account?
284 * \param[in] inc indicator
285 * - 0 when solving an increment
286 * - 1 otherwise
287 * \param[in] eqp pointer to a cs_equation_param_t structure which
288 * contains variable calculation options
289 * \param[in] pvar solved variable (current time step)
290 * may be NULL if pvara != NULL
291 * \param[in] pvara solved variable (previous time step)
292 * may be NULL if pvar != NULL
293 * \param[in] bc_coeffs boundary condition structure for the variable
294 * \param[in] val_f boundary face value for gradient
295 * \param[in] flux boundary flux
296 * \param[in] i_massflux mass flux at interior faces
297 * \param[in] b_massflux mass flux at boundary faces
298 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
299 * at interior faces for the r.h.s.
300 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
301 * at boundary faces for the r.h.s.
302 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
303 * \param[in] xcpp array of specific heat (Cp)
304 * \param[in] weighf internal face weight between cells i j in case
305 * of tensor diffusion
306 * \param[in] weighb boundary face weight for cells i in case
307 * of tensor diffusion
308 * \param[in] icvflb global indicator of boundary convection flux
309 * - 0 upwind scheme at all boundary faces
310 * - 1 imposed flux at some boundary faces
311 * \param[in] icvfli boundary face indicator array of convection flux
312 * - 0 upwind scheme
313 * - 1 imposed flux
314 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
315 * \param[in,out] i_flux interior flux (or nullptr)
316 * \param[in,out] b_flux boundary flux (or nullptr)
317 */
318/*----------------------------------------------------------------------------*/
319
320void
322 int f_id,
323 int imucpp,
324 int imasac,
325 int inc,
327 cs_real_t pvar[],
328 const cs_real_t pvara[],
329 const cs_field_bc_coeffs_t *bc_coeffs,
330 const cs_real_t val_f[],
331 const cs_real_t flux[],
332 const cs_real_t i_massflux[],
333 const cs_real_t b_massflux[],
334 const cs_real_t i_visc[],
335 const cs_real_t b_visc[],
336 cs_real_6_t viscel[],
337 const cs_real_t xcpp[],
338 const cs_real_2_t weighf[],
339 const cs_real_t weighb[],
340 int icvflb,
341 const int icvfli[],
342 cs_real_t rhs[],
343 cs_real_2_t i_flux[],
344 cs_real_t b_flux[]);
345
346#endif // defined(__cplusplus)
347
349
350/*----------------------------------------------------------------------------*/
351/*
352 * \brief Wrapper to the function which adds the explicit part of the
353 * convection/diffusion
354 * terms of a transport equation of a vector field \f$ \vect{\varia} \f$.
355 *
356 * More precisely, the right hand side \f$ \vect{Rhs} \f$ is updated as
357 * follows:
358 * \f[
359 * \vect{Rhs} = \vect{Rhs} - \sum_{\fij \in \Facei{\celli}} \left(
360 * \dot{m}_\ij \left( \vect{\varia}_\fij - \vect{\varia}_\celli \right)
361 * - \mu_\fij \gradt_\fij \vect{\varia} \cdot \vect{S}_\ij \right)
362 * \f]
363 *
364 * Remark:
365 * if ivisep = 1, then we also take \f$ \mu \transpose{\gradt\vect{\varia}}
366 * + \lambda \trace{\gradt\vect{\varia}} \f$, where \f$ \lambda \f$ is
367 * the secondary viscosity, i.e. usually \f$ -\frac{2}{3} \mu \f$.
368 *
369 * Warning:
370 * - \f$ \vect{Rhs} \f$ has already been initialized
371 * before calling cs_balance_vector!
372 * - mind the sign minus
373 *
374 * Options for the convective scheme:
375 * - blencp = 0: upwind scheme for the advection
376 * - blencp = 1: no upwind scheme except in the slope test
377 * - ischcp = 0: SOLU
378 * - ischcp = 1: centered
379 * - ischcp = 2: SOLU with upwind gradient reconstruction
380 * - ischcp = 3: blending SOLU centered
381 * - ischcp = 4: NVD-TVD
382 *
383 * \param[in] idtvar indicator of the temporal scheme
384 * \param[in] f_id field id (or -1)
385 * \param[in] imasac take mass accumulation into account?
386 * \param[in] inc indicator
387 * - 0 when solving an increment
388 * - 1 otherwise
389 * \param[in] ivisep indicator to take \f$ \divv
390 * \left(\mu \gradt \transpose{\vect{a}} \right)
391 * -2/3 \grad\left( \mu \dive \vect{a} \right)\f$
392 * - 1 take into account,
393 * - 0 otherwise
394 * \param[in] eqp pointer to a cs_equation_param_t structure which
395 * contains variable calculation options
396 * \param[in] pvar solved velocity (current time step)
397 * \param[in] pvara solved velocity (previous time step)
398 * \param[in] bc_coeffs_v boundary condition structure for the variable
399 * \param[in] val_f boundary face value for gradient
400 * \param[in] flux boundary flux
401 * \param[in] i_massflux mass flux at interior faces
402 * \param[in] b_massflux mass flux at boundary faces
403 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
404 * at interior faces for the r.h.s.
405 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
406 * at boundary faces for the r.h.s.
407 * \param[in] secvif secondary viscosity at interior faces
408 * \param[in] secvib secondary viscosity at boundary faces
409 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
410 * \param[in] weighf internal face weight between cells i j in case
411 * of tensor diffusion
412 * \param[in] weighb boundary face weight for cells i in case
413 * of tensor diffusion
414 * \param[in] icvflb global indicator of boundary convection flux
415 * - 0 upwind scheme at all boundary faces
416 * - 1 imposed flux at some boundary faces
417 * \param[in] icvfli boundary face indicator array of convection flux
418 * - 0 upwind scheme
419 * - 1 imposed flux
420 * \param[in,out] smbr right hand side \f$ \vect{Rhs} \f$
421 */
422/*----------------------------------------------------------------------------*/
423
424void
426 int f_id,
427 int imasac,
428 int inc,
429 int ivisep,
431 cs_real_3_t pvar[],
432 const cs_real_3_t pvara[],
433 const cs_field_bc_coeffs_t *bc_coeffs_v,
434 const cs_real_t val_f[][3],
435 const cs_real_t flux[][3],
436 const cs_real_t i_massflux[],
437 const cs_real_t b_massflux[],
438 const cs_real_t i_visc[],
439 const cs_real_t b_visc[],
440 const cs_real_t secvif[],
441 const cs_real_t secvib[],
442 cs_real_6_t viscel[],
443 const cs_real_2_t weighf[],
444 const cs_real_t weighb[],
445 int icvflb,
446 const int icvfli[],
447 cs_real_3_t i_pvar[],
448 cs_real_3_t b_pvar[],
449 cs_real_3_t smbr[]);
450
451/*----------------------------------------------------------------------------*/
452/*
453 * \brief Wrapper to the function which adds the explicit part of the
454 * convection/diffusion
455 * terms of a transport equation of a tensor field \f$ \tens{\varia} \f$.
456 *
457 * More precisely, the right hand side \f$ \vect{Rhs} \f$ is updated as
458 * follows:
459 * \f[
460 * \tens{Rhs} = \tens{Rhs} - \sum_{\fij \in \Facei{\celli}} \left(
461 * \dot{m}_\ij \left( \tens{\varia}_\fij - \tens{\varia}_\celli \right)
462 * - \mu_\fij \gradt_\fij \tens{\varia} \cdot \tens{S}_\ij \right)
463 * \f]
464 *
465 * Warning:
466 * - \f$ \tens{Rhs} \f$ has already been initialized before calling bilscts!
467 * - mind the sign minus
468 *
469 * Options for the convective scheme:
470 * - blencp = 0: upwind scheme for the advection
471 * - blencp = 1: no upwind scheme except in the slope test
472 * - ischcp = 0: SOLU
473 * - ischcp = 1: centered
474 * - ischcp = 2: SOLU with upwind gradient reconstruction
475 * - ischcp = 3: blending SOLU centered
476 * - ischcp = 4: NVD-TVD
477 *
478 * \param[in] idtvar indicator of the temporal scheme
479 * \param[in] f_id field id (or -1)
480 * \param[in] imasac take mass accumulation into account?
481 * \param[in] inc indicator
482 * \param[in] eqp pointer to a cs_equation_param_t structure which
483 * contains variable calculation options
484 * \param[in] pvar solved velocity (current time step)
485 * \param[in] pvara solved velocity (previous time step)
486 * \param[in] bc_coeffs_ts boundary condition structure for the variable
487 * \param[in] val_f boundary face value for gradient
488 * \param[in] flux boundary flux
489 * \param[in] i_massflux mass flux at interior faces
490 * \param[in] b_massflux mass flux at boundary faces
491 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
492 * at interior faces for the r.h.s.
493 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
494 * at boundary faces for the r.h.s.
495 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
496 * \param[in] weighf internal face weight between cells i j in case
497 * of tensor diffusion
498 * \param[in] weighb boundary face weight for cells i in case
499 * of tensor diffusion
500 * \param[in] icvflb global indicator of boundary convection flux
501 * - 0 upwind scheme at all boundary faces
502 * - 1 imposed flux at some boundary faces
503 * \param[in] icvfli boundary face indicator array of convection flux
504 * - 0 upwind scheme
505 * - 1 imposed flux
506 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
507 */
508/*----------------------------------------------------------------------------*/
509
510void
512 int f_id,
513 int imasac,
514 int inc,
516 cs_real_6_t pvar[],
517 const cs_real_6_t pvara[],
518 const cs_field_bc_coeffs_t *bc_coeffs_ts,
519 const cs_real_t val_f[][6],
520 const cs_real_t flux[][6],
521 const cs_real_t i_massflux[],
522 const cs_real_t b_massflux[],
523 const cs_real_t i_visc[],
524 const cs_real_t b_visc[],
525 cs_real_6_t viscel[],
526 const cs_real_2_t weighf[],
527 const cs_real_t weighb[],
528 int icvflb,
529 const int icvfli[],
530 cs_real_6_t rhs[]);
531
532/*----------------------------------------------------------------------------*/
533
535
536#endif /* __CS_BALANCE_H__ */
void cs_balance_tensor(int idtvar, int f_id, int imasac, int inc, cs_equation_param_t *eqp, cs_real_6_t pvar[], const cs_real_6_t pvara[], const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_t val_f[][6], const cs_real_t flux[][6], const cs_real_t i_massflux[], const cs_real_t b_massflux[], 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[], cs_real_6_t rhs[])
Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport...
Definition: cs_balance.cpp:917
void cs_balance_initialize(void)
Initialize balance timers.
Definition: cs_balance.cpp:116
void cs_balance_scalar(int idtvar, int f_id, int imucpp, int imasac, int inc, cs_equation_param_t *eqp, cs_real_t pvar[], const cs_real_t pvara[], 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_visc[], const cs_real_t b_visc[], cs_real_6_t viscel[], const cs_real_t xcpp[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_t rhs[])
Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport...
Definition: cs_balance.cpp:198
void cs_balance_vector(int idtvar, int f_id, int imasac, int inc, int ivisep, cs_equation_param_t *eqp, cs_real_3_t pvar[], const cs_real_3_t pvara[], const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_t val_f[][3], const cs_real_t flux[][3], const cs_real_t i_massflux[], const cs_real_t b_massflux[], 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_3_t i_pvar[], cs_real_3_t b_pvar[], cs_real_3_t smbr[])
Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport...
Definition: cs_balance.cpp:706
#define BEGIN_C_DECLS
Definition: cs_defs.h:554
double cs_real_t
Floating-point value.
Definition: cs_defs.h:357
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:374
cs_real_t cs_real_2_t[2]
vector of 2 floating-point values
Definition: cs_defs.h:373
cs_real_t cs_real_6_t[6]
vector of 6 floating-point values
Definition: cs_defs.h:376
#define END_C_DECLS
Definition: cs_defs.h:555
integer(c_int), pointer, save idtvar
option for a variable time step
Definition: optcal.f90:85
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition: cs_equation_param.h:194
Field boundary condition descriptor (for variables)
Definition: cs_field.h:120