7.0
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-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_base.h"
42 #include "cs_cdo_connect.h"
43 #include "cs_cdo_quantities.h"
44 #include "cs_cdofb_priv.h"
45 #include "cs_equation_assemble.h"
46 #include "cs_equation_common.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 typedef struct _cs_cdofb_t cs_cdofb_vecteq_t;
70 
71 /*============================================================================
72  * Public function prototypes
73  *============================================================================*/
74 
75 /*----------------------------------------------------------------------------*/
82 /*----------------------------------------------------------------------------*/
83 
84 bool
86 
87 /*----------------------------------------------------------------------------*/
98 /*----------------------------------------------------------------------------*/
99 
100 void
102  const cs_cdo_connect_t *connect,
103  const cs_time_step_t *time_step,
104  const cs_matrix_structure_t *ms);
105 
106 /*----------------------------------------------------------------------------*/
112 /*----------------------------------------------------------------------------*/
113 
114 const cs_matrix_structure_t *
116 
117 /*----------------------------------------------------------------------------*/
124 /*----------------------------------------------------------------------------*/
125 
126 void
128  cs_cell_builder_t **cb);
129 
130 /*----------------------------------------------------------------------------*/
135 /*----------------------------------------------------------------------------*/
136 
137 void
139 
140 /*----------------------------------------------------------------------------*/
152 /*----------------------------------------------------------------------------*/
153 
154 void *
156  int var_id,
157  int bflux_id,
158  cs_equation_builder_t *eqb);
159 
160 /*----------------------------------------------------------------------------*/
168 /*----------------------------------------------------------------------------*/
169 
170 void *
171 cs_cdofb_vecteq_free_context(void *data);
172 
173 /*----------------------------------------------------------------------------*/
186 /*----------------------------------------------------------------------------*/
187 
188 void
190  const int field_id,
191  const cs_mesh_t *mesh,
192  const cs_equation_param_t *eqp,
194  void *context);
195 
196 /*----------------------------------------------------------------------------*/
216 /*----------------------------------------------------------------------------*/
217 
218 void
220  const cs_equation_param_t *eqp,
221  const cs_equation_builder_t *eqb,
222  const cs_real_t dir_values[],
223  const cs_lnum_t forced_ids[],
224  const cs_real_t val_f_n[],
225  const cs_real_t val_c_n[],
226  const cs_real_t val_f_nm1[],
227  const cs_real_t val_c_nm1[],
228  cs_cell_sys_t *csys,
229  cs_cell_builder_t *cb);
230 
231 /*----------------------------------------------------------------------------*/
244 /*----------------------------------------------------------------------------*/
245 
246 void
248  const cs_mesh_t *mesh,
249  const cs_equation_param_t *eqp,
251  cs_real_t *p_dir_values[],
252  cs_lnum_t *p_enforced_ids[]);
253 
254 /*----------------------------------------------------------------------------*/
267 /*----------------------------------------------------------------------------*/
268 
269 void
271  const cs_equation_builder_t *eqb,
272  const cs_cdofb_vecteq_t *eqc,
273  const cs_cell_mesh_t *cm,
274  cs_hodge_t *diff_hodge,
275  cs_cell_sys_t *csys,
276  cs_cell_builder_t *cb);
277 
278 /*----------------------------------------------------------------------------*/
294 /*----------------------------------------------------------------------------*/
295 
296 void
298  const cs_equation_builder_t *eqb,
299  const cs_cdofb_vecteq_t *eqc,
300  const cs_cell_mesh_t *cm,
302  cs_hodge_t *diff_hodge,
303  cs_cell_sys_t *csys,
304  cs_cell_builder_t *cb);
305 
306 /*----------------------------------------------------------------------------*/
322 /*----------------------------------------------------------------------------*/
323 
324 static inline void
325 cs_cdofb_vecteq_sourceterm(const cs_cell_mesh_t *cm,
326  const cs_equation_param_t *eqp,
327  const cs_real_t t_eval,
328  const cs_real_t coef,
330  cs_cell_builder_t *cb,
332  cs_cell_sys_t *csys)
333 {
334  /* Reset the local contribution */
335  memset(csys->source, 0, csys->n_dofs*sizeof(cs_real_t));
336 
338  (cs_xdef_t *const *)eqp->source_terms,
339  cm,
340  eqb->source_mask,
341  eqb->compute_source,
342  t_eval,
343  mass_hodge, /* Mass matrix context */
344  cb,
345  csys->source);
346 
347  /* Only cell-DoFs are involved */
348  const short int _off = 3*cm->n_fc;
349  csys->rhs[_off ] += coef * csys->source[_off ];
350  csys->rhs[_off + 1] += coef * csys->source[_off + 1];
351  csys->rhs[_off + 2] += coef * csys->source[_off + 2];
352 }
353 
354 /*----------------------------------------------------------------------------*/
368 /*----------------------------------------------------------------------------*/
369 
370 void
372  const cs_range_set_t *rs,
373  const cs_cell_mesh_t *cm,
374  const bool has_sourceterm,
375  cs_cdofb_vecteq_t *eqc,
378  cs_real_t rhs[]);
379 
380 /*----------------------------------------------------------------------------*/
390 /*----------------------------------------------------------------------------*/
391 
392 void
394  cs_field_t *fld,
395  cs_cdofb_vecteq_t *eqc,
396  bool cur2prev);
397 
398 /*----------------------------------------------------------------------------*/
411 /*----------------------------------------------------------------------------*/
412 
413 void
415  const cs_mesh_t *mesh,
416  const int field_id,
417  const cs_equation_param_t *eqp,
419  void *context);
420 
421 /*----------------------------------------------------------------------------*/
434 /*----------------------------------------------------------------------------*/
435 
436 void
437 cs_cdofb_vecteq_solve_implicit(bool cur2prev,
438  const cs_mesh_t *mesh,
439  const int field_id,
440  const cs_equation_param_t *eqp,
442  void *context);
443 
444 /*----------------------------------------------------------------------------*/
457 /*----------------------------------------------------------------------------*/
458 
459 void
460 cs_cdofb_vecteq_solve_theta(bool cur2prev,
461  const cs_mesh_t *mesh,
462  const int field_id,
463  const cs_equation_param_t *eqp,
465  void *context);
466 
467 /*----------------------------------------------------------------------------*/
476 /*----------------------------------------------------------------------------*/
477 
478 void
481  void *context);
482 
483 /*----------------------------------------------------------------------------*/
491 /*----------------------------------------------------------------------------*/
492 
493 void
496  void *context);
497 
498 /*----------------------------------------------------------------------------*/
511 /*----------------------------------------------------------------------------*/
512 
513 cs_real_t *
514 cs_cdofb_vecteq_get_cell_values(void *context,
515  bool previous);
516 
517 /*----------------------------------------------------------------------------*/
528 /*----------------------------------------------------------------------------*/
529 
530 cs_real_t *
531 cs_cdofb_vecteq_get_face_values(void *context,
532  bool previous);
533 
534 /*----------------------------------------------------------------------------*/
543 /*----------------------------------------------------------------------------*/
544 
545 void
547  const char *eqname,
548  void *scheme_context);
549 
550 /*----------------------------------------------------------------------------*/
559 /*----------------------------------------------------------------------------*/
560 
561 void
563  const char *eqname,
564  void *scheme_context);
565 
566 /*----------------------------------------------------------------------------*/
567 
569 
570 #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:1615
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:962
Store common elements used when building an algebraic system related to an equation.
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 dir_values[], const cs_lnum_t forced_ids[], 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:332
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:2098
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:202
struct _cs_matrix_assembler_values_t cs_matrix_assembler_values_t
Definition: cs_matrix_assembler.h:65
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:2237
double * source
Definition: cs_cdo_local.h:117
#define BEGIN_C_DECLS
Definition: cs_defs.h:495
struct _cs_equation_assemble_t cs_equation_assemble_t
Definition: cs_equation_assemble.h:52
void * cs_cdofb_vecteq_free_context(void *data)
Destroy a cs_cdofb_vecteq_t structure.
Definition: cs_cdofb_vecteq.c:1853
Set of local quantities and connectivities related to a mesh cell This is a key structure for all cel...
Definition: cs_cdo_local.h:159
const cs_matrix_structure_t * cs_cdofb_vecteq_matrix_structure(void)
Get the pointer to the related cs_matrix_structure_t.
Definition: cs_cdofb_vecteq.c:1542
Definition: cs_cdo_connect.h:76
Structure associated to a discrete Hodge operator.
Definition: cs_hodge.h:189
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:2128
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:1460
cs_hodge_t ** mass_hodge
Definition: cs_cdofb_priv.h:100
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:2013
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:526
double cs_real_t
Floating-point value.
Definition: cs_defs.h:307
Definition: cs_cdo_quantities.h:124
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:1193
void cs_cdofb_vecteq_finalize_common(void)
Free work buffer and general structure related to CDO face-based schemes.
Definition: cs_cdofb_vecteq.c:1579
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:107
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:2156
int n_dofs
Definition: cs_cdo_local.h:111
Definition: cs_mesh.h:84
double * rhs
Definition: cs_cdo_local.h:116
void cs_cdofb_vecteq_assembly(const cs_cell_sys_t *csys, const cs_range_set_t *rs, const cs_cell_mesh_t *cm, const bool has_sourceterm, cs_cdofb_vecteq_t *eqc, cs_equation_assemble_t *eqa, cs_matrix_assembler_values_t *mav, cs_real_t rhs[])
Perform the assembly stage for a vector-valued system obtained with CDO-Fb scheme.
Definition: cs_cdofb_vecteq.c:677
cs_xdef_t ** source_terms
Definition: cs_equation_param.h:705
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:461
short int n_fc
Definition: cs_cdo_local.h:189
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, cs_real_t *p_dir_values[], cs_lnum_t *p_enforced_ids[])
Set the boundary conditions known from the settings Define an indirection array for the enforcement o...
Definition: cs_cdofb_vecteq.c:274
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:1891
struct _cs_restart_t cs_restart_t
Definition: cs_restart.h:93
Definition: cs_range_set.h:57
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:722
void cs_cdofb_vecteq_init_common(const cs_cdo_quantities_t *quant, const cs_cdo_connect_t *connect, const cs_time_step_t *time_step, const cs_matrix_structure_t *ms)
Allocate work buffer and general structures related to CDO vector-valued face-based schemes...
Definition: cs_cdofb_vecteq.c:1482
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:759
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:2044
Set of local and temporary buffers useful for building the algebraic system with a cellwise process...
Definition: cs_cdo_local.h:69
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:704
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:301
#define END_C_DECLS
Definition: cs_defs.h:496
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:1557
Definition: cs_cdofb_priv.h:52
Structure and routines handling the specific settings related to a cs_equation_t structure.
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:846
struct _cs_matrix_structure_t cs_matrix_structure_t
Definition: cs_matrix.h:90
Definition: mesh.f90:26
Definition: cs_timer.h:57