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

ElasticityDispErrResid_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 "Intrepid_FunctionSpaceTools.hpp"
00011 
00012 namespace LCM {
00013 
00014 //**********************************************************************
00015 template<typename EvalT, typename Traits>
00016 ElasticityDispErrResid<EvalT, Traits>::
00017 ElasticityDispErrResid(const Teuchos::ParameterList& p) :
00018   ErrorStress  (p.get<std::string>                   ("Error Stress Name"),
00019                 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00020   wGradBF      (p.get<std::string>                   ("Weighted Gradient BF Name"),
00021                 p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout") ),
00022   DispResid    (p.get<std::string>                   ("Displacement Residual Name"),
00023                 p.get<Teuchos::RCP<PHX::DataLayout> >("Node Vector Data Layout") ),
00024   ExResidual   (p.get<std::string>                   ("Residual Name"),
00025                 p.get<Teuchos::RCP<PHX::DataLayout> >("Node Vector Data Layout") )
00026 {
00027   this->addDependentField(ErrorStress);
00028   this->addDependentField(wGradBF);
00029   this->addDependentField(DispResid);
00030 
00031   this->addEvaluatedField(ExResidual);
00032 
00033   this->setName("ElasticityDispErrResid"+PHX::TypeString<EvalT>::value);
00034 
00035   std::vector<PHX::DataLayout::size_type> dims;
00036   wGradBF.fieldTag().dataLayout().dimensions(dims);
00037   numNodes = dims[1];
00038   numQPs   = dims[2];
00039   numDims  = dims[3];
00040 }
00041 
00042 //**********************************************************************
00043 template<typename EvalT, typename Traits>
00044 void ElasticityDispErrResid<EvalT, Traits>::
00045 postRegistrationSetup(typename Traits::SetupData d,
00046                       PHX::FieldManager<Traits>& fm)
00047 {
00048   this->utils.setFieldData(ErrorStress,fm);
00049   this->utils.setFieldData(DispResid,fm);
00050   this->utils.setFieldData(wGradBF,fm);
00051   
00052   this->utils.setFieldData(ExResidual,fm);
00053 }
00054 
00055 //**********************************************************************
00056 template<typename EvalT, typename Traits>
00057 void ElasticityDispErrResid<EvalT, Traits>::
00058 evaluateFields(typename Traits::EvalData workset)
00059 {
00060   for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00061     for (std::size_t node=0; node < numNodes; ++node) {
00062       for (std::size_t dim=0; dim<numDims; dim++)  ExResidual(cell,node,dim)=0.0;
00063       for (std::size_t qp=0; qp < numQPs; ++qp) {
00064         for (std::size_t i=0; i<numDims; i++) {
00065           ExResidual(cell,node,i) += DispResid(cell,node,i);
00066           for (std::size_t dim=0; dim<numDims; dim++) {
00067             ExResidual(cell,node,i) += ErrorStress(cell, qp, i, dim) * wGradBF(cell, node, qp, dim);
00068     } } } } }
00069 }
00070 
00071 //**********************************************************************
00072 }
00073 

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