8.3
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-2024 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
70typedef struct _cs_cdofb_t cs_cdofb_vecteq_t;
71
72/*============================================================================
73 * Static inline public function prototypes
74 *============================================================================*/
75
76/*----------------------------------------------------------------------------*/
91/*----------------------------------------------------------------------------*/
92
93static inline void
95 const cs_equation_param_t *eqp,
96 const cs_real_t t_eval,
97 const cs_real_t coef,
101 cs_cell_sys_t *csys)
102{
103 /* Reset the local contribution */
104
105 memset(csys->source, 0, csys->n_dofs*sizeof(cs_real_t));
106
108 (cs_xdef_t *const *)eqp->source_terms,
109 cm,
110 eqb->source_mask,
111 eqb->compute_source,
112 t_eval,
113 mass_hodge, /* Mass matrix context */
114 cb,
115 csys->source);
116
117 /* Only cell-DoFs are involved */
118
119 const short int _off = 3*cm->n_fc;
120 csys->rhs[_off ] += coef * csys->source[_off ];
121 csys->rhs[_off + 1] += coef * csys->source[_off + 1];
122 csys->rhs[_off + 2] += coef * csys->source[_off + 2];
123}
124
125/*============================================================================
126 * Public function prototypes
127 *============================================================================*/
128
129/*----------------------------------------------------------------------------*/
136/*----------------------------------------------------------------------------*/
137
138bool
140
141/*----------------------------------------------------------------------------*/
151/*----------------------------------------------------------------------------*/
152
153void
155 const cs_cdo_connect_t *connect,
156 const cs_time_step_t *time_step);
157
158/*----------------------------------------------------------------------------*/
165/*----------------------------------------------------------------------------*/
166
167void
169 cs_cell_builder_t **cb);
170
171/*----------------------------------------------------------------------------*/
176/*----------------------------------------------------------------------------*/
177
178void
180
181/*----------------------------------------------------------------------------*/
193/*----------------------------------------------------------------------------*/
194
195void *
197 int var_id,
198 int bflux_id,
200
201/*----------------------------------------------------------------------------*/
209/*----------------------------------------------------------------------------*/
210
211void *
213
214/*----------------------------------------------------------------------------*/
227/*----------------------------------------------------------------------------*/
228
229void
231 const int field_id,
232 const cs_mesh_t *mesh,
233 const cs_equation_param_t *eqp,
235 void *context);
236
237/*----------------------------------------------------------------------------*/
255/*----------------------------------------------------------------------------*/
256
257void
259 const cs_equation_param_t *eqp,
260 const cs_equation_builder_t *eqb,
261 const cs_real_t val_f_n[],
262 const cs_real_t val_c_n[],
263 const cs_real_t val_f_nm1[],
264 const cs_real_t val_c_nm1[],
265 cs_cell_sys_t *csys,
267
268/*----------------------------------------------------------------------------*/
279/*----------------------------------------------------------------------------*/
280
281void
283 const cs_mesh_t *mesh,
284 const cs_equation_param_t *eqp,
286
287/*----------------------------------------------------------------------------*/
300/*----------------------------------------------------------------------------*/
301
302void
304 const cs_equation_builder_t *eqb,
305 const cs_cdofb_vecteq_t *eqc,
306 const cs_cell_mesh_t *cm,
307 cs_hodge_t *diff_hodge,
308 cs_cell_sys_t *csys,
310
311/*----------------------------------------------------------------------------*/
327/*----------------------------------------------------------------------------*/
328
329void
331 const cs_equation_builder_t *eqb,
332 const cs_cdofb_vecteq_t *eqc,
333 const cs_cell_mesh_t *cm,
335 cs_hodge_t *diff_hodge,
336 cs_cell_sys_t *csys,
338
339/*----------------------------------------------------------------------------*/
350/*----------------------------------------------------------------------------*/
351
352void
355 cs_real_t *rhs,
356 cs_cdofb_vecteq_t *eqc,
357 cs_cdo_assembly_t *asb);
358
359/*----------------------------------------------------------------------------*/
369/*----------------------------------------------------------------------------*/
370
371void
373 cs_field_t *fld,
374 cs_cdofb_vecteq_t *eqc,
375 bool cur2prev);
376
377/*----------------------------------------------------------------------------*/
390/*----------------------------------------------------------------------------*/
391
392void
394 const cs_mesh_t *mesh,
395 const int field_id,
396 const cs_equation_param_t *eqp,
398 void *context);
399
400/*----------------------------------------------------------------------------*/
413/*----------------------------------------------------------------------------*/
414
415void
417 const cs_mesh_t *mesh,
418 const int field_id,
419 const cs_equation_param_t *eqp,
421 void *context);
422
423/*----------------------------------------------------------------------------*/
436/*----------------------------------------------------------------------------*/
437
438void
439cs_cdofb_vecteq_solve_theta(bool cur2prev,
440 const cs_mesh_t *mesh,
441 const int field_id,
442 const cs_equation_param_t *eqp,
444 void *context);
445
446/*----------------------------------------------------------------------------*/
455/*----------------------------------------------------------------------------*/
456
457void
460 void *context);
461
462/*----------------------------------------------------------------------------*/
470/*----------------------------------------------------------------------------*/
471
472void
475 void *context);
476
477/*----------------------------------------------------------------------------*/
490/*----------------------------------------------------------------------------*/
491
492cs_real_t *
494 bool previous);
495
496/*----------------------------------------------------------------------------*/
507/*----------------------------------------------------------------------------*/
508
509cs_real_t *
511 bool previous);
512
513/*----------------------------------------------------------------------------*/
522/*----------------------------------------------------------------------------*/
523
524void
526 const char *eqname,
527 void *scheme_context);
528
529/*----------------------------------------------------------------------------*/
538/*----------------------------------------------------------------------------*/
539
540void
542 const char *eqname,
543 void *scheme_context);
544
545/*----------------------------------------------------------------------------*/
546
548
549#endif /* __CS_CDOFB_VECTEQ_H__ */
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.cpp:1494
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.cpp:2166
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.cpp:654
void * cs_cdofb_vecteq_free_context(void *data)
Destroy a cs_cdofb_vecteq_t structure.
Definition: cs_cdofb_vecteq.cpp:1904
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.cpp:703
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:94
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.cpp:337
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.cpp:500
void cs_cdofb_vecteq_finalize_sharing(void)
Free work buffer and general structure related to CDO face-based schemes.
Definition: cs_cdofb_vecteq.cpp:1572
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.cpp:2224
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.cpp:2196
bool cs_cdofb_vecteq_is_initialized(void)
Check if the generic structures for building a CDO-Fb scheme are allocated.
Definition: cs_cdofb_vecteq.cpp:1473
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.cpp:744
void * cs_cdofb_vecteq_init_context(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.cpp:1608
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.cpp:2313
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.cpp:950
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.cpp:2105
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....
Definition: cs_cdofb_vecteq.cpp:1943
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.cpp:433
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.cpp:289
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.cpp:2074
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.cpp:1555
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.cpp:1192
#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
struct _cs_restart_t cs_restart_t
Definition: cs_restart.h:95
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.cpp:1002
Definition: mesh.f90:26
Definition: cs_cdofb_priv.h:51
cs_hodge_t ** mass_hodge
Definition: cs_cdofb_priv.h:104
Definition: cs_cdo_connect.h:61
Definition: cs_cdo_quantities.h:139
Definition: cs_cdo_system.h:335
Set of local and temporary buffers.
Definition: cs_cdo_local.h:60
Set of local quantities and connectivities related to a mesh cell.
Definition: cs_cdo_local.h:203
short int n_fc
Definition: cs_cdo_local.h:238
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
double * source
Definition: cs_cdo_local.h:157
double * rhs
Definition: cs_cdo_local.h:156
int n_dofs
Definition: cs_cdo_local.h:151
Store common elements used when building an algebraic system related to an equation.
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition: cs_equation_param.h:192
cs_xdef_t ** source_terms
Definition: cs_equation_param.h:730
int n_source_terms
Definition: cs_equation_param.h:729
Field descriptor.
Definition: cs_field.h:131
Structure associated to a discrete Hodge operator *.
Definition: cs_hodge.h:183
Definition: cs_mesh.h:85
time step descriptor
Definition: cs_time_step.h:64
Definition: cs_timer.h:55
Structure storing medata for defining a quantity in a very flexible way.
Definition: cs_xdef.h:160