8.1
general documentation
cs_cdovb_scaleq.h File Reference
#include "cs_base.h"
#include "cs_cdo_connect.h"
#include "cs_cdo_local.h"
#include "cs_cdo_quantities.h"
#include "cs_cdo_toolbox.h"
#include "cs_equation_builder.h"
#include "cs_equation_param.h"
#include "cs_field.h"
#include "cs_matrix.h"
#include "cs_mesh.h"
#include "cs_source_term.h"
#include "cs_time_step.h"
+ Include dependency graph for cs_cdovb_scaleq.h:

Go to the source code of this file.

Functions

bool cs_cdovb_scaleq_is_initialized (void)
 Check if the generic structures for building a CDO-Vb scheme are allocated. More...
 
void cs_cdovb_scaleq_init_sharing (const cs_cdo_quantities_t *quant, const cs_cdo_connect_t *connect, const cs_time_step_t *time_step)
 Allocate work buffer and general structures related to CDO vertex-based schemes. Set shared pointers. More...
 
void cs_cdovb_scaleq_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_scaleq_finalize_sharing (void)
 Free work buffer and general structure related to CDO vertex-based schemes. More...
 
void * cs_cdovb_scaleq_init_context (const cs_equation_param_t *eqp, int var_id, int bflux_id, cs_equation_builder_t *eqb)
 Initialize a cs_cdovb_scaleq_t structure storing data useful for building and managing such a scheme. More...
 
void * cs_cdovb_scaleq_free_context (void *scheme_context)
 Destroy a cs_cdovb_scaleq_t structure. More...
 
void cs_cdovb_scaleq_setup (cs_real_t t_eval, const cs_mesh_t *mesh, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, cs_flag_t vtx_bc_flag[])
 Set the boundary conditions known from the settings Define an indirection array for the enforcement of internal DoFs only if needed. Case of scalar-valued CDO-Vb schemes. More...
 
