Quadrature.h DG++. More...
#include <Quadrature.h>
Public Member Functions | |
Quadrature (const double *const xqdat, const double *const wqdat, const size_t NC, const size_t NQ) | |
Quadrature (const double *const xqdatmap, const double *const xqdatshape, const double *const wqdat, const size_t NCmap, const size_t NCshape, const size_t NQ) | |
virtual | ~Quadrature () |
Quadrature (const Quadrature &) | |
Quadrature * | clone () const |
size_t | getNumQuadraturePoints () const |
size_t | getNumCoordinates () const |
Returns the number of map coordinates. | |
size_t | getNumShapeCoordinates () const |
Returns the number of shape coordinates. | |
const double * | getQuadraturePoint (size_t q) const |
Return map coordinates of quadrature point q. | |
const double * | getQuadraturePointShape (size_t q) const |
Return shape coordinates of quadrature point q. | |
const double | getQuadratureWeights (size_t q) const |
Returns weight of quadrature point q. | |
Private Attributes | |
double * | xqmap |
double * | xqshape |
double * | wq |
size_t | numMapCoordinates |
size_t | numShapeCoordinates |
size_t | numQuadraturePoints |
Quadrature.h DG++.
Created by Adrian Lew on 9/4/06.
Copyright (c) 2006 Adrian Lew
Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
base class for any quadrature rule A quadrature rule provides and approximation of the integral over a domain , and it has:
1) A number of quadrature point coordinates in . (map) 2) Weights at each quadrature point 3) A second set of coordinates for the same quadrature points. (Shape)
The number of coordinates in (1) of each quadrature point should be unisolvent, i.e., each coordinate can be varied independently. The reason for this choice is that otherwise we need a way to convey the constraint between coordinates when functions of these coordinates and their derivatives are considered. For example, the mapping from the parametric configuration to the real space in finite elements.
In addition to the unisolvent set of coordinates above, the class also allows for an alternative set of coordinates for each gauss point. These do not need to be unisolvent, and can be used to localize the quadrature points when embedded in a higher dimensional space. For example, if is a triangle in 3D, then the second set of coordinates provide the coordinates of these points in 3D.
A specific example is the following: consider a triangle in 2D, and the integration over a segment on its boundary. The quadrature object for a segment should have only one (map) coordinate in (1) to be unisolvent. This coordinate is one of the barycentric coordinates of the quadrature point over the segment. However, if we have a function over the triangle and would like to restrict its values to the segment, we need the coordinates of the quadrature points in the triangle. This is the case when integrating a shape function in the triangle over the segment. This second set of coordinates is given in (3).
Based on this example we shall call the coordinates in (1) map coordinates, and the coordinates in (3) Shape coordinates.
Although the coordinates in (3) are not a natural inclusion in a quadrature rule, they provide the simplest implementation for the scenario just described.
A number of specific Quadrature objects are defined in Quadrature.cpp, and declared here.
Quadrature::Quadrature | ( | const double *const | xqdat, | |
const double *const | wqdat, | |||
const size_t | NC, | |||
const size_t | NQ | |||
) |
NQ | number of quadrature points | |
NC | number of coordinates for each point | |
xqdat | vector with coordinates of the quadrature points xqdat[a*NC+i] gives the i-th coordinate of quadrature point a | |
wqdat | vector with quadrature point weights wq[a] is the weight of quadrature point a This constructor sets the two sets of coordinates for the quadrature points to be the same. |
Quadrature::Quadrature | ( | const double *const | xqdatmap, | |
const double *const | xqdatshape, | |||
const double *const | wqdat, | |||
const size_t | NCmap, | |||
const size_t | NCshape, | |||
const size_t | NQ | |||
) |
NQ | number of quadrature points | |
NCmap | number of map coordinates for each point | |
NCshape | number of shape coordinates for each point | |
xqdatmap | vector with coordinates of the quadrature points xqdatmap[a*NC+i] gives the i-th coordinate of quadrature point a | |
xqdatshape | vector with coordinates of the quadrature points xqdatshape[a*NC+i] gives the i-th coordinate of quadrature point a | |
wqdat | vector with quadrature point weights wq[a] is the weight of quadrature point a |
virtual Quadrature::~Quadrature | ( | ) | [inline, virtual] |
Quadrature::Quadrature | ( | const Quadrature & | SQ | ) |
Quadrature * Quadrature::clone | ( | ) | const |
size_t Quadrature::getNumCoordinates | ( | ) | const [inline] |
Returns the number of map coordinates.
size_t Quadrature::getNumQuadraturePoints | ( | ) | const [inline] |
size_t Quadrature::getNumShapeCoordinates | ( | ) | const [inline] |
Returns the number of shape coordinates.
const double* Quadrature::getQuadraturePoint | ( | size_t | q | ) | const [inline] |
Return map coordinates of quadrature point q.
const double* Quadrature::getQuadraturePointShape | ( | size_t | q | ) | const [inline] |
Return shape coordinates of quadrature point q.
const double Quadrature::getQuadratureWeights | ( | size_t | q | ) | const [inline] |
Returns weight of quadrature point q.
size_t Quadrature::numMapCoordinates [private] |
size_t Quadrature::numQuadraturePoints [private] |
size_t Quadrature::numShapeCoordinates [private] |
double* Quadrature::wq [private] |
double* Quadrature::xqmap [private] |
double* Quadrature::xqshape [private] |