AdH Kraken 0.0.0
Next generation Adaptive Hydraulics
Loading...
Searching...
No Matches
fe_sw2_boundary_integrals.c File Reference
#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...
 

Detailed Description

This file collections functions responsible for the 2D shallow water body contributions to the elemental residual.

Function Documentation

◆ fe_sw2_1d_density_pressure()

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.

Author
Charlie Berger, Ph.D.
Gaurav Savant, Ph.D.
Gary Brown, Ph.D.
Corey Trahan, Ph.D.
Bug:
none
Warning
none
Parameters
[in,out]elem_rhsthe 1D elemental residual array
[in]iethe 1D element ID
[in]nrmlthe 1D element normal
[in]djacthe 1D element jacobian
[in]dtthe current model time-step
[in]hthe water depth
[in]dwater density
[in]ggravity

\begin{eqnarray*}
                \residDA{i}{}{mx}  &=& dt * \bcPressureDD{}{e}{\phidd{i}}{(\rho^h \, \depth{h})}{x} \\
                \residDA{i}{}{my}  &=& dt * \bcPressureDD{}{e}{\phidd{i}}{(\rho^h \, \depth{h})}{y}
\end{eqnarray*}

◆ fe_sw2_1d_explicit_flow()

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.

Author
Gaurav Savant, Ph.D.
Corey Trahan, Ph.D.
Bug:
none
Warning
none
Parameters
[in,out]elem_rhsthe 1D elemental residual array
[in]iethe 1D element ID
[in]nrmlthe 1D element normal
[in]djacthe 1D element jacobian
[in]dtthe current model time-step
[in]fluxa user-supplied flux
[in]h_avgthe elementally averaged depth
[in]ux-velocity
[in]vy-velocity

\begin{eqnarray*}  \residDA{i}{}{c}   &=& dt * \intbe{\phidd{i} \, q_n}{e}    \end{eqnarray*}

For the momentum contributions:
– if the flux is into the domain, the normal flux is split into x,y velocities, so that

\begin{eqnarray*}
    u^* = (q_n / \overline{\depth{h}}) \, n_x \\
    v^* = (q_n / \overline{\depth{h}}) \, n_y
 \end{eqnarray*}

and added to the residual as

\begin{eqnarray*}
  \residDA{i}{}{mx}  &=& dt * \bcConv{}{e}{\phidd{i}}{(\velb{h} \, \ub{h} \, \depth{h})} = dt * \intbe{\phidd{i} \, q_n \, u^*}{e} \\
  \residDA{i}{}{my}  &=& dt * \bcConv{}{e}{\phidd{i}}{(\velb{h} \, \vb{h} \, \depth{h})} = dt * \intbe{\phidd{i} \, q_n \, v^*}{e}
 \end{eqnarray*}

– if the flux is out of the domain, an implicity u,v is used so that

\begin{eqnarray*}
  \residDA{i}{}{mx}  &=& dt * \bcConv{}{e}{\phidd{i}}{(\velb{h} \, \ub{h} \, \depth{h})} = dt * \intbe{\phidd{i} \, q_n \, u}{e} \\
  \residDA{i}{}{my}  &=& dt 8 \bcConv{}{e}{\phidd{i}}{(\velb{h} \, \vb{h} \, \depth{h})} = dt * \intbe{\phidd{i} \, q_n \, v}{e}
 \end{eqnarray*}

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

◆ fe_sw2_1d_hvel()

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.

Author
Charlie Berger, Ph.D.
Gaurav Savant, Ph.D.
Corey Trahan, Ph.D.
Bug:
none
Warning
none
Parameters
[in,out]elem_rhsthe 1D elemental residual array
[in]iethe 1D element ID
[in]nrmlthe 1D element normal
[in]djacthe 1D element jacobian
[in]dtthe current model time-step
[in]vthe water velocity
[in]user_velocitythe water velocity prescribed by user
[in]hthe water depth

\begin{eqnarray*}  \residDA{i}{}{c}   &=& dt * \intbe{\phidd{i} \, (\velb{h} \cdot \nrml) \, \depth{h}}{e} \\
                \residDA{i}{}{mx}  &=& dt * \intbe{\phidd{i} \, q_n \, \ub{}}{e}\\
                \residDA{i}{}{my}  &=& dt * \intbe{\phidd{i} \, q_n \, \vb{}}{e}
\end{eqnarray*}

◆ fe_sw2_1d_implicit_flow()

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.

Author
Charlie Berger, Ph.D.
Gaurav Savant, Ph.D.
Gary Brown, Ph.D.
Corey Trahan, Ph.D.
Bug:
none
Warning
none
Parameters
[in,out]elem_rhsthe 1D elemental residual array
[in]iethe 1D element ID
[in]nrmlthe 1D element normal
[in]djacthe 1D element jacobian
[in]dtthe current model time-step
[in]vthe hydrodynamic velocity
[in]hthe water depth

\begin{eqnarray*}  \residDA{i}{}{c}   &=& dt * \intbe{\phidd{i} \, q_n}{e} \\
                \residDA{i}{}{mx}  &=& dt * \bcConv{}{e}{\phidd{i}}{(\velb{h} \, \ub{h} \, \depth{h})} \\
                \residDA{i}{}{my}  &=& dt * \bcConv{}{e}{\phidd{i}}{(\velb{h} \, \vb{h} \, \depth{h})}
\end{eqnarray*}

Note
CJT:: Gaurav *thinks the "r" variable is a prox for v so he can eliminate eddies forming on the boundary

◆ fe_sw2_1d_pressure()

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.

Author
Charlie Berger, Ph.D.
Gaurav Savant, Ph.D.
Corey Trahan, Ph.D.
Bug:
none
Warning
none
Parameters
[in,out]elem_rhsthe 1D elemental residual array
[in]iethe 1D element ID
[in]nrmlthe 1D element normal
[in]djacthe 1D element jacobian
[in]dtthe current model time-step
[in]hthe water depth
[in]ggravity

\begin{eqnarray*}
                \residDA{i}{}{mx}  &=& dt * \bcPressureDD{}{e}{\phidd{i}}{(\depth{h})}{x} \\
                \residDA{i}{}{my}  &=& dt * \bcPressureDD{}{e}{\phidd{i}}{(\depth{h})}{y}
\end{eqnarray*}

◆ fe_sw2_1d_wall_friction()

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.

Author
Charlie Berger, Ph.D.
Gaurav Savant, Ph.D.
Gary Brown, Ph.D.
Corey Trahan, Ph.D.
Bug:
none
Warning
none
Parameters
[in,out]elem_rhsthe 1D elemental residual array
[in]iethe 1D element ID
[in]nrmlthe 1D element normal
[in]djacthe 1D element jacobian
[in]dtthe current model time-step
[in]ux-velocity
[in]vy-velocity
[in]resistanceresistance

Integrates the weak, discrete outflow boundary terms:

\begin{eqnarray*}
  \residDA{i}{}{mx}  &=& dt * \bcFriction{}{e}{\phidd{i}}{\ub{h}}{\velb{h}} \\
  \residDA{i}{}{my}  &=& dt * \bcFriction{}{e}{\phidd{i}}{\vb{h}}{\velb{h}}
 \end{eqnarray*}