00001
00002
00003
00004
00005
00006
00007 #include "Teuchos_TestForException.hpp"
00008 #include "Phalanx_DataLayout.hpp"
00009
00010 namespace PHAL {
00011
00012
00013 template<typename EvalT, typename Traits>
00014 LaplaceBeltramiResid<EvalT, Traits>::
00015 LaplaceBeltramiResid(const Teuchos::ParameterList& p, const Teuchos::RCP<Albany::Layouts>& dl) :
00016 solnVec(p.get<std::string> ("Solution Vector Name"), dl->node_vector),
00017 Gc (p.get<std::string> ("Contravariant Metric Tensor Name"), dl->qp_tensor),
00018 cubature(p.get<Teuchos::RCP <Intrepid::Cubature<RealType> > >("Cubature")),
00019 cellType(p.get<Teuchos::RCP <shards::CellTopology> > ("Cell Type")),
00020 intrepidBasis(p.get<Teuchos::RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > > > ("Intrepid Basis")),
00021 solnResidual(p.get<std::string> ("Residual Name"), dl->node_vector) {
00022
00023 this->addDependentField(Gc);
00024 this->addDependentField(solnVec);
00025 this->addEvaluatedField(solnResidual);
00026
00027 std::vector<PHX::DataLayout::size_type> dims;
00028 dl->node_qp_vector->dimensions(dims);
00029 worksetSize = dims[0];
00030 numNodes = dims[1];
00031 numQPs = dims[2];
00032 numDims = dims[3];
00033
00034
00035 grad_at_cub_points.resize(numNodes, numQPs, numDims);
00036 refPoints.resize(numQPs, numDims);
00037 refWeights.resize(numQPs);
00038 jacobian.resize(worksetSize, numQPs, numDims, numDims);
00039 jacobian_det.resize(worksetSize, numQPs);
00040
00041
00042 cubature->getCubature(refPoints, refWeights);
00043 intrepidBasis->getValues(grad_at_cub_points, refPoints, Intrepid::OPERATOR_GRAD);
00044
00045 this->setName("LaplaceBeltramiResid" + PHX::TypeString<EvalT>::value);
00046
00047 }
00048
00049
00050 template<typename EvalT, typename Traits>
00051 void LaplaceBeltramiResid<EvalT, Traits>::
00052 postRegistrationSetup(typename Traits::SetupData d,
00053 PHX::FieldManager<Traits>& fm) {
00054
00055 this->utils.setFieldData(solnVec, fm);
00056 this->utils.setFieldData(Gc, fm);
00057
00058 this->utils.setFieldData(solnResidual, fm);
00059 }
00060
00061
00062 template<typename EvalT, typename Traits>
00063 void LaplaceBeltramiResid<EvalT, Traits>::
00064 evaluateFields(typename Traits::EvalData workset) {
00065
00066
00067 Intrepid::CellTools<ScalarT>::setJacobian(jacobian, refPoints, solnVec, *cellType);
00068 Intrepid::CellTools<ScalarT>::setJacobianDet(jacobian_det, jacobian);
00069
00070 for(std::size_t cell = 0; cell < workset.numCells; ++cell) {
00071 for(std::size_t node_a = 0; node_a < numNodes; ++node_a) {
00072
00073 for(std::size_t eq = 0; eq < numDims; eq++)
00074 solnResidual(cell, node_a, eq) = 0.0;
00075
00076 for(std::size_t qp = 0; qp < numQPs; ++qp) {
00077 for(std::size_t node_b = 0; node_b < numNodes; ++node_b) {
00078
00079 ScalarT kk = 0.0;
00080
00081 for(std::size_t i = 0; i < numDims; i++) {
00082 for(std::size_t j = 0; j < numDims; j++) {
00083
00084 kk += Gc(cell, qp, i, j) * grad_at_cub_points(node_a, qp, i) * grad_at_cub_points(node_b, qp, j);
00085
00086 }
00087 }
00088
00089 for(std::size_t eq = 0; eq < numDims; eq++) {
00090
00091 solnResidual(cell, node_a, eq) +=
00092 kk * solnVec(cell, node_b, eq) * jacobian_det(cell, qp) * refWeights(qp);
00093
00094 }
00095 }
00096 }
00097 }
00098 }
00099
00100 }
00101
00102
00103 }
00104