Panzer  Version of the Day
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 PHX_EVALUATOR_CTOR(ScatterCellQuantity,p) :
63  mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
64 {
65  using panzer::Cell;
66  using panzer::Point;
67 
68  std::string scatterName = p.get<std::string>("Scatter Name");
69  int worksetSize = p.get<int>("Workset Size");
70 
71  const std::vector<std::string> & names =
72  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
73 
74  Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(new PHX::MDALayout<Cell>(worksetSize));
75 
76  // build dependent fields
77  scatterFields_.resize(names.size());
78  for (std::size_t fd = 0; fd < names.size(); ++fd) {
79  scatterFields_[fd] = PHX::MDField<const ScalarT,Cell>(names[fd],dl_cell);
80  this->addDependentField(scatterFields_[fd]);
81  }
82 
83  // setup a dummy field to evaluate
84  PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
85  this->addEvaluatedField(scatterHolder);
86 
87  this->setName(scatterName+": STK-Scatter Cell Quantity Fields");
88 }
89 
90 PHX_POST_REGISTRATION_SETUP(ScatterCellQuantity,d,fm)
91 {
92  for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd) {
93  std::string fieldName = scatterFields_[fd].fieldTag().name();
94 
95  // setup the field data object
96  this->utils.setFieldData(scatterFields_[fd],fm);
97  }
98 }
99 
100 PHX_EVALUATE_FIELDS(ScatterCellQuantity,workset)
101 {
102  panzer::MDFieldArrayFactory af("",true);
103 
104  // for convenience pull out some objects from workset
105  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
106  std::string blockId = this->wda(workset).block_id;
107 
108  for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
109  PHX::MDField<const ScalarT,panzer::Cell> & field = scatterFields_[fieldIndex];
110  // std::vector<double> value(field.dimension(0),0.0);
111  PHX::MDField<double,panzer::Cell,panzer::NODE> value = af.buildStaticArray<double,panzer::Cell,panzer::NODE>("",field.dimension(0),1);
112 
113  // write to double field
114  for(unsigned i=0; i<field.dimension(0);i++)
115  value(i,0) = Sacado::ScalarValue<ScalarT>::eval(field(i));
116 
117  mesh_->setCellFieldData(field.fieldTag().name(),blockId,localCellIds,value);
118  }
119 }
120 
121 } // end panzer_stk
122 
123 #endif
PHX_POST_REGISTRATION_SETUP(ScatterCellAvgQuantity, d, fm)
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::Point > > scatterFields_
PHX_EVALUATE_FIELDS(ScatterCellAvgQuantity, workset)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< STK_Interface > mesh_
PHX_EVALUATOR_CTOR(ScatterCellAvgQuantity, p)
PHX::MDField< const ScalarT, Cell, IP > field