Code Documentation

class ArtificialTsunami2d : public tsunami_lab::setups::Setup

Two-dimensional artificial tsunami setup.

Public Functions

virtual t_real getHeight(t_real, t_real) const

Gets the water height at a given point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real, t_real) const

Gets the momentum in x-direction.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real i_x, t_real i_y) const

Gets the bathymetry.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

bathymetry

Private Functions

t_real getDisplacement(t_real i_x, t_real i_y) const

Gets displacement.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

t_real getF(t_real i_x, t_real) const

Gets the F.

Parameters:

i_x – x-coordinate of the queried point.

t_real getG(t_real, t_real i_y) const

Gets the G.

Parameters:

i_y – y-coordinate of the queried point.

Private Members

t_real m_delta = 20

delta for heights

class Csv : public tsunami_lab::io::IoWriter

Public Functions

Csv(t_real i_dxy, t_idx i_nx, t_idx i_ny, t_idx i_stride, t_idx i_ghostCellsX, t_idx i_ghostCellsY, t_real i_offsetX, t_real i_offsetY, t_real, t_real const *i_b, bool i_useCheckpoint)

Initialize the CSV File.

Parameters:
  • i_dxy – cell size.

  • i_nx – number of cells in x-direction.

  • i_ny – number of cells in y-direction.

  • i_stride – stride of the data.

  • i_ghostCellsX – number of ghost cells in x-direction.

  • i_ghostCellsY – number of ghost cells in y-direction.

  • i_offsetX – offset in x-direction.

  • i_offsetY – offset in y-direction.

  • i_b – bathymetry.

  • i_useCheckpoint – flag if checkpoint is used.

virtual void write(t_real const *i_h, t_real const *i_hu, t_real const *i_hv, t_real, t_idx i_nOut)

Writes the data to the output.

Parameters:
  • i_h – water height.

  • i_hu – momentum in x-direction.

  • i_hv – momentum in y-direction.

  • i_nOut – number of the output.

Public Static Functions

static void openCSV(const std::string &i_filePath, rapidcsv::Document &o_doc, bool header)

gets rapidcsv::Document and row count from CSV file

Parameters:
  • i_filePath – path to CSV file

  • o_doc – csv file as rapidcsv::Document

Private Members

const t_real m_dxy

cell width in x- and y-direction.

const t_idx m_nx

number of cells in x-direction.

const t_idx m_ny

number of cells in y-direction.

const t_idx m_stride

stride of the data arrays.

const t_idx m_ghostCellsX

number of ghost cells in x-direction.

const t_idx m_ghostCellsY

number of ghost cells in y-direction.

const t_real m_offsetX

offset in x-direction.

const t_real m_offsetY

offset in y-direction.

char *m_outputPath

file stream

t_real const *m_b

bathymetry.

class Custom1d : public tsunami_lab::setups::Setup

Public Functions

Custom1d(t_real i_heightLeft, t_real i_heightRight, t_real i_momentumLeft, t_real i_momentumRight, t_real i_middle)

Constructor.

Parameters:
  • i_heightLeft – water height on the left side of the dam.

  • i_heightRight – water height on the right side of the dam.

  • i_momentumLeft – momentum on the left side of the dam.

  • i_momentumRight – momentum on the right side of the dam.

  • i_middle – location (x-coordinate) of the middle.

virtual t_real getHeight(t_real i_x, t_real) const

Gets the water height at a given point.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real, t_real) const

Gets the momentum in x-direction.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real, t_real) const

Gets the bathymetry.

Returns:

bathymetry

Private Members

t_real m_heightLeft = 0

height on the left side

t_real m_heightRight = 0

height on the right side

t_real m_momentumLeft = 0

momentum on the left side

t_real m_momentumRight = 0

momentum on the right side

t_real m_middlePoint = 0

middlepoint

class DamBreak1d : public tsunami_lab::setups::Setup
#include <DamBreak1d.h>

1d dam break setup.

Public Functions

DamBreak1d(t_real i_heightLeft, t_real i_heightRight, t_real i_locationDam)

Constructor.

Parameters:
  • i_heightLeft – water height on the left side of the dam.

  • i_heightRight – water height on the right side of the dam.

  • i_locationDam – location (x-coordinate) of the dam.

virtual t_real getHeight(t_real i_x, t_real) const

Gets the water height at a given point.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real, t_real) const

Gets the momentum in x-direction.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real, t_real) const

Gets the bathymetry.

Returns:

bathymetry

Private Members

t_real m_heightLeft = 0

height on the left side

t_real m_heightRight = 0

height on the right side

t_real m_locationDam = 0

location of the dam

class DamBreak2d : public tsunami_lab::setups::Setup
#include <DamBreak2d.h>

2d dam break setup.

Public Functions

virtual t_real getHeight(t_real i_x, t_real i_y) const

Gets the water height at a given point.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real, t_real) const

Gets the momentum in x-direction.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real i_x, t_real i_y) const

Gets the bathymetry.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

bathymetry

class FWave

Public Static Functions

static void netUpdates(t_real i_hL, t_real i_hR, t_real i_huL, t_real i_huR, t_real i_bL, t_real i_bR, t_real o_netUpdateL[2], t_real o_netUpdateR[2])

Computes the net-updates.

Parameters:
  • i_hL – height of the left side.

  • i_hR – height of the right side.

  • i_huL – momentum of the left side.

  • i_huR – momentum of the right side.

  • i_bL – left bathymetry

  • i_bR – right bathymetry

  • o_netUpdateL – will be set to the net-updates for the left side; 0: height, 1: momentum.

  • o_netUpdateR – will be set to the net-updates for the right side; 0: height, 1: momentum.

