2#include "Sacado_Fad_DFad.hpp"
7 template<
class Scalar,
class ArrayScalar>
18 template<
class Scalar,
class ArrayScalar>
52 template<
class Scalar,
class ArrayScalar>
54 const ArrayScalar& inputPoints,
55 const EOperator operatorType)
const{
56 TEUCHOS_TEST_FOR_EXCEPTION (
true, std::logic_error,
57 ">>>ERROR (Basis_HGRAD_POLY_C1_FEM): Polygonal basis calling FEM member function");
61 template<
class Scalar,
class ArrayScalar>
63 const ArrayScalar& inputPoints,
64 const ArrayScalar& cellVertices,
65 const EOperator operatorType)
const{
69 switch (operatorType) {
78 for (
int i=0;i<inputPoints.dimension(0);i++){
79 for (
int j=0;j<2;j++){
80 dInput(i,j) = Sacado::Fad::DFad<Scalar>( inputPoints(i,j));
81 dInput(i,j).diff(j,2);
85 for (
int i=0;i<cellVertices.dimension(0);i++){
86 for (
int j=0;j<cellVertices.dimension(1);j++){
87 cellVerticesFad(i,j) = Sacado::Fad::DFad<Scalar>( cellVertices(i,j) );
95 for (
int i=0;i<outputValues.dimension(0);i++){
96 for (
int j=0;j<outputValues.dimension(1);j++){
97 for (
int k=0;k<outputValues.dimension(2);k++){
98 outputValues(i,j,k) = dOutput(i,j).dx(k);
115 TEUCHOS_TEST_FOR_EXCEPTION (
true, std::invalid_argument,
116 ">>> ERROR (Basis_HGRAD_POLY_C1_FEM): operator not implemented yet");
121 ">>> ERROR (Basis_HGRAD_POLY_C1_FEM): Invalid operator type");
128 template<
class Scalar,
class ArrayScalar>
129 template<
class Scalar1,
class ArrayScalar1>
131 const ArrayScalar1& inputPoints,
132 const ArrayScalar1& cellVertices)
const{
133 int numPoints = inputPoints.dimension(0);
136 for (
int pointIndex = 0;pointIndex<outputValues.dimension(1);pointIndex++){
137 Scalar1 denominator=0;
139 for (
int k=0;k<weightFunctions.
dimension(0);k++){
140 denominator += weightFunctions(k,pointIndex);
142 for (
int k=0;k<outputValues.dimension(0);k++){
143 outputValues(k,pointIndex) = weightFunctions(k,pointIndex)/denominator;
150 template<
class Scalar,
class ArrayScalar>
151 template<
class Scalar1,
class ArrayScalar1>
153 const ArrayScalar1& p2,
154 const ArrayScalar1& p3)
const{
155 Scalar1 area = 0.5*(p1(0)*(p2(1)-p3(1))
157 +p3(0)*(p1(1)-p2(1)));
162 template<
class Scalar,
class ArrayScalar>
163 template<
class Scalar1,
class ArrayScalar1>
165 const ArrayScalar1& inputValues,
166 const ArrayScalar1& cellVertices)
const{
170 for (
int k=0;k<outputValues.dimension(0);k++){
174 int adjIndex1 = -1, adjIndex2 = -1;
181 TEUCHOS_TEST_FOR_EXCEPTION( (adjIndex1 == -1 || adjIndex2 == -1), std::invalid_argument,
182 ">>> ERROR (Intrepid_HGRAD_POLY_C1_FEM): cannot find adjacent nodes when evaluating Wachspress weight function.");
186 for (
int i=0;i<spaceDim;i++){
187 p1(i) = cellVertices(adjIndex1,i);
188 p2(i) = cellVertices(k,i);
189 p3(i) = cellVertices(adjIndex2,i);
193 for (
int pointIndex = 0;pointIndex < inputValues.dimension(0);pointIndex++){
194 Scalar1 product = a_k;
195 for (
int edgeIndex = 0;edgeIndex < this->
basisCellTopology_.getEdgeCount();edgeIndex++){
198 if ( edgeNode1 != k && edgeNode2 != k ){
199 for (
int i=0;i<spaceDim;i++){
200 p1(i) = inputValues(pointIndex,i);
201 p2(i) = cellVertices(edgeNode1,i);
202 p3(i) = cellVertices(edgeNode2,i);
207 outputValues(k,pointIndex) = product;
215#if defined(Intrepid_SHOW_DEPRECATED_WARNINGS)
217#warning "The Intrepid package is deprecated"
Header file for utility class to provide multidimensional containers.
int isValidOperator(const EOperator operatorType)
Verifies validity of an operator enum.
void setOrdinalTagData(std::vector< std::vector< std::vector< int > > > &tagToOrdinal, std::vector< std::vector< int > > &ordinalToTag, const int *tags, const int basisCard, const int tagSize, const int posScDim, const int posScOrd, const int posDfOrd)
Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data.
Scalar1 computeArea(const ArrayScalar1 &p1, const ArrayScalar1 &p2, const ArrayScalar1 &p3) const
Helper function to compute area of triangle formed by 3 points.
void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
FEM reference basis evaluation: invocation of this method throws an exception.
void evaluateWeightFunctions(ArrayScalar1 &outputValues, const ArrayScalar1 &inputValues, const ArrayScalar1 &cellVertices) const
Evaluation of the Wachspress weight functions.
void shapeFunctions(ArrayScalar1 &outputValues, const ArrayScalar1 &inputValues, const ArrayScalar1 &cellVertices) const
Evaluation of Wachspress shape functions.
Basis_HGRAD_POLY_C1_FEM(const shards::CellTopology &cellTopology)
Constructor.
bool basisTagsAreSet_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
std::vector< std::vector< std::vector< int > > > tagToOrdinal_
DoF tag to ordinal lookup table.
int basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom.
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
EBasis basisType_
Type of the basis.
int basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
std::vector< std::vector< int > > ordinalToTag_
DoF ordinal to tag lookup table.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package http://trilinos....
virtual int getCardinality() const
Returns cardinality of the basis.
Implementation of a templated lexicographical container for a multi-indexed scalar quantity....
int dimension(const int whichDim) const
Returns the specified dimension.