8.3
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#include "cs_defs.h"
31#include "cs_math.h"
32#include "cs_parameters.h" // for BC types
33
34/*----------------------------------------------------------------------------
35 * Local headers
36 *----------------------------------------------------------------------------*/
37
38#include "cs_base.h"
39#include "cs_halo.h"
40
41/*----------------------------------------------------------------------------*/
42
43/*=============================================================================
44 * Local Macro definitions
45 *============================================================================*/
46
47/*============================================================================
48 * Type definition
49 *============================================================================*/
50
51/*============================================================================
52 * Global variables
53 *============================================================================*/
54
55/*=============================================================================
56 * Public function prototypes
57 *============================================================================*/
58
59#if defined(__cplusplus)
60
61/*----------------------------------------------------------------------------*/
62/*
63 * \brief Build the diagonal of the advection/diffusion matrix
64 * for determining the variable time step, flow, Fourier.
65 *
66 * \param[in, out] a pointer to matrix structure
67 * \param[in] f pointer to field, or null
68 * \param[in] iconvp indicator
69 * - 1 advection
70 * - 0 otherwise
71 * \param[in] idiffp indicator
72 * - 1 diffusion
73 * - 0 otherwise
74 * \param[in] ndircp number of Dirichlet BCs
75 * \param[in] thetap time scheme parameter
76 * \param[in] relaxp relaxation coefficient (if < 1)
77 * \param[in] imucp 1 for temperature (with Cp), 0 otherwise
78 * \param[in] bc_coeffs boundary condition structure
79 * \param[in] rovsdt implicit terms (rho / dt)
80 * \param[in] i_massflux mass flux at interior faces
81 * \param[in] b_massflux mass flux at border faces
82 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
83 * at interior faces for the matrix
84 * \param[in] b_visc \f$ S_\fib \f$
85 * at border faces for the matrix
86 * \param[in] xcpp Cp per cell, or null
87 */
88/*----------------------------------------------------------------------------*/
89
90void
92 const cs_field_t *f,
93 int iconvp,
94 int idiffp,
95 int ndircp,
96 double thetap,
97 double relaxp,
98 int imucpp,
99 const cs_field_bc_coeffs_t *bc_coeffs,
100 const cs_real_t rovsdt[],
101 const cs_real_t i_massflux[],
102 const cs_real_t b_massflux[],
103 const cs_real_t i_visc[],
104 const cs_real_t b_visc[],
105 const cs_real_t xcpp[]);
106
107/*----------------------------------------------------------------------------*/
137/*----------------------------------------------------------------------------*/
138
139template <cs_lnum_t stride>
140void
142(
143 cs_matrix_t *a,
144 const cs_field_t *f,
145 int iconvp,
146 int idiffp,
147 int tensorial_diffusion,
148 int ndircp,
149 cs_lnum_t eb_size,
150 double thetap,
151 double relaxp,
152 const cs_field_bc_coeffs_t *bc_coeffs,
153 const cs_real_t fimp[][stride][stride],
154 const cs_real_t i_massflux[],
155 const cs_real_t b_massflux[],
156 const cs_real_t i_visc[],
157 const cs_real_t b_visc[]
158);
159
160/*----------------------------------------------------------------------------
161 * Compute legacy matrix coefficients
162 *----------------------------------------------------------------------------*/
163
164void
165cs_matrix_wrapper(int iconvp,
166 int idiffp,
167 int ndircp,
168 int isym,
169 double thetap,
170 int imucpp,
171 const cs_field_bc_coeffs_t *bc_coeffs,
172 const cs_real_t rovsdt[],
173 const cs_real_t i_massflux[],
174 const cs_real_t b_massflux[],
175 const cs_real_t i_visc[],
176 const cs_real_t b_visc[],
177 const cs_real_t xcpp[],
178 cs_real_t da[],
179 cs_real_t xa[]);
180
181/*----------------------------------------------------------------------------
182 * Compute legacy matrix coefficients
183 *----------------------------------------------------------------------------*/
184
185template <cs_lnum_t stride>
186void
187cs_matrix_wrapper(int iconvp,
188 int idiffp,
189 int tensorial_diffusion,
190 int ndircp,
191 int isym,
192 cs_lnum_t eb_size,
193 double thetap,
194 const cs_field_bc_coeffs_t *bc_coeffs_v,
195 const cs_real_t fimp[][stride][stride],
196 const cs_real_t i_massflux[],
197 const cs_real_t b_massflux[],
198 const cs_real_t i_visc[],
199 const cs_real_t b_visc[],
200 cs_real_t da[][stride][stride],
201 cs_real_t xa[]);
202
203/*----------------------------------------------------------------------------*/
204
205#endif //defined(__cplusplus)
206
208
209/*----------------------------------------------------------------------------*/
210/*
211 * \brief Build the diagonal of the advection/diffusion matrix
212 * for determining the variable time step, flow, Fourier.
213 *
214 * \param[in] m pointer to mesh structure
215 * \param[in] iconvp indicator
216 * - 1 advection
217 * - 0 otherwise
218 * \param[in] idiffp indicator
219 * - 1 diffusion
220 * - 0 otherwise
221 * \param[in] isym indicator
222 * - 1 symmetric matrix
223 * - 2 non symmmetric matrix
224 * \param[in] bc_coeffs boundary condition structure for the variable
225 * \param[in] i_massflux mass flux at interior faces
226 * \param[in] b_massflux mass flux at border faces
227 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
228 * at interior faces for the matrix
229 * \param[in] b_visc \f$ S_\fib \f$
230 * at border faces for the matrix
231 * \param[out] da diagonal part of the matrix
232 */
233/*----------------------------------------------------------------------------*/
234
235void
237 int iconvp,
238 int idiffp,
239 int isym,
240 const cs_field_bc_coeffs_t *bc_coeffs,
241 const cs_real_t i_massflux[],
242 const cs_real_t b_massflux[],
243 const cs_real_t i_visc[],
244 const cs_real_t b_visc[],
245 cs_real_t *da);
246
247/*----------------------------------------------------------------------------*/
248
250
251#endif /* __CS_MATRIX_BUILDING_H__ */
#define BEGIN_C_DECLS
Definition: cs_defs.h:542
double cs_real_t
Floating-point value.
Definition: cs_defs.h:342
#define END_C_DECLS
Definition: cs_defs.h:543
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:335
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:110
void cs_matrix_wrapper(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.cpp:2668
void cs_matrix_compute_coeffs(cs_matrix_t *a, const cs_field_t *f, int iconvp, int idiffp, int ndircp, double thetap, double relaxp, 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[])
Build the diagonal of the advection/diffusion matrix for determining the variable time step,...
Definition: cs_matrix_building.cpp:1967
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 *da)
Field boundary condition descriptor (for variables)
Definition: cs_field.h:104
Field descriptor.
Definition: cs_field.h:131
Definition: cs_mesh.h:85