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 FELIX {
00013
00014
00015 template<typename EvalT, typename Traits>
00016 StokesContravarientMetricTensor<EvalT, Traits>::
00017 StokesContravarientMetricTensor(const Teuchos::ParameterList& p,
00018 const Teuchos::RCP<Albany::Layouts>& dl) :
00019 coordVec (p.get<std::string> ("Coordinate Vector Name"), dl->vertices_vector),
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"), dl->qp_tensor)
00023 {
00024 this->addDependentField(coordVec);
00025 this->addEvaluatedField(Gc);
00026
00027
00028 std::vector<PHX::DataLayout::size_type> dim;
00029 dl->qp_gradient->dimensions(dim);
00030 int containerSize = dim[0];
00031 numQPs = dim[1];
00032 numDims = dim[2];
00033
00034
00035 refPoints.resize(numQPs, numDims);
00036 refWeights.resize(numQPs);
00037 jacobian.resize(containerSize, numQPs, numDims, numDims);
00038 jacobian_inv.resize(containerSize, numQPs, numDims, numDims);
00039
00040
00041 cubature->getCubature(refPoints, refWeights);
00042
00043 this->setName("StokesContravarientMetricTensor"+PHX::TypeString<EvalT>::value);
00044 }
00045
00046
00047 template<typename EvalT, typename Traits>
00048 void StokesContravarientMetricTensor<EvalT, Traits>::
00049 postRegistrationSetup(typename Traits::SetupData d,
00050 PHX::FieldManager<Traits>& fm)
00051 {
00052 this->utils.setFieldData(coordVec,fm);
00053 this->utils.setFieldData(Gc,fm);
00054 }
00055
00056
00057 template<typename EvalT, typename Traits>
00058 void StokesContravarientMetricTensor<EvalT, Traits>::
00059 evaluateFields(typename Traits::EvalData workset)
00060 {
00061
00070
00071
00072
00073 Intrepid::CellTools<RealType>::setJacobian(jacobian, refPoints, coordVec, *cellType);
00074 Intrepid::CellTools<MeshScalarT>::setJacobianInv(jacobian_inv, jacobian);
00075
00076 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00077 for (std::size_t qp=0; qp < numQPs; ++qp) {
00078 for (std::size_t i=0; i < numDims; ++i) {
00079 for (std::size_t j=0; j < numDims; ++j) {
00080 Gc(cell,qp,i,j) = 0.0;
00081 for (std::size_t alpha=0; alpha < numDims; ++alpha) {
00082 Gc(cell,qp,i,j) += jacobian_inv(cell,qp,alpha,i)*jacobian_inv(cell,qp,alpha,j);
00083 }
00084 }
00085 }
00086 }
00087 }
00088
00089 }
00090
00091
00092 }