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

CalcInstantaneousCoords_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 CalcInstantaneousCoords<EvalT, Traits>::
00017 CalcInstantaneousCoords(const Teuchos::ParameterList& p, const Teuchos::RCP<Albany::Layouts>& dl) :
00018 
00019   coordVec      (p.get<std::string> ("Coordinate Vector Name"), dl->vertices_vector),
00020   dispVec       (p.get<std::string> ("Solution Vector Name"), dl->node_vector),
00021   instCoords    (p.get<std::string> ("Instantaneous Coordinates Name"), dl->node_vector)
00022 
00023 {
00024 
00025   this->addDependentField(dispVec);
00026   this->addDependentField(coordVec);
00027   this->addEvaluatedField(instCoords);
00028 
00029   // Get Dimensions
00030   std::vector<PHX::DataLayout::size_type> dim;
00031   dl->node_qp_vector->dimensions(dim);
00032   int containerSize = dim[0];
00033   numNodes = dim[1];
00034   numQPs = dim[2];
00035   numDims = dim[3];
00036 
00037   this->setName("CalcInstantaneousCoords"+PHX::TypeString<EvalT>::value);
00038 }
00039 
00040 //**********************************************************************
00041 template<typename EvalT, typename Traits>
00042 void CalcInstantaneousCoords<EvalT, Traits>::
00043 postRegistrationSetup(typename Traits::SetupData d,
00044                       PHX::FieldManager<Traits>& fm)
00045 {
00046   this->utils.setFieldData(dispVec, fm);
00047   this->utils.setFieldData(coordVec, fm);
00048   this->utils.setFieldData(instCoords, fm);
00049 }
00050 
00051 //**********************************************************************
00052 template<typename EvalT, typename Traits>
00053 void CalcInstantaneousCoords<EvalT, Traits>::
00054 evaluateFields(typename Traits::EvalData workset)
00055 {
00056 
00057   for(std::size_t cell = 0; cell < workset.numCells; ++cell) {
00058     for(std::size_t node = 0; node < numNodes; ++node) {
00059       for(std::size_t eq = 0; eq < numDims; eq++)  {
00060 
00061           instCoords(cell, node, eq) = coordVec(cell, node, eq) + dispVec(cell, node, eq);
00062 
00063       }
00064     }
00065   }
00066 }
00067 
00068 //**********************************************************************
00069 }

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