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 ContravariantTargetMetricTensor<EvalT, Traits>::
00017 ContravariantTargetMetricTensor(const Teuchos::ParameterList& p, const Teuchos::RCP<Albany::Layouts>& dl) :
00018
00019 solnVec (p.get<std::string> ("Solution Vector Name"), dl->node_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> ("Contravariant Metric Tensor Name"), dl->qp_tensor)
00023
00024 {
00025
00026 this->addDependentField(solnVec);
00027 this->addEvaluatedField(Gc);
00028
00029
00030 std::vector<PHX::DataLayout::size_type> dim;
00031 dl->node_vector->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("ContravariantTargetMetricTensor"+PHX::TypeString<EvalT>::value);
00046 }
00047
00048
00049 template<typename EvalT, typename Traits>
00050 void ContravariantTargetMetricTensor<EvalT, Traits>::
00051 postRegistrationSetup(typename Traits::SetupData d,
00052 PHX::FieldManager<Traits>& fm)
00053 {
00054 this->utils.setFieldData(solnVec, fm);
00055 this->utils.setFieldData(Gc, fm);
00056 }
00057
00058
00059 template<typename EvalT, typename Traits>
00060 void ContravariantTargetMetricTensor<EvalT, Traits>::
00061 evaluateFields(typename Traits::EvalData workset)
00062 {
00063
00064
00065 Intrepid::CellTools<ScalarT>::setJacobian(jacobian, refPoints, solnVec, *cellType);
00066 Intrepid::CellTools<ScalarT>::setJacobianInv(jacobian_inv, jacobian);
00067
00068 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00069 for (std::size_t qp=0; qp < numQPs; ++qp) {
00070 for (std::size_t i=0; i < numDims; ++i) {
00071 for (std::size_t j=0; j < numDims; ++j) {
00072 Gc(cell, qp, i, j) = 0.0;
00073 for (std::size_t alpha=0; alpha < numDims; ++alpha) {
00074 Gc(cell, qp, i, j) += jacobian_inv(cell, qp, alpha, i) * jacobian_inv(cell, qp, alpha, j);
00075 }
00076 }
00077 }
00078 }
00079 }
00080
00081 }
00082
00083
00084 }