Private Static Functions

static void waveSpeeds(t_real i_hL, t_real i_hR, t_real i_uL, t_real i_uR, t_real &o_waveSpeedL, t_real &o_waveSpeedR)

Computes the wave speeds.

Parameters:
  • i_hL – height of the left side.

  • i_hR – height of the right side.

  • i_uL – particle velocity of the leftside.

  • i_uR – particles velocity of the right side.

  • o_waveSpeedL – will be set to the speed of the wave propagating to the left.

  • o_waveSpeedR – will be set to the speed of the wave propagating to the right.

static void waveStrengths(t_real i_hL, t_real i_hR, t_real i_huL, t_real i_huR, t_real i_bL, t_real i_bR, t_real i_waveSpeedL, t_real i_waveSpeedR, t_real &o_strengthL, t_real &o_strengthR)

Computes the wave strengths.

Parameters:
  • i_hL – height of the left side.

  • i_hR – height of the right side.

  • i_huL – momentum of the left side.

  • i_huR – momentum of the right side.

  • i_bl – bathymetry of the left side.

  • i_br – bathymetry of the right side.

  • i_waveSpeedL – speed of the wave propagating to the left.

  • i_waveSpeedR – speed of the wave propagating to the right.

  • o_strengthL – will be set to the strength of the wave propagating to the left.

  • o_strengthR – will be set to the strength of the wave propagating to the right.

static void flux(t_real i_h, t_real i_hu, t_real &o_flux0, t_real &o_flux1)

Computes the Flux .

Parameters:
  • i_h – height

  • i_hu – momentum

  • o_flux0 – fluxfunction[0]

  • o_flux1 – fluxfunction[1]

static void deltaXPsi(t_real i_bL, t_real i_bR, t_real i_hL, t_real i_hR, t_real &o_deltaXPsi)

calculates the difference in the source term

Parameters:
  • i_bL – left bathymetry

  • i_bR – right bathymetry

  • i_hL – left height

  • i_hR – right height

  • o_deltaXPsi – difference in the source term

Private Static Attributes

static constexpr t_real m_g = 9.80665

gravity

static constexpr t_real m_gSqrt = 3.131557121

square root of gravity

class IoWriter

Subclassed by tsunami_lab::io::Csv, tsunami_lab::io::NetCdf

Public Functions

inline virtual ~IoWriter()

Virtual destructor for base class.

virtual void write(t_real const *i_h, t_real const *i_hu, t_real const *i_hv, t_real i_time, t_idx i_nOut) = 0

Writes the data to the output.

Parameters:
  • i_h – water height.

  • i_hu – momentum in x-direction.

  • i_hv – momentum in y-direction.

  • i_time – time in simulation.

  • i_nOut – number of the output.

class NetCdf : public tsunami_lab::io::IoWriter

Public Functions

NetCdf(t_real i_dxy, t_idx i_nx, t_idx i_ny, t_idx i_stride, t_idx i_ghostCellsX, t_idx i_ghostCellsY, t_real i_offsetX, t_real i_offsetY, t_real i_k, t_real const *i_b, bool i_useCheckpoint)

Initialize the netCdf File.

Parameters:
  • i_dxy – cell size.

  • i_nx – number of cells in x-direction.

  • i_ny – number of cells in y-direction.

  • i_stride – stride of the data.

  • i_ghostCellsX – number of ghost cells in x-direction.

  • i_ghostCellsY – number of ghost cells in y-direction.

  • i_offsetX – offset in x-direction.

  • i_offsetY – offset in y-direction.

  • i_k – cell size to be averaged.

  • i_b – bathymetry.

  • i_useCheckpoint – flag if checkpoint is used.

virtual void write(t_real const *i_h, t_real const *i_hu, t_real const *i_hv, t_real i_time, t_idx i_nOut)

Writes the data to the output.

Parameters:
  • i_h – water height.

  • i_hu – momentum in x-direction.

  • i_hv – momentum in y-direction.

  • i_time – time in simulation.

Public Static Functions

static void read(char *i_filePath, t_idx *o_nx, t_idx *o_ny, t_real **o_x, t_real **o_y, t_real **o_z)

read 3D data from netCdf file.

Parameters:
  • i_filePath – path to the netCdf file.

  • o_nx – number of cells in x-direction.

  • o_ny – number of cells in y-direction.

  • o_nz – number of cells in z-direction.

  • o_x – pointer to array of x-coordinates (Important: Gets a new dynamically allocated array written on it).

  • o_y – pointer to array of y-coordinates (Important: Gets a new dynamically allocated array written on it).

  • o_z – pointer to array of z-coordinates (Important: Gets a new dynamically allocated array written on it).

static void ncCheck(int i_status, char const *i_file, int i_line)

checks if netCdf operation was successful and prints error on failure

Parameters:
  • i_status – status code of the operation

  • i_file – file name of the operation

  • i_line – line number of the operation

static void readCheckpoint(char *i_filePath, t_idx *o_nx, t_idx *o_ny, bool *o_useFWave, bool *o_useCuda, tsunami_lab::t_boundary *o_boundaryL, tsunami_lab::t_boundary *o_boundaryR, tsunami_lab::t_boundary *o_boundaryB, tsunami_lab::t_boundary *o_boundaryT, tsunami_lab::t_real *o_endTime, tsunami_lab::t_real *o_width, tsunami_lab::t_real *o_xOffset, tsunami_lab::t_real *o_yOffset, tsunami_lab::t_real *o_hMax, std::string *o_stationFilePath, tsunami_lab::t_idx *o_nFrames, tsunami_lab::t_idx *o_k, tsunami_lab::t_idx *o_timeStep, tsunami_lab::t_idx *o_nOut, int *o_maxHours, t_real **o_b, t_real **o_h, t_real **o_hu, t_real **o_hv)

