• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

TPSLaplaceResid_Def.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 #include "Teuchos_TestForException.hpp"
00008 #include "Phalanx_DataLayout.hpp"
00009 
00010 #include "Sacado.hpp"
00011 #include "Sacado_Traits.hpp"
00012 
00013 namespace PHAL {
00014 
00015 //**********************************************************************
00016 template<typename EvalT, typename Traits>
00017 TPSLaplaceResid<EvalT, Traits>::
00018 TPSLaplaceResid(const Teuchos::ParameterList& p, const Teuchos::RCP<Albany::Layouts>& dl) :
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   intrepidBasis(p.get<Teuchos::RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > > > ("Intrepid Basis")),
00023   solnResidual(p.get<std::string> ("Residual Name"), dl->node_vector) {
00024 
00025 
00026   this->addDependentField(solnVec);
00027   this->addEvaluatedField(solnResidual);
00028 
00029   std::vector<PHX::DataLayout::size_type> dims;
00030   dl->node_qp_vector->dimensions(dims);
00031   worksetSize = dims[0];
00032   numNodes = dims[1];
00033   numQPs  = dims[2];
00034   numDims = dims[3];
00035 
00036   // Allocate Temporary FieldContainers
00037   grad_at_cub_points.resize(numNodes, numQPs, numDims);
00038   refPoints.resize(numQPs, numDims);
00039   refWeights.resize(numQPs);
00040   jacobian.resize(worksetSize, numQPs, numDims, numDims);
00041   jacobian_det.resize(worksetSize, numQPs);
00042 
00043   // Pre-Calculate reference element quantitites
00044   cubature->getCubature(refPoints, refWeights);
00045   intrepidBasis->getValues(grad_at_cub_points, refPoints, Intrepid::OPERATOR_GRAD);
00046 
00047   this->setName("LaplaceResid" + PHX::TypeString<EvalT>::value);
00048 
00049 }
00050 
00051 //**********************************************************************
00052 template<typename EvalT, typename Traits>
00053 void TPSLaplaceResid<EvalT, Traits>::
00054 postRegistrationSetup(typename Traits::SetupData d,
00055                       PHX::FieldManager<Traits>& fm) {
00056   this->utils.setFieldData(solnVec, fm);
00057   this->utils.setFieldData(solnResidual, fm);
00058 }
00059 
00060 //**********************************************************************
00061 template<typename EvalT, typename Traits>
00062 void TPSLaplaceResid<EvalT, Traits>::
00063 evaluateFields(typename Traits::EvalData workset) {
00064 
00065   // Note that all the integration operations are ScalarT! We need the partials of the Jacobian to show up in the
00066   // system Jacobian as the solution is a function of coordinates (it IS the coordinates!).
00067 
00068   // This adds significant time to the compile
00069 
00070   Intrepid::CellTools<ScalarT>::setJacobian(jacobian, refPoints, solnVec, *cellType);
00071   Intrepid::CellTools<ScalarT>::setJacobianDet(jacobian_det, jacobian);
00072 
00073    // Straight Laplace's equation evaluation for the nodal coord solution
00074 
00075     for(std::size_t cell = 0; cell < workset.numCells; ++cell) {
00076       for(std::size_t node_a = 0; node_a < numNodes; ++node_a) {
00077 
00078         for(std::size_t eq = 0; eq < numDims; eq++)  {
00079 
00080           solnResidual(cell, node_a, eq) = 0.0;
00081 
00082         }
00083 
00084         for(std::size_t qp = 0; qp < numQPs; ++qp) {
00085           for(std::size_t node_b = 0; node_b < numNodes; ++node_b) {
00086 
00087             ScalarT kk = 0.0;
00088 
00089             for(std::size_t i = 0; i < numDims; i++) {
00090 
00091               kk += grad_at_cub_points(node_a, qp, i) * grad_at_cub_points(node_b, qp, i);
00092 
00093             }
00094 
00095             for(std::size_t eq = 0; eq < numDims; eq++) {
00096 
00097               solnResidual(cell, node_a, eq) +=
00098                 kk * solnVec(cell, node_b, eq) * jacobian_det(cell, qp) * refWeights(qp);
00099 
00100             }
00101           }
00102         }
00103       }
00104     }
00105 }
00106 
00107 //**********************************************************************
00108 }
00109 

Generated on Wed Mar 26 2014 18:36:45 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1