WARPXM v1.10.0
Loading...
Searching...
No Matches
wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts Class Reference

This class solves: \(\boldsymbol{v} + \nabla\cdot \bar{\bar{F}} = \boldsymbol{S}\) in weak form for auxiliary variable \(\boldsymbol{v}\) by parts in same manner as with DG. More...

#include <divergence_integral_by_parts.h>

Inheritance diagram for wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts:
wxm::dfem::variable_adjuster::variable_adjuster_t WmPatchProcess

Detailed Description

This class solves: \(\boldsymbol{v} + \nabla\cdot \bar{\bar{F}} = \boldsymbol{S}\) in weak form for auxiliary variable \(\boldsymbol{v}\) by parts in same manner as with DG.

Useful for solving terms of Ohm's Law for electric field, which is a state equation, as opposed to a time advancing pde.

For warpy class see warpy.variable_adjusters.auxiliary_variables.auxiliary_variables.divergence_integral_by_parts

Explicitly, it calculates:

\begin{align*} v_{ij}^{\lambda} = {J_{lm}^{\lambda}}^{-1} \Upsilon_{jlk} F_{imk}^{\lambda} - \sum_{\gamma \in \bar{\Gamma}_{\lambda}} G_{\lambda\gamma}^{-1} \Xi_{jk}^{\lambda\gamma} F_{ik}^{\lambda\gamma} + S_{ij} \end{align*}

For the internal flux part \({J_{lm}^{\lambda}}^{-1} \Upsilon_{jlk} F_{imk}^{\lambda}\):

  • \({J_{lm}^{\lambda}}^{-1} = \frac{\partial \xi_{l}}{\partial x_{m}}\) is the inverse jacobian for the element.
  • \(\Upsilon_{jlk}\) is the internal flux matrix constructed as \(\Upsilon_{jlk} = M_{jp}^{-1}A_{klp}\).
  • \(M_{jp}\) is mass matrix \(\int\limits_{\lambda} \psi_{j} \psi_{p}d\boldsymbol{\xi}\).
  • \(A_{klp}\) is the advection matrix \(\int\limits_{\lambda} \psi_{k} \frac{\partial}{\partial\xi_{l}} \psi_{p} d\boldsymbol{\xi}\).
  • \(F_{imk}^{\lambda}\) is the analytical flux.

For the numerical flux part \(\sum_{\gamma \in \bar{\Gamma}_{\lambda}} G_{\lambda\gamma}^{-1}\Xi_{jk}^{\lambda\gamma} F_{ik}^{\lambda\gamma}\):

  • \(G_{\lambda\gamma}^{-1} = \left|\left|{J_{ji}^{\lambda}}^{-1} m_{j}^{\lambda\gamma}\right|\right|\) where \(m\) is the normal of the isoparametric element for the face of element \(\lambda\) facing element \(\gamma\).
  • \(\Xi_{jk}^{\lambda\gamma}\) is the numerical flux matrix constructed as \(\Xi_{jk}^{\lambda\gamma} = M_{jl}^{-1}\Gamma_{lk}^{\lambda\gamma}\).
  • \(\Gamma_{lk}^{\lambda\gamma}\) is the area mass matrix \(\int\limits_{\partial\lambda^{\gamma}} \psi_{l} \psi_{k} d\boldsymbol{\xi}\) integrated along the simplex face between elements \(\lambda\) and \(\gamma\).
  • \(F_{ik}^{\lambda\gamma}\) is the numerical flux on the the simplex face between elements \(\lambda\) and \(\gamma\).

For the source part \(S_{ij}\):

  • \(S_{ij}\) is calulated and accumulated into \(v_{ij}\) for LGL node source calculation.
  • For GQ node source calculation, we need \(^{3}\Upsilon_{jk} S_{ik}\) where
  • \(^{3}\Upsilon_{jk} = M_{jl}^{-1}\Theta_{lk}\) where
  • \(\Theta_{lk}\) is the source mass matrix \(w_{l}\psi_{k}(\boldsymbol{\xi}_{l})\) where
  • \(w_{l}\) is the gaussian quadrature basis function weight for its node located at \(\boldsymbol{\xi}_{l}\) and
  • \(\psi_{k}\) is the LGL basis function used for the DG scheme.

Note that in general:

  • \(\lambda\) = "inside" element
  • \(\gamma\) = "outside" element
  • i = auxiliary variable component
  • j = local node
  • k = cross local node
  • l = isoparametric space dimension index
  • m = real space dimension index

Public Member Functions

 DivergenceIntegralByParts ()=default
 
void setup (const WxCryptSet &wxc) override
 Setup the spatial solver using the cryptset.
 
void solve (real time, variables_type &input) override
 Solves the spatial system and puts result in _rhs.
 
- Public Member Functions inherited from wxm::dfem::variable_adjuster::variable_adjuster_t
 variable_adjuster_t ()
 Constructor.
 
void setup (const WxCryptSet &wxc) override
 Setup the spatial solver using the cryptset.
 
virtual void solve (real time, variables_type &input)=0
 applies the variable adjuster to input
 
