File Reimann.h
Go to the source code of this file
#include "General.h"#include "Param.h"#include "Arrays.h"#include "Forcing.h"#include "MemManagement.h"#include "Util_CPU.h"
Public Functions
| Type | Name |
|---|---|
| __host__ void | UpdateButtingerXCPU (Param XParam, BlockP< T > XBlock, EvolvingP< T > XEv, GradientsP< T > XGrad, FluxP< T > XFlux, T * dtmax, T * zb) "Adaptive" second-order hydrostatic reconstruction. CPU version for the X-axis |
| __global__ void | UpdateButtingerXGPU (Param XParam, BlockP< T > XBlock, EvolvingP< T > XEv, GradientsP< T > XGrad, FluxP< T > XFlux, T * dtmax, T * zb) "Adaptive" second-order hydrostatic reconstruction. GPU version for the X-axis |
| __host__ void | UpdateButtingerYCPU (Param XParam, BlockP< T > XBlock, EvolvingP< T > XEv, GradientsP< T > XGrad, FluxP< T > XFlux, T * dtmax, T * zb) "Adaptive" second-order hydrostatic reconstruction. CPU version for the Y-axis |
| __global__ void | UpdateButtingerYGPU (Param XParam, BlockP< T > XBlock, EvolvingP< T > XEv, GradientsP< T > XGrad, FluxP< T > XFlux, T * dtmax, T * zb) "Adaptive" second-order hydrostatic reconstruction. GPU version for the Y-axis |
| __host__ __device__ T | hllc (T g, T delta, T epsi, T CFL, T cm, T fm, T hm, T hp, T um, T up, T & fh, T & fq) Calculate the Harten-Lax-van Leer-contact (HLLC) flux. |
Public Functions Documentation
function UpdateButtingerXCPU
"Adaptive" second-order hydrostatic reconstruction. CPU version for the X-axis
template<class T>
__host__ void UpdateButtingerXCPU (
Param XParam,
BlockP < T > XBlock,
EvolvingP < T > XEv,
GradientsP < T > XGrad,
FluxP < T > XFlux,
T * dtmax,
T * zb
)
Description
This function computes the flux term at the cell interface using the hydrostatic reconstruction from Buttinger et al (2019). This reconstruction is safe for steep slopes with thin water depth and is well-balanced, conserving "lake-at-rest" states.
For optimizing the code on CPU and GPU, there are 4 versions of this function: X or Y and CPU or GPU.
Where does this come from:
This scheme was adapted/modified from the Basilisk / B-Flood source code. I (CypB) changed the zr and zl term back to the Audusse type reconstruction http://basilisk.fr/sandbox/b-flood/saint-venant-topo.h
Reference: * Buttinger-Kreuzhuber, A., Horvath, Z., Noelle, S., Bloschl, G., and Waser, J.: A fast second-order shallow water scheme on two-dimensional structured grids over abrupt topography, Advances in water resources, 127, 89-108, 2019. * Kirstetter, G., Delestre, O., Lagree, P.-Y., Popinet, S., and Josserand, C.: B-flood 1.0: an open-source Saint-Venant model for flash flood simulation using adaptive refinement, Geosci. Model Dev. Discuss. [preprint], https://doi.org/10.5194/gmd-2021-15, in review, 2021.
Parameters:
XParamStructure containing the simulation parametersXBlockStructure containing the block informationXEvStructure containing the evolving variablesXGradStructure containing the gradients of the evolving variablesXFluxStructure containing the fluxes to be updateddtmaxArray to store the maximum allowable time step for each cellzbArray containing the bed elevation
Template parameters:
TData type, either float or double
Note:
This function is designed to be run on the CPU and should be called within a loop over all blocks.
function UpdateButtingerXGPU
"Adaptive" second-order hydrostatic reconstruction. GPU version for the X-axis
template<class T>
__global__ void UpdateButtingerXGPU (
Param XParam,
BlockP < T > XBlock,
EvolvingP < T > XEv,
GradientsP < T > XGrad,
FluxP < T > XFlux,
T * dtmax,
T * zb
)
Description
This function computes the flux term at the cell interface using the hydrostatic reconstruction from Buttinger et al (2019). This reconstruction is safe for steep slopes with thin water depth and is well-balanced, conserving "lake-at-rest" states.
For optimising the code on CPU and GPU there are 4 versions of this function: X or Y and CPU or GPU
Where does this come from:
This scheme was adapted/modified from the Basilisk / B-Flood source code. I (CypB) changed the zr and zl term back to the Audusse type reconstruction http://basilisk.fr/sandbox/b-flood/saint-venant-topo.h
Reference:
Kirstetter, G., Delestre, O., Lagree, P.-Y., Popinet, S., and Josserand, C.: B-flood 1.0: an open-source Saint-Venant model for flash flood simulation using adaptive refinement, Geosci. Model Dev. Discuss. [preprint], https://doi.org/10.5194/gmd-2021-15, in review, 2021.* Buttinger-Kreuzhuber, A., Horvath, Z., Noelle, S., Bloschl, G., and Waser, J.: A fast second-order shallow water scheme on two-dimensional structured grids over abrupt topography, Advances in water resources, 127, 89-108, 2019.
Parameters:
XParamStructure containing the simulation parametersXBlockStructure containing the block informationXEvStructure containing the evolving variablesXGradStructure containing the gradients of the evolving variablesXFluxStructure containing the fluxes to be updateddtmaxArray to store the maximum allowable time step for each cellzbArray containing the bed elevationTData type, either float or double
Note:
This function is designed to be run on the GPU and should be launched with a grid and block configuration that matches the number of blocks and block size.
function UpdateButtingerYCPU
"Adaptive" second-order hydrostatic reconstruction. CPU version for the Y-axis
template<class T>
__host__ void UpdateButtingerYCPU (
Param XParam,
BlockP < T > XBlock,
EvolvingP < T > XEv,
GradientsP < T > XGrad,
FluxP < T > XFlux,
T * dtmax,
T * zb
)
Description
This function computes the flux term at the cell interface using the hydrostatic reconstruction from Buttinger et al (2019). This reconstruction is safe for steep slope with thin water depth and is well-balanced meaning that it conserve the "lake-at-rest" states.
For optimising the code on CPU and GPU there are 4 versions of this function: X or Y and CPU or GPU
Where does this come from:
This scheme was adapted/modified from the Basilisk / B-Flood source code. I (CypB) changed the zr and zl term back to the Audusse type reconstruction http://basilisk.fr/sandbox/b-flood/saint-venant-topo.h
Reference: Kirstetter, G., Delestre, O., Lagree, P.-Y., Popinet, S., and Josserand, C.: B-flood 1.0: an open-source Saint-Venant model for flash flood simulation using adaptive refinement, Geosci. Model Dev. Discuss. [preprint], https://doi.org/10.5194/gmd-2021-15, in review, 2021.* Buttinger-Kreuzhuber, A., Horvath, Z., Noelle, S., Bloschl, G., and Waser, J.: A fast second-order shallow water scheme on two-dimensional structured grids over abrupt topography, Advances in water resources, 127, 89-108, 2019.
Parameters:
XParamStructure containing the simulation parametersXBlockStructure containing the block informationXEvStructure containing the evolving variablesXGradStructure containing the gradients of the evolving variablesXFluxStructure containing the fluxes to be updateddtmaxArray to store the maximum allowable time step for each cellzbArray containing the bed elevationTData type, either float or double
Note:
This function is designed to be run on the CPU and should be called within a loop over all blocks.
function UpdateButtingerYGPU
"Adaptive" second-order hydrostatic reconstruction. GPU version for the Y-axis
template<class T>
__global__ void UpdateButtingerYGPU (
Param XParam,
BlockP < T > XBlock,
EvolvingP < T > XEv,
GradientsP < T > XGrad,
FluxP < T > XFlux,
T * dtmax,
T * zb
)
Description
This function computes the flux term at the cell interface using the hydrostatic reconstruction from Buttinger et al (2019). This reconstruction is safe for steep slope with thin water depth and is well-balanced meaning that it conserve the "lake-at-rest" states.
For optimising the code on CPU and GPU there are 4 versions of this function: X or Y and CPU or GPU
Where does this come from:
This scheme was adapted/modified from the Basilisk / B-Flood source code. I (CypB) changed the zr and zl term back to the Audusse type reconstruction http://basilisk.fr/sandbox/b-flood/saint-venant-topo.h
Reference: Kirstetter, G., Delestre, O., Lagree, P.-Y., Popinet, S., and Josserand, C.: B-flood 1.0: an open-source Saint-Venant model for flash flood simulation using adaptive refinement, Geosci. Model Dev. Discuss. [preprint], https://doi.org/10.5194/gmd-2021-15, in review, 2021.* Buttinger-Kreuzhuber, A., Horvath, Z., Noelle, S., Bloschl, G., and Waser, J.: A fast second-order shallow water scheme on two-dimensional structured grids over abrupt topography, Advances in water resources, 127, 89-108, 2019.
Parameters:
XParamStructure containing the simulation parametersXBlockStructure containing the block informationXEvStructure containing the evolving variablesXGradStructure containing the gradients of the evolving variablesXFluxStructure containing the fluxes to be updateddtmaxArray to store the maximum allowable time step for each cellzbArray containing the bed elevation
Template parameters:
TData type, either float or double
Note:
This function is designed to be run on the GPU and should be launched with a grid and block configuration that matches the number of blocks and block size.
function hllc
Calculate the Harten-Lax-van Leer-contact (HLLC) flux.
template<class T>
__host__ __device__ T hllc (
T g,
T delta,
T epsi,
T CFL,
T cm,
T fm,
T hm,
T hp,
T um,
T up,
T & fh,
T & fq
)
Description
This an implementation of the HLLC solver.
Where does this come from:
This scheme was adapted/modified from the Basilisk source code. http://basilisk.fr/src/riemann.h
Reference: (Basilisk reference the scheme from Kurganov reference below) Kurganov, A., & Levy, D. (2002). Central-upwind schemes for the Saint-Venant system. Mathematical Modelling and Numerical Analysis, 36(3), 397-425.
Parameters:
gGravitational accelerationdeltaGrid resolution at the current levelepsiSmall number to prevent division by zeroCFLCourant-Friedrichs-Lewy number for stability conditioncmMetric term for spherical coordinates (1.0 for Cartesian)fmMetric term for spherical coordinates (1.0 for Cartesian)hmWater depth on the left side of the interfacehpWater depth on the right side of the interfaceumVelocity in the x-direction on the left side of the interfaceupVelocity in the x-direction on the right side of the interfacefhReference to store the computed flux for water depthfqReference to store the computed flux for momentum
Returns:
The maximum allowable time step based on the wave speeds and CFL condition
Template parameters:
TData type, either float or double
Note:
This function is designed to be run on both CPU and GPU.
The documentation for this class was generated from the following file src/Reimann.h