Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007 #include "Teuchos_TestForException.hpp"
00008 #include "Phalanx_DataLayout.hpp"
00009
00010 #include "Intrepid_FunctionSpaceTools.hpp"
00011
00012 namespace PHAL {
00013
00014
00015 template<typename EvalT, typename Traits>
00016 NSContravarientMetricTensor<EvalT, Traits>::
00017 NSContravarientMetricTensor(const Teuchos::ParameterList& p) :
00018 coordVec (p.get<std::string> ("Coordinate Vector Name"),
00019 p.get<Teuchos::RCP<PHX::DataLayout> >("Coordinate Data Layout") ),
00020 cubature (p.get<Teuchos::RCP <Intrepid::Cubature<RealType> > >("Cubature")),
00021 cellType (p.get<Teuchos::RCP <shards::CellTopology> > ("Cell Type")),
00022 Gc (p.get<std::string> ("Contravarient Metric Tensor Name"),
00023 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") )
00024 {
00025 this->addDependentField(coordVec);
00026 this->addEvaluatedField(Gc);
00027
00028
00029 Teuchos::RCP<PHX::DataLayout> vector_dl = p.get< Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout");
00030 std::vector<PHX::DataLayout::size_type> dim;
00031 vector_dl->dimensions(dim);
00032 int containerSize = dim[0];
00033 numQPs = dim[1];
00034 numDims = dim[2];
00035
00036
00037 refPoints.resize(numQPs, numDims);
00038 refWeights.resize(numQPs);
00039 jacobian.resize(containerSize, numQPs, numDims, numDims);
00040 jacobian_inv.resize(containerSize, numQPs, numDims, numDims);
00041
00042
00043 cubature->getCubature(refPoints, refWeights);
00044
00045 this->setName("NSContravarientMetricTensor"+PHX::TypeString<EvalT>::value);
00046 }
00047
00048
00049 template<typename EvalT, typename Traits>
00050 void NSContravarientMetricTensor<EvalT, Traits>::
00051 postRegistrationSetup(typename Traits::SetupData d,
00052 PHX::FieldManager<Traits>& fm)
00053 {
00054 this->utils.setFieldData(coordVec,fm);
00055 this->utils.setFieldData(Gc,fm);
00056 }
00057
00058
00059 template<typename EvalT, typename Traits>
00060 void NSContravarientMetricTensor<EvalT, Traits>::
00061 evaluateFields(typename Traits::EvalData workset)
00062 {
00063
00072
00073
00074
00075 Intrepid::CellTools<RealType>::setJacobian(jacobian, refPoints, coordVec, *cellType);
00076 Intrepid::CellTools<MeshScalarT>::setJacobianInv(jacobian_inv, jacobian);
00077
00078 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00079 for (std::size_t qp=0; qp < numQPs; ++qp) {
00080 for (std::size_t i=0; i < numDims; ++i) {
00081 for (std::size_t j=0; j < numDims; ++j) {
00082 Gc(cell,qp,i,j) = 0.0;
00083 for (std::size_t alpha=0; alpha < numDims; ++alpha) {
00084 Gc(cell,qp,i,j) += jacobian_inv(cell,qp,alpha,i)*jacobian_inv(cell,qp,alpha,j);
00085 }
00086 }
00087 }
00088 }
00089 }
00090
00091 }
00092
00093
00094 }