My Project
programmer's documentation
Loading...
Searching...
No Matches
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 */
67typedef struct _cs_cdofb_t cs_cdofb_vecteq_t;
68
69/*============================================================================
70 * Public function prototypes
71 *============================================================================*/
72
73/*----------------------------------------------------------------------------*/
80/*----------------------------------------------------------------------------*/
81
82bool
84
85/*----------------------------------------------------------------------------*/
96/*----------------------------------------------------------------------------*/
97
98void
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
114
115/*----------------------------------------------------------------------------*/
122/*----------------------------------------------------------------------------*/
123
124void
126 cs_cell_builder_t **cb);
127
128/*----------------------------------------------------------------------------*/
133/*----------------------------------------------------------------------------*/
134
135void
137
138/*----------------------------------------------------------------------------*/
150/*----------------------------------------------------------------------------*/
151
152void *
154 int var_id,
155 int bflux_id,
157
158/*----------------------------------------------------------------------------*/
166/*----------------------------------------------------------------------------*/
167
168void *
170
171/*----------------------------------------------------------------------------*/
184/*----------------------------------------------------------------------------*/
185
186void
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
211void
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,
222
223/*----------------------------------------------------------------------------*/
233/*----------------------------------------------------------------------------*/
234
235void
237 const cs_mesh_t *mesh,
238 const cs_equation_param_t *eqp,
240 cs_real_t *p_dir_values[]);
241
242/*----------------------------------------------------------------------------*/
256/*----------------------------------------------------------------------------*/
257
258void
259cs_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,
267
268/*----------------------------------------------------------------------------*/
280/*----------------------------------------------------------------------------*/
281
282void
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,
289
290/*----------------------------------------------------------------------------*/
304/*----------------------------------------------------------------------------*/
305
306static inline void
308 const cs_equation_param_t *eqp,
309 const cs_real_t t_eval,
310 const cs_real_t coef,
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
349int
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
372static 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
419void
421 const int field_id,
422 const cs_equation_param_t *eqp,
424 void *context);
425
426/*----------------------------------------------------------------------------*/
438/*----------------------------------------------------------------------------*/
439
440void
442 const int field_id,
443 const cs_equation_param_t *eqp,
445 void *context);
446
447/*----------------------------------------------------------------------------*/
459/*----------------------------------------------------------------------------*/
460
461void
463 const int field_id,
464 const cs_equation_param_t *eqp,
466 void *context);
467
468/*----------------------------------------------------------------------------*/
480/*----------------------------------------------------------------------------*/
481
482void
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
502void
503cs_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
523cs_real_t *
525
526/*----------------------------------------------------------------------------*/
536/*----------------------------------------------------------------------------*/
537
538cs_real_t *
540
541/*----------------------------------------------------------------------------*/
550/*----------------------------------------------------------------------------*/
551
552void
554 const char *eqname,
555 void *scheme_context);
556
557/*----------------------------------------------------------------------------*/
566/*----------------------------------------------------------------------------*/
567
568void
570 const char *eqname,
571 void *scheme_context);
572
573/*----------------------------------------------------------------------------*/
574
576
577#endif /* __CS_CDOFB_VECTEQ_H__ */
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
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
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
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
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
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
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
void * cs_cdofb_vecteq_free_context(void *data)
Destroy a cs_cdofb_vecteq_t structure.
Definition cs_cdofb_vecteq.c:1677
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
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
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
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
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
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 (...
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
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
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
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_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
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_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
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
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
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
#define BEGIN_C_DECLS
Definition cs_defs.h:467
double cs_real_t
Floating-point value.
Definition cs_defs.h:302
#define END_C_DECLS
Definition cs_defs.h:468
unsigned short int cs_flag_t
Definition cs_defs.h:304
struct _cs_equation_assemble_t cs_equation_assemble_t
Definition cs_equation_assemble.h:52
Structure and routines handling the specific settings related to a cs_equation_t structure.
@ k
Definition cs_field_pointer.h:70
@ fm
Definition cs_field_pointer.h:134
struct _cs_matrix_structure_t cs_matrix_structure_t
Definition cs_matrix.h:86
struct _cs_matrix_t cs_matrix_t
Definition cs_matrix.h:90
struct _cs_matrix_assembler_values_t cs_matrix_assembler_values_t
Definition cs_matrix_assembler.h:65
void matrix(const int *iconvp, const int *idiffp, const int *ndircp, const int *isym, const cs_real_t *thetap, const int *imucpp, const cs_real_t coefbp[], const cs_real_t cofbfp[], const cs_real_t rovsdt[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t xcpp[], cs_real_t da[], cs_real_t xa[])
Definition cs_matrix_building.c:111
struct _cs_restart_t cs_restart_t
Definition cs_restart.h:87
struct _cs_sles_t cs_sles_t
Definition cs_sles.h:68
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
Definition field.f90:27
Definition mesh.f90:26
Definition cs_cdofb_priv.h:51
Definition cs_cdo_connect.h:74
Definition cs_cdo_quantities.h:94
Set of local and temporary buffers useful for building the algebraic system with a cellwise process....
Definition cs_cdo_local.h:56
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_lnum_t c_id
Definition cs_cdo_local.h:157
short int n_fc
Definition cs_cdo_local.h:175
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
double * source
Definition cs_cdo_local.h:105
cs_lnum_t * dof_ids
Definition cs_cdo_local.h:100
double * rhs
Definition cs_cdo_local.h:104
int n_dofs
Definition cs_cdo_local.h:98
Store common elements used when building an algebraic system related to an equation.
Definition cs_equation_common.h:61
cs_mask_t * source_mask
Definition cs_equation_common.h:98
cs_source_term_cellwise_t * compute_source[CS_N_MAX_SOURCE_TERMS]
Definition cs_equation_common.h:107
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition cs_equation_param.h:148
cs_xdef_t ** source_terms
Definition cs_equation_param.h:353
int n_source_terms
Definition cs_equation_param.h:352
Set of local quantities and connectivities related to a mesh face Structure used to get a better memo...
Definition cs_cdo_local.h:212
Field descriptor.
Definition cs_field.h:124
Definition cs_mesh.h:63
Definition cs_range_set.h:57
time step descriptor
Definition cs_time_step.h:51
Structure storing medata for defining a quantity in a very flexible way.
Definition cs_xdef.h:126