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

DamageResid_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 DamageResid<EvalT, Traits>::
00017 DamageResid(const Teuchos::ParameterList& p) :
00018   wBF         (p.get<std::string>                   ("Weighted BF Name"),
00019          p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Scalar Data Layout") ),
00020   damage      (p.get<std::string>                   ("QP Variable Name"),
00021          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00022   damage_dot  (p.get<std::string>                   ("QP Time Derivative Variable Name"),
00023          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00024   damageLS    (p.get<std::string>                   ("Damage Length Scale Name"),
00025          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00026   wGradBF     (p.get<std::string>                   ("Weighted Gradient BF Name"),
00027          p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout") ),
00028   damage_grad (p.get<std::string>                   ("Gradient QP Variable Name"),
00029          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout") ),
00030   source      (p.get<std::string>                   ("Damage Source Name"),
00031          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00032   gc          (p.get<double>("gc Name")),
00033   dResidual   (p.get<std::string>                   ("Residual Name"),
00034          p.get<Teuchos::RCP<PHX::DataLayout> >("Node Scalar Data Layout") )
00035 {
00036   if (p.isType<bool>("Disable Transient"))
00037     enableTransient = !p.get<bool>("Disable Transient");
00038   else enableTransient = true;
00039 
00040   this->addDependentField(wBF);
00041   this->addDependentField(damage);
00042   this->addDependentField(damageLS);
00043   if (enableTransient) this->addDependentField(damage_dot);
00044   this->addDependentField(damage_grad);
00045   this->addDependentField(wGradBF);
00046   this->addDependentField(source);
00047 
00048   this->addEvaluatedField(dResidual);
00049 
00050   Teuchos::RCP<PHX::DataLayout> vector_dl =
00051     p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00052   std::vector<PHX::DataLayout::size_type> dims;
00053   vector_dl->dimensions(dims);
00054   numQPs  = dims[1];
00055   numDims = dims[2];
00056 
00057   // Allocate workspace
00058   flux.resize(dims[0], numQPs, numDims);
00059 
00060   this->setName("DamageResid"+PHX::TypeString<EvalT>::value);
00061 }
00062 
00063 //**********************************************************************
00064 template<typename EvalT, typename Traits>
00065 void DamageResid<EvalT, Traits>::
00066 postRegistrationSetup(typename Traits::SetupData d,
00067                       PHX::FieldManager<Traits>& fm)
00068 {
00069   this->utils.setFieldData(wBF,fm);
00070   this->utils.setFieldData(damage,fm);
00071   this->utils.setFieldData(damageLS,fm);
00072   this->utils.setFieldData(damage_grad,fm);
00073   this->utils.setFieldData(wGradBF,fm);
00074   this->utils.setFieldData(source,fm);
00075   if (enableTransient) this->utils.setFieldData(damage_dot,fm);
00076 
00077   this->utils.setFieldData(dResidual,fm);
00078 }
00079 
00080 //**********************************************************************
00081 template<typename EvalT, typename Traits>
00082 void DamageResid<EvalT, Traits>::
00083 evaluateFields(typename Traits::EvalData workset)
00084 {
00085   typedef Intrepid::FunctionSpaceTools FST;
00086   typedef Intrepid::RealSpaceTools<ScalarT> RST;
00087 
00088   FST::scalarMultiplyDataData<ScalarT> (flux, damageLS, damage_grad);
00089   RST::scale(flux,-gc);
00090 
00091   FST::integrate<ScalarT>(dResidual, flux, wGradBF, Intrepid::COMP_CPP, false); // "false" overwrites
00092 
00093   //for (int i=0; i < source.size(); i++) source[i] *= -1.0;
00094   FST::integrate<ScalarT>(dResidual, source, wBF, Intrepid::COMP_CPP, true); // "true" sums into
00095   
00096   if (workset.transientTerms && enableTransient) 
00097     FST::integrate<ScalarT>(dResidual, damage_dot, wBF, Intrepid::COMP_CPP, true); // "true" sums into
00098 }
00099 
00100 //**********************************************************************
00101 }
00102 

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