![]() |
AdH Kraken 0.0.0
Next generation Adaptive Hydraulics
|
#include "adh.h"
Functions | |
void | fe_sw2_1d_explicit_flow (double *elem_rhs, SELEM_1D elem1d, int ie, SVECT2D nrml, double djac, double dt, double flux, double h_avg, double *u, double *v, int DEBUG, int DEBUG_LOCAL) |
Returns the 1D element explicit flow contributions to 2D-SW equations. More... | |
void | fe_sw2_1d_implicit_flow (double *elem_rhs, SELEM_1D elem1d, int ie, SVECT2D nrml, double djac, double dt, SVECT2D *v, double *h, int DEBUG, int DEBUG_LOCAL) |
Returns the 1D element implicit flow contributions to 2D-SW equations. More... | |
void | fe_sw2_1d_wall_friction (double *elem_rhs, SELEM_1D elem1d, int ie, SVECT2D nrml, double djac, double dt, double *u, double *v, double *resistance, int DEBUG, int DEBUG_LOCAL) |
Returns the 1D element wall friction contributions to 2D-SW equations. More... | |
void | fe_sw2_1d_pressure (double *elem_rhs, SELEM_1D elem1d, int ie, SVECT2D nrml, double djac, double dt, double *h, double g, int DEBUG, int DEBUG_LOCAL) |
Returns the 1D element pressure contributions to 2D-SW equations. More... | |
void | fe_sw2_1d_density_pressure (double *elem_rhs, SELEM_1D elem1d, int ie, SVECT2D nrml, double djac, double dt, double *h, double *d, double g, int DEBUG, int DEBUG_LOCAL) |
Returns the 1D element density pressure contributions to 2D-SW equations. More... | |
void | fe_sw2_1d_hvel (double *elem_rhs, SELEM_1D elem1d, int ie, SVECT2D nrml, double djac, double dt, SVECT2D *v, SVECT2D user_velocity, double *h, int DEBUG, int DEBUG_LOCAL) |
Returns the 1D element explicit flow and pressure contributions to 2D-SW equations. More... | |
This file collections functions responsible for the 2D shallow water body contributions to the elemental residual.
void fe_sw2_1d_density_pressure | ( | double * | elem_rhs, |
SELEM_1D | elem1d, | ||
int | ie, | ||
SVECT2D | nrml, | ||
double | djac, | ||
double | dt, | ||
double * | h, | ||
double * | d, | ||
double | g, | ||
int | DEBUG, | ||
int | DEBUG_LOCAL | ||
) |
Returns the 1D element density pressure contributions to 2D-SW equations.
[in,out] | elem_rhs | the 1D elemental residual array |
[in] | ie | the 1D element ID |
[in] | nrml | the 1D element normal |
[in] | djac | the 1D element jacobian |
[in] | dt | the current model time-step |
[in] | h | the water depth |
[in] | d | water density |
[in] | g | gravity |
void fe_sw2_1d_explicit_flow | ( | double * | elem_rhs, |
SELEM_1D | elem1d, | ||
int | ie, | ||
SVECT2D | nrml, | ||
double | djac, | ||
double | dt, | ||
double | flux, | ||
double | h_avg, | ||
double * | u, | ||
double * | v, | ||
int | DEBUG, | ||
int | DEBUG_LOCAL | ||
) |
Returns the 1D element explicit flow contributions to 2D-SW equations.
[in,out] | elem_rhs | the 1D elemental residual array |
[in] | ie | the 1D element ID |
[in] | nrml | the 1D element normal |
[in] | djac | the 1D element jacobian |
[in] | dt | the current model time-step |
[in] | flux | a user-supplied flux |
[in] | h_avg | the elementally averaged depth |
[in] | u | x-velocity |
[in] | v | y-velocity |
For the momentum contributions:
– if the flux is into the domain, the normal flux is split into x,y velocities, so that
and added to the residual as
– if the flux is out of the domain, an implicity u,v is used so that
Forms the 1d element matrix for the flux terms in the shallow water equations.
This is the unit discharge out of the boundary. Positive is out, Negative is in. NOTE: CJT :: all units here must be equivalent to Int [d(uh)/dt * dxdy] = m^4/s^2
void fe_sw2_1d_hvel | ( | double * | elem_rhs, |
SELEM_1D | elem1d, | ||
int | ie, | ||
SVECT2D | nrml, | ||
double | djac, | ||
double | dt, | ||
SVECT2D * | v, | ||
SVECT2D | user_velocity, | ||
double * | h, | ||
int | DEBUG, | ||
int | DEBUG_LOCAL | ||
) |
Returns the 1D element explicit flow and pressure contributions to 2D-SW equations.
[in,out] | elem_rhs | the 1D elemental residual array |
[in] | ie | the 1D element ID |
[in] | nrml | the 1D element normal |
[in] | djac | the 1D element jacobian |
[in] | dt | the current model time-step |
[in] | v | the water velocity |
[in] | user_velocity | the water velocity prescribed by user |
[in] | h | the water depth |
void fe_sw2_1d_implicit_flow | ( | double * | elem_rhs, |
SELEM_1D | elem1d, | ||
int | ie, | ||
SVECT2D | nrml, | ||
double | djac, | ||
double | dt, | ||
SVECT2D * | v, | ||
double * | h, | ||
int | DEBUG, | ||
int | DEBUG_LOCAL | ||
) |
Returns the 1D element implicit flow contributions to 2D-SW equations.
[in,out] | elem_rhs | the 1D elemental residual array |
[in] | ie | the 1D element ID |
[in] | nrml | the 1D element normal |
[in] | djac | the 1D element jacobian |
[in] | dt | the current model time-step |
[in] | v | the hydrodynamic velocity |
[in] | h | the water depth |
void fe_sw2_1d_pressure | ( | double * | elem_rhs, |
SELEM_1D | elem1d, | ||
int | ie, | ||
SVECT2D | nrml, | ||
double | djac, | ||
double | dt, | ||
double * | h, | ||
double | g, | ||
int | DEBUG, | ||
int | DEBUG_LOCAL | ||
) |
Returns the 1D element pressure contributions to 2D-SW equations.
[in,out] | elem_rhs | the 1D elemental residual array |
[in] | ie | the 1D element ID |
[in] | nrml | the 1D element normal |
[in] | djac | the 1D element jacobian |
[in] | dt | the current model time-step |
[in] | h | the water depth |
[in] | g | gravity |
void fe_sw2_1d_wall_friction | ( | double * | elem_rhs, |
SELEM_1D | elem1d, | ||
int | ie, | ||
SVECT2D | nrml, | ||
double | djac, | ||
double | dt, | ||
double * | u, | ||
double * | v, | ||
double * | resistance, | ||
int | DEBUG, | ||
int | DEBUG_LOCAL | ||
) |
Returns the 1D element wall friction contributions to 2D-SW equations.
[in,out] | elem_rhs | the 1D elemental residual array |
[in] | ie | the 1D element ID |
[in] | nrml | the 1D element normal |
[in] | djac | the 1D element jacobian |
[in] | dt | the current model time-step |
[in] | u | x-velocity |
[in] | v | y-velocity |
[in] | resistance | resistance |
Integrates the weak, discrete outflow boundary terms: