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

PHAL_DOFInterpolation_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 PHAL {
00013 
00014 //**********************************************************************
00015 template<typename EvalT, typename Traits>
00016 DOFInterpolation<EvalT, Traits>::
00017 DOFInterpolation(const Teuchos::ParameterList& p,
00018                               const Teuchos::RCP<Albany::Layouts>& dl) :
00019   val_node    (p.get<std::string>   ("Variable Name"), dl->node_scalar),
00020   BF          (p.get<std::string>   ("BF Name"), dl->node_qp_scalar),
00021   val_qp      (p.get<std::string>   ("Variable Name"), dl->qp_scalar )
00022 {
00023   this->addDependentField(val_node);
00024   this->addDependentField(BF);
00025   this->addEvaluatedField(val_qp);
00026 
00027   this->setName("DOFInterpolation"+PHX::TypeString<EvalT>::value);
00028 
00029   std::vector<PHX::DataLayout::size_type> dims;
00030   BF.fieldTag().dataLayout().dimensions(dims);
00031   numNodes = dims[1];
00032   numQPs   = dims[2];
00033 }
00034 
00035 //**********************************************************************
00036 template<typename EvalT, typename Traits>
00037 void DOFInterpolation<EvalT, Traits>::
00038 postRegistrationSetup(typename Traits::SetupData d,
00039                       PHX::FieldManager<Traits>& fm)
00040 {
00041   this->utils.setFieldData(val_node,fm);
00042   this->utils.setFieldData(BF,fm);
00043   this->utils.setFieldData(val_qp,fm);
00044 }
00045 
00046 //**********************************************************************
00047 template<typename EvalT, typename Traits>
00048 void DOFInterpolation<EvalT, Traits>::
00049 evaluateFields(typename Traits::EvalData workset)
00050 {
00051   //Intrepid version:
00052   // for (int i=0; i < val_qp.size() ; i++) val_qp[i] = 0.0;
00053   // Intrepid::FunctionSpaceTools:: evaluate<ScalarT>(val_qp, val_node, BF);
00054 
00055   for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00056     for (std::size_t qp=0; qp < numQPs; ++qp) {
00057       ScalarT& vqp = val_qp(cell,qp);
00058       vqp = val_node(cell, 0) * BF(cell, 0, qp);
00059       for (std::size_t node=1; node < numNodes; ++node) {
00060         vqp += val_node(cell, node) * BF(cell, node, qp);
00061       }
00062     }
00063   }
00064 }
00065 
00066 //**********************************************************************
00067 template<typename Traits>
00068 DOFInterpolation<PHAL::AlbanyTraits::Jacobian, Traits>::
00069 DOFInterpolation(const Teuchos::ParameterList& p,
00070                               const Teuchos::RCP<Albany::Layouts>& dl) :
00071   val_node    (p.get<std::string>   ("Variable Name"), dl->node_scalar),
00072   BF          (p.get<std::string>   ("BF Name"), dl->node_qp_scalar),
00073   val_qp      (p.get<std::string>   ("Variable Name"), dl->qp_scalar )
00074 {
00075   this->addDependentField(val_node);
00076   this->addDependentField(BF);
00077   this->addEvaluatedField(val_qp);
00078 
00079   this->setName("DOFInterpolation"+PHX::TypeString<PHAL::AlbanyTraits::Jacobian>::value);
00080 
00081   std::vector<PHX::DataLayout::size_type> dims;
00082   BF.fieldTag().dataLayout().dimensions(dims);
00083   numNodes = dims[1];
00084   numQPs   = dims[2];
00085 
00086   offset = p.get<int>("Offset of First DOF");
00087 }
00088 
00089 //**********************************************************************
00090 template<typename Traits>
00091 void DOFInterpolation<PHAL::AlbanyTraits::Jacobian, Traits>::
00092 postRegistrationSetup(typename Traits::SetupData d,
00093                       PHX::FieldManager<Traits>& fm)
00094 {
00095   this->utils.setFieldData(val_node,fm);
00096   this->utils.setFieldData(BF,fm);
00097   this->utils.setFieldData(val_qp,fm);
00098 }
00099 
00100 //**********************************************************************
00101 template<typename Traits>
00102 void DOFInterpolation<PHAL::AlbanyTraits::Jacobian, Traits>::
00103 evaluateFields(typename Traits::EvalData workset)
00104 {
00105   //Intrepid version:
00106   // for (int i=0; i < val_qp.size() ; i++) val_qp[i] = 0.0;
00107   // Intrepid::FunctionSpaceTools:: evaluate<ScalarT>(val_qp, val_node, BF);
00108 
00109   int num_dof = val_node(0,0).size();
00110   int neq = num_dof / numNodes;
00111 
00112   for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00113     for (std::size_t qp=0; qp < numQPs; ++qp) {
00114       ScalarT& vqp = val_qp(cell,qp);
00115       vqp = FadType(num_dof, val_node(cell, 0).val() * BF(cell, 0, qp));
00116       if (num_dof) vqp.fastAccessDx(offset) = val_node(cell, 0).fastAccessDx(offset) * BF(cell, 0, qp);
00117       for (std::size_t node=1; node < numNodes; ++node) {
00118         vqp.val() += val_node(cell, node).val() * BF(cell, node, qp);
00119         if (num_dof) vqp.fastAccessDx(neq*node+offset) += val_node(cell, node).fastAccessDx(neq*node+offset) * BF(cell, node, qp);
00120       }
00121     }
00122   }
00123 }
00124 
00125 //**********************************************************************
00126 
00127 }
00128 

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