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

QCAD_PoissonResid_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 
00013 //**********************************************************************
00014 template<typename EvalT, typename Traits>
00015 QCAD::PoissonResid<EvalT, Traits>::
00016 PoissonResid(const Teuchos::ParameterList& p,
00017                  const Teuchos::RCP<Albany::Layouts>& dl) :
00018   wBF         (p.get<std::string>  ("Weighted BF Name"), dl->node_qp_scalar),
00019   Potential   (p.get<std::string>  ("QP Variable Name"), dl->qp_scalar),
00020   Permittivity (p.get<std::string>  ("Permittivity Name"), dl->qp_scalar),
00021   wGradBF     (p.get<std::string>  ("Weighted Gradient BF Name"), dl->node_qp_gradient),
00022   PhiGrad     (p.get<std::string>  ("Gradient QP Variable Name"), dl->qp_gradient),
00023   PhiFlux     (p.get<std::string>  ("Flux QP Variable Name"), dl->qp_gradient),
00024   Source      (p.get<std::string>  ("Source Name"), dl->qp_scalar ),
00025   haveSource  (p.get<bool>("Have Source")),
00026   PhiResidual (p.get<std::string>  ("Residual Name"),  dl->node_scalar )
00027 {
00028   this->addDependentField(wBF);
00029   //this->addDependentField(Potential);
00030   this->addDependentField(Permittivity);
00031   this->addDependentField(PhiGrad);
00032   this->addDependentField(wGradBF);
00033   if (haveSource) this->addDependentField(Source);
00034 
00035   this->addEvaluatedField(PhiResidual);
00036   this->addEvaluatedField(PhiFlux);
00037 
00038   this->setName("PoissonResid"+PHX::TypeString<EvalT>::value);
00039 }
00040 
00041 //**********************************************************************
00042 template<typename EvalT, typename Traits>
00043 void QCAD::PoissonResid<EvalT, Traits>::
00044 postRegistrationSetup(typename Traits::SetupData d,
00045                       PHX::FieldManager<Traits>& fm)
00046 {
00047   this->utils.setFieldData(wBF,fm);
00048   //this->utils.setFieldData(Potential,fm);
00049   this->utils.setFieldData(Permittivity,fm);
00050   this->utils.setFieldData(PhiGrad,fm);
00051   this->utils.setFieldData(PhiFlux,fm);
00052   this->utils.setFieldData(wGradBF,fm);
00053   if (haveSource)  this->utils.setFieldData(Source,fm);
00054 
00055   this->utils.setFieldData(PhiResidual,fm);
00056 }
00057 
00058 //**********************************************************************
00059 template<typename EvalT, typename Traits>
00060 void QCAD::PoissonResid<EvalT, Traits>::
00061 evaluateFields(typename Traits::EvalData workset)
00062 {
00063   typedef Intrepid::FunctionSpaceTools FST;
00064 
00065   // Scale gradient into a flux, reusing same memory
00066   FST::scalarMultiplyDataData<ScalarT> (PhiFlux, Permittivity, PhiGrad);
00067 
00068   FST::integrate<ScalarT>(PhiResidual, PhiFlux, wGradBF, Intrepid::COMP_CPP, false); // "false" overwrites
00069 
00070   if (haveSource) {
00071     for (int i=0; i<Source.size(); i++) Source[i] *= -1.0;
00072     FST::integrate<ScalarT>(PhiResidual, Source, wBF, Intrepid::COMP_CPP, true); // "true" sums into
00073   }
00074 }
00075 //**********************************************************************
00076 

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