My Project
programmer's 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-2019 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_equation_common.h"
45 #include "cs_equation_param.h"
46 #include "cs_field.h"
47 #include "cs_matrix.h"
48 #include "cs_mesh.h"
49 #include "cs_restart.h"
50 #include "cs_sles.h"
51 #include "cs_source_term.h"
52 #include "cs_time_step.h"
53 
54 /*----------------------------------------------------------------------------*/
55 
57 
58 /*============================================================================
59  * Macro definitions
60  *============================================================================*/
61 
62 /*============================================================================
63  * Type definitions
64  *============================================================================*/
65 
66 /* Algebraic system for CDO face-based discretization */
67 typedef struct _cs_cdofb_t cs_cdofb_vecteq_t;
68 
69 /*============================================================================
70  * Public function prototypes
71  *============================================================================*/
72 
73 /*----------------------------------------------------------------------------*/
80 /*----------------------------------------------------------------------------*/
81 
82 bool
84 
85 /*----------------------------------------------------------------------------*/
96 /*----------------------------------------------------------------------------*/
97 
98 void
100  const cs_cdo_connect_t *connect,
101  const cs_time_step_t *time_step,
102  const cs_matrix_structure_t *ms);
103 
104 /*----------------------------------------------------------------------------*/
110 /*----------------------------------------------------------------------------*/
111 
112 const cs_matrix_structure_t *
114 
115 /*----------------------------------------------------------------------------*/
122 /*----------------------------------------------------------------------------*/
123 
124 void
126  cs_cell_builder_t **cb);
127 
128 /*----------------------------------------------------------------------------*/
133 /*----------------------------------------------------------------------------*/
134 
135 void
137 
138 /*----------------------------------------------------------------------------*/
150 /*----------------------------------------------------------------------------*/
151 
152 void *
154  int var_id,
155  int bflux_id,
156  cs_equation_builder_t *eqb);
157 
158 /*----------------------------------------------------------------------------*/
166 /*----------------------------------------------------------------------------*/
167 
168 void *
169 cs_cdofb_vecteq_free_context(void *data);
170 
171 /*----------------------------------------------------------------------------*/
184 /*----------------------------------------------------------------------------*/
185 
186 void
188  const int field_id,
189  const cs_mesh_t *mesh,
190  const cs_equation_param_t *eqp,
192  void *context);
193 
194 /*----------------------------------------------------------------------------*/
209 /*----------------------------------------------------------------------------*/
210 
211 void
213  const cs_cell_mesh_t *cm,
214  const cs_equation_param_t *eqp,
215  const cs_equation_builder_t *eqb,
216  const cs_cdofb_vecteq_t *eqc,
217  const cs_real_t dir_values[],
218  const cs_real_t field_tn[],
219  cs_real_t t_eval,
220  cs_cell_sys_t *csys,
221  cs_cell_builder_t *cb);
222 
223 /*----------------------------------------------------------------------------*/
233 /*----------------------------------------------------------------------------*/
234 
235 void
237  const cs_mesh_t *mesh,
238  const cs_equation_param_t *eqp,
240  cs_real_t *p_dir_values[]);
241 
242 /*----------------------------------------------------------------------------*/
256 /*----------------------------------------------------------------------------*/
257 
258 void
259 cs_cdofb_vecteq_diffusion(double time_eval,
260  const cs_equation_param_t *eqp,
261  const cs_equation_builder_t *eqb,
262  const cs_cdofb_vecteq_t *eqc,
263  const cs_cell_mesh_t *cm,
265  cs_cell_sys_t *csys,
266  cs_cell_builder_t *cb);
267 
268 /*----------------------------------------------------------------------------*/
280 /*----------------------------------------------------------------------------*/
281 
282 void
283 cs_cdofb_vecteq_advection_diffusion(double time_eval,
284  const cs_equation_param_t *eqp,
285  const cs_cdofb_vecteq_t *eqc,
286  const cs_cell_mesh_t *cm,
287  cs_cell_sys_t *csys,
288  cs_cell_builder_t *cb);
289 
290 /*----------------------------------------------------------------------------*/
304 /*----------------------------------------------------------------------------*/
305 
306 static inline void
308  const cs_equation_param_t *eqp,
309  const cs_real_t t_eval,
310  const cs_real_t coef,
311  cs_cell_builder_t *cb,
313  cs_cell_sys_t *csys)
314 {
315  /* Reset the local contribution */
316  memset(csys->source, 0, csys->n_dofs*sizeof(cs_real_t));
317 
319  (cs_xdef_t *const *)eqp->source_terms,
320  cm,
321  eqb->source_mask,
322  eqb->compute_source,
323  t_eval,
324  NULL, /* No input structure */
325  cb, /* mass matrix is cb->hdg */
326  csys->source);
327 
328  /* Only cell-DoFs are involved */
329  const short int _off = 3*cm->n_fc;
330  csys->rhs[_off ] += coef * csys->source[_off ];
331  csys->rhs[_off + 1] += coef * csys->source[_off + 1];
332  csys->rhs[_off + 2] += coef * csys->source[_off + 2];
333 }
334 
335 /*----------------------------------------------------------------------------*/
347 /*----------------------------------------------------------------------------*/
348 
349 int
351  const cs_matrix_t *matrix,
352  const cs_equation_param_t *eqp,
353  cs_real_t *x,
354  cs_real_t *b);
355 
356 /*----------------------------------------------------------------------------*/
370 /*----------------------------------------------------------------------------*/
371 
372 static inline void
374  const cs_range_set_t *rs,
375  const cs_cell_mesh_t *cm,
376  const bool has_sourceterm,
377  cs_cdofb_vecteq_t *eqc,
380  cs_real_t rhs[])
381 {
382  assert(eqa != NULL); /* Sanity check */
383 
384  const short int n_f_dofs = 3*cm->n_fc;
385 
386  /* Matrix assembly */
387  eqc->assemble(csys, rs, eqa, mav);
388 
389  /* RHS assembly */
390 # pragma omp critical
391  {
392  for (short int f = 0; f < n_f_dofs; f++)
393  rhs[csys->dof_ids[f]] += csys->rhs[f];
394  }
395 
396  /* Reset the value of the source term for the cell DoF
397  Source term is only hold by the cell DoF in face-based schemes */
398  if (has_sourceterm) {
399  cs_real_t *st = eqc->source_terms + 3*cm->c_id;
400  for (int k = 0; k < 3; k++)
401  st[k] = csys->source[n_f_dofs + k];
402  }
403 }
404 
405 /*----------------------------------------------------------------------------*/
417 /*----------------------------------------------------------------------------*/
418 
419 void
421  const int field_id,
422  const cs_equation_param_t *eqp,
424  void *context);
425 
426 /*----------------------------------------------------------------------------*/
438 /*----------------------------------------------------------------------------*/
439 
440 void
442  const int field_id,
443  const cs_equation_param_t *eqp,
445  void *context);
446 
447 /*----------------------------------------------------------------------------*/
459 /*----------------------------------------------------------------------------*/
460 
461 void
463  const int field_id,
464  const cs_equation_param_t *eqp,
466  void *context);
467 
468 /*----------------------------------------------------------------------------*/
480 /*----------------------------------------------------------------------------*/
481 
482 void
484  const cs_real_t *rhs,
485  const cs_equation_param_t *eqp,
487  void *data,
488  cs_real_t *field_val);
489 
490 /*----------------------------------------------------------------------------*/
500 /*----------------------------------------------------------------------------*/
501 
502 void
503 cs_cdofb_vecteq_extra_op(const char *eqname,
504  const cs_field_t *field,
505  const cs_equation_param_t *eqp,
507  void *data);
508 
509 /*----------------------------------------------------------------------------*/
521 /*----------------------------------------------------------------------------*/
522 
523 cs_real_t *
524 cs_cdofb_vecteq_get_cell_values(void *context);
525 
526 /*----------------------------------------------------------------------------*/
536 /*----------------------------------------------------------------------------*/
537 
538 cs_real_t *
539 cs_cdofb_vecteq_get_face_values(void *context);
540 
541 /*----------------------------------------------------------------------------*/
550 /*----------------------------------------------------------------------------*/
551 
552 void
554  const char *eqname,
555  void *scheme_context);
556 
557 /*----------------------------------------------------------------------------*/
566 /*----------------------------------------------------------------------------*/
567 
568 void
570  const char *eqname,
571  void *scheme_context);
572 
573 /*----------------------------------------------------------------------------*/
574 
576 
577 #endif /* __CS_CDOFB_VECTEQ_H__ */
cs_cdofb_vecteq_solve_theta
void cs_cdofb_vecteq_solve_theta(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:1045
cs_equation_builder_t::source_mask
cs_mask_t * source_mask
Definition: cs_equation_common.h:98
cs_defs.h
cs_cdofb_vecteq_get_face_values
cs_real_t * cs_cdofb_vecteq_get_face_values(void *context)
Retrieve an array of values at mesh faces for the current context. The lifecycle of this array is man...
Definition: cs_cdofb_vecteq.c:1888
cs_cdofb_vecteq_assembly
static 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.h:373
cs_cdofb_vecteq_solve_steady_state
void cs_cdofb_vecteq_solve_steady_state(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:638
cs_xdef_t
Structure storing medata for defining a quantity in a very flexible way.
Definition: cs_xdef.h:126
cs_equation_common.h
fm
Definition: cs_field_pointer.h:134
field
Definition: field.f90:27
cs_cdofb_vecteq_finalize_common
void cs_cdofb_vecteq_finalize_common(void)
Free work buffer and general structure related to CDO face-based schemes.
Definition: cs_cdofb_vecteq.c:1415
cs_matrix_structure_t
struct _cs_matrix_structure_t cs_matrix_structure_t
Definition: cs_matrix.h:86
cs_cdofb_vecteq_get
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:1393
cs_source_term.h
cs_fuel_incl::b
double precision, save b
Definition: cs_fuel_incl.f90:146
cs_cdofb_vecteq_sourceterm
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_cell_builder_t *cb, cs_equation_builder_t *eqb, cs_cell_sys_t *csys)
Compute the term source for a vector-valued CDO-Fb scheme and add it to the local rhs.
Definition: cs_cdofb_vecteq.h:307
cs_cdofb_vecteq_diffusion
void cs_cdofb_vecteq_diffusion(double time_eval, 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_face_mesh_t *fm, cs_cell_sys_t *csys, cs_cell_builder_t *cb)
Build the local matrices arising from the diffusion, advection, reaction terms in CDO-Fb schemes.
Definition: cs_cdofb_vecteq.c:391
cs_matrix_assembler_values_t
struct _cs_matrix_assembler_values_t cs_matrix_assembler_values_t
Definition: cs_matrix_assembler.h:65
cs_face_mesh_t
Set of local quantities and connectivities related to a mesh face Structure used to get a better memo...
Definition: cs_cdo_local.h:212
END_C_DECLS
#define END_C_DECLS
Definition: cs_defs.h:468
cs_equation_builder_t::compute_source
cs_source_term_cellwise_t * compute_source[CS_N_MAX_SOURCE_TERMS]
Definition: cs_equation_common.h:107
cs_cell_mesh_t::c_id
cs_lnum_t c_id
Definition: cs_cdo_local.h:157
cs_cell_sys_t::dof_ids
cs_lnum_t * dof_ids
Definition: cs_cdo_local.h:100
cs_real_t
double cs_real_t
Floating-point value.
Definition: cs_defs.h:302
cs_equation_assemble_t
struct _cs_equation_assemble_t cs_equation_assemble_t
Definition: cs_equation_assemble.h:52
BEGIN_C_DECLS
#define BEGIN_C_DECLS
Definition: cs_defs.h:467
cs_cdofb_vecteq_init_cell_system
void cs_cdofb_vecteq_init_cell_system(const cs_flag_t cell_flag, const cs_cell_mesh_t *cm, const cs_equation_param_t *eqp, const cs_equation_builder_t *eqb, const cs_cdofb_vecteq_t *eqc, const cs_real_t dir_values[], const cs_real_t field_tn[], cs_real_t t_eval, cs_cell_sys_t *csys, cs_cell_builder_t *cb)
Initialize the local structure for the current cell.
Definition: cs_cdofb_vecteq.c:303
cs_matrix_t
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:90
cs_cdofb_vecteq_update_field
void cs_cdofb_vecteq_update_field(const cs_real_t *solu, const cs_real_t *rhs, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *data, cs_real_t *field_val)
Store solution(s) of the linear system into a field structure Update extra-field values if required (...
cs_cdofb_vecteq_write_restart
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:1991
cs_cdo_quantities.h
mesh
Definition: mesh.f90:26
cs_equation_param_t::source_terms
cs_xdef_t ** source_terms
Definition: cs_equation_param.h:353
cs_mesh.h
cs_equation_param_t
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition: cs_equation_param.h:148
cs_cdofb_vecteq_solve_system
int cs_cdofb_vecteq_solve_system(cs_sles_t *sles, const cs_matrix_t *matrix, const cs_equation_param_t *eqp, cs_real_t *x, cs_real_t *b)
Solve a linear system arising from a scalar-valued CDO-Fb scheme.
Definition: cs_cdofb_vecteq.c:545
cs_cdofb_vecteq_solve_implicit
void cs_cdofb_vecteq_solve_implicit(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:829
cs_cdofb_vecteq_init_values
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.c:1711
cs_time_step_t
time step descriptor
Definition: cs_time_step.h:51
cs_cell_sys_t::n_dofs
int n_dofs
Definition: cs_cdo_local.h:98
cs_field.h
cs_cdo_quantities_t
Definition: cs_cdo_quantities.h:94
cs_cdofb_vecteq_init_context
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:1451
cs_cdofb_vecteq_get_cell_values
cs_real_t * cs_cdofb_vecteq_get_cell_values(void *context)
Get the computed values at mesh cells from the inverse operation w.r.t. the static condensation (DoF ...
Definition: cs_cdofb_vecteq.c:1867
cs_cdofb_vecteq_free_context
void * cs_cdofb_vecteq_free_context(void *data)
Destroy a cs_cdofb_vecteq_t structure.
Definition: cs_cdofb_vecteq.c:1677
cs_restart_t
struct _cs_restart_t cs_restart_t
Definition: cs_restart.h:87
cs_cdo_connect_t
Definition: cs_cdo_connect.h:74
cs_range_set_t
Definition: cs_range_set.h:57
cs_cell_sys_t::source
double * source
Definition: cs_cdo_local.h:105
cs_cell_builder_t
Set of local and temporary buffers useful for building the algebraic system with a cellwise process....
Definition: cs_cdo_local.h:56
cs_equation_param_t::n_source_terms
int n_source_terms
Definition: cs_equation_param.h:352
cs_cell_mesh_t::n_fc
short int n_fc
Definition: cs_cdo_local.h:175
cs_cdo_connect.h
cs_time_step.h
cs_sles_t
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:68
cs_cdofb_vecteq_setup_bc
void cs_cdofb_vecteq_setup_bc(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[])
Set the boundary conditions known from the settings.
Definition: cs_cdofb_vecteq.c:259
cs_cdofb_vecteq_is_initialized
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:1296
cs_equation_builder_t
Store common elements used when building an algebraic system related to an equation.
Definition: cs_equation_common.h:61
cs_cdofb_vecteq_init_common
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:1318
cs_cdofb_vecteq_read_restart
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:1910
cs_flag_t
unsigned short int cs_flag_t
Definition: cs_defs.h:304
cs_sles.h
atincl::solu
double precision, dimension(:,:), allocatable solu
Definition: atincl.f90:265
cs_cdofb_vecteq_extra_op
void cs_cdofb_vecteq_extra_op(const char *eqname, const cs_field_t *field, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *data)
Predefined extra-operations related to this equation.
Definition: cs_cdofb_vecteq.c:1812
cs_source_term_compute_cellwise
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:799
cs_matrix.h
cs_equation_param.h
Structure and routines handling the specific settings related to a cs_equation_t structure.
_cs_cdofb_t
Definition: cs_cdofb_priv.h:51
cs_cell_sys_t
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:93
cs_cdofb_vecteq_matrix_structure
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:1378
cs_restart.h
cs_field_t
Field descriptor.
Definition: cs_field.h:124
cs_cdofb_vecteq_advection_diffusion
void cs_cdofb_vecteq_advection_diffusion(double time_eval, const cs_equation_param_t *eqp, const cs_cdofb_vecteq_t *eqc, const cs_cell_mesh_t *cm, cs_cell_sys_t *csys, cs_cell_builder_t *cb)
Build the local matrices arising from the diffusion, advection, reaction terms in CDO-Fb schemes.
Definition: cs_cdofb_vecteq.c:455
cs_cell_mesh_t
Set of local quantities and connectivities related to a mesh cell This is a key structure for all cel...
Definition: cs_cdo_local.h:146
cs_mesh_t
Definition: cs_mesh.h:63
cs_base.h
cs_cell_sys_t::rhs
double * rhs
Definition: cs_cdo_local.h:104
k
Definition: cs_field_pointer.h:70