reads a checkpoint file

Parameters:
  • i_filePath – path to the checkpoint file

  • o_nx – number of cells in x-direction

  • o_ny – number of cells in y-direction

  • o_useFWave – flag if f-wave is used

  • o_useCuda – flag if cuda is used

  • o_boundaryL – left boundary

  • o_boundaryR – right boundary

  • o_boundaryB – bottom boundary

  • o_boundaryT – top boundary

  • o_endTime – end time of simulation

  • o_width – width of the domain

  • o_xOffset – x-offset of the domain

  • o_yOffset – y-offset of the domain

  • o_hMax – maximum water height

  • o_stationFilePath – path to the station file

  • o_nFrames – number of frames

  • o_k – cell size to be averaged

  • o_timeStep – time step

  • o_nOut – number of outputs

  • o_maxHours – maximum hours of simulation

  • o_b – bathymetry

  • o_h – water height

  • o_hu – momentum in x-direction

  • o_hv – momentum in y-direction

static void writeCheckpoint(t_idx i_nx, t_idx i_ny, t_idx i_stride, t_idx i_ghostCellsX, t_idx i_ghostCellsY, bool i_useFWave, bool i_useCuda, tsunami_lab::t_boundary i_boundaryL, tsunami_lab::t_boundary i_boundaryR, tsunami_lab::t_boundary i_boundaryB, tsunami_lab::t_boundary i_boundaryT, tsunami_lab::t_real i_endTime, tsunami_lab::t_real i_width, tsunami_lab::t_real i_xOffset, tsunami_lab::t_real i_yOffset, tsunami_lab::t_real i_hMax, std::string i_stationFilePath, tsunami_lab::t_idx i_nFrames, tsunami_lab::t_idx i_k, tsunami_lab::t_idx i_timeStep, tsunami_lab::t_idx i_nOut, int i_maxHours, const t_real *i_b, const t_real *i_h, const t_real *i_hu, const t_real *i_hv)

writes a checkpoint file

Parameters:
  • i_nx – number of cells in x-direction

  • i_ny – number of cells in y-direction

  • i_stride – stride of the data

  • i_ghostCellsX – number of ghost cells in x-direction

  • i_ghostCellsY – number of ghost cells in y-direction

  • i_useFWave – flag if f-wave is used

  • i_useCuda – flag if cuda is used

  • i_boundaryL – left boundary

  • i_boundaryR – right boundary

  • i_boundaryB – bottom boundary

  • i_boundaryT – top boundary

  • i_endTime – end time of simulation

  • i_width – width of the domain

  • i_xOffset – x-offset of the domain

  • i_yOffset – y-offset of the domain

  • i_hMax – maximum water height

  • i_stationFilePath – path to the station file

  • i_nFrames – number of frames

  • i_k – cell size to be averaged

  • i_timeStep – time step

  • i_nOut – number of outputs

  • i_maxHours – maximum hours of simulation

  • i_b – bathymetry

  • i_h – water height

  • i_hu – momentum in x-direction

  • i_hv – momentum in y-direction

Private Functions

void putVaraWithGhostcells(t_real const *i_data, int l_ncidp, int i_var, t_idx i_nOut, bool i_hasTime)

Prune Ghost Cells of Data.

Parameters:
  • i_data – input data

  • i_var – variable id to input

  • i_nOut – output time step

  • i_hasTime – true if data has time

Private Members

const t_real m_dxy

cell width in x- and y-direction.

const t_idx m_nx

number of cells in x-direction.

const t_idx m_ny

number of cells in y-direction.

const t_idx m_stride

stride of the data arrays.

const t_idx m_ghostCellsX

number of ghost cells in x-direction.

const t_idx m_ghostCellsY

number of ghost cells in y-direction.

const t_real m_offsetX

offset in x-direction.

const t_real m_offsetY

offset in y-direction.

const t_idx m_k

cell size to be averaged.

Private Static Functions

static tsunami_lab::t_boundary intToTBoundary(int i_boundary)

Converts saved boundary int to enum.

Parameters:

i_boundary – int representing t_boundary

static int tBoundaryToInt(tsunami_lab::t_boundary i_boundary)

Converts t_boundary enum to int.

Parameters:

i_boundary – t_boundary enum

class RareRare1d : public tsunami_lab::setups::Setup
#include <RareRare1d.h>

1d rarefaction rarefaction wave setup.

Public Functions

RareRare1d(t_real i_height, t_real i_momentum, t_real i_middlePoint)

Constructor.

Parameters:
  • i_height – water height on both sides.

  • i_momentum – momentum in x-direction.

  • i_middlePoint – location of the middle point.

virtual t_real getHeight(t_real, t_real) const

Gets the water height at a given point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real i_x, t_real) const

Gets the momentum in x-direction.

Parameters:

i_x – momentum of the queried point in x-direction.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real, t_real) const

Gets the bathymetry.

Returns:

bathymetry

Private Members

t_real m_height = 0

height of the water, same for both sides

t_real m_momentum = 0

momentum in x-direction

t_real m_middlePoint = 0

location of the middle point

class Roe

Public Static Functions

static void netUpdates(t_real i_hL, t_real i_hR, t_real i_huL, t_real i_huR, t_real o_netUpdateL[2], t_real o_netUpdateR[2])

Computes the net-updates.

Parameters:
  • i_hL – height of the left side.

  • i_hR – height of the right side.

  • i_huL – momentum of the left side.

  • i_huR – momentum of the right side.

  • o_netUpdateL – will be set to the net-updates for the left side; 0: height, 1: momentum.

  • o_netUpdateR – will be set to the net-updates for the right side; 0: height, 1: momentum.

