AdH Kraken 0.0.0
Next generation Adaptive Hydraulics
Loading...
Searching...
No Matches
sw2_body_resid.c File Reference
#include "adh.h"

Functions

void fe_sw2_temporal (int ie, SELEM_2D *elem2d, int nnodes, SVECT *elem_nds, double djac, double drying_lower_limit, double *elem_head, SVECT2D *elem_vel, double wd_factor, double dt_factor, double *elem_rhs, char *string, int DEBUG, int DEBUG_LOCAL, int wd_flag)
 Returns the 2D body temporal contributions to the shallow water equations. More...
 
int fe_sw2_body_resid (SMODEL_SUPER *mod, double *elem_rhs, int ie, double perturbation, int perturb_node, int perturb_var, int perturb_sign, int DEBUG)
 Returns the 2D shallow water elemental residual. More...
 

Detailed Description

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

Function Documentation

◆ fe_sw2_body_resid()

int fe_sw2_body_resid ( SMODEL_SUPER mod,
double *  elem_rhs,
int  ie,
double  perturbation,
int  perturb_node,
int  perturb_var,
int  perturb_sign,
int  DEBUG 
)

Returns the 2D shallow water elemental residual.

Author
Charlie Berger, Ph.D.
Gaurav Savant, Ph.D.
Gary Brown
Corey Trahan, Ph.D.
Mark Loveland, Ph.D.
Bug:
none
Warning
none
Parameters
[in]mod(SMODEL_SUPER *) - the super model, contains pointer to grid, dependent vars, maps, and such
[in,out]elem_rhs(double *) - the 2D elemental residual array
[in]ie(int) - the elemental id
[in]pertubation(double) - the F-D approximation size aka Newton pertubation
[in]perturb_node(int) - the index of node to be perturbed (local to element)
[in]perturb_var(int) - the index of the variable to be perturbed
[in]perturb_sign(int) - the direction of Newton perturbation (either -1 or 1)
[in]DEBUG(int) - a debug flag
Returns
a possible error code

Solves the body integals of the following weak, discrete body terms of the 2D shallow water equation:

\begin{eqnarray*}
 \weakSWDAcont{e}{i}{h} \\
 \weakSWMxDD{e}{i}{h} \\
 \weakSWMyDD{e}{i}{h}
 \end{eqnarray*}

+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

FINITE ELEMENT INTEGRATIONS

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

SHOCK CAPTURING CONTRIBUTION

Calculates residual contributions from momentum and mass shock capturing.

Note
CJT :: wd_flag = 0 (element fully wet), 1 (element partially dry), 2 (element full dry)
CJT :: shock capturing terms only contribute for partially or fully dry elements
CJT :: formulation?

+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

ADVECTION CONTRIBUTION

Calculates the advection addition to the SW body terms.

\begin{eqnarray*}
   \residDA{i}{}{c}   &=& -  dt * \bodyConv{\,2d}{e}{\phidd{i}}{(\velb{h} \, \depth{h})} \\
   \residDA{i}{}{mx}  &=& -  dt * \bodyConv{\,2d}{e}{\phidd{i}}{(\velb{h} \, \ub{h}\depth{h})} \\
   \residDA{i}{}{my}  &=& -  dt * \bodyConv{\,2d}{e}{\phidd{i}}{(\velb{h} \, \vb{h}\depth{h})}
  \end{eqnarray*}

Note
GS :: Wrap in wet/dry integration
CJT :: Until we include quadrilaterals in wetting and drying, this is only for triangles

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

TEMPORAL CONTRIBUTION

Calculates the temporal addition to the elemental residual.

\begin{eqnarray*}
   \residDA{i}{}{c}   &=& dt * \bodyTime{\,2d}{e}{\phidd{i}}{\depth{h}} \\
   \residDA{i}{}{mx}  &=& dt * \bodyTime{\,2d}{e}{\phidd{i}}{\,\ub{h}\depth{h}} \\
   \residDA{i}{}{my}  &=& dt * \bodyTime{\,2d}{e}{\phidd{i}}{\,\vb{h}\depth{h}}
  \end{eqnarray*}

Note
GLB :: consistent mass integration of continuity so only wet portion is itegrated, momentum temporal is lumped
CJT :: consistent mass temporal integral: $ \bodyTime{\,2d}{e}{\phidd{i}}{\depth{h}_1 \, \phidd{1} + \depth{h}_2 \, \phidd{2} + \depth{h}_3 \, \phidd{3}} $
CJT :: lumped & grouped momentum temporal integral: $ \bodyTime{\,2d}{e}{\phidd{i}}{\depth{h}_i \, \ub{h}_i} $

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

DIFFUSION CONTRIBUTION

Calculates the diffusion addition to the elemental residual.

\begin{eqnarray*}
   \residDA{i}{}{mx}  &=& dt * \bodyDiffusion{\,2d}{e}{\phidd{i}}{\depth{h}\,\diffTensorDA{x}{h}{\diffTensorSW}{}}  \\
   \residDA{i}{}{my}  &=& dt * \bodyDiffusion{\,2d}{e}{\phidd{i}}{\depth{h}\,\diffTensorDA{y}{h}{\diffTensorSW}{}}
  \end{eqnarray*}

Note
CJT :: currently uses elementally averaged depth

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

FRICTION CONTRIBUTION

Calculates the friction addition to the elemental residual.

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

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

VORTICITY CONTRIBUTION

Calculates the vorticity addition to the elemental residual.

\begin{eqnarray*}
   \residDA{i}{}{mx}  &=& dt * \bodySource{\,2d}{e}{\phidd{i}}{stream_x}  \\
   \residDA{i}{}{my}  &=& dt * \bodySource{\,2d}{e}{\phidd{i}}{stream_y}
  \end{eqnarray*}

