1 #ifndef INTREPID_HGRAD_WEDGE_C1_FEMDEF_HPP 
    2 #define INTREPID_HGRAD_WEDGE_C1_FEMDEF_HPP 
   53   template<
class Scalar, 
class ArrayScalar>
 
   56     this -> basisCardinality_  = 6;
 
   57     this -> basisDegree_       = 1;    
 
   58     this -> basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Wedge<6> >() );
 
   59     this -> basisType_         = BASIS_FEM_DEFAULT;
 
   60     this -> basisCoordinates_  = COORDINATES_CARTESIAN;
 
   61     this -> basisTagsAreSet_   = 
false;
 
   65 template<
class Scalar, 
class ArrayScalar>
 
   75   int tags[]  = { 0, 0, 0, 1,
 
   83   Intrepid::setOrdinalTagData(
this -> tagToOrdinal_,
 
   84                               this -> ordinalToTag_,
 
   86                               this -> basisCardinality_,
 
   95 template<
class Scalar, 
class ArrayScalar>
 
   97                                                              const ArrayScalar &  inputPoints,
 
   98                                                              const EOperator      operatorType)
 const {
 
  101 #ifdef HAVE_INTREPID_DEBUG 
  102   Intrepid::getValues_HGRAD_Args<Scalar, ArrayScalar>(outputValues,
 
  105                                                       this -> getBaseCellTopology(),
 
  106                                                       this -> getCardinality() );
 
  110   int dim0 = inputPoints.dimension(0);  
 
  117   switch (operatorType) {
 
  120       for (
int i0 = 0; i0 < dim0; i0++) {
 
  121         x = inputPoints(i0, 0);
 
  122         y = inputPoints(i0, 1);
 
  123         z = inputPoints(i0, 2);
 
  126         outputValues(0, i0) = (1.0 - x - y)*(1.0 - z)/2.0;
 
  127         outputValues(1, i0) = x*(1.0 - z)/2.0;
 
  128         outputValues(2, i0) = y*(1.0 - z)/2.0;
 
  129         outputValues(3, i0) = (1.0 - x - y)*(1.0 + z)/2.0;
 
  130         outputValues(4, i0) = x*(1.0 + z)/2.0;
 
  131         outputValues(5, i0) = y*(1.0 + z)/2.0;
 
  137       for (
int i0 = 0; i0 < dim0; i0++) {
 
  138         x = inputPoints(i0,0);
 
  139         y = inputPoints(i0,1);
 
  140         z = inputPoints(i0, 2);
 
  143         outputValues(0, i0, 0) = -(1.0 - z)/2.0;
 
  144         outputValues(0, i0, 1) = -(1.0 - z)/2.0;
 
  145         outputValues(0, i0, 2) = -(1.0 - x - y)/2.0;
 
  147         outputValues(1, i0, 0) =  (1.0 - z)/2.0;
 
  148         outputValues(1, i0, 1) =  0.0;
 
  149         outputValues(1, i0, 2) = -x/2.0;
 
  151         outputValues(2, i0, 0) =  0.0;
 
  152         outputValues(2, i0, 1) =  (1.0 - z)/2.0;
 
  153         outputValues(2, i0, 2) = -y/2.0;
 
  156         outputValues(3, i0, 0) = -(1.0 + z)/2.0;
 
  157         outputValues(3, i0, 1) = -(1.0 + z)/2.0;
 
  158         outputValues(3, i0, 2) =  (1.0 - x - y)/2.0;
 
  160         outputValues(4, i0, 0) =  (1.0 + z)/2.0;
 
  161         outputValues(4, i0, 1) =  0.0;
 
  162         outputValues(4, i0, 2) =  x/2.0;
 
  164         outputValues(5, i0, 0) =  0.0;
 
  165         outputValues(5, i0, 1) =  (1.0 + z)/2.0;
 
  166         outputValues(5, i0, 2) =  y/2.0;
 
  171       TEUCHOS_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_CURL), std::invalid_argument,
 
  172                           ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): CURL is invalid operator for rank-0 (scalar) functions in 3D");
 
  176       TEUCHOS_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_DIV), std::invalid_argument,
 
  177                           ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): DIV is invalid operator for rank-0 (scalar) functions in 3D");
 
  181       for (
int i0 = 0; i0 < dim0; i0++) {
 
  182         outputValues(0, i0, 0) = 0.0;     outputValues(3, i0, 0) = 0.0; 
 
  183         outputValues(0, i0, 1) = 0.0;     outputValues(3, i0, 1) = 0.0;
 
  184         outputValues(0, i0, 2) = 0.5;     outputValues(3, i0, 2) =-0.5; 
 
  185         outputValues(0, i0, 3) = 0.0;     outputValues(3, i0, 3) = 0.0;
 
  186         outputValues(0, i0, 4) = 0.5;     outputValues(3, i0, 4) =-0.5;
 
  187         outputValues(0, i0, 5) = 0.0;     outputValues(3, i0, 5) = 0.0;
 
  189         outputValues(1, i0, 0) = 0.0;     outputValues(4, i0, 0) = 0.0; 
 
  190         outputValues(1, i0, 1) = 0.0;     outputValues(4, i0, 1) = 0.0; 
 
  191         outputValues(1, i0, 2) =-0.5;     outputValues(4, i0, 2) = 0.5;
 
  192         outputValues(1, i0, 3) = 0.0;     outputValues(4, i0, 3) = 0.0;
 
  193         outputValues(1, i0, 4) = 0.0;     outputValues(4, i0, 4) = 0.0;
 
  194         outputValues(1, i0, 5) = 0.0;     outputValues(4, i0, 5) = 0.0;
 
  196         outputValues(2, i0, 0) = 0.0;     outputValues(5, i0, 0) = 0.0;
 
  197         outputValues(2, i0, 1) = 0.0;     outputValues(5, i0, 1) = 0.0;
 
  198         outputValues(2, i0, 2) = 0.0;     outputValues(5, i0, 2) = 0.0;
 
  199         outputValues(2, i0, 3) = 0.0;     outputValues(5, i0, 3) = 0.0;
 
  200         outputValues(2, i0, 4) =-0.5;     outputValues(5, i0, 4) = 0.5;
 
  201         outputValues(2, i0, 5) = 0.0;     outputValues(5, i0, 5) = 0.0;
 
  215         int DkCardinality = Intrepid::getDkCardinality(operatorType, 
 
  216                                                        this -> basisCellTopology_.getDimension() );
 
  217         for(
int dofOrd = 0; dofOrd < 
this -> basisCardinality_; dofOrd++) {
 
  218           for (
int i0 = 0; i0 < dim0; i0++) {
 
  219             for(
int dkOrd = 0; dkOrd < DkCardinality; dkOrd++){
 
  220               outputValues(dofOrd, i0, dkOrd) = 0.0;
 
  228       TEUCHOS_TEST_FOR_EXCEPTION( !( Intrepid::isValidOperator(operatorType) ), std::invalid_argument,
 
  229                           ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): Invalid operator type");
 
  235 template<
class Scalar, 
class ArrayScalar>
 
  237                                                              const ArrayScalar &    inputPoints,
 
  238                                                              const ArrayScalar &    cellVertices,
 
  239                                                              const EOperator        operatorType)
 const {
 
  240   TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::logic_error,
 
  241                       ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): FEM Basis calling an FVD member function");
 
void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays. 
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const 
FEM basis evaluation on a reference Wedge cell. 
Basis_HGRAD_WEDGE_C1_FEM()
Constructor.