Private Static Functions

static void waveSpeeds(t_real i_hL, t_real i_hR, t_real i_uL, t_real i_uR, t_real &o_waveSpeedL, t_real &o_waveSpeedR)

Computes the wave speeds.

Parameters:
  • i_hL – height of the left side.

  • i_hR – height of the right side.

  • i_uL – particle velocity of the leftside.

  • i_uR – particles velocity of the right side.

  • o_waveSpeedL – will be set to the speed of the wave propagating to the left.

  • o_waveSpeedR – will be set to the speed of the wave propagating to the right.

static void waveStrengths(t_real i_hL, t_real i_hR, t_real i_huL, t_real i_huR, t_real i_waveSpeedL, t_real i_waveSpeedR, t_real &o_strengthL, t_real &o_strengthR)

Computes the wave strengths.

Parameters:
  • i_hL – height of the left side.

  • i_hR – height of the right side.

  • i_huL – momentum of the left side.

  • i_huR – momentum of the right side.

  • i_waveSpeedL – speed of the wave propagating to the left.

  • i_waveSpeedR – speed of the wave propagating to the right.

  • o_strengthL – will be set to the strength of the wave propagating to the left.

  • o_strengthR – will be set to the strength of the wave propagating to the right.

Private Static Attributes

static constexpr t_real m_gSqrt = 3.131557121

square root of gravity

class Setup
#include <Setup.h>

Base setup.

Subclassed by tsunami_lab::setups::ArtificialTsunami2d, tsunami_lab::setups::Custom1d, tsunami_lab::setups::DamBreak1d, tsunami_lab::setups::DamBreak2d, tsunami_lab::setups::RareRare1d, tsunami_lab::setups::ShockShock1d, tsunami_lab::setups::Subcritical1d, tsunami_lab::setups::Supercritical1d, tsunami_lab::setups::TsunamiEvent1d, tsunami_lab::setups::TsunamiEvent2d

Public Functions

inline virtual ~Setup()

Virtual destructor for base class.

virtual t_real getHeight(t_real i_x, t_real i_y) const = 0

Gets the water height at a given point.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

water height at the given point.

virtual t_real getMomentumX(t_real i_x, t_real i_y) const = 0

Gets the momentum in x-direction.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real i_x, t_real i_y) const = 0

Gets the momentum in y-direction.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real i_x, t_real i_y) const = 0

Gets the bathymetry.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

bathymetry.

class ShockShock1d : public tsunami_lab::setups::Setup
#include <ShockShock1d.h>

1d shock shock setup.

Public Functions

ShockShock1d(t_real i_height, t_real i_momentum, t_real i_middlePoint)

Constructor.

Parameters:
  • i_height – water height.

  • i_momentum – water momentum.

  • i_middlePoint – location (x-coordinate) of the middlepoint.

virtual t_real getHeight(t_real, t_real) const

Gets the water height at a given point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real i_x, t_real) const

Gets the momentum in x-direction.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real, t_real) const

Gets the bathymetry.

Returns:

bathymetry

Private Members

t_real m_height = 0

height

t_real m_momentum = 0

momentum on the left side (and right side)

t_real m_middlePoint = 0

location of the middle

class Stations

Public Functions

Stations(const std::string path, t_real i_dxy, t_idx i_nx, t_idx i_ny, t_idx i_stride, t_idx i_ghostCellsX, t_idx i_ghostCellsY, t_real i_offsetX, t_real i_offsetY, t_real const *i_b, bool i_useCheckpoint)

Constructs the stations controller.

Parameters:
  • path – path to the stations file

  • i_dxy – cell width in x- and y-direction.

  • i_nx – number of cells in x-direction.

  • i_ny – number of cells in y-direction.

  • i_stride – stride of the data arrays in y-direction (x is assumed to be stride-1).

  • i_ghostCellsX – number of ghost cells in x-direction.

  • i_ghostCellsY – number of ghost cells in y-direction.

  • i_offsetX – offset in x-direction.

  • i_offsetY – offset in y-direction.

  • i_b – bathymetry of the cells; optional: use nullptr if not required.

  • i_useCheckpoint – flag if checkpoint is used.

t_real getT() const

Returns the output period.

Returns:

output period

std::vector<t_station> getStations() const

Returns the station data vector.

Returns:

station data vector

bool hasStations() const

Returns the number of stations.

std::string getPath() const

Returns the path to the stations file.

void write(t_real i_simTime, t_real const *i_h, t_real const *i_hu, t_real const *i_hv)

Writes the data as CSV to the given stream.

Parameters:
  • i_simTime – simulation time.

  • i_h – water height of the cells; optional: use nullptr if not required.

  • i_hu – momentum in x-direction of the cells; optional: use nullptr if not required.

  • i_hv – momentum in y-direction of the cells; optional: use nullptr if not required.

void init()

Initializes output files.

Private Members

t_real m_T = 0

output period

std::vector<t_station> m_stations

station data vector

bool m_hasStations = false

has stations

const std::string m_path

stations file path

const t_real m_dxy

cell width in x- and y-direction.

const t_idx m_nx

number of cells in x-direction.

const t_idx m_ny

number of cells in y-direction.

const t_idx m_stride

stride of the data arrays.

const t_idx m_ghostCellsX

number of ghost cells in x-direction.

const t_idx m_ghostCellsY

number of ghost cells in y-direction.

const t_real m_offsetX

offset in x-direction.

const t_real m_offsetY

offset in y-direction.

const t_real *m_b

bathymetry.

