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

PHAL_ScatterScalarResponse_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 // **********************************************************************
00011 // Base Class Generic Implemtation
00012 // **********************************************************************
00013 namespace PHAL {
00014 
00015 template<typename EvalT, typename Traits>
00016 ScatterScalarResponseBase<EvalT, Traits>::
00017 ScatterScalarResponseBase(const Teuchos::ParameterList& p,
00018         const Teuchos::RCP<Albany::Layouts>& dl)
00019 {
00020   setup(p, dl);
00021 }
00022 
00023 template<typename EvalT, typename Traits>
00024 void ScatterScalarResponseBase<EvalT, Traits>::
00025 postRegistrationSetup(typename Traits::SetupData d,
00026                       PHX::FieldManager<Traits>& fm)
00027 {
00028   this->utils.setFieldData(global_response,fm);
00029 }
00030 
00031 template<typename EvalT, typename Traits>
00032 void
00033 ScatterScalarResponseBase<EvalT, Traits>::
00034 setup(const Teuchos::ParameterList& p, const Teuchos::RCP<Albany::Layouts>& dl)
00035 {
00036   bool stand_alone = p.get<bool>("Stand-alone Evaluator");
00037 
00038   // Setup fields we require
00039   PHX::Tag<ScalarT> global_response_tag =
00040     p.get<PHX::Tag<ScalarT> >("Global Response Field Tag");
00041   global_response = PHX::MDField<ScalarT>(global_response_tag);
00042   std::cout << "global_response layout = " << std::endl;
00043   global_response_tag.dataLayout().print(std::cout);
00044   if (stand_alone)
00045     this->addDependentField(global_response);
00046   else
00047     this->addEvaluatedField(global_response);
00048 
00049   // Setup field we evaluate
00050   std::string fieldName = global_response_tag.name() + " Scatter Response";
00051   scatter_operation = Teuchos::rcp(new PHX::Tag<ScalarT>(fieldName, dl->dummy));
00052   this->addEvaluatedField(*scatter_operation);
00053 
00055   Teuchos::ParameterList* plist = 
00056     p.get<Teuchos::ParameterList*>("Parameter List");
00057   if (stand_alone) {
00058     Teuchos::RCP<const Teuchos::ParameterList> reflist = 
00059       this->getValidResponseParameters();
00060     plist->validateParameters(*reflist,0);
00061   }
00062 
00063   if (stand_alone)
00064     this->setName(fieldName+" Scatter Response" + 
00065       PHX::TypeString<EvalT>::value);
00066 }
00067 
00068 template<typename EvalT,typename Traits>
00069 Teuchos::RCP<const Teuchos::ParameterList>
00070 ScatterScalarResponseBase<EvalT, Traits>::
00071 getValidResponseParameters() const
00072 {
00073   Teuchos::RCP<Teuchos::ParameterList> validPL =
00074     rcp(new Teuchos::ParameterList("Valid ScatterScalarResponse Params"));
00075   return validPL;
00076 }
00077 
00078 // **********************************************************************
00079 // Specialization: Residual
00080 // **********************************************************************
00081 template<typename Traits>
00082 ScatterScalarResponse<PHAL::AlbanyTraits::Residual,Traits>::
00083 ScatterScalarResponse(const Teuchos::ParameterList& p,
00084     const Teuchos::RCP<Albany::Layouts>& dl) 
00085 {
00086   this->setup(p,dl);
00087 }
00088 
00089 template<typename Traits>
00090 void ScatterScalarResponse<PHAL::AlbanyTraits::Residual, Traits>::
00091 postEvaluate(typename Traits::PostEvalData workset)
00092 {
00093   // Here we scatter the *global* response
00094   Teuchos::RCP<Epetra_Vector> g = workset.g;
00095   for (std::size_t res = 0; res < this->global_response.size(); res++) {
00096     (*g)[res] = this->global_response[res];
00097   }
00098 }
00099 
00100 // **********************************************************************
00101 // Specialization: Tangent
00102 // **********************************************************************
00103 
00104 template<typename Traits>
00105 ScatterScalarResponse<PHAL::AlbanyTraits::Tangent, Traits>::
00106 ScatterScalarResponse(const Teuchos::ParameterList& p,
00107     const Teuchos::RCP<Albany::Layouts>& dl) 
00108 {
00109   this->setup(p,dl);
00110 }
00111 
00112 template<typename Traits>
00113 void ScatterScalarResponse<PHAL::AlbanyTraits::Tangent, Traits>::
00114 postEvaluate(typename Traits::PostEvalData workset)
00115 {
00116   // Here we scatter the *global* response and tangent
00117   Teuchos::RCP<Epetra_Vector> g = workset.g;
00118   Teuchos::RCP<Epetra_MultiVector> gx = workset.dgdx;
00119   Teuchos::RCP<Epetra_MultiVector> gp = workset.dgdp;
00120   for (std::size_t res = 0; res < this->global_response.size(); res++) {
00121     ScalarT& val = this->global_response[res];
00122     if (g != Teuchos::null)
00123       (*g)[res] = val.val();
00124     if (gx != Teuchos::null)
00125       for (int col=0; col<workset.num_cols_x; col++)
00126   gx->ReplaceMyValue(res, col, val.dx(col));
00127     if (gp != Teuchos::null)
00128       for (int col=0; col<workset.num_cols_p; col++)
00129   gp->ReplaceMyValue(res, col, val.dx(col+workset.param_offset));
00130   }
00131 }
00132 
00133 // **********************************************************************
00134 // Specialization: Stochastic Galerkin Residual
00135 // **********************************************************************
00136 
00137 #ifdef ALBANY_SG_MP
00138 template<typename Traits>
00139 ScatterScalarResponse<PHAL::AlbanyTraits::SGResidual, Traits>::
00140 ScatterScalarResponse(const Teuchos::ParameterList& p,
00141     const Teuchos::RCP<Albany::Layouts>& dl) 
00142 {
00143   this->setup(p,dl);
00144 }
00145 
00146 template<typename Traits>
00147 void ScatterScalarResponse<PHAL::AlbanyTraits::SGResidual, Traits>::
00148 postEvaluate(typename Traits::PostEvalData workset)
00149 {
00150   // Here we scatter the *global* SG response
00151   Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > g_sg = workset.sg_g;
00152   for (std::size_t res = 0; res < this->global_response.size(); res++) {
00153     ScalarT& val = this->global_response[res];
00154     for (int block=0; block<g_sg->size(); block++)
00155       (*g_sg)[block][res] = val.coeff(block);
00156   }
00157 }
00158 
00159 // **********************************************************************
00160 // Specialization: Stochastic Galerkin Tangent
00161 // **********************************************************************
00162 
00163 template<typename Traits>
00164 ScatterScalarResponse<PHAL::AlbanyTraits::SGTangent, Traits>::
00165 ScatterScalarResponse(const Teuchos::ParameterList& p,
00166     const Teuchos::RCP<Albany::Layouts>& dl) 
00167 {
00168   this->setup(p,dl);
00169 }
00170 
00171 template<typename Traits>
00172 void ScatterScalarResponse<PHAL::AlbanyTraits::SGTangent, Traits>::
00173 postEvaluate(typename Traits::PostEvalData workset)
00174 {
00175   // Here we scatter the *global* SG response and tangent
00176   Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly> g_sg = workset.sg_g;
00177   Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly> gx_sg = workset.sg_dgdx;
00178   Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly> gp_sg = workset.sg_dgdp;
00179   for (std::size_t res = 0; res < this->global_response.size(); res++) {
00180     ScalarT& val = this->global_response[res];
00181     if (g_sg != Teuchos::null)
00182       for (int block=0; block<g_sg->size(); block++)
00183   (*g_sg)[block][res] = val.val().coeff(block);
00184     if (gx_sg != Teuchos::null)
00185       for (int col=0; col<workset.num_cols_x; col++)
00186   for (int block=0; block<gx_sg->size(); block++)
00187     (*gx_sg)[block].ReplaceMyValue(res, col, val.dx(col).coeff(block));
00188     if (gp_sg != Teuchos::null)
00189       for (int col=0; col<workset.num_cols_p; col++)
00190   for (int block=0; block<gp_sg->size(); block++)
00191     (*gp_sg)[block].ReplaceMyValue(
00192       res, col, val.dx(col+workset.param_offset).coeff(block));
00193   }
00194 }
00195 
00196 // **********************************************************************
00197 // Specialization: Multi-point Residual
00198 // **********************************************************************
00199 
00200 template<typename Traits>
00201 ScatterScalarResponse<PHAL::AlbanyTraits::MPResidual, Traits>::
00202 ScatterScalarResponse(const Teuchos::ParameterList& p,
00203     const Teuchos::RCP<Albany::Layouts>& dl) 
00204 {
00205   this->setup(p,dl);
00206 }
00207 
00208 template<typename Traits>
00209 void ScatterScalarResponse<PHAL::AlbanyTraits::MPResidual, Traits>::
00210 postEvaluate(typename Traits::PostEvalData workset)
00211 {
00212   // Here we scatter the *global* MP response
00213   Teuchos::RCP<Stokhos::ProductEpetraVector> g_mp = workset.mp_g;
00214   for (std::size_t res = 0; res < this->global_response.size(); res++) {
00215     ScalarT& val = this->global_response[res];
00216     for (int block=0; block<g_mp->size(); block++)
00217       (*g_mp)[block][res] = val.coeff(block);
00218   }
00219 }
00220 
00221 // **********************************************************************
00222 // Specialization: Multi-point Tangent
00223 // **********************************************************************
00224 
00225 template<typename Traits>
00226 ScatterScalarResponse<PHAL::AlbanyTraits::MPTangent, Traits>::
00227 ScatterScalarResponse(const Teuchos::ParameterList& p,
00228     const Teuchos::RCP<Albany::Layouts>& dl) 
00229 {
00230   this->setup(p,dl);
00231 }
00232 
00233 template<typename Traits>
00234 void ScatterScalarResponse<PHAL::AlbanyTraits::MPTangent, Traits>::
00235 postEvaluate(typename Traits::PostEvalData workset)
00236 {
00237   // Here we scatter the *global* MP response and tangent
00238   Teuchos::RCP<Stokhos::ProductEpetraVector> g_mp = workset.mp_g;
00239   Teuchos::RCP<Stokhos::ProductEpetraMultiVector> gx_mp = workset.mp_dgdx;
00240   Teuchos::RCP<Stokhos::ProductEpetraMultiVector> gp_mp = workset.mp_dgdp;
00241   for (std::size_t res = 0; res < this->global_response.size(); res++) {
00242     ScalarT& val = this->global_response[res];
00243     if (g_mp != Teuchos::null)
00244       for (int block=0; block<g_mp->size(); block++)
00245   (*g_mp)[block][res] = val.val().coeff(block);
00246     if (gx_mp != Teuchos::null)
00247       for (int col=0; col<workset.num_cols_x; col++)
00248   for (int block=0; block<gx_mp->size(); block++)
00249     (*gx_mp)[block].ReplaceMyValue(res, col, val.dx(col).coeff(block));
00250     if (gp_mp != Teuchos::null)
00251       for (int col=0; col<workset.num_cols_p; col++)
00252   for (int block=0; block<gp_mp->size(); block++)
00253     (*gp_mp)[block].ReplaceMyValue(
00254       res, col, val.dx(col+workset.param_offset).coeff(block));
00255   }
00256 }
00257 #endif //ALBANY_SG_MP
00258 
00259 }
00260 

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