7.1
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 
70 typedef struct _cs_cdofb_t cs_cdofb_vecteq_t;
71 
72 /*============================================================================
73  * Public function prototypes
74  *============================================================================*/
75 
76 /*----------------------------------------------------------------------------*/
83 /*----------------------------------------------------------------------------*/
84 
85 bool
87 
88 /*----------------------------------------------------------------------------*/
99 /*----------------------------------------------------------------------------*/
100 
101 void
103  const cs_cdo_connect_t *connect,
104  const cs_time_step_t *time_step,
105  const cs_matrix_structure_t *ms);
106 
107 /*----------------------------------------------------------------------------*/
113 /*----------------------------------------------------------------------------*/
114 
115 const cs_matrix_structure_t *
117 
118 /*----------------------------------------------------------------------------*/
125 /*----------------------------------------------------------------------------*/
126 
127 void
129  cs_cell_builder_t **cb);
130 
131 /*----------------------------------------------------------------------------*/
136 /*----------------------------------------------------------------------------*/
137 
138 void
140 
141 /*----------------------------------------------------------------------------*/
153 /*----------------------------------------------------------------------------*/
154 
155 void *
157  int var_id,
158  int bflux_id,
159  cs_equation_builder_t *eqb);
160 
161 /*----------------------------------------------------------------------------*/
169 /*----------------------------------------------------------------------------*/
170 
171 void *
172 cs_cdofb_vecteq_free_context(void *data);
173 
174 /*----------------------------------------------------------------------------*/
187 /*----------------------------------------------------------------------------*/
188 
189 void
191  const int field_id,
192  const cs_mesh_t *mesh,
193  const cs_equation_param_t *eqp,
195  void *context);
196 
197 /*----------------------------------------------------------------------------*/
215 /*----------------------------------------------------------------------------*/
216 
217 void
219  const cs_equation_param_t *eqp,
220  const cs_equation_builder_t *eqb,
221  const cs_real_t val_f_n[],
222  const cs_real_t val_c_n[],
223  const cs_real_t val_f_nm1[],
224  const cs_real_t val_c_nm1[],
225  cs_cell_sys_t *csys,
226  cs_cell_builder_t *cb);
227 
228 /*----------------------------------------------------------------------------*/
239 /*----------------------------------------------------------------------------*/
240 
241 void
243  const cs_mesh_t *mesh,
244  const cs_equation_param_t *eqp,
245  cs_equation_builder_t *eqb);
246 
247 /*----------------------------------------------------------------------------*/
260 /*----------------------------------------------------------------------------*/
261 
262 void
264  const cs_equation_builder_t *eqb,
265  const cs_cdofb_vecteq_t *eqc,
266  const cs_cell_mesh_t *cm,
267  cs_hodge_t *diff_hodge,
268  cs_cell_sys_t *csys,
269  cs_cell_builder_t *cb);
270 
271 /*----------------------------------------------------------------------------*/
287 /*----------------------------------------------------------------------------*/
288 
289 void
291  const cs_equation_builder_t *eqb,
292  const cs_cdofb_vecteq_t *eqc,
293  const cs_cell_mesh_t *cm,
295  cs_hodge_t *diff_hodge,
296  cs_cell_sys_t *csys,
297  cs_cell_builder_t *cb);
298 
299 /*----------------------------------------------------------------------------*/
315 /*----------------------------------------------------------------------------*/
316 
317 static inline void
318 cs_cdofb_vecteq_sourceterm(const cs_cell_mesh_t *cm,
319  const cs_equation_param_t *eqp,
320  const cs_real_t t_eval,
321  const cs_real_t coef,
323  cs_cell_builder_t *cb,
325  cs_cell_sys_t *csys)
326 {
327  /* Reset the local contribution */
328 
329  memset(csys->source, 0, csys->n_dofs*sizeof(cs_real_t));
330 
332  (cs_xdef_t *const *)eqp->source_terms,
333  cm,
334  eqb->source_mask,
335  eqb->compute_source,
336  t_eval,
337  mass_hodge, /* Mass matrix context */
338  cb,
339  csys->source);
340 
341  /* Only cell-DoFs are involved */
342 
343  const short int _off = 3*cm->n_fc;
344  csys->rhs[_off ] += coef * csys->source[_off ];
345  csys->rhs[_off + 1] += coef * csys->source[_off + 1];
346  csys->rhs[_off + 2] += coef * csys->source[_off + 2];
347 }
348 
349 /*----------------------------------------------------------------------------*/
363 /*----------------------------------------------------------------------------*/
364 
365 void
367  const cs_range_set_t *rs,
368  const cs_cell_mesh_t *cm,
369  const bool has_sourceterm,
370  cs_cdofb_vecteq_t *eqc,
373  cs_real_t rhs[]);
374 
375 /*----------------------------------------------------------------------------*/
385 /*----------------------------------------------------------------------------*/
386 
387 void
389  cs_field_t *fld,
390  cs_cdofb_vecteq_t *eqc,
391  bool cur2prev);
392 
393 /*----------------------------------------------------------------------------*/
406 /*----------------------------------------------------------------------------*/
407 
408 void
410  const cs_mesh_t *mesh,
411  const int field_id,
412  const cs_equation_param_t *eqp,
414  void *context);
415 
416 /*----------------------------------------------------------------------------*/
429 /*----------------------------------------------------------------------------*/
430 
431 void
432 cs_cdofb_vecteq_solve_implicit(bool cur2prev,
433  const cs_mesh_t *mesh,
434  const int field_id,
435  const cs_equation_param_t *eqp,
437  void *context);
438 
439 /*----------------------------------------------------------------------------*/
452 /*----------------------------------------------------------------------------*/
453 
454 void
455 cs_cdofb_vecteq_solve_theta(bool cur2prev,
456  const cs_mesh_t *mesh,
457  const int field_id,
458  const cs_equation_param_t *eqp,
460  void *context);
461 
462 /*----------------------------------------------------------------------------*/
471 /*----------------------------------------------------------------------------*/
472 
473 void
476  void *context);
477 
478 /*----------------------------------------------------------------------------*/
486 /*----------------------------------------------------------------------------*/
487 
488 void
491  void *context);
492 
493 /*----------------------------------------------------------------------------*/
506 /*----------------------------------------------------------------------------*/
507 
508 cs_real_t *
509 cs_cdofb_vecteq_get_cell_values(void *context,
510  bool previous);
511 
512 /*----------------------------------------------------------------------------*/
523 /*----------------------------------------------------------------------------*/
524 
525 cs_real_t *
526 cs_cdofb_vecteq_get_face_values(void *context,
527  bool previous);
528 
529 /*----------------------------------------------------------------------------*/
538 /*----------------------------------------------------------------------------*/
539 
540 void
542  const char *eqname,
543  void *scheme_context);
544 
545 /*----------------------------------------------------------------------------*/
554 /*----------------------------------------------------------------------------*/
555 
556 void
558  const char *eqname,
559  void *scheme_context);
560 
561 /*----------------------------------------------------------------------------*/
562 
564 
565 #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:1661
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:971
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:2169
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:177
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:2315
double * source
Definition: cs_cdo_local.h:156
#define BEGIN_C_DECLS
Definition: cs_defs.h:510
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:1913
Set of local quantities and connectivities related to a mesh cell.
Definition: cs_cdo_local.h:202
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:1588
Definition: cs_cdo_connect.h:79
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:2199
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:1504
cs_hodge_t ** mass_hodge
Definition: cs_cdofb_priv.h:110
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:2076
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:507
double cs_real_t
Floating-point value.
Definition: cs_defs.h:322
Definition: cs_cdo_quantities.h:129
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:1219
void cs_cdofb_vecteq_finalize_common(void)
Free work buffer and general structure related to CDO face-based schemes.
Definition: cs_cdofb_vecteq.c:1625
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:146
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:2227
int n_dofs
Definition: cs_cdo_local.h:150
Definition: cs_mesh.h:84
double * rhs
Definition: cs_cdo_local.h:155
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:669
cs_xdef_t ** source_terms
Definition: cs_equation_param.h:689
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:337
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:439
short int n_fc
Definition: cs_cdo_local.h:232
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:1952
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:717
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:1526
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:756
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:2108
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:688
#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:1603
Definition: cs_cdofb_priv.h:53
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:847
struct _cs_matrix_structure_t cs_matrix_structure_t
Definition: cs_matrix.h:89
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:287