class Subcritical1d : public tsunami_lab::setups::Setup
#include <Subcritical1d.h>

1d subcritical setup.

Public Functions

virtual t_real getHeight(t_real, t_real) const

Gets the water height at a given point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real i_x, t_real) const

Gets the momentum in x-direction.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real i_x, t_real) const

Gets the bathymetry.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

bathymetry

class Supercritical1d : public tsunami_lab::setups::Setup
#include <Supercritical1d.h>

1d supercritical setup.

Public Functions

virtual t_real getHeight(t_real, t_real) const

Gets the water height at a given point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real i_x, t_real) const

Gets the momentum in x-direction.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real i_x, t_real) const

Gets the bathymetry.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

bathymetry

struct t_station
#include <constants.h>

station data

Public Members

std::string name
t_real x
t_real y
class TsunamiEvent1d : public tsunami_lab::setups::Setup
#include <TsunamiEvent1d.h>

1d TSUNAMIEVENT setup.

Public Functions

TsunamiEvent1d(rapidcsv::Document i_doc)

Construct a new TsunamiEvent1d object.

Parameters:

i_doc – csv file as rapidcsv::Document

virtual t_real getHeight(t_real i_x, t_real) const

Gets the water height at a given point.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real, t_real) const

Gets the momentum in x-direction.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real i_x, t_real) const

Gets the bathymetry.

Parameters:

i_x – x-coordinate of the queried point.

Returns:

bathymetry

Private Functions

t_real getDisplacement(t_real i_x) const

Get the initial displacement.

Parameters:

i_x

Returns:

displacement

t_real getBathymetryBin(t_real i_x) const

Get the Bathymetry of the CSV file.

Parameters:

i_x – x-coordinate of the queried point. (length)

Returns:

bathymetry

Private Members

rapidcsv::Document m_doc

csv file as rapidcsv::Document

t_real m_delta = 20

row count of csv file

delta for heights

class TsunamiEvent2d : public tsunami_lab::setups::Setup
#include <TsunamiEvent2d.h>

Two-dimensional articial tsunami setup.

Public Functions

TsunamiEvent2d(char *i_displacement, char *i_bathymetry, t_real *o_width, t_real *o_height, t_real *o_offsetX, t_real *o_offsetY)

Constructor.

Parameters:
  • i_displacement – path to the displacement file

  • i_bathymetry – path to the bathymetry file

  • o_width – width of the domain

  • o_height – height of the domain

  • o_offsetX – offset in x-direction

  • o_offsetY – offset in y-direction

~TsunamiEvent2d()

Destructor.

virtual t_real getHeight(t_real i_x, t_real i_y) const

Gets the water height at a given point.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

height at the given point.

virtual t_real getMomentumX(t_real, t_real) const

Gets the momentum in x-direction.

Returns:

momentum in x-direction.

virtual t_real getMomentumY(t_real, t_real) const

Gets the momentum in y-direction.

Returns:

momentum in y-direction.

virtual t_real getBathymetry(t_real i_x, t_real i_y) const

Gets the bathymetry.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

Returns:

bathymetry

Private Functions

t_real getDisplacement(t_real i_x, t_real i_y) const

Gets displacement.

Parameters:
  • i_x – x-coordinate of the queried point.

  • i_y – y-coordinate of the queried point.

t_real getBathymetryBin(t_real i_x, t_real i_y) const

Get the Bathymetry of the NC file.

Parameters:
  • i_x – x-coordinate of the queried point

  • i_y – y-coordinate of the queried point

Returns:

bathymetry

Private Members

t_real m_delta = 20

delta for heights

t_idx m_nbX

length of bathymetry in x direction.

t_idx m_nbY

length of bathymetry in y direction.

t_real *m_bathymetryX

x-value array of the bathymetry.

t_real *m_bathymetryY

y-value array of the bathymetry.

t_real *m_bathymetry

bathymetry array.

t_real m_stepBathX

cell size in x direction for bathymetry

t_real m_stepBathY

cell size in y direction for bathymetry

t_real m_offsetBathX

offset in x direction for bathymetry

t_real m_offsetBathY

offset in y direction for bathymetry

t_idx m_ndX

length of displacement in x direction.

t_idx m_ndY

length of displacement in y direction.

t_real *m_displacementX

x-value array of the displacement.

t_real *m_displacementY

y-value array of the displacement.

t_real *m_displacement

displacement array.

t_real m_stepDisplX

cell size in x direction for displacement

t_real m_stepDisplY

cell size in y direction for displacement

t_real m_offsetDisplX

offset in x direction for displacement

t_real m_offsetDisplY

offset in y direction for displacement

class WavePropagation

Subclassed by tsunami_lab::patches::WavePropagation1d, tsunami_lab::patches::WavePropagation2d, tsunami_lab::patches::WavePropagationCUDA

Public Functions

inline virtual ~WavePropagation()

Virtual destructor for base class.

virtual void timeStep(t_real i_scaling) = 0

Performs a time step.

Parameters:

i_scaling – scaling of the time step.

virtual t_idx getStride() = 0

Gets the stride in y-direction. x-direction is stride-1.

Returns:

stride in y-direction.

virtual t_idx getGhostCellsX() = 0

Gets number of ghost cells in x-direction.

Returns:

number of ghost cells in x-direction.

virtual t_idx getGhostCellsY() = 0

Gets number of ghost cells in y-direction.

Returns:

number of ghost cells in y-direction.

virtual t_real const *getHeight() = 0

Gets cells’ water heights.

Returns:

water heights.

virtual t_real const *getMomentumX() = 0

Gets the cells’ momenta in x-direction.

Returns:

momenta in x-direction.

