00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00036 #ifndef INTREPID_CELLTOOLS_HPP
00037 #define INTREPID_CELLTOOLS_HPP
00038
00039
00040 #include "Intrepid_FieldContainer.hpp"
00041 #include "Intrepid_RealSpaceTools.hpp"
00042 #include "Intrepid_ConfigDefs.hpp"
00043 #include "Intrepid_Types.hpp"
00044 #include "Intrepid_Utils.hpp"
00045 #include "Intrepid_Basis.hpp"
00046 #include "Intrepid_HGRAD_TRI_C1_FEM.hpp"
00047 #include "Intrepid_HGRAD_QUAD_C1_FEM.hpp"
00048 #include "Intrepid_HGRAD_TET_C1_FEM.hpp"
00049 #include "Intrepid_HGRAD_WEDGE_C1_FEM.hpp"
00050 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
00051 #include "Intrepid_HGRAD_LINE_C1_FEM.hpp"
00052
00053 #include "Intrepid_HGRAD_TRI_C2_FEM.hpp"
00054 #include "Intrepid_HGRAD_QUAD_C2_FEM.hpp"
00055 #include "Intrepid_HGRAD_TET_C2_FEM.hpp"
00056 #include "Intrepid_HGRAD_WEDGE_C2_FEM.hpp"
00057 #include "Intrepid_HGRAD_HEX_C2_FEM.hpp"
00058
00059 #include "Shards_CellTopology.hpp"
00060 #include "Shards_BasicTopologies.hpp"
00061
00062 #include "Teuchos_TestForException.hpp"
00063 #include "Teuchos_RCP.hpp"
00064
00065 namespace Intrepid {
00066
00067
00068
00069
00070
00071
00072
00073
00082 template<class Scalar>
00083 class CellTools {
00084 private:
00085
00086
00087
00088
00089
00090
00091
00092
00105 static const FieldContainer<double>& getSubcellParametrization(const int subcellDim,
00106 const shards::CellTopology& parentCell);
00107
00108
00109
00149 static void setSubcellParametrization(FieldContainer<double>& subcellParametrization,
00150 const int subcellDim,
00151 const shards::CellTopology& parentCell);
00152
00153
00154
00155
00156
00157
00158
00166 template<class ArrayScalar>
00167 static void validateArguments_setJacobian(const ArrayScalar & jacobian,
00168 const ArrayScalar & points,
00169 const ArrayScalar & cellWorkset,
00170 const int & whichCell,
00171 const shards::CellTopology & cellTopo);
00172
00173
00174
00179 template<class ArrayScalar>
00180 static void validateArguments_setJacobianInv(const ArrayScalar & jacobianInv,
00181 const ArrayScalar & jacobian);
00182
00183
00184
00189 template<class ArrayScalar>
00190 static void validateArguments_setJacobianDetArgs(const ArrayScalar & jacobianDet,
00191 const ArrayScalar & jacobian);
00192
00193
00194
00202 template<class ArrayScalar>
00203 static void validateArguments_mapToPhysicalFrame(const ArrayScalar & physPoints,
00204 const ArrayScalar & refPoints,
00205 const ArrayScalar & cellWorkset,
00206 const shards::CellTopology & cellTopo,
00207 const int& whichCell);
00208
00209
00210
00218 template<class ArrayScalar>
00219 static void validateArguments_mapToReferenceFrame(const ArrayScalar & refPoints,
00220 const ArrayScalar & physPoints,
00221 const ArrayScalar & cellWorkset,
00222 const shards::CellTopology & cellTopo,
00223 const int& whichCell);
00224
00225
00226
00235 template<class ArrayType1, class ArrayType2>
00236 static void validateArguments_mapToReferenceFrame(const ArrayType1 & refPoints,
00237 const ArrayType2 & initGuess,
00238 const ArrayType1 & physPoints,
00239 const ArrayType1 & cellWorkset,
00240 const shards::CellTopology & cellTopo,
00241 const int& whichCell);
00242
00243
00244
00252 template<class ArrayInt, class ArrayPoint, class ArrayScalar>
00253 static void validateArguments_checkPointwiseInclusion(ArrayInt & inCell,
00254 const ArrayPoint & physPoints,
00255 const ArrayScalar & cellWorkset,
00256 const int & whichCell,
00257 const shards::CellTopology & cell);
00258 public:
00259
00262 CellTools(){ };
00263
00264
00267 ~CellTools(){ };
00268
00269
00270
00271
00272
00273
00274
00322 template<class ArrayScalar>
00323 static void setJacobian(ArrayScalar & jacobian,
00324 const ArrayScalar & points,
00325 const ArrayScalar & cellWorkset,
00326 const shards::CellTopology & cellTopo,
00327 const int & whichCell = -1);
00328
00329
00330
00343 template<class ArrayScalar>
00344 static void setJacobianInv(ArrayScalar & jacobianInv,
00345 const ArrayScalar & jacobian);
00346
00347
00348
00361 template<class ArrayScalar>
00362 static void setJacobianDet(ArrayScalar & jacobianDet,
00363 const ArrayScalar & jacobian);
00364
00365
00366
00367
00368
00369
00370
00426 template<class ArrayScalar>
00427 static void mapToPhysicalFrame(ArrayScalar & physPoints,
00428 const ArrayScalar & refPoints,
00429 const ArrayScalar & cellWorkset,
00430 const shards::CellTopology & cellTopo,
00431 const int & whichCell = -1);
00432
00433
00434
00493 template<class ArrayScalar>
00494 static void mapToReferenceFrame(ArrayScalar & refPoints,
00495 const ArrayScalar & physPoints,
00496 const ArrayScalar & cellWorkset,
00497 const shards::CellTopology & cellTopo,
00498 const int & whichCell = -1);
00499
00500
00501
00547 template<class ArrayType1, class ArrayType2>
00548 static void mapToReferenceFrame(ArrayType1 & refPoints,
00549 const ArrayType2 & initGuess,
00550 const ArrayType1 & physPoints,
00551 const ArrayType1 & cellWorkset,
00552 const shards::CellTopology & cellTopo,
00553 const int & whichCell = -1);
00554
00555
00556
00607 template<class ArrayTypeOut, class ArrayTypeIn>
00608 static void mapToReferenceSubcell(ArrayTypeOut & refSubcellPoints,
00609 const ArrayTypeIn & paramPoints,
00610 const int subcellDim,
00611 const int subcellOrd,
00612 const shards::CellTopology & parentCell);
00613
00614
00615
00641 template<class ArrayTypeOut>
00642 static void getReferenceEdgeTangent(ArrayTypeOut & refEdgeTangent,
00643 const int & edgeOrd,
00644 const shards::CellTopology & parentCell);
00645
00646
00647
00684 template<class ArrayTypeOut>
00685 static void getReferenceFaceTangents(ArrayTypeOut & refFaceTanU,
00686 ArrayTypeOut & refFaceTanV,
00687 const int & faceOrd,
00688 const shards::CellTopology & parentCell);
00689
00690
00691
00754 template<class ArrayTypeOut>
00755 static void getReferenceSideNormal(ArrayTypeOut & refSideNormal,
00756 const int & sideOrd,
00757 const shards::CellTopology & parentCell);
00758
00759
00760
00799 template<class ArrayTypeOut>
00800 static void getReferenceFaceNormal(ArrayTypeOut & refFaceNormal,
00801 const int & faceOrd,
00802 const shards::CellTopology & parentCell);
00803
00804
00805
00835 template<class ArrayTypeOut, class ArrayTypeIn>
00836 static void getPhysicalEdgeTangents(ArrayTypeOut & edgeTangents,
00837 const ArrayTypeIn & worksetJacobians,
00838 const int & worksetEdgeOrd,
00839 const shards::CellTopology & parentCell);
00840
00841
00842
00882 template<class ArrayTypeOut, class ArrayTypeIn>
00883 static void getPhysicalFaceTangents(ArrayTypeOut & faceTanU,
00884 ArrayTypeOut & faceTanV,
00885 const ArrayTypeIn & worksetJacobians,
00886 const int & worksetFaceOrd,
00887 const shards::CellTopology & parentCell);
00888
00889
00890
00951 template<class ArrayTypeOut, class ArrayTypeIn>
00952 static void getPhysicalSideNormals(ArrayTypeOut & sideNormals,
00953 const ArrayTypeIn & worksetJacobians,
00954 const int & worksetSideOrd,
00955 const shards::CellTopology & parentCell);
00956
00957
00958
00997 template<class ArrayTypeOut, class ArrayTypeIn>
00998 static void getPhysicalFaceNormals(ArrayTypeOut & faceNormals,
00999 const ArrayTypeIn & worksetJacobians,
01000 const int & worksetFaceOrd,
01001 const shards::CellTopology & parentCell);
01002
01003
01004
01005
01006
01007
01008
01009
01020 static int checkPointInclusion(const Scalar* point,
01021 const int pointDim,
01022 const shards::CellTopology & cellTopo,
01023 const double & threshold = INTREPID_THRESHOLD);
01024
01025
01026
01039 template<class ArrayPoint>
01040 static int checkPointsetInclusion(const ArrayPoint& points,
01041 const shards::CellTopology & cellTopo,
01042 const double & threshold = INTREPID_THRESHOLD);
01043
01044
01045
01072 template<class ArrayInt, class ArrayPoint>
01073 static void checkPointwiseInclusion(ArrayInt & inRefCell,
01074 const ArrayPoint & points,
01075 const shards::CellTopology & cellTopo,
01076 const double & threshold = INTREPID_THRESHOLD);
01077
01078
01079
01115 template<class ArrayInt, class ArrayPoint, class ArrayScalar>
01116 static void checkPointwiseInclusion(ArrayInt & inCell,
01117 const ArrayPoint & points,
01118 const ArrayScalar & cellWorkset,
01119 const shards::CellTopology & cell,
01120 const int & whichCell = -1,
01121 const double & threshold = INTREPID_THRESHOLD);
01122
01123
01124
01135 static const double* getReferenceVertex(const shards::CellTopology& cell,
01136 const int vertexOrd);
01137
01138
01139
01154 template<class ArrayOut>
01155 static void getReferenceSubcellVertices(ArrayOut& subcellVertices,
01156 const int subcellDim,
01157 const int subcellOrd,
01158 const shards::CellTopology& parentCell);
01159
01160
01161
01177 static const double* getReferenceNode(const shards::CellTopology& cell,
01178 const int nodeOrd);
01179
01180
01181
01195 template<class ArrayOut>
01196 static void getReferenceSubcellNodes(ArrayOut& subcellNodes,
01197 const int subcellDim,
01198 const int subcellOrd,
01199 const shards::CellTopology& parentCell);
01200
01201
01202
01208 static int hasReferenceCell(const shards::CellTopology & cellTopo);
01209
01210
01211
01212
01213
01214
01215
01216
01217
01218
01224 static void printSubcellVertices(const int subcellDim,
01225 const int subcellOrd,
01226 const shards::CellTopology & parentCell);
01227
01228
01229
01233 template<class ArrayTypeIn>
01234 static void printWorksetSubcell(const ArrayTypeIn& cellWorkset,
01235 const shards::CellTopology& parentCell,
01236 const int& pCellOrd,
01237 const int& subcellDim,
01238 const int& subcellOrd,
01239 const int& fieldWidth = 3);
01240
01241 };
01242
01243 }
01244
01245
01246 #include "Intrepid_CellToolsDef.hpp"
01247
01248 #endif
01249
01250
01251
01252
01253
01254
01255