8.2
general documentation
cs_matrix_building.h
Go to the documentation of this file.
1 #ifndef __CS_MATRIX_BUILDING_H__
2 #define __CS_MATRIX_BUILDING_H__
3 
4 /*============================================================================
5  * Matrix building
6  *============================================================================*/
7 
8 /*
9  This file is part of code_saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2024 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_base.h"
35 #include "cs_halo.h"
36 
37 /*----------------------------------------------------------------------------*/
38 
40 
41 /*=============================================================================
42  * Local Macro definitions
43  *============================================================================*/
44 
45 /*============================================================================
46  * Type definition
47  *============================================================================*/
48 
49 /*============================================================================
50  * Global variables
51  *============================================================================*/
52 
53 /*============================================================================
54  * Public function prototypes for Fortran API
55  *============================================================================*/
56 
57 /*----------------------------------------------------------------------------
58  * Fortran wrapper to cs_matrix_scalar (or its counterpart for
59  * symmetric matrices)
60  *----------------------------------------------------------------------------*/
61 
62 void CS_PROCF (matrix, MATRIX)
63 (
64  const int *iconvp,
65  const int *idiffp,
66  const int *ndircp,
67  const int *isym,
68  const cs_real_t *thetap,
69  const int *imucpp,
70  const cs_real_t coefbp[],
71  const cs_real_t cofbfp[],
72  const cs_real_t rovsdt[],
73  const cs_real_t i_massflux[],
74  const cs_real_t b_massflux[],
75  const cs_real_t i_visc[],
76  const cs_real_t b_visc[],
77  const cs_real_t xcpp[],
78  cs_real_t da[],
79  cs_real_t xa[]
80 );
81 
82 /*=============================================================================
83  * Public function prototypes
84  *============================================================================*/
85 
86 /*----------------------------------------------------------------------------
87  * Wrapper to cs_matrix_scalar (or its counterpart for
88  * symmetric matrices)
89  *----------------------------------------------------------------------------*/
90 
91 void
92 cs_matrix_wrapper_scalar(int iconvp,
93  int idiffp,
94  int ndircp,
95  int isym,
96  double thetap,
97  int imucpp,
98  const cs_field_bc_coeffs_t *bc_coeffs,
99  const cs_real_t rovsdt[],
100  const cs_real_t i_massflux[],
101  const cs_real_t b_massflux[],
102  const cs_real_t i_visc[],
103  const cs_real_t b_visc[],
104  const cs_real_t xcpp[],
105  cs_real_t da[],
106  cs_real_t xa[]);
107 
108 /*----------------------------------------------------------------------------
109  * Wrapper to cs_matrix_vector (or its counterpart for
110  * symmetric matrices)
111  *----------------------------------------------------------------------------*/
112 
113 void
114 cs_matrix_wrapper_vector(int iconvp,
115  int idiffp,
116  int tensorial_diffusion,
117  int ndircp,
118  int isym,
119  cs_lnum_t eb_size,
120  double thetap,
121  const cs_field_bc_coeffs_t *bc_coeffs_v,
122  const cs_real_33_t fimp[],
123  const cs_real_t i_massflux[],
124  const cs_real_t b_massflux[],
125  const cs_real_t i_visc[],
126  const cs_real_t b_visc[],
127  cs_real_33_t da[],
128  cs_real_t xa[]);
129 
130 /*----------------------------------------------------------------------------
131  * Wrapper to cs_matrix_tensor (or its counterpart for
132  * symmetric matrices)
133  *----------------------------------------------------------------------------*/
134 
135 void
136 cs_matrix_wrapper_tensor(int iconvp,
137  int idiffp,
138  int tensorial_diffusion,
139  int ndircp,
140  int isym,
141  double thetap,
142  const cs_field_bc_coeffs_t *bc_coeffs_ts,
143  const cs_real_66_t fimp[],
144  const cs_real_t i_massflux[],
145  const cs_real_t b_massflux[],
146  const cs_real_t i_visc[],
147  const cs_real_t b_visc[],
148  cs_real_66_t da[],
149  cs_real_t xa[]);
150 
151 /*----------------------------------------------------------------------------*/
181 /*----------------------------------------------------------------------------*/
182 
183 void
185  int idiffp,
186  double thetap,
187  const cs_field_bc_coeffs_t *bc_coeffs,
188  const cs_real_t rovsdt[],
189  const cs_real_t i_visc[],
190  const cs_real_t b_visc[],
191  cs_real_t *restrict da,
192  cs_real_t *restrict xa);
193 
194 /*----------------------------------------------------------------------------*/
224 /*----------------------------------------------------------------------------*/
225 
226 void
228  int idiffp,
229  double thetap,
230  const cs_field_bc_coeffs_t *bc_coeffs_ts,
231  const cs_real_66_t fimp[],
232  const cs_real_t i_visc[],
233  const cs_real_t b_visc[],
235  cs_real_t *restrict xa);
236 
237 /*----------------------------------------------------------------------------*/
274 /*----------------------------------------------------------------------------*/
275 
276 void
277 cs_matrix_tensor(const cs_mesh_t *m,
278  int iconvp,
279  int idiffp,
280  double thetap,
281  const cs_field_bc_coeffs_t *bc_coeffs_ts,
282  const cs_real_66_t fimp[],
283  const cs_real_t i_massflux[],
284  const cs_real_t b_massflux[],
285  const cs_real_t i_visc[],
286  const cs_real_t b_visc[],
288  cs_real_2_t *restrict xa);
289 
290 /*----------------------------------------------------------------------------*/
327 /*----------------------------------------------------------------------------*/
328 
329 void
330 cs_matrix_scalar(const cs_mesh_t *m,
331  int iconvp,
332  int idiffp,
333  double thetap,
334  int imucpp,
335  const cs_field_bc_coeffs_t *bc_coeffs,
336  const cs_real_t rovsdt[],
337  const cs_real_t i_massflux[],
338  const cs_real_t b_massflux[],
339  const cs_real_t i_visc[],
340  const cs_real_t b_visc[],
341  const cs_real_t xcpp[],
342  cs_real_t *restrict da,
343  cs_real_2_t *restrict xa);
344 
345 /*----------------------------------------------------------------------------*/
375 /*----------------------------------------------------------------------------*/
376 
377 void
379  int idiffp,
380  double thetap,
381  const cs_field_bc_coeffs_t *bc_coeffs_v,
382  const cs_real_33_t fimp[],
383  const cs_real_t i_visc[],
384  const cs_real_t b_visc[],
386  cs_real_t *restrict xa);
387 
388 /*----------------------------------------------------------------------------*/
426 /*----------------------------------------------------------------------------*/
427 
428 void
429 cs_matrix_vector(const cs_mesh_t *m,
430  const cs_mesh_quantities_t *mq,
431  int iconvp,
432  int idiffp,
433  cs_lnum_t eb_size,
434  double thetap,
435  const cs_field_bc_coeffs_t *bc_coeffs_v,
436  const cs_real_33_t fimp[],
437  const cs_real_t i_massflux[],
438  const cs_real_t b_massflux[],
439  const cs_real_t i_visc[],
440  const cs_real_t b_visc[],
442  cs_real_2_t *restrict xa);
443 
444 /*----------------------------------------------------------------------------*/
468 /*----------------------------------------------------------------------------*/
469 
470 void
472  int iconvp,
473  int idiffp,
474  int isym,
475  const cs_field_bc_coeffs_t *bc_coeffs,
476  const cs_real_t i_massflux[],
477  const cs_real_t b_massflux[],
478  const cs_real_t i_visc[],
479  const cs_real_t b_visc[],
480  cs_real_t *restrict da);
481 
482 /*----------------------------------------------------------------------------*/
517 /*----------------------------------------------------------------------------*/
518 
519 void
521  const cs_mesh_quantities_t *mq,
522  int iconvp,
523  int idiffp,
524  double thetap,
525  const cs_field_bc_coeffs_t *bc_coeffs_v,
526  const cs_real_33_t fimp[],
527  const cs_real_t i_massflux[],
528  const cs_real_t b_massflux[],
529  const cs_real_33_t i_visc[],
530  const cs_real_t b_visc[],
532  cs_real_332_t *restrict xa);
533 
534 /*----------------------------------------------------------------------------*/
568 /*----------------------------------------------------------------------------*/
569 
570 void
572  int iconvp,
573  int idiffp,
574  double thetap,
575  const cs_field_bc_coeffs_t *bc_coeffs_ts,
576  const cs_real_66_t fimp[],
577  const cs_real_t i_massflux[],
578  const cs_real_t b_massflux[],
579  const cs_real_66_t i_visc[],
580  const cs_real_t b_visc[],
582  cs_real_662_t *restrict xa);
583 
584 /*----------------------------------------------------------------------------*/
612 /*----------------------------------------------------------------------------*/
613 
614 void
616  int idiffp,
617  double thetap,
618  const cs_field_bc_coeffs_t *bc_coeffs_v,
619  const cs_real_33_t fimp[],
620  const cs_real_33_t i_visc[],
621  const cs_real_t b_visc[],
623  cs_real_33_t *restrict xa);
624 
625 /*----------------------------------------------------------------------------*/
655 /*----------------------------------------------------------------------------*/
656 
657 void
659 (const cs_mesh_t *m,
660  int idiffp,
661  double thetap,
662  const cs_field_bc_coeffs_t *bc_coeffs_ts,
663  const cs_real_66_t fimp[],
664  const cs_real_66_t i_visc[],
665  const cs_real_t b_visc[],
667  cs_real_66_t *restrict xa);
668 
669 /*----------------------------------------------------------------------------*/
670 
672 
673 #endif /* __CS_MATRIX_BUILDING_H__ */
#define restrict
Definition: cs_defs.h:141
#define BEGIN_C_DECLS
Definition: cs_defs.h:528
double cs_real_t
Floating-point value.
Definition: cs_defs.h:332
cs_real_33_t cs_real_332_t[2]
vector of 2 3x3 matrices of floating-point values
Definition: cs_defs.h:368
#define CS_PROCF(x, y)
Definition: cs_defs.h:560
cs_real_t cs_real_66_t[6][6]
6x6 matrix of floating-point values
Definition: cs_defs.h:357
cs_real_66_t cs_real_662_t[2]
Definition: cs_defs.h:370
cs_real_t cs_real_2_t[2]
vector of 2 floating-point values
Definition: cs_defs.h:346
#define END_C_DECLS
Definition: cs_defs.h:529
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition: cs_defs.h:356
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:325
void cs_matrix_anisotropic_diffusion_tensor(const cs_mesh_t *m, int iconvp, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_66_t fimp[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_66_t i_visc[], const cs_real_t b_visc[], cs_real_66_t *restrict da, cs_real_662_t *restrict xa)
Build the advection/diffusion matrix for a tensor field with a tensorial diffusivity.
Definition: cs_matrix_building.c:1875
void cs_sym_matrix_vector(const cs_mesh_t *m, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_33_t fimp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_33_t *restrict da, cs_real_t *restrict xa)
Build the diffusion matrix for a vector field (symmetric matrix).
Definition: cs_matrix_building.c:827
void cs_matrix_anisotropic_diffusion(const cs_mesh_t *m, const cs_mesh_quantities_t *mq, int iconvp, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_33_t fimp[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_33_t i_visc[], const cs_real_t b_visc[], cs_real_33_t *restrict da, cs_real_332_t *restrict xa)
Build the advection/diffusion matrix for a vector field with a tensorial diffusivity.
Definition: cs_matrix_building.c:1664
void cs_matrix_vector(const cs_mesh_t *m, const cs_mesh_quantities_t *mq, int iconvp, int idiffp, cs_lnum_t eb_size, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_33_t fimp[], 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_33_t *restrict da, cs_real_2_t *restrict xa)
Build the advection/diffusion matrix for a vector field (non-symmetric matrix).
Definition: cs_matrix_building.c:1077
void matrix(const int *iconvp, const int *idiffp, const int *ndircp, const int *isym, const cs_real_t *thetap, const int *imucpp, const cs_real_t coefbp[], const cs_real_t cofbfp[], const cs_real_t rovsdt[], 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 xcpp[], cs_real_t da[], cs_real_t xa[])
void cs_matrix_time_step(const cs_mesh_t *m, int iconvp, int idiffp, int isym, 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_t *restrict da)
Build the diagonal of the advection/diffusion matrix for determining the variable time step,...
Definition: cs_matrix_building.c:1506
void cs_sym_matrix_anisotropic_diffusion(const cs_mesh_t *m, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_33_t fimp[], const cs_real_33_t i_visc[], const cs_real_t b_visc[], cs_real_33_t *restrict da, cs_real_33_t *restrict xa)
Build the diffusion matrix for a vector field with a tensorial diffusivity (symmetric matrix).
Definition: cs_matrix_building.c:2034
void cs_matrix_wrapper_vector(int iconvp, int idiffp, int tensorial_diffusion, int ndircp, int isym, cs_lnum_t eb_size, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_33_t fimp[], 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_33_t da[], cs_real_t xa[])
Definition: cs_matrix_building.c:198
void cs_matrix_wrapper_scalar(int iconvp, int idiffp, int ndircp, int isym, double thetap, int imucpp, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t rovsdt[], 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 xcpp[], cs_real_t da[], cs_real_t xa[])
Definition: cs_matrix_building.c:111
void cs_matrix_scalar(const cs_mesh_t *m, int iconvp, int idiffp, double thetap, int imucpp, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t rovsdt[], 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 xcpp[], cs_real_t *restrict da, cs_real_2_t *restrict xa)
Build the advection/diffusion matrix for a scalar field.
Definition: cs_matrix_building.c:606
void cs_sym_matrix_tensor(const cs_mesh_t *m, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_66_t fimp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_66_t *restrict da, cs_real_t *restrict xa)
Build the diffusion matrix for a tensor field (symmetric matrix).
Definition: cs_matrix_building.c:948
void cs_sym_matrix_anisotropic_diffusion_tensor(const cs_mesh_t *m, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_66_t fimp[], const cs_real_66_t i_visc[], const cs_real_t b_visc[], cs_real_66_t *restrict da, cs_real_66_t *restrict xa)
Build the diffusion matrix for a vector field with a tensorial diffusivity (symmetric matrix).
Definition: cs_matrix_building.c:2164
void cs_matrix_tensor(const cs_mesh_t *m, int iconvp, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_66_t fimp[], 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_66_t *restrict da, cs_real_2_t *restrict xa)
Build the advection/diffusion matrix for a tensor field (non-symmetric matrix).
Definition: cs_matrix_building.c:1361
void cs_matrix_wrapper_tensor(int iconvp, int idiffp, int tensorial_diffusion, int ndircp, int isym, double thetap, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_66_t fimp[], 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_66_t da[], cs_real_t xa[])
Definition: cs_matrix_building.c:323
void cs_sym_matrix_scalar(const cs_mesh_t *m, int idiffp, double thetap, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t rovsdt[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t *restrict da, cs_real_t *restrict xa)
Build the diffusion matrix for a scalar field. (symmetric matrix).
Definition: cs_matrix_building.c:469
Field boundary condition descriptor (for variables)
Definition: cs_field.h:104
Definition: cs_mesh_quantities.h:92
Definition: cs_mesh.h:85