virtual t_real const *getMomentumY() = 0

Gets the cells’ momenta in y-direction.

Returns:

momenta in y-direction.

virtual t_real const *getBathymetry() = 0

Gets the cells bathymetry.

Returns:

bathymetry.

virtual void setHeight(t_idx i_ix, t_idx i_iy, t_real i_h) = 0

Sets the height of the cell to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_iy – id of the cell in y-direction.

  • i_h – water height.

virtual void setMomentumX(t_idx i_ix, t_idx i_iy, t_real i_hu) = 0

Sets the momentum in x-direction to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_iy – id of the cell in y-direction.

  • i_hu – momentum in x-direction.

virtual void setMomentumY(t_idx i_ix, t_idx i_iy, t_real i_hv) = 0

Sets the momentum in y-direction to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_iy – id of the cell in y-direction.

  • i_hv – momentum in y-direction.

virtual void setBathymetry(t_idx i_ix, t_idx i_iy, t_real i_b) = 0

Sets the bathymetry to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_iy – id of the cell in y-direction.

  • i_b – bathymetry.

virtual void initGhostCells() = 0

Initializes the ghost cells. You can’t use set Operations afterwards.

virtual void prepareDataAccess() = 0

Prepares Data Access (you can use get Operations afterwards until the next time step)

class WavePropagation1d : public tsunami_lab::patches::WavePropagation

Public Functions

WavePropagation1d(t_idx i_nCells, bool i_useFWave, t_boundary i_boundaryLeft, t_boundary i_boundaryRight)

Constructs the 1d wave propagation solver.

Parameters:
  • i_nCells – number of cells.

  • i_useFWave – bool: true if FWave solver should be used, false if Roe solver should be used.

  • i_boundaryLeft – left boundary condition.

  • i_boundaryRight – right boundary condition.

~WavePropagation1d()

Destructor which frees all allocated memory.

virtual void timeStep(t_real i_scaling)

Performs a time step.

Parameters:

i_scaling – scaling of the time step (dt / dx).

inline virtual t_idx getStride()

Gets the stride in y-direction. x-direction is stride-1.

Returns:

stride in y-direction.

inline virtual t_real const *getHeight()

Gets cells’ water heights.

Returns:

water heights.

inline virtual t_real const *getMomentumX()

Gets the cells’ momenta in x-direction.

Returns:

momenta in x-direction.

inline virtual t_real const *getMomentumY()

Dummy function which returns a nullptr.

inline virtual t_real const *getBathymetry()

Gets the cells bathymetry.

Returns:

bathymetry.

inline virtual t_idx getGhostCellsX()

Gets number of ghost cells in x-direction.

Returns:

number of ghost cells in x-direction.

inline virtual t_idx getGhostCellsY()

Gets number of ghost cells in y-direction.

Returns:

number of ghost cells in y-direction.

inline virtual void setHeight(t_idx i_ix, t_idx, t_real i_h)

Sets the height of the cell to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_h – water height.

inline virtual void setMomentumX(t_idx i_ix, t_idx, t_real i_hu)

Sets the momentum in x-direction to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_hu – momentum in x-direction.

inline virtual void setMomentumY(t_idx, t_idx, t_real)

Dummy function since there is no y-momentum in the 1d solver.

inline virtual void setBathymetry(t_idx i_ix, t_idx, t_real i_b)

Sets the bathymetry of the cell to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_b – bathymetry.

inline virtual void initGhostCells()

Initializes the ghost cells.

inline virtual void prepareDataAccess()

Prepares data access (does nothing on non CUDA)

Private Functions

void setGhostCells()

Sets the values of the ghost cells according to t_boundary set in the class.

Private Members

unsigned short m_step = 0

current step which indicates the active values in the arrays below

t_idx m_nCells = 0

number of cells discretizing the computational domain

bool m_useFWave = true

bool if FWave solver is used

t_boundary m_boundaryLeft = t_boundary::OPEN

left boundary

t_boundary m_boundaryRight = t_boundary::OPEN

right boundary

t_real *m_h[2] = {nullptr, nullptr}

water heights for the current and next time step for all cells

t_real *m_hu[2] = {nullptr, nullptr}

momenta for the current and next time step for all cells

t_real *m_b = nullptr

bathymetry for the current and next time step for all cells

class WavePropagation2d : public tsunami_lab::patches::WavePropagation

Public Functions

WavePropagation2d(t_idx i_nCellsx, t_idx i_nCellsy, bool i_useFWave, t_boundary i_boundaryLeft, t_boundary i_boundaryRight, t_boundary i_boundaryBottom, t_boundary i_boundaryTop)

Constructs the 1d wave propagation solver.

Parameters:
  • i_nCellsx – number of cells in x direction.

  • i_nCellsy – number of cells in y direction.

  • i_useFWave – bool: true if FWave solver should be used, false if Roe solver should be used.

  • i_boundaryLeft – left boundary condition.

  • i_boundaryRight – right boundary condition.

  • i_boundaryBottom – bottom boundary condition.

  • i_boundaryTop – top boundary condition.

~WavePropagation2d()

Destructor which frees all allocated memory.

virtual void timeStep(t_real i_scaling)

Performs a time step.

Parameters:

i_scaling – scaling of the time step (dt / dx).

inline virtual t_idx getStride()

Gets the stride in y-direction. x-direction is stride-1.

Returns:

stride in y-direction.

inline virtual t_idx getGhostCellsX()

Gets number of ghost cells in x-direction.

Returns:

number of ghost cells in x-direction.

inline virtual t_idx getGhostCellsY()

Gets number of ghost cells in y-direction.

Returns:

number of ghost cells in y-direction.

inline virtual t_real const *getHeight()

