Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_ScatterCellQuantity_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef PANZER_STK_SCATTER_CELL_QUANTITY_IMPL_HPP
12 #define PANZER_STK_SCATTER_CELL_QUANTITY_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 
16 #include "Phalanx_config.hpp"
17 #include "Phalanx_Evaluator_Macros.hpp"
18 #include "Phalanx_MDField.hpp"
19 #include "Phalanx_DataLayout.hpp"
20 #include "Phalanx_DataLayout_MDALayout.hpp"
21 
24 
25 #include "Teuchos_FancyOStream.hpp"
26 #include "Teuchos_ArrayRCP.hpp"
27 
28 namespace panzer_stk {
29 
30 template<typename EvalT, typename Traits>
33  const Teuchos::ParameterList& p) :
34  mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
35 {
36  using panzer::Cell;
37  using panzer::Point;
38 
39  std::string scatterName = p.get<std::string>("Scatter Name");
40  int worksetSize = p.get<int>("Workset Size");
41 
42  // if it's there pull the ouptut scaling
43  if (p.isParameter("Variable Scale Factors Map"))
44  varScaleFactors_ = p.get<Teuchos::RCP<std::map<std::string,double>>>("Variable Scale Factors Map");
45 
46  const std::vector<std::string> & names =
47  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
48 
50 
51  // build dependent fields
52  scatterFields_.resize(names.size());
53  for (std::size_t fd = 0; fd < names.size(); ++fd) {
54  scatterFields_[fd] = PHX::MDField<const ScalarT,Cell>(names[fd],dl_cell);
55  this->addDependentField(scatterFields_[fd]);
56  }
57 
58  // setup a dummy field to evaluate
59  PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
60  this->addEvaluatedField(scatterHolder);
61 
62  this->setName(scatterName+": STK-Scatter Cell Quantity Fields");
63 }
64 
65 template<typename EvalT, typename Traits>
66 void
69  typename Traits::SetupData /* d */,
70  PHX::FieldManager<Traits>& /* fm */)
71 {
72  for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd)
73  std::string fieldName = scatterFields_[fd].fieldTag().name();
74 }
75 
76 template<typename EvalT, typename Traits>
77 void
80  typename Traits::EvalData workset)
81 {
82  panzer::MDFieldArrayFactory af("",true);
83 
84  // for convenience pull out some objects from workset
85  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
86  std::string blockId = this->wda(workset).block_id;
87 
88  for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
89  PHX::MDField<const ScalarT,panzer::Cell> & field = scatterFields_[fieldIndex];
90  // std::vector<double> value(field.extent(0),0.0);
92 
93  // write to double field
94  for(unsigned i=0; i<field.extent(0);i++)
95  value(i,0) = Sacado::scalarValue(field(i));
96 
97  std::string varname = field.fieldTag().name();
98  double scalef = 1.0;
99 
100  if (!varScaleFactors_.is_null())
101  {
102  std::map<std::string,double> *tmp_sfs = varScaleFactors_.get();
103  if(tmp_sfs->find(varname) != tmp_sfs->end())
104  scalef = (*tmp_sfs)[varname];
105  }
106 
107  mesh_->setCellFieldData(field.fieldTag().name(),blockId,localCellIds,value.get_view(),scalef);
108  }
109 }
110 
111 } // end panzer_stk
112 
113 #endif
ScatterCellQuantity(const Teuchos::ParameterList &p)
std::vector< PHX::MDField< const ScalarT, panzer::Cell > > scatterFields_
T & get(const std::string &name, T def_value)
Teuchos::RCP< std::map< std::string, double > > varScaleFactors_
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void evaluateFields(typename Traits::EvalData d)
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)