My Project
programmer's documentation
cs_sles_petsc.h
Go to the documentation of this file.
1 #ifndef __CS_SLES_PETSC_H__
2 #define __CS_SLES_PETSC_H__
3 
4 /*============================================================================
5  * Sparse Linear Equation Solvers using PETSc
6  *============================================================================*/
7 
8 /*
9  This file is part of Code_Saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2019 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  * PETSc headers
32  *----------------------------------------------------------------------------*/
33 
34 #include <petscviewer.h>
35 #include <petscksp.h>
36 
37 /*----------------------------------------------------------------------------
38  * Local headers
39  *----------------------------------------------------------------------------*/
40 
41 #include "cs_base.h"
42 #include "cs_halo_perio.h"
43 #include "cs_matrix.h"
44 #include "cs_time_plot.h"
45 
46 /*----------------------------------------------------------------------------*/
47 
49 
50 /*============================================================================
51  * Macro definitions
52  *============================================================================*/
53 
54 /*============================================================================
55  * Type definitions
56  *============================================================================*/
57 
58 /*----------------------------------------------------------------------------
59  * Function pointer for user settings of a PETSc KSP solver setup.
60  *
61  * This function is called the end of the setup stage for a KSP solver.
62  *
63  * Note that using the advanced KSPSetPostSolve and KSPSetPreSolve functions,
64  * this also allows setting furthur function pointers for pre and post-solve
65  * operations (see the PETSc documentation).
66  *
67  * Note: if the context pointer is non-NULL, it must point to valid data
68  * when the selection function is called so that value or structure should
69  * not be temporary (i.e. local);
70  *
71  * parameters:
72  * context <-> pointer to optional (untyped) value or structure
73  * a <-> PETSc matrix context
74  * ksp <-> pointer to PETSc KSP context
75  *----------------------------------------------------------------------------*/
76 
77 typedef void
78 (cs_sles_petsc_setup_hook_t) (void *context,
79  Mat a,
80  KSP ksp);
81 
82 /* Iterative linear solver context (opaque) */
83 
84 typedef struct _cs_sles_petsc_t cs_sles_petsc_t;
85 
86 /*============================================================================
87  * Global variables
88  *============================================================================*/
89 
90 /*=============================================================================
91  * User function prototypes
92  *============================================================================*/
93 
94 /*----------------------------------------------------------------------------
95  * Function pointer for user settings of a PETSc KSP solver setup.
96  *
97  * This function is called the end of the setup stage for a KSP solver.
98  *
99  * Note that using the advanced KSPSetPostSolve and KSPSetPreSolve functions,
100  * this also allows setting furthur function pointers for pre and post-solve
101  * operations (see the PETSc documentation).
102  *
103  * Note: if the context pointer is non-NULL, it must point to valid data
104  * when the selection function is called so that value or structure should
105  * not be temporary (i.e. local);
106  *
107  * parameters:
108  * context <-> pointer to optional (untyped) value or structure
109  * a <-> PETSc matrix context
110  * ksp <-> pointer to PETSc KSP context
111  *----------------------------------------------------------------------------*/
112 
113 void
114 cs_user_sles_petsc_hook(void *context,
115  Mat a,
116  KSP ksp);
117 
118 /*=============================================================================
119  * Static inline public function prototypes
120  *============================================================================*/
121 
122 /*----------------------------------------------------------------------------
123  * Initialized PETSc if needed
124  *----------------------------------------------------------------------------*/
125 
126 static inline void
128 {
129  /* Initialization must be called before setting options;
130  it does not need to be called before calling
131  cs_sles_petsc_define(), as this is handled automatically. */
132 
133  PetscBool is_initialized;
134  PetscInitialized(&is_initialized);
135 
136  if (is_initialized == PETSC_FALSE) {
137 #if defined(HAVE_MPI)
138  PETSC_COMM_WORLD = cs_glob_mpi_comm;
139 #endif
140  PetscInitializeNoArguments();
141  }
142 }
143 
144 /*----------------------------------------------------------------------------
145  * \brief Output the settings of a KSP structure
146  *
147  * \param[in] ksp Krylov SubSpace structure
148  *----------------------------------------------------------------------------*/
149 
150 static inline void
152 {
153  PetscViewer v;
154 
155  PetscViewerCreate(PETSC_COMM_WORLD, &v);
156  PetscViewerSetType(v, PETSCVIEWERASCII);
157  PetscViewerFileSetMode(v, FILE_MODE_APPEND);
158  PetscViewerFileSetName(v, "petsc_setup.log");
159 
160  KSPView(ksp, v);
161  PetscViewerDestroy(&v);
162 }
163 
164 /*=============================================================================
165  * Public function prototypes
166  *============================================================================*/
167 
168 /*----------------------------------------------------------------------------
169  * Define and associate a PETSc linear system solver
170  * for a given field or equation name.
171  *
172  * If this system did not previously exist, it is added to the list of
173  * "known" systems. Otherwise, its definition is replaced by the one
174  * defined here.
175  *
176  * This is a utility function: if finer control is needed, see
177  * cs_sles_define() and cs_sles_petsc_create().
178  *
179  * In case of rotational periodicity for a block (non-scalar) matrix,
180  * the matrix type will be forced to MATSHELL ("shell") regardless
181  * of the option used.
182  *
183  * Note that this function returns a pointer directly to the iterative solver
184  * management structure. This may be used to set further options.
185  * If needed, cs_sles_find() may be used to obtain a pointer to the matching
186  * cs_sles_t container.
187  *
188  * parameters:
189  * f_id <-- associated field id, or < 0
190  * name <-- associated name if f_id < 0, or NULL
191  * matrix_type <-- PETSc matrix type
192  * setup_hook <-- pointer to optional setup epilogue function
193  * context <-> pointer to optional (untyped) value or structure
194  * for setup_hook, or NULL
195  *
196  * returns:
197  * pointer to newly created iterative solver info object.
198  *----------------------------------------------------------------------------*/
199 
202  const char *name,
203  MatType matrix_type,
204  cs_sles_petsc_setup_hook_t *setup_hook,
205  void *context);
206 
207 /*----------------------------------------------------------------------------
208  * Create PETSc linear system solver info and context.
209  *
210  * In case of rotational periodicity for a block (non-scalar) matrix,
211  * the matrix type will be forced to MATSHELL ("shell") regardless
212  * of the option used.
213  *
214  * parameters:
215  * matrix_type <-- PETSc matrix type
216  * setup_hook <-- pointer to optional setup epilogue function
217  * context <-> pointer to optional (untyped) value or structure
218  * for setup_hook, or NULL
219  *
220  * returns:
221  * pointer to newly created solver info object.
222  *----------------------------------------------------------------------------*/
223 
225 cs_sles_petsc_create(MatType matrix_type,
226  cs_sles_petsc_setup_hook_t *setup_hook,
227  void *context);
228 
229 /*----------------------------------------------------------------------------
230  * Create PETSc linear system solver info and context
231  * based on existing info and context.
232  *
233  * parameters:
234  * context <-- pointer to reference info and context
235  * (actual type: cs_sles_petsc_t *)
236  *
237  * returns:
238  * pointer to newly created solver info object
239  * (actual type: cs_sles_petsc_t *)
240  *----------------------------------------------------------------------------*/
241 
242 void *
243 cs_sles_petsc_copy(const void *context);
244 
245 /*----------------------------------------------------------------------------
246  * Destroy PETSc linear system solver info and context.
247  *
248  * parameters:
249  * context <-> pointer to PETSc linear solver info
250  * (actual type: cs_sles_petsc_t **)
251  *----------------------------------------------------------------------------*/
252 
253 void
254 cs_sles_petsc_destroy(void **context);
255 
256 /*----------------------------------------------------------------------------
257  * Setup PETSc linear equation solver.
258  *
259  * parameters:
260  * context <-> pointer to PETSc linear solver info
261  * (actual type: cs_sles_petsc_t *)
262  * name <-- pointer to system name
263  * a <-- associated matrix
264  * verbosity <-- verbosity level
265  *----------------------------------------------------------------------------*/
266 
267 void
268 cs_sles_petsc_setup(void *context,
269  const char *name,
270  const cs_matrix_t *a,
271  int verbosity);
272 
273 /*----------------------------------------------------------------------------
274  * Call PETSc linear equation solver.
275  *
276  * parameters:
277  * context <-> pointer to PETSc linear solver info
278  * (actual type: cs_sles_petsc_t *)
279  * name <-- pointer to system name
280  * a <-- matrix
281  * verbosity <-- verbosity level
282  * rotation_mode <-- halo update option for rotational periodicity
283  * precision <-- solver precision
284  * r_norm <-- residue normalization
285  * n_iter --> number of iterations
286  * residue --> residue
287  * rhs <-- right hand side
288  * vx <-> system solution
289  * aux_size <-- number of elements in aux_vectors (in bytes)
290  * aux_vectors --- optional working area (internal allocation if NULL)
291  *
292  * returns:
293  * convergence state
294  *----------------------------------------------------------------------------*/
295 
297 cs_sles_petsc_solve(void *context,
298  const char *name,
299  const cs_matrix_t *a,
300  int verbosity,
301  cs_halo_rotation_t rotation_mode,
302  double precision,
303  double r_norm,
304  int *n_iter,
305  double *residue,
306  const cs_real_t *rhs,
307  cs_real_t *vx,
308  size_t aux_size,
309  void *aux_vectors);
310 
311 /*----------------------------------------------------------------------------
312  * Free PETSc linear equation solver setup context.
313  *
314  * This function frees resolution-related data, such as
315  * buffers and preconditioning but does not free the whole context,
316  * as info used for logging (especially performance data) is maintained.
317 
318  * parameters:
319  * context <-> pointer to PETSc linear solver info
320  * (actual type: cs_sles_petsc_t *)
321  *----------------------------------------------------------------------------*/
322 
323 void
324 cs_sles_petsc_free(void *context);
325 
326 /*----------------------------------------------------------------------------*/
343 /*----------------------------------------------------------------------------*/
344 
345 bool
348  const cs_matrix_t *a,
349  cs_halo_rotation_t rotation_mode,
350  const cs_real_t *rhs,
351  cs_real_t *vx);
352 
353 /*----------------------------------------------------------------------------
354  * Log sparse linear equation solver info.
355  *
356  * parameters:
357  * context <-> pointer to PETSc linear solver info
358  * (actual type: cs_sles_petsc_t *)
359  * log_type <-- log type
360  *----------------------------------------------------------------------------*/
361 
362 void
363 cs_sles_petsc_log(const void *context,
364  cs_log_t log_type);
365 
366 /*----------------------------------------------------------------------------*/
367 
369 
370 #endif /* __CS_SLES_PETSC_H__ */
f_id
void const int * f_id
Definition: cs_gui.h:292
cs_sles_petsc_destroy
void cs_sles_petsc_destroy(void **context)
Destroy PETSc linear system solver info and context.
Definition: cs_sles_petsc.c:683
cs_fuel_incl::a
double precision, save a
Definition: cs_fuel_incl.f90:146
cs_sles_petsc_free
void cs_sles_petsc_free(void *context)
Free PETSc linear equation solver setup context.
Definition: cs_sles_petsc.c:1477
cs_sles_petsc_init
static void cs_sles_petsc_init(void)
Definition: cs_sles_petsc.h:127
cs_sles_petsc_t
struct _cs_sles_petsc_t cs_sles_petsc_t
Definition: cs_sles_petsc.h:84
cs_sles_petsc_log
void cs_sles_petsc_log(const void *context, cs_log_t log_type)
Log sparse linear equation solver info.
Definition: cs_sles_petsc.c:1566
cs_sles_petsc_solve
cs_sles_convergence_state_t cs_sles_petsc_solve(void *context, const char *name, const cs_matrix_t *a, int verbosity, cs_halo_rotation_t rotation_mode, double precision, double r_norm, int *n_iter, double *residue, const cs_real_t *rhs, cs_real_t *vx, size_t aux_size, void *aux_vectors)
Call PETSc linear equation solver.
Definition: cs_sles_petsc.c:1301
END_C_DECLS
#define END_C_DECLS
Definition: cs_defs.h:468
cs_sles_petsc_setup
void cs_sles_petsc_setup(void *context, const char *name, const cs_matrix_t *a, int verbosity)
Setup PETSc linear equation solver.
Definition: cs_sles_petsc.c:729
cs_time_plot.h
cs_sles_petsc_setup_hook_t
void() cs_sles_petsc_setup_hook_t(void *context, Mat a, KSP ksp)
Function pointer for user settings of a PETSc KSP solver setup.
Definition: cs_sles_petsc.h:78
atimbr::v
double precision, dimension(:,:,:), allocatable v
Definition: atimbr.f90:114
cs_real_t
double cs_real_t
Floating-point value.
Definition: cs_defs.h:302
BEGIN_C_DECLS
#define BEGIN_C_DECLS
Definition: cs_defs.h:467
cs_matrix_t
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:90
cs_sles_petsc_log_setup
static void cs_sles_petsc_log_setup(KSP ksp)
Definition: cs_sles_petsc.h:151
cs_sles_convergence_state_t
cs_sles_convergence_state_t
Convergence status indicator.
Definition: cs_sles.h:56
cs_glob_mpi_comm
MPI_Comm cs_glob_mpi_comm
Definition: cs_defs.c:181
cs_sles_petsc_copy
void * cs_sles_petsc_copy(const void *context)
Create PETSc linear system solver info and context based on existing info and context.
Definition: cs_sles_petsc.c:659
cs_halo_perio.h
cs_user_sles_petsc_hook
void cs_user_sles_petsc_hook(void *context, Mat a, KSP ksp)
Definition: cs_sles_petsc.c:490
cs_log_t
cs_log_t
Definition: cs_log.h:48
cs_sles_t
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:68
cs_matrix.h
cs_halo_rotation_t
cs_halo_rotation_t
Definition: cs_halo.h:60
cs_sles_petsc_error_post_and_abort
bool cs_sles_petsc_error_post_and_abort(cs_sles_t *sles, cs_sles_convergence_state_t state, const cs_matrix_t *a, cs_halo_rotation_t rotation_mode, const cs_real_t *rhs, cs_real_t *vx)
Error handler for PETSc solver.
Definition: cs_sles_petsc.c:1524
cs_sles_petsc_define
cs_sles_petsc_t * cs_sles_petsc_define(int f_id, const char *name, MatType matrix_type, cs_sles_petsc_setup_hook_t *setup_hook, void *context)
Define and associate a PETSc linear system solver for a given field or equation name.
Definition: cs_sles_petsc.c:536
cs_base.h
cs_sles_petsc_create
cs_sles_petsc_t * cs_sles_petsc_create(MatType matrix_type, cs_sles_petsc_setup_hook_t *setup_hook, void *context)
Create PETSc linear system solver info and context.
Definition: cs_sles_petsc.c:582