Gets cells’ water heights.

Returns:

water heights.

inline virtual t_real const *getMomentumX()

Gets the cells’ momenta in x-direction.

Returns:

momenta in x-direction.

inline virtual t_real const *getMomentumY()

Dummy function which returns a nullptr.

inline virtual t_real const *getBathymetry()

Gets the cells bathymetry.

Returns:

bathymetry.

inline virtual void setHeight(t_idx i_ix, t_idx i_iy, t_real i_h)

Sets the height of the cell to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_h – water height.

inline virtual void setMomentumX(t_idx i_ix, t_idx i_iy, t_real i_hu)

Sets the momentum in x-direction to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_hu – momentum in x-direction.

inline virtual void setMomentumY(t_idx i_ix, t_idx i_iy, t_real i_hv)

Dummy function since there is no y-momentum in the 1d solver.

inline virtual void setBathymetry(t_idx i_ix, t_idx i_iy, t_real i_b)

Sets the bathymetry of the cell to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_b – bathymetry.

virtual void initGhostCells()

Initializes Bathymetry Ghost Cells.

inline virtual void prepareDataAccess()

Prepares data access (does nothing on non CUDA)

Private Functions

t_idx getCoord(t_idx i_x, t_idx i_y)

get Coordinate from x and y index

Parameters:
  • i_x – x index

  • i_y – y index

Returns:

coordinate

void setGhostCellsX()

Sets the values of the ghost cells according to t_boundary set in the class for the X Sweep.

void setGhostCellsY()

Sets the values of the ghost cells according to t_boundary set in the class for the Y Sweep.

Private Members

const t_idx m_nCellsx = 0

current step which indicates the active values in the arrays below

number of cells discretizing the computational domain in x direction

const t_idx m_nCellsy = 0

number of cells discretizing the computational domain in y direction

const bool m_useFWave = true

bool if FWave solver is used

const t_boundary m_boundaryLeft = t_boundary::OPEN

left boundary

const t_boundary m_boundaryRight = t_boundary::OPEN

right boundary

const t_boundary m_boundaryBottom = t_boundary::OPEN

bottom boundary

const t_boundary m_boundaryTop = t_boundary::OPEN

top boundary

t_real *m_h = nullptr

water heights for the current and next time step for all cells access array like m_h[i_x + i_y * (m_nCells + 2)]

t_real *m_hu = nullptr

momenta for the current and next time step for all cells in x-direction access array like m_hu[i_x + i_y * (m_nCells + 2)]

t_real *m_hv = nullptr

momenta for the current and next time step for all cells in y-direction access array like m_hv[i_x + i_y * (m_nCells + 2)]

t_real *m_hTemp = nullptr

temp array for height

t_real *m_huvTemp = nullptr

temp array for momentum

t_real *m_b = nullptr

bathymetry for the current and next time step for all cells access array like m_h[i_x + i_y * (m_nCells + 2)]

class WavePropagationCUDA : public tsunami_lab::patches::WavePropagation

Public Functions

WavePropagationCUDA(t_idx i_nCellsx, t_idx i_nCellsy, t_boundary i_boundaryLeft, t_boundary i_boundaryRight, t_boundary i_boundaryBottom, t_boundary i_boundaryTop)

Constructs the 1d wave propagation solver.

Parameters:
  • i_nCellsx – number of cells in x direction.

  • i_nCellsy – number of cells in y direction.

  • i_boundaryLeft – left boundary condition.

  • i_boundaryRight – right boundary condition.

  • i_boundaryBottom – bottom boundary condition.

  • i_boundaryTop – top boundary condition.

~WavePropagationCUDA()

Destructor which frees all allocated memory.

virtual void timeStep(t_real i_scaling)

Performs a time step.

Parameters:

i_scaling – scaling of the time step (dt / dx).

inline virtual t_idx getStride()

Gets the stride in y-direction. x-direction is stride-1.

Returns:

stride in y-direction.

inline virtual t_idx getGhostCellsX()

Gets number of ghost cells in x-direction.

Returns:

number of ghost cells in x-direction.

inline virtual t_idx getGhostCellsY()

Gets number of ghost cells in y-direction.

Returns:

number of ghost cells in y-direction.

inline virtual t_real const *getHeight()

Gets cells’ water heights.

Returns:

water heights.

inline virtual t_real const *getMomentumX()

Gets the cells’ momenta in x-direction.

Returns:

momenta in x-direction.

inline virtual t_real const *getMomentumY()

Dummy function which returns a nullptr.

inline virtual t_real const *getBathymetry()

Gets the cells bathymetry.

Returns:

bathymetry.

inline virtual void setHeight(t_idx i_ix, t_idx i_iy, t_real i_h)

Sets the height of the cell to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_h – water height.

inline virtual void setMomentumX(t_idx i_ix, t_idx i_iy, t_real i_hu)

Sets the momentum in x-direction to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_hu – momentum in x-direction.

inline virtual void setMomentumY(t_idx i_ix, t_idx i_iy, t_real i_hv)

Dummy function since there is no y-momentum in the 1d solver.

inline virtual void setBathymetry(t_idx i_ix, t_idx i_iy, t_real i_b)

Sets the bathymetry of the cell to the given value.

Parameters:
  • i_ix – id of the cell in x-direction.

  • i_b – bathymetry.

virtual void prepareDataAccess()

Prepares Data Access.

virtual void initGhostCells()

Sets the bathymetry in Ghost Cells.

Private Functions

t_idx getCoord(t_idx i_x, t_idx i_y)

get Coordinate from x and y index

Parameters:
  • i_x – x index

  • i_y – y index

Returns:

coordinate

Private Members

