Plugins
|
Finite Volume Element Discretization for the incompressible Navier-Stokes Equation. More...
#include <navier_stokes_fv1.h>
Public Member Functions | |
template<typename TElem , typename TFVGeom > | |
void | add_def_A_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
adds the stiffness part to the local defect | |
template<typename TElem , typename TFVGeom > | |
void | add_def_M_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
adds the mass part to the local defect | |
template<typename TElem , typename TFVGeom > | |
void | add_jac_A_elem (LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
adds the stiffness part to the local jacobian | |
template<typename TElem , typename TFVGeom > | |
void | add_jac_M_elem (LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
adds the mass part to the local jacobian | |
template<typename TElem , typename TFVGeom > | |
void | add_rhs_elem (LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
adds the source part to the local defect | |
SmartPtr< CplUserData< number, dim > > | density () |
returns density | |
virtual std::string | disc_type () const |
returns string identifying disc type | |
template<typename TElem , typename TFVGeom > | |
void | ex_nodal_velocity (MathVector< dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< TFVGeom::dim > vLocIP[], const size_t nip, bool bDeriv, std::vector< std::vector< MathVector< dim > > > vvvDeriv[]) |
export value of the velocity | |
template<typename TElem , typename TFVGeom > | |
void | ex_velocity_grad (MathMatrix< dim, dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< TFVGeom::dim > vLocIP[], const size_t nip, bool bDeriv, std::vector< std::vector< MathMatrix< dim, dim > > > vvvDeriv[]) |
computes the value of the gradient of the velocity | |
template<typename TElem , typename TFVGeom > | |
void | fsh_elem_loop () |
finishes the element loop | |
SmartPtr< CplUserData< number, dim > > | kinematic_viscosity () |
returns kinematic viscosity | |
template<typename TFVGeom > | |
number | peclet_blend (MathVector< dim > &UpwindVel, const TFVGeom &geo, size_t ip, const MathVector< dim > &StdVel, number kinVisco) |
computes the pecled blended Upwind veloctity | |
template<typename TElem , typename TFVGeom > | |
void | prep_elem (const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[]) |
prepares the element for evaluation | |
template<typename TElem , typename TFVGeom > | |
void | prep_elem_loop (const ReferenceObjectID roid, const int si) |
prepares the element loop | |
virtual void | prepare_setting (const std::vector< LFEID > &vLfeID, bool bNonRegularGrid) |
type of trial space for each function used | |
void | set_density (SmartPtr< CplUserData< number, dim > > user) |
sets the density | |
void | set_kinematic_viscosity (SmartPtr< CplUserData< number, dim > > user) |
sets the kinematic viscosity | |
void | set_pac_upwind (bool bPac) |
void | set_source (SmartPtr< CplUserData< MathVector< dim >, dim > > user) |
sets the source function | |
void | set_stabilization (const std::string &name) |
sets stabilization based on string identifier | |
void | set_stabilization (const std::string &name, const std::string &diffLength) |
sets stabilization and diff length method based on string identifier (only for Schneider-Raw stabilizations) | |
void | set_stabilization (SmartPtr< INavierStokesFV1Stabilization< dim > > spStab) |
sets the stabilization used to compute the stabilized velocities | |
void | set_upwind (const std::string &name) |
sets the upwind based on a string identifier | |
void | set_upwind (SmartPtr< INavierStokesFV1Stabilization< dim > > spStab) |
sets a stabilization for upwinding (Physical Advection Correction) | |
void | set_upwind (SmartPtr< INavierStokesUpwind< dim > > spUpwind) |
sets an upwinding for the convective term of momentum equation | |
DataImport< MathVector< dim >, dim > | sourceSCVF () |
returns scvf source | |
SmartPtr< INavierStokesFV1Stabilization< dim > > | stabilization () |
returns stabilization | |
NavierStokesFV1 (const char *functions, const char *subsets) | |
NavierStokesFV1 (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset) | |
![]() | |
IncompressibleNavierStokesBase (const char *functions, const char *subsets) | |
IncompressibleNavierStokesBase (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset) | |
void | set_density (number val) |
sets the kinematic viscosity | |
void | set_stokes (bool Stokes) |
switches the convective terms off (to solve the Stokes equation) | |
bool | stokes () |
flag if using Peclet Blending | |
void | set_laplace (bool bLaplace) |
sets assembling of diffusive term to laplace | |
bool | laplace () |
flag if using Peclet Blending | |
void | set_peclet_blend (bool pecletBlend) |
sets if peclet blending is used in momentum equation | |
void | set_grad_div (number factor) |
flag if using Peclet Blending | |
virtual bool | requests_local_time_series () |
returns if local time series is needed | |
SmartPtr< CplUserData< MathVector< dim >, dim > > | velocity () |
returns the export of the velocity | |
SmartPtr< CplUserData< MathMatrix< dim, dim >, dim > > | velocity_grad () |
returns the export of the velocity gradient | |
![]() | |
void | set_exact_jacobian (bool bExactJacobian) |
sets if the exact jacobian is computed (fixpoint approximation else) | |
void | set_exact_jacobian (number fullNewtonFactor) |
NavierStokesBase (const char *functions, const char *subsets) | |
NavierStokesBase (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset) | |
void | set_kinematic_viscosity (number val) |
sets the kinematic viscosity | |
void | set_source (const std::vector< number > &vSource) |
sets the source function | |
![]() | |
void | add_elem_modifier (SmartPtr< IElemDiscModifier< TDomain > > elemModifier) |
std::vector< SmartPtr< IElemDiscModifier< TDomain > > > & | get_elem_modifier () |
IElemDisc (const char *functions, const char *subsets) | |
IElemDisc (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset) | |
![]() | |
virtual void | add_def_A_expl_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
TLeaf & | asLeaf () |
void | check_roid (ReferenceObjectID roid, int discType) |
void | do_add_def_A_elem (LocalVector &d, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_add_def_A_expl_elem (LocalVector &d, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_add_def_M_elem (LocalVector &d, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_add_jac_A_elem (LocalMatrix &J, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_add_jac_M_elem (LocalMatrix &J, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_add_rhs_elem (LocalVector &rhs, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_fsh_elem_loop () |
void | do_fsh_timestep (const number time, VectorProxyBase *u, size_t algebra_id) |
void | do_fsh_timestep_elem (const number time, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_prep_elem (LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[]) |
void | do_prep_elem_loop (const ReferenceObjectID roid, const int si) |
void | do_prep_timestep (number future_time, const number time, VectorProxyBase *u, size_t algebra_id) |
void | do_prep_timestep_elem (const number time, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
virtual void | fsh_timestep (number time, VectorProxyBase *u) |
virtual void | fsh_timestep_elem (const number time, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
IElemAssembleFuncs () | |
virtual void | post_assemble_loop () |
virtual void | prep_assemble_loop () |
virtual void | prep_timestep (number future_time, number time, VectorProxyBase *u) |
virtual void | prep_timestep_elem (const number time, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | set_roid (ReferenceObjectID id, int discType) |
virtual | ~IElemAssembleFuncs () |
![]() | |
IElemError (const char *functions, const char *subsets) | |
IElemError (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset) | |
![]() | |
SmartPtr< ApproximationSpace< TDomain > > | approx_space () |
ConstSmartPtr< ApproximationSpace< TDomain > > | approx_space () const |
void | check_setup (bool bNonRegularGrid) |
void | clear_imports () |
SmartPtr< TDomain > | domain () |
ConstSmartPtr< TDomain > | domain () const |
const FunctionGroup & | function_group () const |
ConstSmartPtr< FunctionPattern > | function_pattern () const |
IDataImport< dim > & | get_import (size_t i) |
IElemDiscBase (const char *functions="", const char *subsets="") | |
IElemDiscBase (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset) | |
bool | is_time_dependent () const |
bool | local_time_series_needed () |
const LocalVectorTimeSeries * | local_time_solutions () const |
const FunctionIndexMapping & | map () const |
number | mass_scale () const |
number | mass_scale () const |
number | mass_scale (const size_t timePoint) const |
number | mass_scale (const size_t timePoint) const |
const std::vector< number > & | mass_scales () const |
const std::vector< number > & | mass_scales () const |
size_t | num_fct () const |
size_t | num_imports () const |
size_t | num_subsets () const |
void | register_import (IDataImport< dim > &Imp) |
void | set_approximation_space (SmartPtr< ApproximationSpace< TDomain > > approxSpace) |
void | set_functions (const std::string &functions) |
void | set_functions (const std::vector< std::string > &functions) |
void | set_stationary () |
void | set_stationary (bool bStationaryForced=true) |
void | set_subsets (const std::string &subsets) |
void | set_subsets (const std::vector< std::string > &subsets) |
void | set_time_dependent (LocalVectorTimeSeries &locTimeSeries, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff) |
void | set_time_independent () |
void | set_time_point (const size_t timePoint) |
number | stiff_scale () const |
number | stiff_scale () const |
number | stiff_scale (const size_t timePoint) const |
number | stiff_scale (const size_t timePoint) const |
const std::vector< number > & | stiff_scales () const |
const std::vector< number > & | stiff_scales () const |
TDomain::subset_handler_type & | subset_handler () |
const TDomain::subset_handler_type & | subset_handler () const |
const std::vector< std::string > & | symb_fcts () const |
const std::vector< std::string > & | symb_subsets () const |
number | time () const |
size_t | time_point () const |
virtual int | type () const |
virtual bool | use_hanging () const |
virtual | ~IElemDiscBase () |
![]() | |
TLeaf & | asLeaf () |
void | check_roid (ReferenceObjectID roid, int discType) |
virtual void | compute_err_est_A_elem (const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale) |
virtual void | compute_err_est_M_elem (const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale) |
virtual void | compute_err_est_rhs_elem (GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale) |
void | do_compute_err_est_A_elem (LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale) |
void | do_compute_err_est_M_elem (LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale) |
void | do_compute_err_est_rhs_elem (GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale) |
void | do_fsh_err_est_elem_loop () |
void | do_prep_err_est_elem (LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
void | do_prep_err_est_elem_loop (const ReferenceObjectID roid, const int si) |
virtual SmartPtr< IErrEstData< TDomain > > | err_est_data () |
bool | err_est_enabled () const |
virtual void | fsh_err_est_elem_loop () |
IElemEstimatorFuncs () | |
virtual void | prep_err_est_elem (const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[]) |
virtual void | prep_err_est_elem_loop (const ReferenceObjectID roid, const int si) |
void | set_error_estimator (SmartPtr< IErrEstData< TDomain > > ee) |
void | set_roid (ReferenceObjectID id, int discType) |
virtual | ~IElemEstimatorFuncs () |
Static Public Attributes | |
static const int | dim = base_type::dim |
World dimension. | |
![]() | |
static const int | dim = base_type::dim |
World dimension. | |
![]() | |
static const int | dim = base_type::dim |
World dimension. | |
![]() | |
static const int | dim |
![]() | |
static const int | dim |
![]() | |
static const int | dim |
![]() | |
static const int | dim |
![]() | |
static const int | dim |
Protected Types | |
typedef IncompressibleNavierStokesBase< TDomain > | base_type |
Base class type. | |
typedef NavierStokesFV1< TDomain > | this_type |
own type | |
![]() | |
typedef NavierStokesBase< TDomain > | base_type |
Base class type. | |
typedef IncompressibleNavierStokesBase< TDomain > | this_type |
own type | |
![]() | |
typedef IElemDisc< TDomain > | base_type |
Base class type. | |
typedef NavierStokesBase< TDomain > | this_type |
own type | |
Static Protected Attributes | |
static const size_t | _P_ = dim |
abbreviation for pressure | |
Additional Inherited Members | |
![]() | |
typedef IElemAssembleFuncs< IElemDisc< TDomain >, TDomain > | assemble_base_type |
typedef IElemError< TDomain > | base_type |
typedef TDomain | domain_type |
typedef IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain > | estimator_base_type |
![]() | |
typedef TDomain | domain_type |
typedef TLeaf | leaf_type |
![]() | |
typedef TDomain | domain_type |
typedef IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain > | estimator_base_type |
![]() | |
typedef TDomain | domain_type |
typedef TDomain::position_type | position_type |
![]() | |
typedef TDomain | domain_type |
typedef TLeaf | leaf_type |
Finite Volume Element Discretization for the incompressible Navier-Stokes Equation.
This class implements the IElemDisc interface to provide element local assemblings for the incompressible Navier-Stokes equation.
The unknowns of the equation are:
The equation takes the form
\begin{align*} \frac{\partial \rho \vec{u}}{\partial t} - \nabla \left( \rho \nu (\nabla \vec{u} + (\nabla \vec{u})^T) \right) + \nabla \cdot \left( \rho \vec{u} \vec{u}^T \right) + \nabla p &= \vec{f}\\ \nabla \cdot (\rho \vec{u}) &= 0 \end{align*}
with
The first equation models the conservation of momentum and is therefore referred to as the Momentum equation. The second equation models the conservation of mass and is known as the Mass continuity equation or simply Continuity equation.
In component-wise notation, the equation reads
\begin{align*} \frac{\partial \rho u_{d_1}}{\partial t} - \sum_{d_2 = 1}^{\text{dim}} \frac{\partial}{\partial x_{d_2}} \left( \rho \nu \left( \frac{\partial u_{d_1}}{\partial x_{d_2}} + \frac{\partial u_{d_2}}{\partial x_{d_1}} \right)\right) + \sum_{d_2 = 1}^{\text{dim}} \frac{\partial}{\partial x_{d_2}} \left( \rho u_{d_2} u_{d_1} \right) + \frac{\partial p}{\partial x_{d_1}} &= f_{d_1} \qquad \forall \quad d_1 = 1, \dots, \text{dim}\\ \sum_{d = 1}^{\text{dim}} \frac{\partial \rho u_d}{\partial x_{d}} &= 0 \end{align*}
Since the density is assumed to be constant, it set to \(\rho \equiv 1 \).
The finite volume element discretization uses a dual mesh consisting of Control Volumes \(\mathcal{B}_h\), that cover the domain. For each control volume \(B \in \mathcal{B}_h\) the following equation is solved
\begin{align*} \frac{\partial}{\partial t} \int_{B} \begin{pmatrix} \vec{u} \\ 0 \end{pmatrix} dV + \int_{\partial B} \begin{pmatrix} - \rho \nu \left(\nabla \vec{u} + (\nabla \vec{u})^T \right) \vec{n} + \vec{u} \vec{u}^T \vec{n} + p \vec{n} \\ \vec{u} \vec{n} \end{pmatrix} dS = \int_B \begin{pmatrix} \vec{f} \\ 0 \end{pmatrix} dV \end{align*}
where \( \vec{n}\) is the outer normal on the boundary of the control volume and \(\int_B \cdot \; dV \) and \(\int_{\partial B} \cdot \; dS \) are the integration over the control volume and the integration over the boundary of the control volume. By $f B_{sh} $f is denoted the control volume associated to a shape function (i.e. vertices, since we use P1 Lagrange ansatz functions), T usually denotes the Element.
In order to number the local unknowns, we use the following notation.
The access to local unknowns of local vectors and matrices is now given by
\begin{align*} &\mathbf{d}(d_1, sh), \mathbf{d}(\text{\_P\_}, sh),\\ &\mathcal{J}(d_1, sh_1, d_2, sh_2), \mathcal{J}(d_1, sh_1, \text{\_P\_}, sh_2), \end{align*}
TDomain | Domain |
TAlgebra | Algebra |
|
protected |
Base class type.
|
protected |
own type
ug::NavierStokes::NavierStokesFV1< TDomain >::NavierStokesFV1 | ( | const char * | functions, |
const char * | subsets | ||
) |
Constructor (setting default values)
References ug::NavierStokes::NavierStokesFV1< TDomain >::init().
ug::NavierStokes::NavierStokesFV1< TDomain >::NavierStokesFV1 | ( | const std::vector< std::string > & | vFct, |
const std::vector< std::string > & | vSubset | ||
) |
Constructor (setting default values)
References ug::NavierStokes::NavierStokesFV1< TDomain >::init().
|
virtual |
adds the stiffness part to the local defect
This function adds the contribution of one element of the stiffness part to the local defect.
The local solution \( \mathbf{u} \) passes the p1 conform lagrange dofs of the element, i.e. the unknowns in the corners of the element.
Define the functional matrix of the velocity at some integration point by
\begin{align*} \nabla \vec{u}|_{ip} := \left(\left.\frac{\partial u_i}{\partial x_j} \right|_{ip}\right)_{i,j = 1,\dots,dim} \end{align*}
and the diffusive flux to
\begin{align*} \vec{f}^{\text{diff}}|_{ip} := - \nu \left( \nabla \vec{u}|_{ip} + (\nabla \vec{u}|_{ip})^T) \right) \cdot \vec{n}_{ip} \end{align*}
Suppose that a procedure producing some upwind velocity \( \vec{u}^{up}\) is available for each ip. Then, we define the convective flux to
\begin{align*} \vec{f}^{\text{conv}}|_{ip} := \left( \vec{u}^{up}|_{ip} \vec{u}^{up}|_{ip}^T \right) \cdot \vec{n}_{ip} \end{align*}
Define the pressure at the ip to be interpolated as
\begin{align*} p|_{ip} := \sum_{sh=0}^{n_{sh}-1} \mathbf{u}(\text{\_P\_}, sh) \phi_{sh} \end{align*}
The contribution added for the momentum equation is:
\begin{align*} \forall d=0,\dots, \text{dim}-1, \quad \forall sh=0,\dots, n_{sh}-1: \qquad \mathbf{d}(d, sh) &\text{ +=} \int_{\partial B_{sh} \cap T} \vec{f}^{\text{diff}}_d + \vec{f}^{\text{conv}}_d + p \vec{n}_d \;dS \approx \sum_{i=1}^{n_{\partial B_{sh}}} |\partial B^i_{sh} \cap T| \left( \left.\vec{f}^{\text{diff}}_d\right|_{ip} + \left.\vec{f}^{\text{conv}}_d\right|_{ip} + p|_{ip} \vec{n}|_{ip} \right)\\ \end{align*}
Assume, that some stabilized velocity \( \vec{u}^{stab}\) is available for each ip. Then, the contribution added for the continuity equation is:
\begin{align*} \forall sh=0,\dots, n_{sh}-1: \qquad \mathbf{d}(\text{\_P\_}, sh) &\text{ +=} \int_{\partial B_{sh} \cap T} \vec{u}^{stab} \cdot \vec{n} \;dS \approx \sum_{i=1}^{n_{\partial B_{sh}}} |\partial B^i_{sh} \cap T| \; \vec{u}^{stab}|_{ip} \cdot \vec{n}|_{ip}\\ \end{align*}
[in,out] | d | local defect with values added on output |
[in] | u | local vector of current solution |
[in] | time | current time step |
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
References ug::GeomProvider< typename TGeom >::get(), ug::MatVecMult(), ug::LocalVectorTimeSeries::size(), ug::LocalVectorTimeSeries::solution(), ug::NavierStokes::INavierStokesFV1Stabilization< dim >::stab_vel(), ug::LocalVectorTimeSeries::time(), ug::TransposedMatVecMultAdd(), UG_ASSERT, UG_THROW, ug::NavierStokes::INavierStokesUpwind< dim >::upwind_vel(), ug::VecProd(), ug::VecScale(), and ug::VecSet().
|
virtual |
adds the mass part to the local defect
This function adds the contribution of one element of the mass part to the local defect.
The local solution \( \mathbf{u} \) passes the p1 conform lagrange dofs of the element, i.e. the unknowns in the corners of the element.
The contribution added is:
\begin{align*} \forall d=0,\dots, \text{dim}-1, \quad \forall sh=0,\dots, n_{sh}-1: \qquad \mathbf{d}(d, sh) &\text{ +=} \int_{B_{sh} \cap T} \vec{u}_d dV \approx |B_{sh} \cap T| \mathbf{u}(d, sh)\\ \mathbf{d}(\text{\_P\_, sh}) &\text{ +=} 0 \end{align*}
[in,out] | d | local defect with values added on output |
[in] | u | local vector of current solution |
[in] | time | current time step |
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
References ug::GeomProvider< typename TGeom >::get(), and UG_ASSERT.
|
virtual |
adds the stiffness part to the local jacobian
This function adds the local contributions of the stiffness part to the local jacobian.
For the definition of \( \vec{f}^{\text{diff}}|_{ip}\) and \( \vec{f}^{\text{conv}}|_{ip}\) see add_def_A_elem.
The derivative of the diffusive flux is given by
\begin{align*} \frac{\partial \vec{f}^{\text{diff}}_{d_1}|_{ip}}{\partial \vec{u}_{d_2}^{sh}} &= - \nu \left( \delta_{d_1, d_2} \sum_{k=1}^{\text{dim}} \left. \frac{\partial \phi^{sh}}{\partial x_k}\right|_{ip} \vec{n_k}|_{ip} + \left. \frac{\partial \phi^{sh}}{\partial x_{d_1}}\right|_{ip} \vec{n_{d_2}}|_{ip} \right)\\ \frac{\partial \vec{f}^{\text{diff}}_{d_1}|_{ip}}{\partial p^{sh}} &= 0 \end{align*}
For the derivative of the convective term, we can use a fixpoint linearization if we only differentiate the first factor of \( \vec{f}^{\text{conv}}_{d_1} = (\vec{u}^{\text{blend}} (\vec{u}^{\text{blend}})^T)_{d_1 k} \vec{n}_k = \vec{u}^{\text{blend}}_{d_1} \sum_{k=1}^{\text{dim}} \vec{u}^{\text{blend}}_k \vec{n}_k = \vec{u}^{\text{blend}}_{d_1} \langle \vec{u}^{\text{blend}}, \vec{n} \rangle \)
This gives
\begin{align*} \frac{\partial \vec{f}^{\text{conv}}_{d_1}}{\partial \vec{u}_{d_2}^{sh}} &= \langle \vec{u}^{\text{blend}}, \vec{n} \rangle \frac{\partial \vec{u}^{\text{blend}}_{d_1}}{\partial \vec{u}_{d_2}^{sh}} &= \langle \vec{u}^{\text{blend}}, \vec{n} \rangle \left( (1-\omega) \delta_{d_1 d_2} \phi^{sh} + \omega \frac{\partial \vec{u}^{\text{up}}_{d_1}}{\partial \vec{u}_{d_2}^{sh}} \right) \end{align*}
The derivative of the pressure term is given by
\begin{align*} \frac{\partial p \vec{n}|_{ip}}{\partial \vec{u}_{d_2}^{sh}} &= 0\\ \frac{\partial p \vec{n}|_{ip}}{\partial p^{sh}} &= \phi^{sh}|_{ip} \vec{n} \end{align*}
The derivative of the continuity term is given by
\begin{align*} \frac{\partial \vec{u}^{\text{stab}}|_{ip} \vec{n}|_{ip}}{\partial \vec{u}_{d_1}^{sh}} &= \sum_{d_2=1}^{dim} \left. \frac{\partial u_{d_2}^{\text{stab}}} {\partial \vec{u}_{d_1}^{sh}} \right|_{ip} \vec{n}_{d_2}\\ \frac{\partial \vec{u}^{\text{stab}}|_{ip} \vec{n}|_{ip}}{\partial p^{sh}} &=\sum_{d_2=1}^{dim} \left. \frac{\partial u_{d_2}^{\text{stab}}} {\partial p^{sh}} \right|_{ip} \vec{n}_{d_2}\\ \end{align*}
[in,out] | J | local jacobian with values added on output |
[in] | u | local vector of current solution |
[in] | time | current time step |
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
References ug::GeomProvider< typename TGeom >::get(), ug::LocalVectorTimeSeries::size(), ug::LocalVectorTimeSeries::solution(), ug::NavierStokes::INavierStokesFV1Stabilization< dim >::stab_shape_p(), ug::NavierStokes::INavierStokesFV1Stabilization< dim >::stab_shape_vel(), ug::NavierStokes::INavierStokesFV1Stabilization< dim >::stab_vel(), ug::LocalVectorTimeSeries::time(), UG_ASSERT, UG_THROW, ug::NavierStokes::INavierStokesUpwind< dim >::upwind_shape_ip(), ug::NavierStokes::INavierStokesUpwind< dim >::upwind_shape_sh(), ug::NavierStokes::INavierStokesUpwind< dim >::upwind_vel(), ug::VecDot(), ug::VecProd(), ug::VecSet(), and ug::NavierStokes::INavierStokesFV1Stabilization< dim >::vel_comp_connected().
|
virtual |
adds the mass part to the local jacobian
This function adds the contribution of one element of the mass part to the local jacobian.
The local solution \( \mathbf{u} \) passes the p1 conform lagrange dofs of the element, i.e. the unknowns in the corners of the element.
The contribution added is:
\begin{align*} \forall d=0,\dots, \text{dim}-1, \quad \forall sh=0,\dots, n_{sh}-1: \qquad \mathcal{J}(d, sh, d, sh) &\text{ +=} \int_{B_{sh} \cap T} dV = |B_{sh} \cap T|\\ \end{align*}
[in,out] | J | local jacobian with values added on output |
[in] | u | local vector of current solution |
[in] | time | current time step |
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
References ug::GeomProvider< typename TGeom >::get(), and UG_ASSERT.
|
virtual |
adds the source part to the local defect
This function adds the contribution of one element of the source part to the local defect, if a source function \( \vec{f} \) has been specified. Otherwise a zero source function is assumed.
The contribution added is:
\begin{align*} \forall d=0,\dots, \text{dim}-1, \quad \forall sh=0,\dots, n_{sh}-1: \qquad \mathbf{d}(d, sh) &\text{ +=} \int_{B_{sh} \cap T} \vec{f}_d \; dV \approx |B_{sh} \cap T| \; \vec{f}_d|_{ip} \end{align*}
[in,out] | d | local defect with values added on output |
[in] | time | current time step |
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
References ug::GeomProvider< typename TGeom >::get(), and UG_ASSERT.
|
inlinevirtual |
returns density
Implements ug::NavierStokes::IncompressibleNavierStokesBase< TDomain >.
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_imDensitySCVF, and ug::DataImport< typename TData, int dim >::user_data().
|
inlinevirtual |
returns string identifying disc type
Implements ug::NavierStokes::IncompressibleNavierStokesBase< TDomain >.
void ug::NavierStokes::NavierStokesFV1< TDomain >::ex_nodal_velocity | ( | MathVector< dim > | vValue[], |
const MathVector< dim > | vGlobIP[], | ||
number | time, | ||
int | si, | ||
const LocalVector & | u, | ||
GridObject * | elem, | ||
const MathVector< dim > | vCornerCoords[], | ||
const MathVector< TFVGeom::dim > | vLocIP[], | ||
const size_t | nip, | ||
bool | bDeriv, | ||
std::vector< std::vector< MathVector< dim > > > | vvvDeriv[] | ||
) |
export value of the velocity
References ug::GeomProvider< typename TGeom >::get(), ug::BaseLSFS< typename TImpl, int TDim, typename TShape, typename TGrad >::shapes(), and ug::VecSet().
void ug::NavierStokes::NavierStokesFV1< TDomain >::ex_velocity_grad | ( | MathMatrix< dim, dim > | vValue[], |
const MathVector< dim > | vGlobIP[], | ||
number | time, | ||
int | si, | ||
const LocalVector & | u, | ||
GridObject * | elem, | ||
const MathVector< dim > | vCornerCoords[], | ||
const MathVector< TFVGeom::dim > | vLocIP[], | ||
const size_t | nip, | ||
bool | bDeriv, | ||
std::vector< std::vector< MathMatrix< dim, dim > > > | vvvDeriv[] | ||
) |
computes the value of the gradient of the velocity
References ug::GeomProvider< typename TGeom >::get(), ug::BaseLSFS< typename TImpl, int TDim, typename TShape, typename TGrad >::grads(), and ug::MatSet().
|
virtual |
finishes the element loop
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
|
protectedvirtual |
References UG_THROW.
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::NavierStokesFV1(), and ug::NavierStokes::NavierStokesFV1< TDomain >::NavierStokesFV1().
|
inlinevirtual |
returns kinematic viscosity
Implements ug::NavierStokes::IncompressibleNavierStokesBase< TDomain >.
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_imKinViscosity, and ug::DataImport< typename TData, int dim >::user_data().
|
inline |
computes the pecled blended Upwind veloctity
This function compute the peclet blended velocity. First, the standard interpolated (no upwind) velocity is computed at the integration point of the scvf:
\begin{align*} \vec{u}^{\text{std}}_d|_{ip} := \sum_{sh=0}^{n_{sh}-1} \mathbf{u}(d, sh) \; \phi_{sh} \end{align*}
Defining the Peclet Number
\begin{align*} Pe := \frac{\vec{u}^{\text{std}}|_{ip} \cdot \vec{n}|_{ip} \; L}{\nu} \end{align*}
with the Viscosity \( \nu \) and the diffusion length \( L \), a weighting factor is computed:
\begin{align*} \omega := \frac{Pe^2}{5 + Pe^2} \end{align*}
and the blended velocity is computed to
\begin{align*} \vec{u}^{\text{blend}}|_{ip} := \omega \vec{u}^{\text{up}} + (1-\omega) \vec{u}^{\text{std}} \end{align*}
[in,out] | UpwindVel | upwind vel on entry, blended vel on exit |
[in] | scvf | SubControlVolumeFace of the Velocity |
[in] | StdVel | standard interpolation vel |
[in] | kinVisco | kinematic Viscosity at scvf |
References ug::VecDistance(), ug::VecProd(), ug::VecScaleAdd(), and ug::VecTwoNormSq().
|
virtual |
prepares the element for evaluation
This function prepare a specific element for assembling. This function is always called before any of the assemble_... functions is performed. Here, the Finite Volume Geometry is prepared and global positions of the evaluation points of the DataImports are requested (and local positions in case of hanging node assembling)
[in] | elem | element to prepare the data for |
[in] | u | current local solution vector |
[in] | glob_ind | global indices of the local vector components |
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
References ug::GeomProvider< typename TGeom >::get(), subset_handler(), and UG_CATCH_THROW.
|
virtual |
prepares the element loop
This function is used to prepare the element loop. It gets called, before the loop over all elements starts. Here, general checks are performed:
In addition, local coordinates of evaluation points are requested by the DataImports in case of element-fixed points.
Reimplemented from ug::IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >.
References ug::GeomProvider< typename TGeom >::get(), and UG_THROW.
|
virtual |
type of trial space for each function used
Implements ug::IElemDiscBase< typename TDomain >.
References ug::LFEID::LAGRANGE, type(), and UG_THROW.
|
protectedvirtual |
register utils
|
protected |
register utils
|
virtual |
sets the density
Implements ug::NavierStokes::IncompressibleNavierStokesBase< TDomain >.
|
virtual |
sets the kinematic viscosity
Implements ug::NavierStokes::IncompressibleNavierStokesBase< TDomain >.
|
inline |
|
virtual |
sets the source function
Implements ug::NavierStokes::IncompressibleNavierStokesBase< TDomain >.
|
inline |
sets stabilization based on string identifier
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvUpwind, and ug::NavierStokes::NavierStokesFV1< TDomain >::m_spStab.
|
inline |
sets stabilization and diff length method based on string identifier (only for Schneider-Raw stabilizations)
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvUpwind, and ug::NavierStokes::NavierStokesFV1< TDomain >::m_spStab.
|
inline |
sets the stabilization used to compute the stabilized velocities
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_spStab.
|
inline |
sets the upwind based on a string identifier
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvStab, ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvUpwind, ug::NavierStokes::NavierStokesFV1< TDomain >::m_spStab, and SPNULL.
|
inline |
sets a stabilization for upwinding (Physical Advection Correction)
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvStab, ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvUpwind, and SPNULL.
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::set_pac_upwind().
|
inline |
sets an upwinding for the convective term of momentum equation
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvStab, ug::NavierStokes::NavierStokesFV1< TDomain >::m_spConvUpwind, and SPNULL.
|
inline |
returns scvf source
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_imSourceSCVF.
|
inline |
returns stabilization
References ug::NavierStokes::NavierStokesFV1< TDomain >::m_spStab.
|
staticprotected |
abbreviation for pressure
|
static |
World dimension.
|
protected |
flag if using only laplace term
|
protected |
flag if using Peclet Blending
|
protected |
flag if solving the Stokes equation
|
protected |
Export for the velocity.
|
protected |
Export for the velocity gradient.
|
protected |
|
protected |
Data import for density.
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::density().
|
protected |
Data import for kinematic viscosity.
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::kinematic_viscosity().
|
protected |
Data import for source.
|
protected |
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::sourceSCVF().
|
protected |
Stabilization for velocity in convective term of momentum equation Here, the stabilization is used as an upwinding
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::set_upwind(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_upwind(), and ug::NavierStokes::NavierStokesFV1< TDomain >::set_upwind().
|
protected |
Upwinding for velocity in convective term of momentum equation.
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::set_pac_upwind(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_stabilization(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_stabilization(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_upwind(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_upwind(), and ug::NavierStokes::NavierStokesFV1< TDomain >::set_upwind().
|
protected |
Stabilization for velocity in continuity equation.
Referenced by ug::NavierStokes::NavierStokesFV1< TDomain >::set_pac_upwind(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_stabilization(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_stabilization(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_stabilization(), ug::NavierStokes::NavierStokesFV1< TDomain >::set_upwind(), and ug::NavierStokes::NavierStokesFV1< TDomain >::stabilization().