Note
CJT :: currently only for triangles, but can use elementally averaged gradients and include quads

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

CORIOLIS CONTRIBUTION

Calculates the coriolis addition to the elemental residual.

\begin{eqnarray*}
   \residDA{i}{}{mx}  &=& - dt * \bodyCoriolis{\,2d}{e}{\phidd{i}}{\depth{h}\,\vb{h}}   \\
   \residDA{i}{}{my}  &=& + dt * \bodyCoriolis{\,2d}{e}{\phidd{i}}{\depth{h}\,\ub{h}}
  \end{eqnarray*}

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

WIND AND WAVE CONTRIBUTION

Calculates the wind and wave addition to the elemental residual.

\begin{eqnarray*}
   \residDA{i}{}{mx}  &=& -dt * \bodyStress{\,2d}{e}{\phidd{i}}{\lrpb{\tau_{winds,x}^{\,h} + \tau_{waves,x}^{\,h}}}   \\
   \residDA{i}{}{my}  &=& -dt * \bodyStress{\,2d}{e}{\phidd{i}}{\lrpb{\tau_{winds,y}^{\,h} + \tau_{waves,y}^{\,h}}}
  \end{eqnarray*}

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

RAINFALL CONTRIBUTION

Calculates the rainfall addition to the elemental residual.
$ \residDA{i}{}{c} = dt * \bodySource{\,2d}{e}{\phidd{i}}{R} $

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

PRESSURE CONTRIBUTION

Calculates the pressure additions to the elemental residual.

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Calculates the body pressure on element.

\begin{eqnarray*}
  \residDA{i}{}{mx}  &=& -dt * \bodyPressureDD{\,2d}{e}{\phidd{i}}{x}{(\depth{h})}    \\
  \residDA{i}{}{my}  &=& -dt * \bodyPressureDD{\,2d}{e}{\phidd{i}}{y}{(\depth{h})}
 \end{eqnarray*}

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Calculates the bathymetry gradient integral on element.

\begin{eqnarray*}
  \residDA{i}{}{mx}  &=& -dt * \bodySource{\,2d}{e}{\phidd{i}}{ g \, \depth{h} \frac{ \partial z_\bed{}}{\partial x} }    \\
  \residDA{i}{}{my}  &=& -dt * \bodySource{\,2d}{e}{\phidd{i}}{ g \, \depth{h} \frac{ \partial z_\bed{}}{\partial y} }
 \end{eqnarray*}

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Calculates the boundary integral on element edges. This is to auto-find no-flow boundaries.

\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*}

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

BAROCLINIC PRESSURE CONTRIBUTION

Calculates the density pressure additions to the elemental residual.

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Calculates the density driven body pressure on element.

\begin{eqnarray*}
  \residDA{i}{}{mx}  &=& -dt * \bodyDensityPressureDD{\,2d}{e}{\phidd{i}}{x}{(\depth{h})}    \\
  \residDA{i}{}{my}  &=& -dt * \bodyDensityPressureDD{\,2d}{e}{\phidd{i}}{y}{(\depth{h})}
 \end{eqnarray*}

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Calculates the density bathymetry gradient integral on element.

\begin{eqnarray*}
  \residDA{i}{}{mx}  &=& -dt * \bodySource{\,2d}{e}{\phidd{i}}{ g \, \depth{h} \frac{ \partial z_\bed{}}{\partial x} \, \frac{\bigtriangleup \rho}{\rho_o}}    \\
  \residDA{i}{}{my}  &=& -dt * \bodySource{\,2d}{e}{\phidd{i}}{ g \, \depth{h} \frac{ \partial z_\bed{}}{\partial y} \, \frac{\bigtriangleup \rho}{\rho_o}}
 \end{eqnarray*}

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Calculates the density driven pressure boundary integral on element edges. This is to auto-find no-flow boundaries.

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

++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

SUPG CONTRIBUTION

Calculates the SUPG addition to the elemental residual.

\begin{eqnarray*}
  \weakSwCsupgDD{2d}{e}{i} {\ub{h}} {\depth{h}} {\vb{h}} {\depth{h}} {c} \\
  \weakSwMXsupgDD{2d}{e}{i}{g\,\depth{h}}{\ub{h} \depth{h}}{}{\vb{h} \depth{h}}{mx} \\
  \weakSwMYsupgDD{2d}{e}{i}{g\,\depth{h}}{\vb{h} \depth{h}}{}{\ub{h} \depth{h}}{my}
  \end{eqnarray*}

Note
CJT :: For tau, use method flag == 2, the elemental area to get elemental length estimate
CJT :: Not sure what g_factor is in the tau calculation and why area is X factor_old

+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

SET DIRICHLET BCS

Sets the Dirichlet boundary conditions.

◆ fe_sw2_temporal()

void fe_sw2_temporal ( int  ie,
SELEM_2D elem2d,
int  nnodes,
SVECT elem_nds,
double  djac,
double  drying_lower_limit,
double *  elem_head,
SVECT2D elem_vel,
double  wd_factor,
double  dt_factor,
double *  elem_rhs,
char *  string,
int  DEBUG,
int  DEBUG_LOCAL,
int  wd_flag 
)

Returns the 2D body temporal contributions to the shallow water equations.

Author
Charlie Berger, Ph.D.
Gaurav Savant, Ph.D.
Gary Brown
Corey Trahan, Ph.D.
Bug:
none
Warning
none
Note
CJT:: this calculation groups hu,v in the momentum temporal terms