const t_idx m_nCellsx = 0

current step which indicates the active values in the arrays below

number of cells discretizing the computational domain in x direction

const t_idx m_nCellsy = 0

number of cells discretizing the computational domain in y direction

const t_boundary m_boundaryLeft = t_boundary::OPEN

left boundary

const t_boundary m_boundaryRight = t_boundary::OPEN

right boundary

const t_boundary m_boundaryBottom = t_boundary::OPEN

bottom boundary

const t_boundary m_boundaryTop = t_boundary::OPEN

top boundary

t_real *m_h = nullptr

water heights for the current and next time step for all cells access array like m_h[i_x + i_y * (m_nCells + 2)]

t_real *m_h_host = nullptr
t_real *m_hu = nullptr

momenta for the current and next time step for all cells in x-direction access array like m_hu[i_x + i_y * (m_nCells + 2)]

t_real *m_hu_host = nullptr
t_real *m_hv = nullptr

momenta for the current and next time step for all cells in y-direction access array like m_hv[i_x + i_y * (m_nCells + 2)]

t_real *m_hv_host = nullptr
t_real *m_hTemp = nullptr

temp array for height

t_real *m_huvTemp = nullptr

temp array for momentum

t_real *m_b = nullptr

bathymetry for the current and next time step for all cells access array like m_h[i_x + i_y * (m_nCells + 2)]

t_real *m_b_host = nullptr
namespace tsunami_lab

Author

Alexander Breuer (alex.breuer AT uni-jena.de)

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Constants / typedefs used throughout the code.

Author

Alexander Breuer (alex.breuer AT uni-jena.de)

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

IO-routines for writing a snapshot as Comma Separated Values (CSV).

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Abstract base class for IO-writes.

IO-routines for reading / writing netCdf data.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Stations Controller Class to get data at specific points in space

Author

Alexander Breuer (alex.breuer AT uni-jena.de)

DESCRIPTION

Base class of the wave propagation patches.

Author

Alexander Breuer (alex.breuer AT uni-jena.de)

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

One-dimensional wave propagation patch.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Two-dimensional wave propagation patch.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Two-dimensional wave propagation patch for CUDA.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Two-dimensional artificial tsunami problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

One-dimensional custom wave problem.

Author

Alexander Breuer (alex.breuer AT uni-jena.de)

DESCRIPTION

One-dimensional dam break problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Two-dimensional dam break problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

One-dimensional rarefaction rarefaction wave problem.

Author

Alexander Breuer (alex.breuer AT uni-jena.de)

DESCRIPTION

Simulation setup.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

One-dimensional shock shock problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

One-dimensional subcritical problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

One-dimensional supercritical problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

One-dimensional TSUNAMIEVENT problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

Two-dimensional tsunami event problem.

Author

Justus Dreßler (justus.dressler AT uni-jena.de)

Author

Thorsten Kröhl (thorsten.kroehl AT uni-jena.de)

Author

Julius Halank (julius.halank AT uni-jena.de)

DESCRIPTION

F-Wave solver for one dimensional shallow water equations.

Author

Alexander Breuer (alex.breuer AT uni-jena.de)

DESCRIPTION

Roe Riemann solver for the one-dimensional shallow water equations.

Typedefs

typedef std::size_t t_idx

integral type for cell-ids, pointer arithmetic, etc.

typedef float t_real

floating point type

Enums

enum class t_boundary

boundary conditions

Values:

enumerator WALL
enumerator OPEN
namespace io
namespace patches
namespace setups
namespace solvers
file constants.h
#include <cstddef>
#include <string>
file Csv.h
#include “../../constants.h
#include “../IoWriter.h
#include <cstring>
#include <string>
#include <iostream>
#include <vector>
#include “rapidcsv.h”
file IoWriter.h
#include “../constants.h
file NetCdf.h
#include “../../constants.h
#include “../IoWriter.h
file Stations.h
#include <string>
#include <vector>
#include “../../constants.h
#include <nlohmann/json.hpp>
file WavePropagation.h
#include “../constants.h
file WavePropagation1d.h
#include “../WavePropagation.h
file WavePropagation2d.h
#include “../WavePropagation.h
file WavePropagationCUDA.h
#include “../WavePropagation.h
file ArtificialTsunami2d.h
#include <cmath>
#include “../Setup.h
file Custom1d.h
#include “../Setup.h
file DamBreak1d.h
#include “../Setup.h
file DamBreak2d.h
#include <cmath>
#include “../Setup.h
file RareRare1d.h
#include “../Setup.h
file Setup.h
#include “../constants.h
file ShockShock1d.h
#include “../Setup.h
file Subcritical1d.h
#include “../Setup.h
file Supercritical1d.h
#include “../Setup.h
file TsunamiEvent1d.h
#include “../Setup.h
#include <cmath>
#include <algorithm>
#include “../../io/csv/Csv.h
#include “../../io/csv/rapidcsv.h”
file TsunamiEvent2d.h
#include <cmath>
#include “../Setup.h
file FWave.h
#include “../../constants.h
file Roe.h
#include “../../constants.h
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/artificialTsunami2d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/io/csv
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/custom1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/damBreak1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/damBreak2d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/solvers/fWave
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/io
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/io/netCdf
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/patches
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/rareRare1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/solvers/roe
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/shockShock1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/solvers
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/io/stations
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/subcritical1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/supercritical1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/tsunamiEvent1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/setups/tsunamiEvent2d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/patches/wavePropagation1d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/patches/wavePropagation2d
dir /home/docs/checkouts/readthedocs.org/user_builds/tsunami-lab/checkouts/latest/src/patches/wavePropagationCUDA