![]() |
AdH Kraken 0.0.0
Next generation Adaptive Hydraulics
|
#include "adh.h"
Functions | |
void | fe_sw2_wd_average (SVECT *elem_nds, double *depth, SVECT2D *v_wd, double *f_wd, double djac, double *f_wd_avg, SVECT2D *v_wd_avg) |
Wet/Dry routine for calculating the wet-dry integral of shallow water a generic function, f, and/or vector, v: ![]() | |
void | fe_sw2_wd_convection_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *elem_rhs) |
Add the 2D shallow water wet-dry wrapped convection terms to the FE residual. More... | |
void | fe_sw2_wd_continuity_temporal_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Routine for wet-dry integrating the SW 2D continuity temporal term on a triangle. More... | |
void | fe_sw2_wd_pressure_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Add the 2D shallow water wet-dry wrapped pressure terms to the FE residual. More... | |
void | fe_sw2_wd_bodyForce_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Add the 2D shallow water wet-dry wrapped body-force terms to the FE residual. More... | |
void | fe_sw2_wd_boundaryPressure_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *elem_rhs) |
Add the 2D shallow water wet-dry wrapped boundary pressure terms to the FE residual. More... | |
void | fe_sw2_wd_densityPressure_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Add the 2D shallow water wet-dry wrapped density pressure terms to the FE residual. More... | |
void | fe_sw2_wd_densityBodyForce_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Add the 2D shallow water wet-dry wrapped density body-force terms to the FE residual. More... | |
void | fe_sw2_wd_densityBoundaryPressure_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Add the 2D shallow water wet-dry wrapped density boundary pressure terms to the FE residual. More... | |
void | fe_sw2_wd_integrate_triangle_f (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Wet/Dry routine for calculating the wet-dry integral of shallow water H, u and v: ![]() | |
void | fe_sw2_wd_gls_convection_triangle (SVECT *elem_nds, double *elem_head, SVECT2D *grad_phi, SVECT2D *v, SVECT2D *v_wet_dry, double *f, double *f_wet_dry, double djac, double *vars, double *rhs) |
Add the 2D shallow water wet-dry wrapped *nonconservative convection terms to the FE residual. More... | |
This file collects the 2D shallow water wet/dry integrations. All function arguements must be consistent with fe_sw2_wet_dry_wrapper.
|
inline |
Wet/Dry routine for calculating the wet-dry integral of shallow water a generic function, f, and/or vector, v: .
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | not used |
[in] | v | not used |
[in] | v_wet_dry | wet-dry velocities (averaged here) |
[in] | f | not used |
[in] | f_wet_dry | the scalar function to be averaged (averaged here) |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
|
inline |
Add the 2D shallow water wet-dry wrapped body-force terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | not used |
[in] | v | the fully wet bathymetry gradient |
[in] | v_wet_dry | not used |
[in] | f | not used |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
|
inline |
Add the 2D shallow water wet-dry wrapped boundary pressure terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | the wet-dry nodal positions |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | not used |
[in] | v | not used |
[in] | v_wet_dry | not used |
[in] | f | not used |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
|
inline |
Routine for wet-dry integrating the SW 2D continuity temporal term on a triangle.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | not used |
[in] | v | not used |
[in] | v_wet_dry | not used |
[in] | f | not used |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | the second order time-derivative weight |
|
inline |
Add the 2D shallow water wet-dry wrapped convection terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | the elemental shape function gradients in configuration space |
[in] | v | fully wet velocites |
[in] | v_wet_dry | the wet-dry velocities |
[in] | f | not used |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
Solves the following weak, wet/dry discrete body terms of the 2D shallow water equation:
|
inline |
Add the 2D shallow water wet-dry wrapped density body-force terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | not used |
[in] | v | the fully wet bathymetry gradient |
[in] | v_wet_dry | not used |
[in] | f | fully we density |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
|
inline |
Add the 2D shallow water wet-dry wrapped density boundary pressure terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | the wet-dry nodal positions |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | not used |
[in] | v | not used |
[in] | v_wet_dry | not used |
[in] | f | fully wet density |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
|
inline |
Add the 2D shallow water wet-dry wrapped density pressure terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | the elemental shape function gradients in configuration space |
[in] | v | not used |
[in] | v_wet_dry | not used |
[in] | f | fully wet density |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
|
inline |
Add the 2D shallow water wet-dry wrapped *nonconservative convection terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | the elemental shape function gradients in configuration space |
[in] | v | not used |
[in] | v_wet_dry | wet-dry velocities |
[in] | f | not used |
[in] | f_wet_dry | wet-dry velocities |
[in] | djac | the 2d element area |
[in] | dt | time-step |
[in] | c | a constant |
Solves the following weak, non-conservative wet/dry discrete equation:
|
inline |
Wet/Dry routine for calculating the wet-dry integral of shallow water H, u and v: .
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | not used |
[in] | v | not used |
[in] | v_wet_dry | wet-dry velocities (integrated here) |
[in] | f | not used |
[in] | f_wet_dry | the scalar function (integrated here) |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |
|
inline |
Add the 2D shallow water wet-dry wrapped pressure terms to the FE residual.
[in,out] | rhs | the 2D elemental residual array contribution with wet/dry inclusion |
[in] | elem_nds | not used |
[in] | elem_head | the 2d element nodal wet-dry depths |
[in] | grad_phi | the elemental shape function gradients in configuration space |
[in] | v | not used |
[in] | v_wet_dry | not used |
[in] | f | not used |
[in] | f_wet_dry | not used |
[in] | djac | the 2d element area |
[in] | dt | not used |
[in] | c | a constant |