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 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_STK_SCATTER_CELL_QUANTITY_IMPL_HPP
44 #define PANZER_STK_SCATTER_CELL_QUANTITY_IMPL_HPP
45 
46 #include "Teuchos_Assert.hpp"
47 
48 #include "Phalanx_config.hpp"
49 #include "Phalanx_Evaluator_Macros.hpp"
50 #include "Phalanx_MDField.hpp"
51 #include "Phalanx_DataLayout.hpp"
52 #include "Phalanx_DataLayout_MDALayout.hpp"
53 
56 
57 #include "Teuchos_FancyOStream.hpp"
58 #include "Teuchos_ArrayRCP.hpp"
59 
60 namespace panzer_stk {
61 
62 template<typename EvalT, typename Traits>
65  const Teuchos::ParameterList& p) :
66  mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
67 {
68  using panzer::Cell;
69  using panzer::Point;
70 
71  std::string scatterName = p.get<std::string>("Scatter Name");
72  int worksetSize = p.get<int>("Workset Size");
73 
74  const std::vector<std::string> & names =
75  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
76 
77  Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(new PHX::MDALayout<Cell>(worksetSize));
78 
79  // build dependent fields
80  scatterFields_.resize(names.size());
81  for (std::size_t fd = 0; fd < names.size(); ++fd) {
82  scatterFields_[fd] = PHX::MDField<const ScalarT,Cell>(names[fd],dl_cell);
83  this->addDependentField(scatterFields_[fd]);
84  }
85 
86  // setup a dummy field to evaluate
87  PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
88  this->addEvaluatedField(scatterHolder);
89 
90  this->setName(scatterName+": STK-Scatter Cell Quantity Fields");
91 }
92 
93 template<typename EvalT, typename Traits>
94 void
97  typename Traits::SetupData /* d */,
98  PHX::FieldManager<Traits>& /* fm */)
99 {
100  for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd)
101  std::string fieldName = scatterFields_[fd].fieldTag().name();
102 }
103 
104 template<typename EvalT, typename Traits>
105 void
108  typename Traits::EvalData workset)
109 {
110  panzer::MDFieldArrayFactory af("",true);
111 
112  // for convenience pull out some objects from workset
113  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
114  std::string blockId = this->wda(workset).block_id;
115 
116  for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
117  PHX::MDField<const ScalarT,panzer::Cell> & field = scatterFields_[fieldIndex];
118  // std::vector<double> value(field.extent(0),0.0);
119  PHX::MDField<double,panzer::Cell,panzer::NODE> value = af.buildStaticArray<double,panzer::Cell,panzer::NODE>("",field.extent(0),1);
120 
121  // write to double field
122  for(unsigned i=0; i<field.extent(0);i++)
123  value(i,0) = Sacado::ScalarValue<ScalarT>::eval(field(i));
124 
125  mesh_->setCellFieldData(field.fieldTag().name(),blockId,localCellIds,value);
126  }
127 }
128 
129 } // end panzer_stk
130 
131 #endif
ScatterCellQuantity(const Teuchos::ParameterList &p)
std::vector< PHX::MDField< const ScalarT, panzer::Cell > > scatterFields_
T & get(const std::string &name, T def_value)
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) 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...
BASIS NODE
Spatial Dimension Tag.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)