7.2
general documentation
cs_cdofb_vecteq.h
Go to the documentation of this file.
1 #ifndef __CS_CDOFB_VECTEQ_H__
2 #define __CS_CDOFB_VECTEQ_H__
3 
4 /*============================================================================
5  * Build an algebraic CDO face-based system for unsteady convection/diffusion
6  * reaction of vector-valued equations with source terms
7  *============================================================================*/
8 
9 /*
10  This file is part of code_saturne, a general-purpose CFD tool.
11 
12  Copyright (C) 1998-2022 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_base.h"
42 #include "cs_cdo_assembly.h"
43 #include "cs_cdo_connect.h"
44 #include "cs_cdo_quantities.h"
45 #include "cs_cdofb_priv.h"
46 #include "cs_equation_builder.h"
47 #include "cs_equation_param.h"
48 #include "cs_field.h"
49 #include "cs_matrix.h"
50 #include "cs_mesh.h"
51 #include "cs_restart.h"
52 #include "cs_sles.h"
53 #include "cs_source_term.h"
54 #include "cs_time_step.h"
55 
56 /*----------------------------------------------------------------------------*/
57 
59 
60 /*============================================================================
61  * Macro definitions
62  *============================================================================*/
63 
64 /*============================================================================
65  * Type definitions
66  *============================================================================*/
67 
68 /* Algebraic system for CDO face-based discretization */
69 
70 typedef struct _cs_cdofb_t cs_cdofb_vecteq_t;
71 
72 /*============================================================================
73  * Static inline public function prototypes
74  *============================================================================*/
75 
76 /*----------------------------------------------------------------------------*/
92 /*----------------------------------------------------------------------------*/
93 
94 static inline void
96  const cs_equation_param_t *eqp,
97  const cs_real_t t_eval,
98  const cs_real_t coef,
100  cs_cell_builder_t *cb,
102  cs_cell_sys_t *csys)
103 {
104  /* Reset the local contribution */
105 
106  memset(csys->source, 0, csys->n_dofs*sizeof(cs_real_t));
107 
109  (cs_xdef_t *const *)eqp->source_terms,
110  cm,
111  eqb->source_mask,
112  eqb->compute_source,
113  t_eval,
114  mass_hodge, /* Mass matrix context */
115  cb,
116  csys->source);
117 
118  /* Only cell-DoFs are involved */
119 
120  const short int _off = 3*cm->n_fc;
121  csys->rhs[_off ] += coef * csys->source[_off ];
122  csys->rhs[_off + 1] += coef * csys->source[_off + 1];
123  csys->rhs[_off + 2] += coef * csys->source[_off + 2];
124 }
125 
126 /*============================================================================
127  * Public function prototypes
128  *============================================================================*/
129 
130 /*----------------------------------------------------------------------------*/
137 /*----------------------------------------------------------------------------*/
138 
139 bool
141 
142 /*----------------------------------------------------------------------------*/
152 /*----------------------------------------------------------------------------*/
153 
154 void
156  const cs_cdo_connect_t *connect,
157  const cs_time_step_t *time_step);
158 
159 /*----------------------------------------------------------------------------*/
166 /*----------------------------------------------------------------------------*/
167 
168 void
170  cs_cell_builder_t **cb);
171 
172 /*----------------------------------------------------------------------------*/
177 /*----------------------------------------------------------------------------*/
178 
179 void
181 
182 /*----------------------------------------------------------------------------*/
194 /*----------------------------------------------------------------------------*/
195 
196 void *
198  int var_id,
199  int bflux_id,
200  cs_equation_builder_t *eqb);
201 
202 /*----------------------------------------------------------------------------*/
210 /*----------------------------------------------------------------------------*/
211 
212 void *
213 cs_cdofb_vecteq_free_context(void *data);
214 
215 /*----------------------------------------------------------------------------*/
228 /*----------------------------------------------------------------------------*/
229 
230 void
232  const int field_id,
233  const cs_mesh_t *mesh,
234  const cs_equation_param_t *eqp,
236  void *context);
237 
238 /*----------------------------------------------------------------------------*/
256 /*----------------------------------------------------------------------------*/
257 
258 void
260  const cs_equation_param_t *eqp,
261  const cs_equation_builder_t *eqb,
262  const cs_real_t val_f_n[],
263  const cs_real_t val_c_n[],
264  const cs_real_t val_f_nm1[],
265  const cs_real_t val_c_nm1[],
266  cs_cell_sys_t *csys,
267  cs_cell_builder_t *cb);
268 
269 /*----------------------------------------------------------------------------*/
280 /*----------------------------------------------------------------------------*/
281 
282 void
284  const cs_mesh_t *mesh,
285  const cs_equation_param_t *eqp,
286  cs_equation_builder_t *eqb);
287 
288 /*----------------------------------------------------------------------------*/
301 /*----------------------------------------------------------------------------*/
302 
303 void
305  const cs_equation_builder_t *eqb,
306  const cs_cdofb_vecteq_t *eqc,
307  const cs_cell_mesh_t *cm,
308  cs_hodge_t *diff_hodge,
309  cs_cell_sys_t *csys,
310  cs_cell_builder_t *cb);
311 
312 /*----------------------------------------------------------------------------*/
328 /*----------------------------------------------------------------------------*/
329 
330 void
332  const cs_equation_builder_t *eqb,
333  const cs_cdofb_vecteq_t *eqc,
334  const cs_cell_mesh_t *cm,
336  cs_hodge_t *diff_hodge,
337  cs_cell_sys_t *csys,
338  cs_cell_builder_t *cb);
339 
340 /*----------------------------------------------------------------------------*/
351 /*----------------------------------------------------------------------------*/
352 
353 void
355  cs_cdo_system_block_t *block,
356  cs_real_t *rhs,
357  cs_cdofb_vecteq_t *eqc,
358  cs_cdo_assembly_t *asb);
359 
360 /*----------------------------------------------------------------------------*/
370 /*----------------------------------------------------------------------------*/
371 
372 void
374  cs_field_t *fld,
375  cs_cdofb_vecteq_t *eqc,
376  bool cur2prev);
377 
378 /*----------------------------------------------------------------------------*/
391 /*----------------------------------------------------------------------------*/
392 
393 void
395  const cs_mesh_t *mesh,
396  const int field_id,
397  const cs_equation_param_t *eqp,
399  void *context);
400 
401 /*----------------------------------------------------------------------------*/
414 /*----------------------------------------------------------------------------*/
415 
416 void
417 cs_cdofb_vecteq_solve_implicit(bool cur2prev,
418  const cs_mesh_t *mesh,
419  const int field_id,
420  const cs_equation_param_t *eqp,
422  void *context);
423 
424 /*----------------------------------------------------------------------------*/
437 /*----------------------------------------------------------------------------*/
438 
439 void
440 cs_cdofb_vecteq_solve_theta(bool cur2prev,
441  const cs_mesh_t *mesh,
442  const int field_id,
443  const cs_equation_param_t *eqp,
445  void *context);
446 
447 /*----------------------------------------------------------------------------*/
456 /*----------------------------------------------------------------------------*/
457 
458 void
461  void *context);
462 
463 /*----------------------------------------------------------------------------*/
471 /*----------------------------------------------------------------------------*/
472 
473 void
476  void *context);
477 
478 /*----------------------------------------------------------------------------*/
491 /*----------------------------------------------------------------------------*/
492 
493 cs_real_t *
494 cs_cdofb_vecteq_get_cell_values(void *context,
495  bool previous);
496 
497 /*----------------------------------------------------------------------------*/
508 /*----------------------------------------------------------------------------*/
509 
510 cs_real_t *
511 cs_cdofb_vecteq_get_face_values(void *context,
512  bool previous);
513 
514 /*----------------------------------------------------------------------------*/
523 /*----------------------------------------------------------------------------*/
524 
525 void
527  const char *eqname,
528  void *scheme_context);
529 
530 /*----------------------------------------------------------------------------*/
539 /*----------------------------------------------------------------------------*/
540 
541 void
543  const char *eqname,
544  void *scheme_context);
545 
546 /*----------------------------------------------------------------------------*/
547 
549 
550 #endif /* __CS_CDOFB_VECTEQ_H__ */
void * cs_cdofb_vecteq_init_context(const cs_equation_param_t *eqp, int var_id, int bflux_id, cs_equation_builder_t *eqb)
Initialize a cs_cdofb_vecteq_t structure storing data useful for building and managing such a scheme...
Definition: cs_cdofb_vecteq.c:1618
void cs_cdofb_vecteq_solve_implicit(bool cur2prev, const cs_mesh_t *mesh, const int field_id, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
Build and solve the linear system arising from a vector diffusion equation with a CDO-Fb scheme and a...
Definition: cs_cdofb_vecteq.c:959
Store common elements used when building an algebraic system related to an equation.
time step descriptor
Definition: cs_time_step.h:64
cs_real_t * cs_cdofb_vecteq_get_cell_values(void *context, bool previous)
Get the computed values at mesh cells from the inverse operation w.r.t. the static condensation (DoF ...
Definition: cs_cdofb_vecteq.c:2176
Field descriptor.
Definition: cs_field.h:125
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources...
Definition: cs_equation_param.h:186
void cs_cdofb_vecteq_write_restart(cs_restart_t *restart, const char *eqname, void *scheme_context)
Write additional arrays (not defined as fields) but useful for the checkpoint/restart process...
Definition: cs_cdofb_vecteq.c:2322
double * source
Definition: cs_cdo_local.h:157
#define BEGIN_C_DECLS
Definition: cs_defs.h:510
void * cs_cdofb_vecteq_free_context(void *data)
Destroy a cs_cdofb_vecteq_t structure.
Definition: cs_cdofb_vecteq.c:1914
Set of local quantities and connectivities related to a mesh cell.
Definition: cs_cdo_local.h:203
Definition: cs_cdo_connect.h:61
Structure associated to a discrete Hodge operator *.
Definition: cs_hodge.h:186
cs_real_t * cs_cdofb_vecteq_get_face_values(void *context, bool previous)
Retrieve an array of values at mesh faces for the current context. The lifecycle of this array is man...
Definition: cs_cdofb_vecteq.c:2206
bool cs_cdofb_vecteq_is_initialized(void)
Check if the generic structures for building a CDO-Fb scheme are allocated.
Definition: cs_cdofb_vecteq.c:1478
void cs_cdofb_vecteq_assembly(const cs_cell_sys_t *csys, cs_cdo_system_block_t *block, cs_real_t *rhs, cs_cdofb_vecteq_t *eqc, cs_cdo_assembly_t *asb)
Perform the assembly stage for a vector-valued system obtained with CDO-Fb scheme.
Definition: cs_cdofb_vecteq.c:667
cs_hodge_t ** mass_hodge
Definition: cs_cdofb_priv.h:104
void cs_cdofb_vecteq_current_to_previous(const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
Operate a current to previous operation for the field associated to this equation and potentially for...
Definition: cs_cdofb_vecteq.c:2083
void cs_cdofb_vecteq_conv_diff_reac(const cs_equation_param_t *eqp, const cs_equation_builder_t *eqb, const cs_cdofb_vecteq_t *eqc, const cs_cell_mesh_t *cm, cs_hodge_t *mass_hodge, cs_hodge_t *diff_hodge, cs_cell_sys_t *csys, cs_cell_builder_t *cb)
Build the local matrices arising from the convection, diffusion, reaction terms in vector-valued CDO-...
Definition: cs_cdofb_vecteq.c:508
double cs_real_t
Floating-point value.
Definition: cs_defs.h:322
Definition: cs_cdo_quantities.h:132
void cs_cdofb_vecteq_solve_theta(bool cur2prev, const cs_mesh_t *mesh, const int field_id, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
Build and solve the linear system arising from a vector diffusion equation with a CDO-Fb scheme and a...
Definition: cs_cdofb_vecteq.c:1199
Set of arrays and local (small) dense matrices related to a mesh cell This is a key structure for bui...
Definition: cs_cdo_local.h:147
void cs_cdofb_vecteq_init_sharing(const cs_cdo_quantities_t *quant, const cs_cdo_connect_t *connect, const cs_time_step_t *time_step)
Allocate work buffer and general structures related to CDO vector-valued face-based schemes...
Definition: cs_cdofb_vecteq.c:1499
void cs_cdofb_vecteq_read_restart(cs_restart_t *restart, const char *eqname, void *scheme_context)
Read additional arrays (not defined as fields) but useful for the checkpoint/restart process...
Definition: cs_cdofb_vecteq.c:2234
int n_dofs
Definition: cs_cdo_local.h:151
Definition: cs_mesh.h:84
double * rhs
Definition: cs_cdo_local.h:156
cs_xdef_t ** source_terms
Definition: cs_equation_param.h:698
void cs_cdofb_vecteq_init_cell_system(const cs_cell_mesh_t *cm, const cs_equation_param_t *eqp, const cs_equation_builder_t *eqb, const cs_real_t val_f_n[], const cs_real_t val_c_n[], const cs_real_t val_f_nm1[], const cs_real_t val_c_nm1[], cs_cell_sys_t *csys, cs_cell_builder_t *cb)
Initialize the local structure for the current cell The algebraic system for time t^{n+1} is going to...
Definition: cs_cdofb_vecteq.c:338
void cs_cdofb_vecteq_diffusion(const cs_equation_param_t *eqp, const cs_equation_builder_t *eqb, const cs_cdofb_vecteq_t *eqc, const cs_cell_mesh_t *cm, cs_hodge_t *diff_hodge, cs_cell_sys_t *csys, cs_cell_builder_t *cb)
Build the local matrices arising from the diffusion term in the vector-valued CDO-Fb schemes...
Definition: cs_cdofb_vecteq.c:440
short int n_fc
Definition: cs_cdo_local.h:238
void cs_cdofb_vecteq_init_values(cs_real_t t_eval, const int field_id, const cs_mesh_t *mesh, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
Set the initial values of the variable field taking into account the boundary conditions. Case of vector-valued CDO-Fb schemes.
Definition: cs_cdofb_vecteq.c:1953
struct _cs_restart_t cs_restart_t
Definition: cs_restart.h:93
void cs_cdofb_vecteq_update_cell_fields(cs_timer_counter_t *tce, cs_field_t *fld, cs_cdofb_vecteq_t *eqc, bool cur2prev)
Update the variables associated to cells in case of a CDO-Fb scheme. This has to be done after a reso...
Definition: cs_cdofb_vecteq.c:716
Definition: cs_cdo_system.h:326
void cs_cdofb_vecteq_solve_steady_state(bool cur2prev, const cs_mesh_t *mesh, const int field_id, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
Build and solve the linear system arising from a vector steady-state diffusion equation with a CDO-Fb...
Definition: cs_cdofb_vecteq.c:755
void cs_cdofb_vecteq_extra_post(const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
Predefined extra-operations related to this equation.
Definition: cs_cdofb_vecteq.c:2115
Set of local and temporary buffers.
Definition: cs_cdo_local.h:60
Structure storing medata for defining a quantity in a very flexible way.
Definition: cs_xdef.h:154
int n_source_terms
Definition: cs_equation_param.h:697
#define END_C_DECLS
Definition: cs_defs.h:511
void cs_cdofb_vecteq_get(cs_cell_sys_t **csys, cs_cell_builder_t **cb)
Retrieve work buffers used for building a CDO system cellwise.
Definition: cs_cdofb_vecteq.c:1560
Definition: cs_cdofb_priv.h:51
static void cs_cdofb_vecteq_sourceterm(const cs_cell_mesh_t *cm, const cs_equation_param_t *eqp, const cs_real_t t_eval, const cs_real_t coef, cs_hodge_t *mass_hodge, cs_cell_builder_t *cb, cs_equation_builder_t *eqb, cs_cell_sys_t *csys)
Compute the source term for a vector-valued CDO-Fb scheme and add it to the local rhs Mass matrix is ...
Definition: cs_cdofb_vecteq.h:95
void cs_cdofb_vecteq_finalize_sharing(void)
Free work buffer and general structure related to CDO face-based schemes.
Definition: cs_cdofb_vecteq.c:1582
void cs_source_term_compute_cellwise(const int n_source_terms, cs_xdef_t *const *source_terms, const cs_cell_mesh_t *cm, const cs_mask_t *source_mask, cs_source_term_cellwise_t *compute_source[], cs_real_t time_eval, void *input, cs_cell_builder_t *cb, cs_real_t *result)
Compute the local contributions of source terms in a cell.
Definition: cs_source_term.c:926
struct _cs_cdo_assembly_t cs_cdo_assembly_t
Definition: cs_cdo_assembly.h:52
Definition: mesh.f90:26
Definition: cs_timer.h:55
void cs_cdofb_vecteq_setup(cs_real_t t_eval, const cs_mesh_t *mesh, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb)
Set the boundary conditions known from the settings Define an indirection array for the enforcement o...
Definition: cs_cdofb_vecteq.c:288