void cs_cdovb_scaleq_init_properties (int t_id, cs_real_t t_eval, const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
 Set the main properties before the main loop on cells. Case of scalar-valued CDO-Vb schemes. More...
 
double cs_cdovb_scaleq_build_block_implicit (int t_id, cs_lnum_t c_id, bool diag_block, const cs_real_t f_val[], const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context, cs_cell_builder_t *cb, cs_cell_sys_t *csys)
 Build the cell system for the given cell id when the build occurs in a coupled system. Case of scalar-valued CDO-Vb schemes. More...
 
void cs_cdovb_scaleq_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 scalar-valued CDO-Vb schemes. More...
 
void cs_cdovb_scaleq_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 scalar steady-state convection/diffusion/reaction equation with a CDO-Vb scheme One works cellwise and then process to the assembly. More...
 
void cs_cdovb_scaleq_solve_steady_state_incr (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 scalar steady-state convection/diffusion/reaction equation with a CDO-Vb scheme One works cellwise and then process to the assembly. More...
 
void cs_cdovb_scaleq_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 scalar unsteady convection/diffusion/reaction equation with a CDO-Vb scheme Implicit time scheme is used to progress in time. One works cellwise and then process to the assembly. More...
 
void cs_cdovb_scaleq_solve_implicit_incr (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 scalar unsteady convection/diffusion/reaction equation with a CDO-Vb scheme Implicit time scheme is used to progress in time. One works cellwise and then process to the assembly. More...
 
void cs_cdovb_scaleq_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 scalar unsteady convection/diffusion/reaction equation with a CDO-Vb scheme Theta time scheme is used to progress in time. One works cellwise and then process to the assembly. More...
 
cs_real_tcs_cdovb_scaleq_get_vertex_values (void *context, bool previous)
 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_tcs_cdovb_scaleq_get_cell_values (void *context, bool previous)
 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...
 
cs_real_tcs_cdovb_scaleq_get_source_term_values (void *context)
 Retrieve the array storing the source term values at dual mesh cells. The lifecycle of this array is managed by the code. So one does not have to free the return pointer. More...
 
cs_cdo_balance_tcs_cdovb_scaleq_balance (const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
 Compute the balance for an equation over the full computational domain between time t_cur and t_cur + dt_cur Case of scalar-valued CDO vertex-based scheme. More...
 
void cs_cdovb_scaleq_apply_stiffness (const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context, const cs_property_t *pty, const cs_real_t *pot, cs_flag_t res_loc, cs_real_t *res)
 Compute the cellwise stiffness matrix associated to the property given as a parameter and multiply it by the "pot" array to define the resulting array associated to entities defined at res_loc. Case of scalar-valued CDO vertex-based scheme. More...
 
void cs_cdovb_scaleq_boundary_diff_flux (const cs_real_t t_eval, const cs_equation_param_t *eqp, const cs_real_t *pdi, cs_equation_builder_t *eqb, void *context, cs_real_t *vf_flux)
 Compute for each vertex of a boundary face, the portion of diffusive flux across the boundary face. The surface attached to each vertex corresponds to the intersection of its dual cell (associated to a vertex of the face) with the face. Case of scalar-valued CDO-Vb schemes. More...
 
void cs_cdovb_scaleq_flux_across_plane (const cs_real_t normal[], const cs_real_t *pdi, const cs_equation_param_t *eqp, int ml_id, cs_equation_builder_t *eqb, void *context, double *d_flux, double *c_flux)
 Compute the diffusive and convective flux across a list of faces Case of scalar-valued CDO-Vb schemes. More...
 
void cs_cdovb_scaleq_diff_flux_in_cells (const cs_real_t *values, const cs_equation_param_t *eqp, const cs_property_t *diff_pty, cs_real_t t_eval, cs_equation_builder_t *eqb, void *eq_context, cs_real_t *diff_flux)
 Cellwise and threaded computation of an approximation of a constant diffusive flux (a vector) in each cell. Case of scalar-valued CDO-Vb schemes. More...
 
void cs_cdovb_scaleq_diff_flux_dfaces (const cs_real_t *values, const cs_equation_param_t *eqp, const cs_property_t *diff_pty, cs_real_t t_eval, cs_equation_builder_t *eqb, void *eq_context, cs_real_t *diff_flux)
 Cellwise and threaded computation of the diffusive flux accross dual faces (a scalar) in each cell. Case of scalar-valued CDO-Vb schemes. More...
 
void cs_cdovb_scaleq_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 related fields/arrays. More...
 
void cs_cdovb_scaleq_extra_post (const cs_equation_param_t *eqp, cs_equation_builder_t *eqb, void *context)
 Predefined extra-operations related to this equation. More...
 

Function Documentation

◆ cs_cdovb_scaleq_apply_stiffness()

void cs_cdovb_scaleq_apply_stiffness ( const cs_equation_param_t eqp,
cs_equation_builder_t eqb,
void *  context,
const cs_property_t pty,
const cs_real_t pot,
cs_flag_t  res_loc,
cs_real_t res 
)

Compute the cellwise stiffness matrix associated to the property given as a parameter and multiply it by the "pot" array to define the resulting array associated to entities defined at res_loc. Case of scalar-valued CDO vertex-based scheme.

Parameters
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to a scheme builder structure
[in]ptypointer to the property related to the stiffness op.
[in]potarray to multiply with the stiffness matrix
[in]res_loclocation of entities in the resulting array
[in,out]resresulting array

◆ cs_cdovb_scaleq_balance()

cs_cdo_balance_t* cs_cdovb_scaleq_balance ( const cs_equation_param_t eqp,
cs_equation_builder_t eqb,
void *  context 
)

Compute the balance for an equation over the full computational domain between time t_cur and t_cur + dt_cur Case of scalar-valued CDO vertex-based scheme.

Parameters
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to a scheme builder structure
Returns
a pointer to a cs_cdo_balance_t structure

◆ cs_cdovb_scaleq_boundary_diff_flux()

void cs_cdovb_scaleq_boundary_diff_flux ( const cs_real_t  t_eval,
const cs_equation_param_t eqp,
const cs_real_t pdi,
cs_equation_builder_t eqb,
void *  context,
cs_real_t vf_flux 
)

Compute for each vertex of a boundary face, the portion of diffusive flux across the boundary face. The surface attached to each vertex corresponds to the intersection of its dual cell (associated to a vertex of the face) with the face. Case of scalar-valued CDO-Vb schemes.

Parameters
[in]t_evaltime at which one performs the evaluation
[in]eqppointer to a cs_equation_param_t structure
[in]pdipointer to an array of field values
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to a scheme builder structure
[in,out]vf_fluxpointer to the values of the diffusive flux

◆ cs_cdovb_scaleq_build_block_implicit()

double cs_cdovb_scaleq_build_block_implicit ( int  t_id,
cs_lnum_t  c_id,
bool  diag_block,
const cs_real_t  f_val[],
const cs_equation_param_t eqp,
cs_equation_builder_t eqb,
void *  context,
cs_cell_builder_t cb,
cs_cell_sys_t csys 
)

Build the cell system for the given cell id when the build occurs in a coupled system. Case of scalar-valued CDO-Vb schemes.

Warning: The treatment of the BCs differs from the "standard" case. Up to now, one assumes a Dirichlet or a Neumann for all equations (i.e. all blocks) and only an algebraic treatment is performed.

This function may be called inside an openMP loop.

Parameters
[in]t_idthread id if openMP is used
[in]c_idcell id
[in]diag_blocktrue if this a diagonal block in the full system
[in]f_valcurrent field values
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to a scheme context structure
[in,out]cbcell builder structure
[in,out]csyscell system structure
Returns
the value of the rhs_norm for the cellwise system

◆ cs_cdovb_scaleq_current_to_previous()

void cs_cdovb_scaleq_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 related fields/arrays.

Parameters
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to cs_cdovb_scaleq_t structure

◆ cs_cdovb_scaleq_diff_flux_dfaces()

void cs_cdovb_scaleq_diff_flux_dfaces ( const cs_real_t values,
const cs_equation_param_t eqp,
const cs_property_t diff_pty,
cs_real_t  t_eval,
cs_equation_builder_t eqb,
void *  eq_context,
cs_real_t diff_flux 
)

Cellwise and threaded computation of the diffusive flux accross dual faces (a scalar) in each cell. Case of scalar-valued CDO-Vb schemes.

Parameters
[in]valuesdiscrete values for the potential
[in]eqppointer to a cs_equation_param_t structure
[in]diff_ptypointer to the diffusion property to use
[in]t_evaltime at which one performs the evaluation
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]eq_contextpointer to cs_cdovb_scaleq_t structure
[in,out]diff_fluxvalues of the diffusive flux

◆ cs_cdovb_scaleq_diff_flux_in_cells()

void cs_cdovb_scaleq_diff_flux_in_cells ( const cs_real_t values,
const cs_equation_param_t eqp,
const cs_property_t diff_pty,
cs_real_t  t_eval,
cs_equation_builder_t eqb,
void *  eq_context,
cs_real_t diff_flux 
)

Cellwise and threaded computation of an approximation of a constant diffusive flux (a vector) in each cell. Case of scalar-valued CDO-Vb schemes.

Parameters
[in]valuesdiscrete values for the potential
[in]eqppointer to a cs_equation_param_t structure
[in]diff_ptypointer to the diffusion property to use
[in]t_evaltime at which one performs the evaluation
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]eq_contextpointer to cs_cdovb_scaleq_t structure
[in,out]diff_fluxvalue of the diffusive flux

◆ cs_cdovb_scaleq_extra_post()

void cs_cdovb_scaleq_extra_post ( const cs_equation_param_t eqp,
cs_equation_builder_t eqb,
void *  context 
)

Predefined extra-operations related to this equation.

Parameters
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to cs_cdovb_scaleq_t structure

◆ cs_cdovb_scaleq_finalize_sharing()

void cs_cdovb_scaleq_finalize_sharing ( void  )

Free work buffer and general structure related to CDO vertex-based schemes.

◆ cs_cdovb_scaleq_flux_across_plane()

void cs_cdovb_scaleq_flux_across_plane ( const cs_real_t  normal[],
const cs_real_t pdi,
const cs_equation_param_t eqp,
int  ml_id,
cs_equation_builder_t eqb,
void *  context,
double *  d_flux,
double *  c_flux 
)

Compute the diffusive and convective flux across a list of faces Case of scalar-valued CDO-Vb schemes.

Parameters
[in]normalindicate in which direction flux is > 0
[in]pdipointer to an array of field values
[in]eqppointer to a cs_equation_param_t structure
[in]ml_idid related to a cs_mesh_location_t struct.
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to data specific for this scheme
[in,out]d_fluxpointer to the value of the diffusive flux
[in,out]c_fluxpointer to the value of the convective flux

◆ cs_cdovb_scaleq_free_context()

void* cs_cdovb_scaleq_free_context ( void *  scheme_context)

Destroy a cs_cdovb_scaleq_t structure.

Parameters
[in,out]scheme_contextpointer to a cs_cdovb_scaleq_t structure
Returns
a NULL pointer

◆ cs_cdovb_scaleq_get()

void cs_cdovb_scaleq_get ( cs_cell_sys_t **  csys,
cs_cell_builder_t **  cb 
)

Retrieve work buffers used for building a CDO system cellwise.

Parameters
[out]csyspointer to a pointer on a cs_cell_sys_t structure
[out]cbpointer to a pointer on a cs_cell_builder_t structure

◆ cs_cdovb_scaleq_get_cell_values()

cs_real_t* cs_cdovb_scaleq_get_cell_values ( void *  context,
bool  previous 
)

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.

Parameters
[in,out]contextpointer to a data structure cast on-the-fly
[in]previousretrieve the previous state (true/false)
Returns
a pointer to an array of cs_real_t (size: n_cells)

◆ cs_cdovb_scaleq_get_source_term_values()

cs_real_t* cs_cdovb_scaleq_get_source_term_values ( void *  context)

Retrieve the array storing the source term values at dual mesh cells. The lifecycle of this array is managed by the code. So one does not have to free the return pointer.

Parameters
[in,out]contextpointer to a data structure cast on-the-fly
Returns
a pointer to an array of cs_real_t (size n_vertices)

◆ cs_cdovb_scaleq_get_vertex_values()

cs_real_t* cs_cdovb_scaleq_get_vertex_values ( void *  context,
bool  previous 
)

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.

Parameters
[in,out]contextpointer to a data structure cast on-the-fly
[in]previousretrieve the previous state (true/false)
Returns
a pointer to an array of cs_real_t (size: n_vertices)

◆ cs_cdovb_scaleq_init_context()

void* cs_cdovb_scaleq_init_context ( const cs_equation_param_t eqp,
int  var_id,
int  bflux_id,
cs_equation_builder_t eqb 
)

Initialize a cs_cdovb_scaleq_t structure storing data useful for building and managing such a scheme.

Parameters
[in]eqppointer to a cs_equation_param_t structure
[in]var_idid of the variable field
[in]bflux_idid of the boundary flux field
[in,out]eqbpointer to a cs_equation_builder_t struct.
Returns
a pointer to a new allocated cs_cdovb_scaleq_t structure

◆ cs_cdovb_scaleq_init_properties()

void cs_cdovb_scaleq_init_properties ( int  t_id,
cs_real_t  t_eval,
const cs_equation_param_t eqp,
cs_equation_builder_t eqb,
void *  context 
)

Set the main properties before the main loop on cells. Case of scalar-valued CDO-Vb schemes.

Parameters
[in]t_idthread id if > 0
[in]t_evaltime at which one evaluates BCs
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to a scheme context structure

◆ cs_cdovb_scaleq_init_sharing()

void cs_cdovb_scaleq_init_sharing ( const cs_cdo_quantities_t quant,
const cs_cdo_connect_t connect,
const cs_time_step_t time_step 
)

Allocate work buffer and general structures related to CDO vertex-based schemes. Set shared pointers.

Parameters
[in]quantadditional mesh quantities struct.
[in]connectpointer to a cs_cdo_connect_t struct.
[in]time_steppointer to a time step structure

◆ cs_cdovb_scaleq_init_values()

void cs_cdovb_scaleq_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 scalar-valued CDO-Vb schemes.

Parameters
[in]t_evaltime at which one evaluates BCs
[in]field_idid related to the variable field of this equation
[in]meshpointer to a cs_mesh_t structure
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to the scheme context (cast on-the-fly)

◆ cs_cdovb_scaleq_is_initialized()

bool cs_cdovb_scaleq_is_initialized ( void  )

Check if the generic structures for building a CDO-Vb scheme are allocated.

Returns
true or false

◆ cs_cdovb_scaleq_setup()

void cs_cdovb_scaleq_setup ( cs_real_t  t_eval,
const cs_mesh_t mesh,
const cs_equation_param_t eqp,
cs_equation_builder_t eqb,
cs_flag_t  vtx_bc_flag[] 
)

Set the boundary conditions known from the settings Define an indirection array for the enforcement of internal DoFs only if needed. Case of scalar-valued CDO-Vb schemes.

Parameters
[in]t_evaltime at which one evaluates BCs
[in]meshpointer to a cs_mesh_t structure
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]vtx_bc_flagpointer to an array of BC flag for each vtx

◆ cs_cdovb_scaleq_solve_implicit()

void cs_cdovb_scaleq_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 scalar unsteady convection/diffusion/reaction equation with a CDO-Vb scheme Implicit time scheme is used to progress in time. One works cellwise and then process to the assembly.

Parameters
[in]cur2prevtrue="current to previous" operation is performed
[in]meshpointer to a cs_mesh_t structure
[in]field_idid of the variable field related to this equation
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to cs_cdovb_scaleq_t structure

◆ cs_cdovb_scaleq_solve_implicit_incr()

void cs_cdovb_scaleq_solve_implicit_incr ( 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 scalar unsteady convection/diffusion/reaction equation with a CDO-Vb scheme Implicit time scheme is used to progress in time. One works cellwise and then process to the assembly.

Variant with an incremental approach. The system is modified to fit the incremental form (unknows are the increments and rhs corresponds to a residual). This is useful when the resolution is embedded into a non-linear process.

Parameters
[in]cur2prevNot used. Should be done before if needed.
[in]meshpointer to a cs_mesh_t structure
[in]field_idid of the variable field related to this equation
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to cs_cdovb_scaleq_t structure

◆ cs_cdovb_scaleq_solve_steady_state()

void cs_cdovb_scaleq_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 scalar steady-state convection/diffusion/reaction equation with a CDO-Vb scheme One works cellwise and then process to the assembly.

Parameters
[in]cur2prevtrue="current to previous" operation is performed
[in]meshpointer to a cs_mesh_t structure
[in]field_idid of the variable field related to this equation
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to cs_cdovb_scaleq_t structure

◆ cs_cdovb_scaleq_solve_steady_state_incr()

void cs_cdovb_scaleq_solve_steady_state_incr ( 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 scalar steady-state convection/diffusion/reaction equation with a CDO-Vb scheme One works cellwise and then process to the assembly.

Variant with an incremental approach. The system is modified to fit the incremental form (unknows are the increments and rhs corresponds to a residual). This is useful when the resolution is embedded into a non-linear process.

Parameters
[in]cur2prevNot used. Should be done before if needed.
[in]meshpointer to a cs_mesh_t structure
[in]field_idid of the variable field related to this equation
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to cs_cdovb_scaleq_t structure

◆ cs_cdovb_scaleq_solve_theta()

void cs_cdovb_scaleq_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 scalar unsteady convection/diffusion/reaction equation with a CDO-Vb scheme Theta time scheme is used to progress in time. One works cellwise and then process to the assembly.

Parameters
[in]cur2prevtrue="current to previous" operation is performed
[in]meshpointer to a cs_mesh_t structure
[in]field_idid of the variable field related to this equation
[in]eqppointer to a cs_equation_param_t structure
[in,out]eqbpointer to a cs_equation_builder_t structure
[in,out]contextpointer to cs_cdovb_scaleq_t structure