9.1
general documentation
cs_runge_kutta_integrator.cpp File Reference
#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_tcs_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 $ \vect{Rhs} $ is updated as follows: More...
 

Variables

cs_runge_kutta_integrator_t ** _rk_lst = nullptr
 
int _n_rk_integrators = 0
 

Function Documentation

◆ _runge_kutta_integrator_free()

static void _runge_kutta_integrator_free ( cs_runge_kutta_integrator_t **  rk)
static

Free memory.

Parameters
[in]rkdouble pointer to a runge-kutta integrator

◆ _runge_kutta_integrator_init_scheme()

static void _runge_kutta_integrator_init_scheme ( cs_runge_kutta_integrator_t rk,
cs_runge_kutta_scheme_t  scheme 
)
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_by_id()

cs_runge_kutta_integrator_t * cs_runge_kutta_integrator_by_id ( int  rk_id)

Return a Runge-Kutta integrator by id.

◆ cs_runge_kutta_integrator_create()

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.

Parameters
[in]schemeRK scheme type (RK_NONE, RK1, RK2, RK3, RK4)
[in]nameassociated equation or field's name
[in]dttime step
[in]dimvariable dimentsion
[in]n_eltsnumber of computational elements

return the RK integrator's id in the RK list

◆ cs_runge_kutta_integrators_destroy()

void cs_runge_kutta_integrators_destroy ( )

Clean all Runge-Kutta integrators.

◆ cs_runge_kutta_integrators_initialize()

void cs_runge_kutta_integrators_initialize ( )

Create RK integrator structures.

◆ cs_runge_kutta_stage_complete_scalar_rhs()

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 $ \vect{Rhs} $ is updated as follows:

\[
\vect{Rhs} = \vect{Rhs} - \sum_{\fij \in \Facei{\celli}}      \left(
       \dot{m}_\ij \left( \vect{\varia}_\fij - \vect{\varia}_\celli \right)
     - \mu_\fij \gradt_\fij \vect{\varia} \cdot \vect{S}_\ij  \right)
\]

Warning:

  • $ \vect{Rhs} $ has already been initialized before calling cs_balance_vector!
  • mind the sign minus
Parameters
[in]ctxReference to dispatch context
[in,out]rkpointer to a Runge-Kutta integrator
[in]idtvarindicator of the temporal scheme
[in]f_idfield id (or -1)
[in]imucppindicator
  • 0 do not multiply the convectiv term by Cp
  • 1 do multiply the convectiv term by Cp
[in]eqppointer to a cs_equation_param_t structure which contains variable calculation options
[in]bc_coeffsboundary condition structure for the variable
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at boundary faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in]pvarsolved variable (current time step)
[in]xcpparray of specific heat (Cp)

Variable Documentation

◆ _n_rk_integrators

int _n_rk_integrators = 0

◆ _rk_lst

cs_runge_kutta_integrator_t** _rk_lst = nullptr