9.1
general documentation
cs_runge_kutta_integrator.h
Go to the documentation of this file.
1#ifndef RK_INTEGRATOR_H
2#define RK_INTEGRATOR_H
3
4/*============================================================================
5 * Explicit Runge-Kutta integrator utilities.
6 *============================================================================*/
7
8/*
9 This file is part of code_saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2025 EDF S.A.
12
13 This program is free software; you can redistribute it and/or modify it under
14 the terms of the GNU General Public License as published by the Free Software
15 Foundation; either version 2 of the License, or (at your option) any later
16 version.
17
18 This program is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21 details.
22
23 You should have received a copy of the GNU General Public License along with
24 this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25 Street, Fifth Floor, Boston, MA 02110-1301, USA.
26*/
27
28/*----------------------------------------------------------------------------*/
29
30#include "base/cs_defs.h"
31
32/*----------------------------------------------------------------------------
33 * Standard C library headers
34 *----------------------------------------------------------------------------*/
35
36#include <assert.h>
37#include <stdarg.h>
38#include <stdio.h>
39#include <stdlib.h>
40#include <string.h>
41
42/*----------------------------------------------------------------------------
43 * Local headers
44 *----------------------------------------------------------------------------*/
45
46#include "alge/cs_balance.h"
47#include "alge/cs_blas.h"
48
49#include "base/cs_array.h"
50#include "base/cs_base_accel.h"
52#include "base/cs_dispatch.h"
53#include "mesh/cs_mesh.h"
54
55/*----------------------------------------------------------------------------
56 * Header for the current file
57 *----------------------------------------------------------------------------*/
58
60
61/*----------------------------------------------------------------------------*/
62
63/*============================================================================
64 * Global variables
65 *============================================================================*/
66
67/*=============================================================================
68 * Public function prototypes
69 *============================================================================*/
70
71/*----------------------------------------------------------------------------*/
83/*----------------------------------------------------------------------------*/
84
85int
87 const char *name,
88 const cs_real_t *dt,
89 int dim,
90 cs_lnum_t n_elts);
91
92/*----------------------------------------------------------------------------*/
96/*----------------------------------------------------------------------------*/
97
98void
100
101/*----------------------------------------------------------------------------*/
112/*----------------------------------------------------------------------------*/
113
114template<cs_lnum_t stride>
115void
118 const cs_real_t *rho,
119 const cs_real_t *vol,
120 cs_real_t *pvara)
121{
122 if (rk == nullptr)
123 return;
124
125 rk->i_stage = 0;
126
127 const cs_lnum_t n_elts = rk->n_elts;
128
129 cs_real_t *mass = rk->mass;
130 cs_real_t *u_old = rk->u_old;
131
132 ctx.parallel_for(n_elts, [=] CS_F_HOST_DEVICE (cs_lnum_t i_elt) {
133 mass[i_elt] = rho[i_elt]*vol[i_elt];
134 for (cs_lnum_t k = 0; k < stride; k++)
135 u_old[stride*i_elt + k] = pvara[stride*i_elt + k];
136 });
137}
138
139/*----------------------------------------------------------------------------*/
147/*----------------------------------------------------------------------------*/
148
149template<cs_lnum_t stride>
150void
153 cs_real_t *pvar_stage)
154{
155 assert(rk != nullptr);
156
157 const int n_elts = rk->n_elts;
158 const cs_real_t *dt = rk->dt;
159 const cs_real_t *mass = rk->mass;
160
161 // get the current stage index
162 const int i_stg = rk->i_stage;
163
164 cs_real_t *u0 = rk->u_old;
165 cs_real_t *u_new = rk->u_new;
166 cs_real_t *rhs_stages = rk->rhs_stages;
167
168 const cs_real_t *a = rk->rk_coeff.a + RK_HIGHEST_ORDER*i_stg;
169
170 ctx.parallel_for (n_elts, [=] CS_F_HOST_DEVICE (cs_lnum_t i_elt) {
171 for (cs_lnum_t k = 0; k < stride; k++)
172 u_new[stride*i_elt + k] = u0[stride*i_elt + k];
173
174 for (int j_stg = 0; j_stg <= i_stg; j_stg++) {
175 cs_real_t *rhs = rhs_stages + j_stg*stride*n_elts;
176 for (cs_lnum_t k = 0; k < stride; k++)
177 u_new[stride*i_elt + k] += dt[i_elt]/mass[i_elt]*a[j_stg]
178 * rhs[i_elt*stride + k];
179 }
180
181 for (cs_lnum_t k = 0; k < stride; k++)
182 pvar_stage[stride*i_elt + k] = u_new[stride*i_elt + k];
183 });
184
185 ctx.wait();
186
187 rk->i_stage++;
188}
189
190/*----------------------------------------------------------------------------*/
203/*----------------------------------------------------------------------------*/
204
205template<int stride>
206void
209 cs_real_t *rhs_pvar)
210{
211 assert(rk != nullptr);
212 const int n_elts = rk->n_elts;
213 // get the current stage index
214 const int i_stg = rk->i_stage;
215
216 cs_real_t *rhs = rk->rhs_stages + i_stg*stride*n_elts;
217 ctx.parallel_for (n_elts, [=] CS_F_HOST_DEVICE (cs_lnum_t i_elt) {
218 for (cs_lnum_t k = 0; k < stride; k++)
219 rhs[stride*i_elt + k] = rhs_pvar[stride*i_elt + k];
220 });
221}
222
223/*----------------------------------------------------------------------------*/
274/*----------------------------------------------------------------------------*/
275
276void
280 int idtvar,
281 int f_id,
282 int imucpp,
284 cs_field_bc_coeffs_t *bc_coeffs,
285 const cs_real_t i_massflux[],
286 const cs_real_t b_massflux[],
287 const cs_real_t i_visc[],
288 const cs_real_t b_visc[],
289 cs_real_t viscel[][6],
290 const cs_real_t weighf[][2],
291 const cs_real_t weighb[],
292 int icvflb,
293 const int icvfli[],
294 cs_real_t pvar[],
295 const cs_real_t xcpp[]);
296
297/*----------------------------------------------------------------------------*/
357/*----------------------------------------------------------------------------*/
358
359template<int stride>
360void
363 int idtvar,
364 int f_id,
365 int ivisep,
367 cs_field_bc_coeffs_t *bc_coeffs,
368 const cs_real_t i_massflux[],
369 const cs_real_t b_massflux[],
370 const cs_real_t i_visc[],
371 const cs_real_t b_visc[],
372 const cs_real_t i_secvis[],
373 const cs_real_t b_secvis[],
374 cs_real_t viscel[][6],
375 const cs_real_t weighf[][2],
376 const cs_real_t weighb[],
377 int icvflb,
378 const int icvfli[],
379 cs_real_t pvar[][stride])
380{
381 assert(rk != nullptr);
382 const int n_elts = rk->n_elts;
383 // get the current stage index
384 const int i_stg = rk->i_stage;
385
386 const cs_lnum_t n_b_faces = cs_glob_mesh->n_b_faces;
387
388 cs_real_t *rhs = rk->rhs_stages + i_stg*stride*n_elts;
389
390 cs_alloc_mode_t amode = ctx.alloc_mode(false);
391
392 /* Allocate non reconstructed face value only if presence of limiter */
393
394 int df_limiter_id = -1;
395 cs_field_t *f = nullptr;
396 if (f_id > -1) {
397 f = cs_field_by_id(f_id);
398 df_limiter_id
399 = cs_field_get_key_int(f, cs_field_key_id("diffusion_limiter_id"));
400 }
401
402 const int ircflp = eqp->ircflu;
403 const int ircflb = (ircflp > 0) ? eqp->b_diff_flux_rc : 0;
404
405 cs_bc_coeffs_solve_t bc_coeffs_solve;
406 cs_init_bc_coeffs_solve(bc_coeffs_solve,
407 n_b_faces,
408 stride,
409 amode,
410 (df_limiter_id > -1 || ircflb != 1));
411
412 using var_t = cs_real_t[stride];
413
414 var_t *val_ip = (var_t *)bc_coeffs_solve.val_ip;
415 var_t *val_f = (var_t *)bc_coeffs_solve.val_f;
416 var_t *flux = (var_t *)bc_coeffs_solve.flux;
417 var_t *flux_lim = (var_t *)bc_coeffs_solve.flux_lim;
418
419 /* We compute the total explicit balance. */
420
421 int inc = 1;
422
423 /* The added convective scalar mass flux is:
424 * (thetex*Y_\face-imasac*Y_\celli)*mf.
425 * When building the explicit part of the rhs, one
426 * has to impose 0 on mass accumulation. */
427 int imasac = 0;
428
429 eqp->theta = 1;
430
431 cs_boundary_conditions_update_bc_coeff_face_values_strided<stride>
432 (ctx, f, bc_coeffs, inc, eqp, pvar,
433 val_ip, val_f, flux, flux_lim);
434
435 if (stride == 3)
437 f_id,
438 imasac,
439 inc, /* inc */
440 ivisep,
441 eqp,
442 nullptr, /* pvar == pvara */
443 (const cs_real_3_t *)pvar,
444 bc_coeffs,
445 (const cs_real_3_t *)val_f,
446 (const cs_real_3_t *)flux_lim,
447 i_massflux,
448 b_massflux,
449 i_visc,
450 b_visc,
451 i_secvis,
452 b_secvis,
453 viscel,
454 weighf,
455 weighb,
456 icvflb,
457 icvfli,
458 nullptr,
459 nullptr,
460 (cs_real_3_t *)rhs);
461
462 else if (stride == 6)
464 f_id,
465 imasac,
466 inc, /* inc */
467 eqp,
468 nullptr, /* pvar == pvara */
469 (const cs_real_6_t *)pvar,
470 bc_coeffs,
471 (const cs_real_6_t *)val_f,
472 (const cs_real_6_t *)flux_lim,
473 i_massflux,
474 b_massflux,
475 i_visc,
476 b_visc,
477 viscel,
478 weighf,
479 weighb,
480 icvflb,
481 icvfli,
482 (cs_real_6_t *)rhs);
483
484 eqp->theta = 0;
485
486 cs_clear_bc_coeffs_solve(bc_coeffs_solve);
487 ctx.wait();
488}
489
490/*----------------------------------------------------------------------------*/
495/*----------------------------------------------------------------------------*/
496
497inline bool
499{
500 if (rk == nullptr)
501 return false;
502 else
503 return (rk->scheme > CS_RK_NONE);
504}
505
506/*----------------------------------------------------------------------------*/
512/*----------------------------------------------------------------------------*/
513
514inline cs_real_t *
516(
519)
520{
521 assert(rk != nullptr);
522 assert(rk->i_stage > 0);
523
524 const cs_lnum_t n_elts = rk->n_elts;
525 // get the current stage index
526 const int i_stg = rk->i_stage;
527 const cs_real_t *dt = rk->dt;
528
529 const cs_real_t scaling_factor = rk->rk_coeff.c[i_stg - 1];
530 cs_real_t *scaled_dt = rk->scaled_dt;
531
532 ctx.parallel_for(n_elts, [=] CS_F_HOST_DEVICE (cs_lnum_t i) {
533 scaled_dt[i] = dt[i] * scaling_factor;
534 });
535 ctx.wait();
536
537 return rk->scaled_dt;
538}
539
540/*----------------------------------------------------------------------------*/
546/*----------------------------------------------------------------------------*/
547
548inline bool
550{
551 if (rk == nullptr)
552 return false;
553
554 return (rk->i_stage < rk->n_stages);
555}
556
557/*----------------------------------------------------------------------------*/
561/*----------------------------------------------------------------------------*/
562
565
566/*----------------------------------------------------------------------------*/
570/*----------------------------------------------------------------------------*/
571
572void
574
575/*----------------------------------------------------------------------------*/
576
577#endif /* __RK_INTEGRATOR_H__ */
auto parallel_for(cs_lnum_t n, F &&f, Args &&... args)
Definition: cs_dispatch.h:1570
void wait(void)
Wait (synchronize) until launched computations have finished.
Definition: cs_dispatch.h:1635
Definition: cs_dispatch.h:1711
void cs_balance_tensor(int idtvar, int f_id, int imasac, int inc, cs_equation_param_t *eqp, cs_real_6_t pvar[], const cs_real_6_t pvara[], const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_t val_f[][6], const cs_real_t flux[][6], 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_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_6_t rhs[])
Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport...
Definition: cs_balance.cpp:917
void cs_balance_vector(int idtvar, int f_id, int imasac, int inc, int ivisep, cs_equation_param_t *eqp, cs_real_3_t pvar[], const cs_real_3_t pvara[], const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_t val_f[][3], const cs_real_t flux[][3], 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 secvif[], const cs_real_t secvib[], cs_real_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_3_t i_pvar[], cs_real_3_t b_pvar[], cs_real_3_t smbr[])
Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport...
Definition: cs_balance.cpp:706
void cs_init_bc_coeffs_solve(cs_bc_coeffs_solve_t &c, cs_lnum_t n_b_faces, cs_lnum_t stride, cs_alloc_mode_t amode, bool limiter)
Initialize boundary condition coefficients solve arrays.
Definition: cs_boundary_conditions_set_coeffs.cpp:5172
void cs_clear_bc_coeffs_solve(cs_bc_coeffs_solve_t &c)
Free boundary condition coefficients solve arrays.
Definition: cs_boundary_conditions_set_coeffs.cpp:5204
#define CS_F_HOST_DEVICE
Definition: cs_defs.h:585
double cs_real_t
Floating-point value.
Definition: cs_defs.h:357
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:374
cs_real_t cs_real_6_t[6]
vector of 6 floating-point values
Definition: cs_defs.h:376
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:350
int cs_field_get_key_int(const cs_field_t *f, int key_id)
Return a integer value for a given key associated with a field.
Definition: cs_field.cpp:3157
cs_field_t * cs_field_by_id(int id)
Return a pointer to a field based on its id.
Definition: cs_field.cpp:2347
int cs_field_key_id(const char *name)
Return an id associated with a given key name.
Definition: cs_field.cpp:2663
@ k
Definition: cs_field_pointer.h:72
@ rho
Definition: cs_field_pointer.h:100
@ dt
Definition: cs_field_pointer.h:65
cs_alloc_mode_t
Definition: cs_mem.h:50
cs_mesh_t * cs_glob_mesh
void cs_runge_kutta_staging(cs_dispatch_context &ctx, cs_runge_kutta_integrator_t *rk, cs_real_t *pvar_stage)
Perform one Runge-Kutta staging.
Definition: cs_runge_kutta_integrator.h:151
cs_real_t * cs_runge_kutta_get_projection_time_scale_by_stage(cs_dispatch_context &ctx, cs_runge_kutta_integrator_t *rk)
Get the scaling factor when conducting a projection per stage.
Definition: cs_runge_kutta_integrator.h:516
void cs_runge_kutta_stage_set_initial_rhs(cs_dispatch_context &ctx, cs_runge_kutta_integrator_t *rk, cs_real_t *rhs_pvar)
Set initial rhs per stage for Runge-Kutta integrator. Align with legacy equations' building sequence....
Definition: cs_runge_kutta_integrator.h:207
void cs_runge_kutta_integrators_initialize()
Create RK integrator structures.
Definition: cs_runge_kutta_integrator.cpp:248
void cs_runge_kutta_integrators_destroy()
Clean all Runge-Kutta integrators.
Definition: cs_runge_kutta_integrator.cpp:301
void cs_runge_kutta_init_state(cs_dispatch_context &ctx, cs_runge_kutta_integrator_t *rk, const cs_real_t *rho, const cs_real_t *vol, cs_real_t *pvara)
Initialize a RK integrator at the begining of a time step.
Definition: cs_runge_kutta_integrator.h:116
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.
Definition: cs_runge_kutta_integrator.cpp:190
void cs_runge_kutta_stage_complete_rhs(cs_dispatch_context &ctx, cs_runge_kutta_integrator_t *rk, int idtvar, int f_id, int ivisep, 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[], const cs_real_t i_secvis[], const cs_real_t b_secvis[], 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[][stride])
prepare and complete rhs per stage for Runge-Kutta integrator. Align with legacy equations' building ...
Definition: cs_runge_kutta_integrator.h:361
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 ...
Definition: cs_runge_kutta_integrator.cpp:365
bool cs_runge_kutta_is_active(cs_runge_kutta_integrator_t *rk)
Indicate if Runge-Kutta integrator is activated.
Definition: cs_runge_kutta_integrator.h:498
cs_runge_kutta_integrator_t * cs_runge_kutta_integrator_by_id(int rk_id)
Return a Runge-Kutta integrator by id.
Definition: cs_runge_kutta_integrator.cpp:285
bool cs_runge_kutta_is_staging(cs_runge_kutta_integrator_t *rk)
indicate if the Runge-Kutta integrator is staging.
Definition: cs_runge_kutta_integrator.h:549
#define RK_HIGHEST_ORDER
Definition: cs_runge_kutta_integrator_priv.h:62
cs_runge_kutta_scheme_t
Definition: cs_runge_kutta_integrator_priv.h:71
@ CS_RK_NONE
Definition: cs_runge_kutta_integrator_priv.h:72
integer(c_int), pointer, save idtvar
option for a variable time step
Definition: optcal.f90:85
Definition: cs_field.h:105
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition: cs_equation_param.h:194
int b_diff_flux_rc
Definition: cs_equation_param.h:514
int ircflu
Definition: cs_equation_param.h:498
double theta
Definition: cs_equation_param.h:501
Field boundary condition descriptor (for variables)
Definition: cs_field.h:120
Field descriptor.
Definition: cs_field.h:156
cs_lnum_t n_b_faces
Definition: cs_mesh.h:99
double * a
Definition: cs_runge_kutta_integrator_priv.h:84
double * c
Definition: cs_runge_kutta_integrator_priv.h:85
Definition: cs_runge_kutta_integrator_priv.h:102
int i_stage
Definition: cs_runge_kutta_integrator_priv.h:107
cs_real_t * scaled_dt
Definition: cs_runge_kutta_integrator_priv.h:110
cs_runge_kutta_coeff_t rk_coeff
Definition: cs_runge_kutta_integrator_priv.h:125
const cs_real_t * dt
Definition: cs_runge_kutta_integrator_priv.h:109
cs_lnum_t n_elts
Definition: cs_runge_kutta_integrator_priv.h:112
cs_real_t * rhs_stages
Definition: cs_runge_kutta_integrator_priv.h:118
cs_real_t * u_old
Definition: cs_runge_kutta_integrator_priv.h:115
cs_runge_kutta_scheme_t scheme
Definition: cs_runge_kutta_integrator_priv.h:104
int n_stages
Definition: cs_runge_kutta_integrator_priv.h:106
cs_real_t * u_new
Definition: cs_runge_kutta_integrator_priv.h:116
cs_real_t * mass
Definition: cs_runge_kutta_integrator_priv.h:121