std::vector< wxm::array::patch_array_t * > get_patch_arrays (const variables_type &in)
 
void process () override
 
virtual const std::vector< size_t > & get_input_indices () const
 
int priority () const
 Getter function for priority.
 
virtual void Barrier (const WxMsgBase &msg, const real time, variables_type &input)
 An opportunity for the variable adjuster to perform any MPI operations it needs to after solve is called.
 
- Public Member Functions inherited from WmPatchProcess
 WmPatchProcess ()
 
virtual ~WmPatchProcess ()=default
 
void step ()
 
virtual void process ()=0
 
void setPatch (const WmUnstructuredPatch *patch)
 
void setParentTaskProcessor (WmPatchProcessor *taskProcessor)
 
virtual std::string name (std::string prefix="patchedProcessor.") const
 
virtual void setup (const WxCryptSet &wxc)
 
void setParentSolver (const WmSolverBase *solver)
 
void setDt (WxStepper::time_t dt)
 
void setTime (WxStepper::time_t time)
 
WxStepper::time_t getSuggestedDt () const
 
const WmUnstructuredPatchget_patch () const
 
const WmSolverget_solver () const
 

Protected Member Functions

void rhs (real time, const WxRange &element_scope_range, std::vector< wxm::array::patch_array_t * > input)
 Actually calculates the right hand side for Auxiliary Variable \(v_{ij}^{\lambda}\).
 

Protected Attributes

std::string _basisSetName
 Spatial basis set name.
 
WmBasisArraySet _basis_set
 Spatial basis set.
 
std::unique_ptr< WmUDGGeometry_dg_geometry
 WmUDGGeometry.
 
wxm::dfem::tools::scope_t _element_flux
 Scope.
 
real _penalty_beta
 Penalties.
 
real _penalty_eta
 
std::string _source_quad_style
 Source quadrature style.
 
bool _has_fluxes
 
bool _has_sources
 
bool _sources_lgl
 
bool _sources_gq
 
std::vector< size_t > _rhs_idcs
 checking of indexes
 
std::vector< size_t > _all_idcs
 
- Protected Attributes inherited from wxm::dfem::variable_adjuster::variable_adjuster_t
std::vector< size_t > _var_idcs
 
std::vector< std::unique_ptr< WmApplication > > _apps
 
bool _grab_all = false
 
int _priority
 
- Protected Attributes inherited from WmPatchProcess
std::string _ppName
 
WmPatchProcessor_parentTaskProcessor
 
const WmUnstructuredPatch_patch
 
std::vector< std::string > _onSubdomains
 
WxStepper::time_t _time
 
WxStepper::time_t _dt
 
WxStepper::time_t _suggested_dt
 
WxLogStream _debStrm = WxLogger::get("warpx-root.console")->getDebugStream()
 
const WmSolverBase_parentSolver = nullptr
 

Additional Inherited Members

- Public Types inherited from wxm::dfem::variable_adjuster::variable_adjuster_t
typedef wxm::temporal_solver::variables_type variables_type
 

Constructor & Destructor Documentation

◆ DivergenceIntegralByParts()

wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::DivergenceIntegralByParts ( )
default

Member Function Documentation

◆ rhs()

void wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::rhs ( real  time,
const WxRange element_scope_range,
std::vector< wxm::array::patch_array_t * >  input 
)
protected

Actually calculates the right hand side for Auxiliary Variable \(v_{ij}^{\lambda}\).

◆ setup()

void wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::setup ( const WxCryptSet wxc)
overridevirtual

Setup the spatial solver using the cryptset.

Parameters
wxcThe cryptset

Reimplemented from WmPatchProcess.

◆ solve()

void wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::solve ( real  time,
variables_type input 
)
overridevirtual

Solves the spatial system and puts result in _rhs.

Implements wxm::dfem::variable_adjuster::variable_adjuster_t.

Member Data Documentation

◆ _all_idcs

std::vector<size_t> wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_all_idcs
protected

◆ _basis_set

WmBasisArraySet wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_basis_set
protected

Spatial basis set.

◆ _basisSetName

std::string wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_basisSetName
protected

Spatial basis set name.

◆ _dg_geometry

std::unique_ptr<WmUDGGeometry> wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_dg_geometry
protected

◆ _element_flux

wxm::dfem::tools::scope_t wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_element_flux
protected

Scope.

◆ _has_fluxes

bool wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_has_fluxes
protected

◆ _has_sources

bool wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_has_sources
protected

◆ _penalty_beta

real wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_penalty_beta
protected

Penalties.

◆ _penalty_eta

real wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_penalty_eta
protected

◆ _rhs_idcs

std::vector<size_t> wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_rhs_idcs
protected

checking of indexes

◆ _source_quad_style

std::string wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_source_quad_style
protected

Source quadrature style.

◆ _sources_gq

bool wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_sources_gq
protected

◆ _sources_lgl

bool wxm::dfem::variable_adjuster::auxiliary_variables::DivergenceIntegralByParts::_sources_lgl
protected

The documentation for this class was generated from the following file: