Class representing the Generalized Duffy quadrature rule in 2d, $ x = u^\beta, y = u^\beta v $, that is for $ \beta=1 $ the standard duffy integration rule. More...

#include <quadRule.hh>

Inheritance diagram for concepts::QuadratureRule2dQuadDuffy:
concepts::QuadratureRule2d concepts::QuadratureRule concepts::OutputOperator

Public Member Functions

virtual const Realabscissas (uint i) const
 Returns the quadrature abcissas in the i-th direction. More...
 
virtual const bool domain () const
 Method delivers if the integration points and weights are computed in [-1,1]^2 or in [0,1]^2. More...
 
const uint n (uint i=0) const
 Returns the number of quadrature points. More...
 
void printData () const
 
virtual bool quadratureData (const uint i, Real3d &q, Real &w) const
 Method delivers the i-th quadrature point. More...
 
 QuadratureRule2dQuadDuffy (Real beta, uint nx, uint ny, uint noVtx=0)
 Constructor. More...
 
virtual const bool tensor () const
 Method delivers information about the quadrature rule structure. More...
 
virtual const Realweights (uint i=0) const
 Returns the i-th weight, belonging to (x_i,y_i) as nontensored. More...
 
virtual ~QuadratureRule2dQuadDuffy ()
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream. More...
 

Protected Attributes

const RealintX_
 Abscissas. More...
 
const RealintY_
 

Private Member Functions

void compute_ (const Real2d S, const Real2d B, const Real2d C, const Real beta, const uint nx, const uint ny, Real *&abscXP, Real *&abscYP, Real *&wP)
 Routine computes generalized duffy quadrature points on the triangle conv(S,B,C), where S is the coordinate of the singularity vertex. More...
 

Private Attributes

Real beta_
 
uint n_
 
const Realweights_
 

Detailed Description

Class representing the Generalized Duffy quadrature rule in 2d, $ x = u^\beta, y = u^\beta v $, that is for $ \beta=1 $ the standard duffy integration rule.

The quadrature rule is not tensored and is most effective to integrate $ 1/r^\alpha $ singularities, for $ 0<\alpha < 2 $.

The weights are given in $ [-1,1]^2 $.

This class is motivated by the paper "Generalized Duffy Transformation for Integrating Vertex Singularites" by Mousavi / Sukumar

The quadrature rule depends on the duffy parameter beta. beta should be chosen such that

\[2\beta -1 -\alpha\beta \geq 0\]

.

Note that for too big beta, the quadrature points in $ [0,1]^2 $ become very close to the origin. This should be taken into account, regarding mashine precision. Due to this fact the singular local vertex by default is in origin.

The abscissas are computed on $ [0,1]^2 $.

The internal compute routine computes generalized duffy points on a triangle. With that, the method can be used for arbitrary reference geometries that can be subdivided into triangles with common singularity vertex.

The quadrature points on the transformed triangles are computed with underlying Gauss-Legrendre points, whose number can be chosen individual in each direction in the transformed quad.

Author
Robert Gruhlke, 2014

Definition at line 481 of file quadRule.hh.

Constructor & Destructor Documentation

◆ QuadratureRule2dQuadDuffy()

concepts::QuadratureRule2dQuadDuffy::QuadratureRule2dQuadDuffy ( Real  beta,
uint  nx,
uint  ny,
uint  noVtx = 0 
)

Constructor.

The singular vertex is by default set to zero. This is ok, since when building the mesh one can place the geometry s.t. the first elements vertex lies on the singularity.

Parameters
betageneralized duffy transformation parameter
nxnumber of quadrature points of underlying Gauss-Legendre points in x-direction on transformed square.
nynumber of quadrature points of underlying Gauss-Legendre points in y-direction on transformed square.
noVtxnumber of local vertex

◆ ~QuadratureRule2dQuadDuffy()

virtual concepts::QuadratureRule2dQuadDuffy::~QuadratureRule2dQuadDuffy ( )
virtual

Member Function Documentation

◆ abscissas()

virtual const Real* concepts::QuadratureRule2d::abscissas ( uint  i) const
inlinevirtualinherited

Returns the quadrature abcissas in the i-th direction.

Parameters
ii = 0 : x-direction, i = 1 : y-direction

Implements concepts::QuadratureRule.

Definition at line 360 of file quadRule.hh.

◆ compute_()

void concepts::QuadratureRule2dQuadDuffy::compute_ ( const Real2d  S,
const Real2d  B,
const Real2d  C,
const Real  beta,
const uint  nx,
const uint  ny,
Real *&  abscXP,
Real *&  abscYP,
Real *&  wP 
)
private

Routine computes generalized duffy quadrature points on the triangle conv(S,B,C), where S is the coordinate of the singularity vertex.

Parameters
S,B,CTriangle points
betaduffy parameter
nxnumber of points of underlying quadrature in x direction
nynumber of points of underlying quadrature in y direction
abscXPabscissa points in x direction
abscYPabscissa points in y direction
wPweights

◆ domain()

virtual const bool concepts::QuadratureRule2dQuadDuffy::domain ( ) const
inlinevirtual

Method delivers if the integration points and weights are computed in [-1,1]^2 or in [0,1]^2.

Return true if integration points are computed in

\[ [-1,1]^2 \]

or false if they are computed in

\[ [0,1]^2 \]

.

Implements concepts::QuadratureRule2d.

Definition at line 520 of file quadRule.hh.

◆ info()

virtual std::ostream& concepts::QuadratureRule2dQuadDuffy::info ( std::ostream &  os) const
protectedvirtual

Returns information in an output stream.

Reimplemented from concepts::OutputOperator.

◆ n()

const uint concepts::QuadratureRule2dQuadDuffy::n ( uint  i = 0) const
inlinevirtual

Returns the number of quadrature points.

Implements concepts::QuadratureRule2d.

Definition at line 511 of file quadRule.hh.

◆ printData()

void concepts::QuadratureRule2dQuadDuffy::printData ( ) const

◆ quadratureData()

virtual bool concepts::QuadratureRule2dQuadDuffy::quadratureData ( const uint  i,
Real3d q,
Real w 
) const
virtual

Method delivers the i-th quadrature point.

\[ q \in [0,1]^n \]

and its belonging i-th weight. If i exceeds the numer of quadrature, false is returned.

Implements concepts::QuadratureRule.

◆ tensor()

virtual const bool concepts::QuadratureRule2dQuadDuffy::tensor ( ) const
inlinevirtual

Method delivers information about the quadrature rule structure.

Returns true for tensor structur and false for non-tensor one.

Implements concepts::QuadratureRule2d.

Definition at line 523 of file quadRule.hh.

◆ weights()

virtual const Real* concepts::QuadratureRule2dQuadDuffy::weights ( uint  i = 0) const
inlinevirtual

Returns the i-th weight, belonging to (x_i,y_i) as nontensored.

Implements concepts::QuadratureRule2d.

Definition at line 503 of file quadRule.hh.

Member Data Documentation

◆ beta_

Real concepts::QuadratureRule2dQuadDuffy::beta_
private

Definition at line 552 of file quadRule.hh.

◆ intX_

const Real* concepts::QuadratureRule2d::intX_
protectedinherited

Abscissas.

Definition at line 382 of file quadRule.hh.

◆ intY_

const Real* concepts::QuadratureRule2d::intY_
protectedinherited

Definition at line 383 of file quadRule.hh.

◆ n_

uint concepts::QuadratureRule2dQuadDuffy::n_
private

Definition at line 554 of file quadRule.hh.

◆ weights_

const Real* concepts::QuadratureRule2dQuadDuffy::weights_
private

Definition at line 550 of file quadRule.hh.


The documentation for this class was generated from the following file:
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich