#include "base/cs_defs.h"#include <assert.h>#include <stdarg.h>#include <stdio.h>#include <stdlib.h>#include <string.h>#include "bft/bft_error.h"#include "bft/bft_printf.h"#include "base/cs_dispatch.h"#include "base/cs_field_default.h"#include "base/cs_field_pointer.h"#include "base/cs_runge_kutta_integrator.h"#include "base/cs_runge_kutta_integrator_priv.h"
Include dependency graph for cs_runge_kutta_integrator.cpp:Functions | |
| static void | _runge_kutta_integrator_init_scheme (cs_runge_kutta_integrator_t *rk, cs_runge_kutta_scheme_t scheme) |
| Fill integrator coeffs for a given scheme. More... | |
| static void | _runge_kutta_integrator_free (cs_runge_kutta_integrator_t **rk) |
| Free memory. More... | |
| int | cs_runge_kutta_integrator_create (cs_runge_kutta_scheme_t scheme, const char *name, const cs_real_t *dt, int dim, cs_lnum_t n_elts) |
| Create a RK integrator. More... | |
| void | cs_runge_kutta_integrators_initialize () |
| Create RK integrator structures. More... | |
| cs_runge_kutta_integrator_t * | cs_runge_kutta_integrator_by_id (int rk_id) |
| Return a Runge-Kutta integrator by id. More... | |
| void | cs_runge_kutta_integrators_destroy () |
| Clean all Runge-Kutta integrators. More... | |
| void | cs_runge_kutta_stage_complete_scalar_rhs (cs_dispatch_context &ctx, cs_runge_kutta_integrator_t *rk, int idtvar, int f_id, int imucpp, cs_equation_param_t *eqp, cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t viscel[][6], const cs_real_t weighf[][2], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_t pvar[], const cs_real_t xcpp[]) |
| prepare and complete rhs per stage for Runge-Kutta integrator. Align with legacy equations' building sequence. 1) Collect initial rhs done in cs_solve_navier_stokes/cs_solve_equation_scalar 2) Complete rhs with adding explicit part of the convection/diffusion balance More precisely, the right hand side | |
Variables | |
| cs_runge_kutta_integrator_t ** | _rk_lst = nullptr |
| int | _n_rk_integrators = 0 |
|
static |
Free memory.
| [in] | rk | double pointer to a runge-kutta integrator |
|
static |
Fill integrator coeffs for a given scheme.
Reference: Accuracy analysis of explicit Runge–Kutta methods applied to the incompressible Navier–Stokes equations. J.Comput.Phys. 231 (2012) 3041–3063
| cs_runge_kutta_integrator_t * cs_runge_kutta_integrator_by_id | ( | int | rk_id | ) |
Return a Runge-Kutta integrator by id.
| int cs_runge_kutta_integrator_create | ( | cs_runge_kutta_scheme_t | scheme, |
| const char * | name, | ||
| const cs_real_t * | dt, | ||
| int | dim, | ||
| cs_lnum_t | n_elts | ||
| ) |
Create a RK integrator.
| [in] | scheme | RK scheme type (RK_NONE, RK1, RK2, RK3, RK4) |
| [in] | name | associated equation or field's name |
| [in] | dt | time step |
| [in] | dim | variable dimentsion |
| [in] | n_elts | number of computational elements |
return the RK integrator's id in the RK list
| void cs_runge_kutta_integrators_destroy | ( | ) |
Clean all Runge-Kutta integrators.
| void cs_runge_kutta_integrators_initialize | ( | ) |
Create RK integrator structures.
| void cs_runge_kutta_stage_complete_scalar_rhs | ( | cs_dispatch_context & | ctx, |
| cs_runge_kutta_integrator_t * | rk, | ||
| int | idtvar, | ||
| int | f_id, | ||
| int | imucpp, | ||
| cs_equation_param_t * | eqp, | ||
| cs_field_bc_coeffs_t * | bc_coeffs, | ||
| const cs_real_t | i_massflux[], | ||
| const cs_real_t | b_massflux[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_t | viscel[][6], | ||
| const cs_real_t | weighf[][2], | ||
| const cs_real_t | weighb[], | ||
| int | icvflb, | ||
| const int | icvfli[], | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | xcpp[] | ||
| ) |
prepare and complete rhs per stage for Runge-Kutta integrator. Align with legacy equations' building sequence. 1) Collect initial rhs done in cs_solve_navier_stokes/cs_solve_equation_scalar 2) Complete rhs with adding explicit part of the convection/diffusion balance More precisely, the right hand side
is updated as follows:
Warning:
| [in] | ctx | Reference to dispatch context |
| [in,out] | rk | pointer to a Runge-Kutta integrator |
| [in] | idtvar | indicator of the temporal scheme |
| [in] | f_id | field id (or -1) |
| [in] | imucpp | indicator
|
| [in] | eqp | pointer to a cs_equation_param_t structure which contains variable calculation options |
| [in] | bc_coeffs | boundary condition structure for the variable |
| [in] | i_massflux | mass flux at interior faces |
| [in] | b_massflux | mass flux at boundary faces |
| [in] | i_visc | |
| [in] | b_visc | |
| [in] | viscel | symmetric cell tensor |
| [in] | weighf | internal face weight between cells i j in case of tensor diffusion |
| [in] | weighb | boundary face weight for cells i in case of tensor diffusion |
| [in] | icvflb | global indicator of boundary convection flux
|
| [in] | icvfli | boundary face indicator array of convection flux
|
| [in] | pvar | solved variable (current time step) |
| [in] | xcpp | array of specific heat (Cp) |
| int _n_rk_integrators = 0 |
| cs_runge_kutta_integrator_t** _rk_lst = nullptr |