7.1
general documentation
cs_cdofb_monolithic_sles.h
Go to the documentation of this file.
1 #ifndef __CS_CDOFB_MONOLITHIC_SLES_H__
2 #define __CS_CDOFB_MONOLITHIC_SLES_H__
3 
4 /*============================================================================
5  * Functions dedicated to the linear algebra settings and operations in case
6  * of CDO face-based schemes with a monolithic velocity-pressure coupling
7  *============================================================================*/
8 
9 /*
10  This file is part of Code_Saturne, a general-purpose CFD tool.
11 
12  Copyright (C) 1998-2021 EDF S.A.
13 
14  This program is free software; you can redistribute it and/or modify it under
15  the terms of the GNU General Public License as published by the Free Software
16  Foundation; either version 2 of the License, or (at your option) any later
17  version.
18 
19  This program is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
21  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
22  details.
23 
24  You should have received a copy of the GNU General Public License along with
25  this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
26  Street, Fifth Floor, Boston, MA 02110-1301, USA.
27 */
28 
29 /*----------------------------------------------------------------------------*/
30 
31 #include "cs_defs.h"
32 
33 /*----------------------------------------------------------------------------
34  * Standard C library headers
35  *----------------------------------------------------------------------------*/
36 
37 /*----------------------------------------------------------------------------
38  * Local headers
39  *----------------------------------------------------------------------------*/
40 
41 #include "cs_navsto_param.h"
42 
43 /*----------------------------------------------------------------------------*/
44 
46 
47 /*============================================================================
48  * Macro definitions
49  *============================================================================*/
50 
51 /*============================================================================
52  * Type definitions
53  *============================================================================*/
54 
55 /* Context related to the resolution of a saddle point problem */
56 typedef struct {
57 
58  /* Block matrices: The gradient operator is the -transpose of div_op */
59 
60  /*
61  * If n_row_blocks = 1 and div_op == NULL then all the velocity-pressure
62  * components are gathered in one block. In this case, the full matrix is a
63  * 1x1 matrix
64  *
65  * If n_row_blocks = 1 and div_op != NULL then all the velocity components
66  * are gathered in one block. In this case, the full matrix is a 2x2 matrix
67 
68  * If n_row_blocks = 3 and div_op == NULL then there is one block dedicated
69  * to each velocity component. In this case, the full matrix is a 4x4 matrix
70  *
71  * Please notice that div_op if not NULL is stored in a non-assembled way.
72  */
73 
75 
76  /* Blocks related to the velocity momentum */
78 
79  /* B*.approx(A^-1).B^t which corresponds to a compatible discretization of
80  the discrete Laplacian on the pressure space */
82 
83  cs_real_t *div_op; /* Block related to the -divergence (block
84  A_{10}) */
85 
86  /* Arrays split according to the block shape. U is interlaced or not
87  * according to the SLES strategy */
88 
89  cs_lnum_t n_faces; /* local number of DoFs for each component
90  * of the velocity */
91  cs_lnum_t n_cells; /* local number of DoFs for the pressure */
92 
93  cs_real_t *u_f; /* velocity values at faces */
94  cs_real_t *p_c; /* pressure values at cells */
95 
96  cs_real_t *b_f; /* RHS for the momentum (size = 3*n_faces) */
97  cs_real_t *b_c; /* RHS for the mass equation (size =
98  n_cells) */
99 
100  cs_sles_t *sles; /* main SLES structure */
101  cs_sles_t *schur_sles; /* auxiliary SLES for the Schur complement
102  * May be NULL */
103 
104  cs_real_t graddiv_coef; /* value of the grad-div coefficient in case
105  * of augmented system */
106 
108 
109 /*============================================================================
110  * Public function prototypes
111  *============================================================================*/
112 
113 /*----------------------------------------------------------------------------*/
119 /*----------------------------------------------------------------------------*/
120 
123 
124 /*----------------------------------------------------------------------------*/
132 /*----------------------------------------------------------------------------*/
133 
134 void
136  cs_lnum_t n_faces,
138 
139 /*----------------------------------------------------------------------------*/
145 /*----------------------------------------------------------------------------*/
146 
147 void
149 
150 /*----------------------------------------------------------------------------*/
156 /*----------------------------------------------------------------------------*/
157 
158 void
160 
161 /*----------------------------------------------------------------------------*/
167 /*----------------------------------------------------------------------------*/
168 
169 void
171 
172 /*----------------------------------------------------------------------------*/
180 /*----------------------------------------------------------------------------*/
181 
182 void
184  const cs_cdo_quantities_t *quant,
185  const cs_range_set_t *rset);
186 
187 /*----------------------------------------------------------------------------*/
198 /*----------------------------------------------------------------------------*/
199 
200 void
202  void *context);
203 
204 /*----------------------------------------------------------------------------*/
217 /*----------------------------------------------------------------------------*/
218 
219 int
221  const cs_equation_param_t *eqp,
223 
224 /*----------------------------------------------------------------------------*/
238 /*----------------------------------------------------------------------------*/
239 
240 int
242  const cs_equation_param_t *eqp,
244 
245 /*----------------------------------------------------------------------------*/
258 /*----------------------------------------------------------------------------*/
259 
260 int
262  const cs_equation_param_t *eqp,
264 
265 /*----------------------------------------------------------------------------*/
277 /*----------------------------------------------------------------------------*/
278 
279 int
281  const cs_equation_param_t *eqp,
283 
284 /*----------------------------------------------------------------------------*/
299 /*----------------------------------------------------------------------------*/
300 
301 int
303  const cs_equation_param_t *eqp,
305 
306 /*----------------------------------------------------------------------------*/
322 /*----------------------------------------------------------------------------*/
323 
324 int
326  const cs_equation_param_t *eqp,
328 
329 /*----------------------------------------------------------------------------*/
341 /*----------------------------------------------------------------------------*/
342 
343 int
345  const cs_equation_param_t *eqp,
347 
348 /*----------------------------------------------------------------------------*/
361 /*----------------------------------------------------------------------------*/
362 
363 int
365  const cs_equation_param_t *eqp,
367 
368 /*----------------------------------------------------------------------------*/
369 
371 
372 #endif /* __CS_CDOFB_MONOLITHIC_SLES_H__ */
int cs_cdofb_monolithic_uzawa_cg_solve(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Use the preconditioned Uzawa-CG algorithm to solve the saddle-point problem arising from CDO-Fb schem...
Definition: cs_cdofb_monolithic_sles.c:4766
void cs_cdofb_monolithic_sles_set_shared(const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_range_set_t *rset)
Set pointers to shared structures.
Definition: cs_cdofb_monolithic_sles.c:3955
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources...
Definition: cs_equation_param.h:177
int cs_cdofb_monolithic_uzawa_al_solve(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Use the Uzawa algorithm with an Augmented Lagrangian technique to solve the saddle-point problem aris...
Definition: cs_cdofb_monolithic_sles.c:5381
cs_real_t * b_c
Definition: cs_cdofb_monolithic_sles.h:97
cs_real_t * div_op
Definition: cs_cdofb_monolithic_sles.h:83
#define BEGIN_C_DECLS
Definition: cs_defs.h:510
int cs_cdofb_monolithic_gkb_solve(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Use the GKB algorithm to solve the saddle-point problem arising from CDO-Fb schemes for Stokes and Na...
Definition: cs_cdofb_monolithic_sles.c:4613
int cs_cdofb_monolithic_by_blocks_solve(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Solve a linear system arising from the discretization of the Navier-Stokes equation with a CDO face-b...
Definition: cs_cdo_connect.h:79
void cs_cdofb_monolithic_sles_reset(cs_cdofb_monolithic_sles_t *msles)
Reset to zero rhs and clean the cs_sles_t structure.
Definition: cs_cdofb_monolithic_sles.c:3871
cs_real_t graddiv_coef
Definition: cs_cdofb_monolithic_sles.h:104
Structure storing the parameters related to the resolution of the Navier-Stokes system.
Definition: cs_navsto_param.h:605
double cs_real_t
Floating-point value.
Definition: cs_defs.h:322
Definition: cs_cdo_quantities.h:129
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:93
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:68
void cs_cdofb_monolithic_sles_init(cs_lnum_t n_cells, cs_lnum_t n_faces, cs_cdofb_monolithic_sles_t *msles)
Allocate and initialize the rhs.
Definition: cs_cdofb_monolithic_sles.c:3838
cs_real_t * u_f
Definition: cs_cdofb_monolithic_sles.h:93
cs_sles_t * sles
Definition: cs_cdofb_monolithic_sles.h:100
Definition: cs_cdofb_monolithic_sles.h:56
cs_real_t * b_f
Definition: cs_cdofb_monolithic_sles.h:96
Definition: cs_range_set.h:57
int cs_cdofb_monolithic_uzawa_n3s_solve(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Use the preconditioned Uzawa-CG algorithm to solve the saddle-point problem arising from CDO-Fb schem...
Definition: cs_cdofb_monolithic_sles.c:5052
cs_matrix_t ** block_matrices
Definition: cs_cdofb_monolithic_sles.h:77
int cs_cdofb_monolithic_krylov_block_precond(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Solve a linear system arising from the discretization of the Navier-Stokes equation with a CDO face-b...
Definition: cs_cdofb_monolithic_sles.c:4432
void cs_cdofb_monolithic_set_sles(cs_navsto_param_t *nsp, void *context)
Start setting-up the Navier-Stokes equations when a monolithic algorithm is used to couple the system...
Definition: cs_cdofb_monolithic_sles.c:3981
cs_real_t * p_c
Definition: cs_cdofb_monolithic_sles.h:94
cs_sles_t * schur_sles
Definition: cs_cdofb_monolithic_sles.h:101
cs_cdofb_monolithic_sles_t * cs_cdofb_monolithic_sles_create(void)
Create an empty cs_cdofb_monolithic_sles_t structure.
Definition: cs_cdofb_monolithic_sles.c:3801
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:316
void cs_cdofb_monolithic_sles_clean(cs_cdofb_monolithic_sles_t *msles)
Free a part of the structure.
Definition: cs_cdofb_monolithic_sles.c:3905
#define END_C_DECLS
Definition: cs_defs.h:511
int n_row_blocks
Definition: cs_cdofb_monolithic_sles.h:74
int cs_cdofb_monolithic_solve(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Solve a linear system arising from the discretization of the Navier-Stokes equation with a CDO face-b...
Definition: cs_cdofb_monolithic_sles.c:4229
cs_lnum_t n_cells
Definition: cs_cdofb_monolithic_sles.h:91
cs_lnum_t n_faces
Definition: cs_cdofb_monolithic_sles.h:89
cs_matrix_t * compatible_laplacian
Definition: cs_cdofb_monolithic_sles.h:81
int cs_cdofb_monolithic_uzawa_al_incr_solve(const cs_navsto_param_t *nsp, const cs_equation_param_t *eqp, cs_cdofb_monolithic_sles_t *msles)
Use the Uzawa algorithm with an Augmented Lagrangian technique in an incremental way to solve the sad...
Definition: cs_cdofb_monolithic_sles.c:5529
void cs_cdofb_monolithic_sles_free(cs_cdofb_monolithic_sles_t **p_msles)
Free memory related to cs_cdofb_monolithic_sles_t structure.
Definition: cs_cdofb_monolithic_sles.c:3929