7.1
general documentation
cs_saddle_itsol.h
Go to the documentation of this file.
1 #ifndef __CS_SADDLE_ITSOL_H__
2 #define __CS_SADDLE_ITSOL_H__
3 
4 /*============================================================================
5  * In-house iterative solvers defined by blocks and associated to CDO
6  * discretizations for saddle-point system
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_iter_algo.h"
42 #include "cs_matrix.h"
43 #include "cs_mesh_adjacencies.h"
44 #include "cs_param_sles.h"
45 #include "cs_range_set.h"
46 #include "cs_sles.h"
47 
48 /*----------------------------------------------------------------------------*/
49 
51 
52 /*============================================================================
53  * Macro definitions
54  *============================================================================*/
55 
56 /*============================================================================
57  * Type definitions
58  *============================================================================*/
59 
60 typedef struct {
61 
62  /*
63  * A hybrid storage is used to represent the system described below
64  *
65  * | M11 | M12 | | x1 | |rhs1 |
66  * M = |-----------| |----| = |-----|
67  * | M21 | 0 | | x2 | |rhs2 |
68  *
69  * One assumes that M12 = M21^T
70  */
71 
72  int n_m11_matrices; /* 1, 3, 6 or 9 is expected */
74  cs_lnum_t x1_size; /* scatter view */
75  cs_lnum_t max_x1_size; /* max(x1_size, n_m11_cols) */
76 
78 
79  cs_lnum_t x2_size; /* scatter view */
80  cs_lnum_t max_x2_size; /* max(x2_size, n_m22_cols) */
82 
85 
86  /* Indexed list used to scan the unassembled m21 operator */
88 
89  /* Structure used for synchronisation (parallel or periodic). Enable to
90  switch from a scatter view (the mesh view) to a gather view (the algebraic
91  view). This structure is shared. */
93 
95 
96 /* Structure handling the block preconditioning */
97 
98 typedef struct {
99 
100  /* Parameters */
103 
104  /* Block 11 settings (mandatory) */
107 
108  /* Schur complement settings (optional) */
113 
114  /* Native arrays for the Schur matrix (optional) */
117 
118  /* Diagonal approximations of block matrices (optional) */
121 
123 
124 /*============================================================================
125  * Public function prototypes
126  *============================================================================*/
127 
128 /*----------------------------------------------------------------------------*/
139 /*----------------------------------------------------------------------------*/
140 
143  cs_param_schur_approx_t schur_type,
144  cs_param_sles_t *m11_slesp,
145  cs_sles_t *m11_sles);
146 
147 /*----------------------------------------------------------------------------*/
153 /*----------------------------------------------------------------------------*/
154 
155 void
157 
158 /*----------------------------------------------------------------------------*/
172 /*----------------------------------------------------------------------------*/
173 
174 void
177  cs_real_t *x1,
178  cs_real_t *x2,
179  cs_iter_algo_t *ia);
180 
181 /*----------------------------------------------------------------------------*/
198 /*----------------------------------------------------------------------------*/
199 
200 void
201 cs_saddle_gcr(int restart,
202  cs_saddle_system_t *ssys,
204  cs_real_t *x1,
205  cs_real_t *x2,
206  cs_iter_algo_t *ia);
207 
208 /*----------------------------------------------------------------------------*/
224 /*----------------------------------------------------------------------------*/
225 
226 void
228  const cs_matrix_t *mat,
229  cs_real_t *vec,
230  cs_real_t *matvec);
231 
232 /*----------------------------------------------------------------------------*/
249 /*----------------------------------------------------------------------------*/
250 
251 void
253  const cs_matrix_t *mat,
254  cs_lnum_t vec_len,
255  cs_real_t *vec,
256  cs_real_t **p_matvec);
257 
258 /*----------------------------------------------------------------------------*/
259 
261 
262 #endif /* __CS_SADDLE_ITSOL_H__ */
cs_lnum_t max_x2_size
Definition: cs_saddle_itsol.h:80
cs_real_t * m11_inv_diag
Definition: cs_saddle_itsol.h:119
cs_param_sles_t * m11_slesp
Definition: cs_saddle_itsol.h:105
cs_real_t * mass22_diag
Definition: cs_saddle_itsol.h:120
cs_matrix_t * schur_matrix
Definition: cs_saddle_itsol.h:109
void cs_saddle_minres(cs_saddle_system_t *ssys, cs_saddle_block_precond_t *sbp, cs_real_t *x1, cs_real_t *x2, cs_iter_algo_t *ia)
Apply the MINRES algorithm to a saddle point problem (the system is stored in a hybrid way)...
Definition: cs_saddle_itsol.c:1567
Definition: cs_mesh_adjacencies.h:68
const cs_range_set_t * rset
Definition: cs_saddle_itsol.h:92
cs_sles_t * m11_sles
Definition: cs_saddle_itsol.h:106
cs_lnum_t max_x1_size
Definition: cs_saddle_itsol.h:75
void cs_saddle_block_precond_free(cs_saddle_block_precond_t **p_sbp)
Free a cs_saddle_block_precond_t structure.
Definition: cs_saddle_itsol.c:1531
#define BEGIN_C_DECLS
Definition: cs_defs.h:510
cs_real_t * schur_xtra
Definition: cs_saddle_itsol.h:116
cs_matrix_t ** m11_matrices
Definition: cs_saddle_itsol.h:73
cs_real_t * rhs2
Definition: cs_saddle_itsol.h:81
cs_param_precond_block_t
Definition: cs_param_types.h:688
cs_lnum_t x1_size
Definition: cs_saddle_itsol.h:74
Structure and routines handling the SLES settings stored inside a cs_param_sles_t structure...
Definition: cs_field_pointer.h:223
double cs_real_t
Floating-point value.
Definition: cs_defs.h:322
cs_sles_t * schur_sles
Definition: cs_saddle_itsol.h:111
cs_param_sles_t * schur_slesp
Definition: cs_saddle_itsol.h:110
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:93
Set of common parameters to manage an iterative algorithm.
Definition: cs_iter_algo.h:101
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:68
cs_param_schur_approx_t
Strategy to build the Schur complement approximation. This appears in block preconditioning or uzawa ...
Definition: cs_param_types.h:621
void cs_saddle_gcr(int restart, cs_saddle_system_t *ssys, cs_saddle_block_precond_t *sbp, cs_real_t *x1, cs_real_t *x2, cs_iter_algo_t *ia)
Apply the GCR algorithm to a saddle point problem (the system is stored in a hybrid way)...
Definition: cs_saddle_itsol.c:1772
Definition: cs_field_pointer.h:222
cs_real_t * schur_diag
Definition: cs_saddle_itsol.h:115
cs_saddle_block_precond_t * cs_saddle_block_precond_create(cs_param_precond_block_t block_type, cs_param_schur_approx_t schur_type, cs_param_sles_t *m11_slesp, cs_sles_t *m11_sles)
Create and initialize a cs_saddle_block_precond_t structure.
Definition: cs_saddle_itsol.c:1484
cs_adjacency_t * m21_adjacency
Definition: cs_saddle_itsol.h:87
double schur_scaling
Definition: cs_saddle_itsol.h:112
Definition: cs_saddle_itsol.h:98
void cs_matrix_vector_multiply_gs(const cs_range_set_t *rset, const cs_matrix_t *mat, cs_lnum_t vec_len, cs_real_t *vec, cs_real_t **p_matvec)
Perform a matrix-vector multiplication in case of scatter-view array as input parameter. Thus, one performs a scatter –> gather (before the multiplication) and a gather –> scatter operation after the multiplication. The output parameter matvec is not allocated. A check on the size is done for the input array.
Definition: cs_saddle_itsol.c:2017
Definition: cs_range_set.h:57
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:316
void cs_matrix_vector_multiply_gs_allocated(const cs_range_set_t *rset, const cs_matrix_t *mat, cs_real_t *vec, cs_real_t *matvec)
Perform a matrix-vector multiplication in case of scatter-view array as input parameter. Thus, one performs a scatter –> gather (before the multiplication) and a gather –> scatter operation after the multiplication. One assumes that matvec is allocated to the right size. No check is done.
Definition: cs_saddle_itsol.c:1964
#define END_C_DECLS
Definition: cs_defs.h:511
cs_real_t * rhs1
Definition: cs_saddle_itsol.h:77
cs_lnum_t x2_size
Definition: cs_saddle_itsol.h:79
cs_param_schur_approx_t schur_type
Definition: cs_saddle_itsol.h:102
int n_m11_matrices
Definition: cs_saddle_itsol.h:72
cs_param_precond_block_t block_type
Definition: cs_saddle_itsol.h:101
Structure storing all metadata related to the resolution of a linear system with an iterative solver...
Definition: cs_param_sles.h:62
cs_real_t * m21_unassembled
Definition: cs_saddle_itsol.h:83
int m21_stride
Definition: cs_saddle_itsol.h:84
Definition: cs_saddle_itsol.h:60