My Project
programmer's documentation
|
Build an algebraic CDO vertex-based system for unsteady convection-diffusion-reaction of vector-valued equations with source terms. More...
#include "cs_defs.h"
#include <assert.h>
#include <float.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <bft_mem.h>
#include "cs_cdo_advection.h"
#include "cs_cdo_bc.h"
#include "cs_cdo_diffusion.h"
#include "cs_cdo_local.h"
#include "cs_cdo_time.h"
#include "cs_cdovb_priv.h"
#include "cs_equation_bc.h"
#include "cs_equation_common.h"
#include "cs_evaluate.h"
#include "cs_hodge.h"
#include "cs_log.h"
#include "cs_math.h"
#include "cs_mesh_location.h"
#include "cs_parall.h"
#include "cs_param.h"
#include "cs_post.h"
#include "cs_quadrature.h"
#include "cs_reco.h"
#include "cs_scheme_geometry.h"
#include "cs_search.h"
#include "cs_sles.h"
#include "cs_source_term.h"
#include "cs_timer.h"
#include "cs_cdovb_vecteq.h"
Functions | |
static cs_cell_builder_t * | _vbv_create_cell_builder (const cs_cdo_connect_t *connect) |
Initialize the local builder structure used for building the system cellwise. More... | |
static void | _vbv_setup (cs_real_t t_eval, const cs_mesh_t *mesh, const cs_equation_param_t *eqp, const cs_equation_builder_t *eqb, cs_flag_t vtx_bc_flag[], 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 of internal DoFs only if needed. More... | |
static void | _vbv_init_cell_system (cs_real_t t_eval, 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_real_t dir_values[], const cs_flag_t vtx_bc_flag[], const cs_lnum_t forced_ids[], const cs_real_t field_tn[], cs_cell_sys_t *csys, cs_cell_builder_t *cb) |
Initialize the local structure for the current cell. Case of vector-valued CDO-Vb schemes. More... | |
static void | _vbv_advection_diffusion_reaction (const cs_equation_param_t *eqp, const cs_equation_builder_t *eqb, const cs_cdovb_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. If asked, a mass matrix is also computed and stored in cb->hdg. Case of vector-valued CDO-Vb schemes. More... | |
static void | _vbv_apply_weak_bc (cs_real_t time_eval, const cs_equation_param_t *eqp, const cs_cdovb_vecteq_t *eqc, const cs_cell_mesh_t *cm, cs_face_mesh_t *fm, cs_cell_sys_t *csys, cs_cell_builder_t *cb) |
First pass to apply boundary conditions enforced weakly. Update the local system before applying the time scheme. Case of vector-valued CDO-Vb schemes. More... | |
static void | _vbv_enforce_values (const cs_equation_param_t *eqp, const cs_cdovb_vecteq_t *eqc, const cs_cell_mesh_t *cm, cs_face_mesh_t *fm, cs_cell_sys_t *csys, cs_cell_builder_t *cb) |
Second pass to apply boundary conditions. Only Dirichlet BCs which are enforced strongly. Apply also the enforcement of internal DoFs. Update the local system after applying the time scheme. Case of vector-valued CDO-Vb schemes. More... | |
static void | _vbv_compute_cw_sles_normalization (const cs_equation_param_t *eqp, const cs_cell_mesh_t *cm, const cs_cell_sys_t *csys, cs_real_t *rhs_norm) |
Cellwise computation of the renormalization coefficient for the the residual norm of the linear system. More... | |
static void | _vbv_sync_sles_normalization (const cs_equation_param_t *eqp, cs_real_t *rhs_norm) |
Last stage to compute of the renormalization coefficient for the the residual norm of the linear system. More... | |
static int | _vbv_solve_system (cs_sles_t *sles, const cs_matrix_t *matrix, const cs_equation_param_t *eqp, const cs_real_t rhs_norm, cs_real_t *x, cs_real_t *b) |
Solve a linear system arising from a vector-valued CDO-Vb scheme. More... | |
bool | cs_cdovb_vecteq_is_initialized (void) |
Check if the generic structures for building a CDO-Vb scheme are allocated. More... | |
void | cs_cdovb_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 vertex-based schemes Set shared pointers. More... | |
void | cs_cdovb_vecteq_get (cs_cell_sys_t **csys, cs_cell_builder_t **cb) |
Retrieve work buffers used for building a CDO system cellwise. More... | |
void | cs_cdovb_vecteq_finalize_common (void) |
Free work buffer and general structure related to CDO vertex-based schemes. More... | |
void * | cs_cdovb_vecteq_init_context (const cs_equation_param_t *eqp, int var_id, int bflux_id, cs_equation_builder_t *eqb) |
Initialize a cs_cdovb_vecteq_t structure storing data useful for building and managing such a scheme. More... | |
void * | cs_cdovb_vecteq_free_context (void *builder) |
Destroy a cs_cdovb_vecteq_t structure. More... | |
void | cs_cdovb_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-Vb schemes. More... | |
void | cs_cdovb_vecteq_set_dir_bc (cs_real_t t_eval, const cs_mesh_t *mesh, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context, cs_real_t field_val[]) |
Set the boundary conditions known from the settings when the fields stem from a vector CDO vertex-based scheme. More... | |
void | cs_cdovb_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 convection/diffusion/reaction equation with a CDO-Vb scheme. One works cellwise and then process to the assembly. More... | |
void | cs_cdovb_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 (for hybrid discretization) More... | |
cs_real_t * | cs_cdovb_vecteq_get_vertex_values (void *context) |
Retrieve an array of values at mesh vertices for the variable field associated to the given context The lifecycle of this array is managed by the code. So one does not have to free the return pointer. More... | |
cs_real_t * | cs_cdovb_vecteq_get_cell_values (void *context) |
Compute an array of values at mesh cells by interpolating the variable field associated to the given context located at mesh vertices The lifecycle of this array is managed by the code. So one does not have to free the return pointer. More... | |
void | cs_cdovb_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. More... | |
Variables | |
static cs_cell_sys_t ** | _vbv_cell_system = NULL |
static cs_cell_builder_t ** | _vbv_cell_builder = NULL |
static const cs_cdo_quantities_t * | cs_shared_quant |
static const cs_cdo_connect_t * | cs_shared_connect |
static const cs_time_step_t * | cs_shared_time_step |
static const cs_matrix_structure_t * | cs_shared_ms |
Build an algebraic CDO vertex-based system for unsteady convection-diffusion-reaction of vector-valued equations with source terms.
|
static |
Build the local matrices arising from the diffusion, advection, reaction terms. If asked, a mass matrix is also computed and stored in cb->hdg. Case of vector-valued CDO-Vb schemes.
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | eqb | pointer to a cs_equation_builder_t structure |
[in] | eqc | context for this kind of discretization |
[in] | cm | pointer to a cellwise view of the mesh |
[in,out] | csys | pointer to a cellwise view of the system |
[in,out] | cb | pointer to a cellwise builder |
|
static |
First pass to apply boundary conditions enforced weakly. Update the local system before applying the time scheme. Case of vector-valued CDO-Vb schemes.
[in] | time_eval | time at which analytical function are evaluated |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | eqc | context for this kind of discretization |
[in] | cm | pointer to a cellwise view of the mesh |
[in,out] | fm | pointer to a facewise view of the mesh |
[in,out] | csys | pointer to a cellwise view of the system |
[in,out] | cb | pointer to a cellwise builder |
|
inlinestatic |
Cellwise computation of the renormalization coefficient for the the residual norm of the linear system.
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | cm | pointer to a cellwise view of the mesh/quantities |
[in] | csys | pointer to a cellwise view of the system |
[in,out] | rhs_norm | quantity used for the RHS normalization |
|
static |
Initialize the local builder structure used for building the system cellwise.
[in] | connect | pointer to a cs_cdo_connect_t structure |
|
static |
Second pass to apply boundary conditions. Only Dirichlet BCs which are enforced strongly. Apply also the enforcement of internal DoFs. Update the local system after applying the time scheme. Case of vector-valued CDO-Vb schemes.
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | eqc | context for this kind of discretization |
[in] | cm | pointer to a cellwise view of the mesh |
[in,out] | fm | pointer to a facewise view of the mesh |
[in,out] | csys | pointer to a cellwise view of the system |
[in,out] | cb | pointer to a cellwise builder |
|
static |
Initialize the local structure for the current cell. Case of vector-valued CDO-Vb schemes.
[in] | t_eval | time at which one performs the evaluation |
[in] | cell_flag | flag related to the current cell |
[in] | cm | pointer to a cellwise view of the mesh |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | eqb | pointer to a cs_equation_builder_t structure |
[in] | dir_values | Dirichlet values associated to each vertex |
[in] | vtx_bc_flag | Flag related to BC associated to each vertex |
[in] | forced_ids | indirection in case of internal enforcement |
[in] | field_tn | values of the field at the last computed time |
[in,out] | csys | pointer to a cellwise view of the system |
[in,out] | cb | pointer to a cellwise builder |
|
static |
Set the boundary conditions known from the settings Define an indirection array for the enforcement of internal DoFs only if needed.
[in] | t_eval | time at which one evaluates BCs |
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | eqb | pointer to a cs_equation_builder_t structure |
[in,out] | vtx_bc_flag | pointer to an array of BC flag for each vtx |
[in,out] | p_dir_values | pointer to the Dirichlet values to set |
[in,out] | p_enforced_ids | pointer to the list of enforced vertices |
|
static |
Solve a linear system arising from a vector-valued CDO-Vb scheme.
[in,out] | sles | pointer to a cs_sles_t structure |
[in] | matrix | pointer to a cs_matrix_t structure |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | rhs_norm | quantity used for the RHS normalization |
[in,out] | x | solution of the linear system (in: initial guess) |
[in,out] | b | right-hand side (scatter/gather if needed) |
|
inlinestatic |
Last stage to compute of the renormalization coefficient for the the residual norm of the linear system.
[in] | eqp | pointer to a cs_equation_param_t structure |
[in,out] | rhs_norm | quantity used for the RHS normalization |
void cs_cdovb_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.
[in] | eqname | name of the equation |
[in] | field | pointer to a field structure |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in,out] | eqb | pointer to a cs_equation_builder_t structure |
[in,out] | data | pointer to cs_cdovb_vecteq_t structure |
void cs_cdovb_vecteq_finalize_common | ( | void | ) |
Free work buffer and general structure related to CDO vertex-based schemes.
void* cs_cdovb_vecteq_free_context | ( | void * | builder | ) |
Destroy a cs_cdovb_vecteq_t structure.
[in,out] | builder | pointer to a cs_cdovb_vecteq_t structure |
void cs_cdovb_vecteq_get | ( | cs_cell_sys_t ** | csys, |
cs_cell_builder_t ** | cb | ||
) |
Retrieve work buffers used for building a CDO system cellwise.
[out] | csys | pointer to a pointer on a cs_cell_sys_t structure |
[out] | cb | pointer to a pointer on a cs_cell_builder_t structure |
cs_real_t* cs_cdovb_vecteq_get_cell_values | ( | void * | context | ) |
Compute an array of values at mesh cells by interpolating the variable field associated to the given context located at mesh vertices The lifecycle of this array is managed by the code. So one does not have to free the return pointer.
[in,out] | context | pointer to a data structure cast on-the-fly |
cs_real_t* cs_cdovb_vecteq_get_vertex_values | ( | void * | context | ) |
Retrieve an array of values at mesh vertices for the variable field associated to the given context The lifecycle of this array is managed by the code. So one does not have to free the return pointer.
[in,out] | context | pointer to a data structure cast on-the-fly |
void cs_cdovb_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 vertex-based schemes Set shared pointers.
[in] | quant | additional mesh quantities struct. |
[in] | connect | pointer to a cs_cdo_connect_t struct. |
[in] | time_step | pointer to a time step structure |
[in] | ms | pointer to a cs_matrix_structure_t structure |
void* cs_cdovb_vecteq_init_context | ( | const cs_equation_param_t * | eqp, |
int | var_id, | ||
int | bflux_id, | ||
cs_equation_builder_t * | eqb | ||
) |
Initialize a cs_cdovb_vecteq_t structure storing data useful for building and managing such a scheme.
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | var_id | id of the variable field |
[in] | bflux_id | id of the boundary flux field |
[in,out] | eqb | pointer to a cs_equation_builder_t structure |
void cs_cdovb_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-Vb schemes.
[in] | t_eval | time at which one evaluates BCs |
[in] | field_id | id related to the variable field of this equation |
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in,out] | eqb | pointer to a cs_equation_builder_t structure |
[in,out] | context | pointer to the scheme context (cast on-the-fly) |
bool cs_cdovb_vecteq_is_initialized | ( | void | ) |
Check if the generic structures for building a CDO-Vb scheme are allocated.
void cs_cdovb_vecteq_set_dir_bc | ( | cs_real_t | t_eval, |
const cs_mesh_t * | mesh, | ||
const cs_equation_param_t * | eqp, | ||
cs_equation_builder_t * | eqb, | ||
void * | context, | ||
cs_real_t | field_val[] | ||
) |
Set the boundary conditions known from the settings when the fields stem from a vector CDO vertex-based scheme.
[in] | t_eval | time at which one evaluates BCs |
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in,out] | eqb | pointer to a cs_equation_builder_t structure |
[in,out] | context | pointer to the scheme context (cast on-the-fly) |
[in,out] | field_val | pointer to the values of the variable field |
void cs_cdovb_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 convection/diffusion/reaction equation with a CDO-Vb scheme. One works cellwise and then process to the assembly.
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | field_id | id of the variable field related to this equation |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in,out] | eqb | pointer to a cs_equation_builder_t structure |
[in,out] | context | pointer to cs_cdovb_vecteq_t structure |
void cs_cdovb_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 (for hybrid discretization)
[in] | solu | solution array |
[in] | rhs | rhs associated to this solution array |
[in] | eqp | pointer to a cs_equation_param_t structure |
[in,out] | eqb | pointer to a cs_equation_builder_t structure |
[in,out] | data | pointer to cs_cdovb_vecteq_t structure |
[in,out] | field_val | pointer to the current value of the field |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |