Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/src/Discretization/Basis/Intrepid_HGRAD_WEDGE_C1_FEMDef.hpp
Go to the documentation of this file.
00001 // @HEADER
00002 // ************************************************************************
00003 //
00004 //                           Intrepid Package
00005 //                 Copyright (2007) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Pavel Bochev (pbboche@sandia.gov) or
00025 //                    Denis Ridzal (dridzal@sandia.gov).
00026 //
00027 // ************************************************************************
00028 // @HEADER
00029 
00035 namespace Intrepid {
00036 
00037   template<class Scalar, class ArrayScalar>
00038   Basis_HGRAD_WEDGE_C1_FEM<Scalar, ArrayScalar>::Basis_HGRAD_WEDGE_C1_FEM()
00039   {
00040     this -> basisCardinality_  = 6;
00041     this -> basisDegree_       = 1;    
00042     this -> basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Wedge<6> >() );
00043     this -> basisType_         = BASIS_FEM_DEFAULT;
00044     this -> basisCoordinates_  = COORDINATES_CARTESIAN;
00045     this -> basisTagsAreSet_   = false;
00046   }
00047   
00048   
00049 template<class Scalar, class ArrayScalar>
00050 void Basis_HGRAD_WEDGE_C1_FEM<Scalar, ArrayScalar>::initializeTags() {
00051   
00052   // Basis-dependent intializations
00053   int tagSize  = 4;        // size of DoF tag
00054   int posScDim = 0;        // position in the tag, counting from 0, of the subcell dim 
00055   int posScOrd = 1;        // position in the tag, counting from 0, of the subcell ordinal
00056   int posDfOrd = 2;        // position in the tag, counting from 0, of DoF ordinal relative to the subcell
00057 
00058   // An array with local DoF tags assigned to basis functions, in the order of their local enumeration 
00059   int tags[]  = { 0, 0, 0, 1,
00060                   0, 1, 0, 1,
00061                   0, 2, 0, 1,
00062                   0, 3, 0, 1,
00063                   0, 4, 0, 1,
00064                   0, 5, 0, 1 };
00065   
00066   // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
00067   Intrepid::setOrdinalTagData(this -> tagToOrdinal_,
00068                               this -> ordinalToTag_,
00069                               tags,
00070                               this -> basisCardinality_,
00071                               tagSize,
00072                               posScDim,
00073                               posScOrd,
00074                               posDfOrd);
00075 }
00076 
00077 
00078 
00079 template<class Scalar, class ArrayScalar>
00080 void Basis_HGRAD_WEDGE_C1_FEM<Scalar, ArrayScalar>::getValues(ArrayScalar &    outputValues,
00081                                                              const ArrayScalar &  inputPoints,
00082                                                              const EOperator      operatorType) const {
00083   
00084   // Verify arguments
00085 #ifdef HAVE_INTREPID_DEBUG
00086   Intrepid::getValues_HGRAD_Args<Scalar, ArrayScalar>(outputValues,
00087                                                       inputPoints,
00088                                                       operatorType,
00089                                                       this -> getBaseCellTopology(),
00090                                                       this -> getCardinality() );
00091 #endif
00092   
00093   // Number of evaluation points = dim 0 of inputPoints
00094   int dim0 = inputPoints.dimension(0);  
00095   
00096   // Temporaries: (x,y,z) coordinates of the evaluation point
00097   Scalar x = 0.0;                                    
00098   Scalar y = 0.0;   
00099   Scalar z = 0.0;
00100   
00101   switch (operatorType) {
00102     
00103     case OPERATOR_VALUE:
00104       for (int i0 = 0; i0 < dim0; i0++) {
00105         x = inputPoints(i0, 0);
00106         y = inputPoints(i0, 1);
00107         z = inputPoints(i0, 2);
00108         
00109         // outputValues is a rank-2 array with dimensions (basisCardinality_, dim0)
00110         outputValues(0, i0) = (1.0 - x - y)*(1.0 - z)/2.0;
00111         outputValues(1, i0) = x*(1.0 - z)/2.0;
00112         outputValues(2, i0) = y*(1.0 - z)/2.0;
00113         outputValues(3, i0) = (1.0 - x - y)*(1.0 + z)/2.0;
00114         outputValues(4, i0) = x*(1.0 + z)/2.0;
00115         outputValues(5, i0) = y*(1.0 + z)/2.0;
00116       }
00117       break;
00118       
00119     case OPERATOR_GRAD:
00120     case OPERATOR_D1:
00121       for (int i0 = 0; i0 < dim0; i0++) {
00122         x = inputPoints(i0,0);
00123         y = inputPoints(i0,1);
00124         z = inputPoints(i0, 2);
00125         
00126         // outputValues is a rank-3 array with dimensions (basisCardinality_, dim0, spaceDim)
00127         outputValues(0, i0, 0) = -(1.0 - z)/2.0;
00128         outputValues(0, i0, 1) = -(1.0 - z)/2.0;
00129         outputValues(0, i0, 2) = -(1.0 - x - y)/2.0;
00130         
00131         outputValues(1, i0, 0) =  (1.0 - z)/2.0;
00132         outputValues(1, i0, 1) =  0.0;
00133         outputValues(1, i0, 2) = -x/2.0;
00134         
00135         outputValues(2, i0, 0) =  0.0;
00136         outputValues(2, i0, 1) =  (1.0 - z)/2.0;
00137         outputValues(2, i0, 2) = -y/2.0;
00138   
00139         
00140         outputValues(3, i0, 0) = -(1.0 + z)/2.0;
00141         outputValues(3, i0, 1) = -(1.0 + z)/2.0;
00142         outputValues(3, i0, 2) =  (1.0 - x - y)/2.0;
00143         
00144         outputValues(4, i0, 0) =  (1.0 + z)/2.0;
00145         outputValues(4, i0, 1) =  0.0;
00146         outputValues(4, i0, 2) =  x/2.0;
00147         
00148         outputValues(5, i0, 0) =  0.0;
00149         outputValues(5, i0, 1) =  (1.0 + z)/2.0;
00150         outputValues(5, i0, 2) =  y/2.0;
00151       }
00152       break;
00153       
00154     case OPERATOR_CURL:
00155       TEST_FOR_EXCEPTION( (operatorType == OPERATOR_CURL), std::invalid_argument,
00156                           ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): CURL is invalid operator for rank-0 (scalar) functions in 3D");
00157       break;
00158       
00159     case OPERATOR_DIV:
00160       TEST_FOR_EXCEPTION( (operatorType == OPERATOR_DIV), std::invalid_argument,
00161                           ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): DIV is invalid operator for rank-0 (scalar) functions in 3D");
00162       break;
00163       
00164     case OPERATOR_D2:
00165       for (int i0 = 0; i0 < dim0; i0++) {
00166         outputValues(0, i0, 0) = 0.0;     outputValues(3, i0, 0) = 0.0; 
00167         outputValues(0, i0, 1) = 0.0;     outputValues(3, i0, 1) = 0.0;
00168         outputValues(0, i0, 2) = 0.5;     outputValues(3, i0, 2) =-0.5; 
00169         outputValues(0, i0, 3) = 0.0;     outputValues(3, i0, 3) = 0.0;
00170         outputValues(0, i0, 4) = 0.5;     outputValues(3, i0, 4) =-0.5;
00171         outputValues(0, i0, 5) = 0.0;     outputValues(3, i0, 5) = 0.0;
00172         
00173         outputValues(1, i0, 0) = 0.0;     outputValues(4, i0, 0) = 0.0; 
00174         outputValues(1, i0, 1) = 0.0;     outputValues(4, i0, 1) = 0.0; 
00175         outputValues(1, i0, 2) =-0.5;     outputValues(4, i0, 2) = 0.5;
00176         outputValues(1, i0, 3) = 0.0;     outputValues(4, i0, 3) = 0.0;
00177         outputValues(1, i0, 4) = 0.0;     outputValues(4, i0, 4) = 0.0;
00178         outputValues(1, i0, 5) = 0.0;     outputValues(4, i0, 5) = 0.0;
00179         
00180         outputValues(2, i0, 0) = 0.0;     outputValues(5, i0, 0) = 0.0;
00181         outputValues(2, i0, 1) = 0.0;     outputValues(5, i0, 1) = 0.0;
00182         outputValues(2, i0, 2) = 0.0;     outputValues(5, i0, 2) = 0.0;
00183         outputValues(2, i0, 3) = 0.0;     outputValues(5, i0, 3) = 0.0;
00184         outputValues(2, i0, 4) =-0.5;     outputValues(5, i0, 4) = 0.5;
00185         outputValues(2, i0, 5) = 0.0;     outputValues(5, i0, 5) = 0.0;
00186       }
00187       break;
00188       
00189     case OPERATOR_D3:
00190     case OPERATOR_D4:
00191     case OPERATOR_D5:
00192     case OPERATOR_D6:
00193     case OPERATOR_D7:
00194     case OPERATOR_D8:
00195     case OPERATOR_D9:
00196     case OPERATOR_D10:
00197       {
00198         // outputValues is a rank-3 array with dimensions (basisCardinality_, dim0, DkCardinality)
00199         int DkCardinality = Intrepid::getDkCardinality(operatorType, 
00200                                                        this -> basisCellTopology_.getDimension() );
00201         for(int dofOrd = 0; dofOrd < this -> basisCardinality_; dofOrd++) {
00202           for (int i0 = 0; i0 < dim0; i0++) {
00203             for(int dkOrd = 0; dkOrd < DkCardinality; dkOrd++){
00204               outputValues(dofOrd, i0, dkOrd) = 0.0;
00205             }
00206           }
00207         }
00208       }
00209       break;
00210       
00211     default:
00212       TEST_FOR_EXCEPTION( !( Intrepid::isValidOperator(operatorType) ), std::invalid_argument,
00213                           ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): Invalid operator type");
00214   }
00215 }
00216 
00217 
00218   
00219 template<class Scalar, class ArrayScalar>
00220 void Basis_HGRAD_WEDGE_C1_FEM<Scalar, ArrayScalar>::getValues(ArrayScalar&       outputValues,
00221                                                              const ArrayScalar &    inputPoints,
00222                                                              const ArrayScalar &    cellVertices,
00223                                                              const EOperator        operatorType) const {
00224   TEST_FOR_EXCEPTION( (true), std::logic_error,
00225                       ">>> ERROR (Basis_HGRAD_WEDGE_C1_FEM): FEM Basis calling an FVD member function");
00226 }
00227 }// namespace Intrepid