Panzer  Version of the Day
Panzer_STK_ScatterCellAvgVector_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_AVG_VECTOR_IMPL_HPP
44 #define PANZER_STK_SCATTER_CELL_AVG_VECTOR_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(ScatterCellAvgVector,p) :
63  mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
64 {
65  using panzer::Cell;
66  using panzer::Point;
67  using panzer::Dim;
68 
69  std::string scatterName = p.get<std::string>("Scatter Name");
70 
71  const std::vector<std::string> & names =
72  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
73 
74  Teuchos::RCP<panzer::IntegrationRule> intRule =
75  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR");
76 
77  // build dependent fields
78  scatterFields_.resize(names.size());
79  stkFields_.resize(names.size());
80  for (std::size_t fd = 0; fd < names.size(); ++fd)
81  {
82  scatterFields_[fd] = PHX::MDField<const ScalarT,Cell,Point,Dim>(names[fd],intRule->dl_vector);
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 Vectors");
91 }
92 
93 
94 PHX_POST_REGISTRATION_SETUP(ScatterCellAvgVector,d,fm)
95 {
96  for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd)
97  {
98  std::string fieldName = scatterFields_[fd].fieldTag().name();
99 
100  stkFields_[fd] = mesh_->getMetaData()->get_field<VariableField>(stk::topology::ELEMENT_RANK, fieldName);
101 
102  // setup the field data object
103  this->utils.setFieldData(scatterFields_[fd],fm);
104  }
105 }
106 
107 
108 PHX_EVALUATE_FIELDS(ScatterCellAvgVector,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  std::string d_mod[3] = {"X","Y","Z"};
116 
117  // loop over the number of vector fields requested for exodus output
118  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++)
119  {
120  PHX::MDField<const ScalarT,panzer::Cell,panzer::Point,panzer::Dim> & field = scatterFields_[fieldIndex];
121  std::string fieldName = field.fieldTag().name();
122  int numCells = field.dimension(0);
123  int numPoints = field.dimension(1);
124  int numDims = field.dimension(2);
125 
126  for (int dim = 0; dim < numDims; dim++)
127  {
128  // std::vector<double> average(numCells,0.0);
129  PHX::MDField<double,panzer::Cell,panzer::NODE> average = af.buildStaticArray<double,panzer::Cell,panzer::NODE>("",numCells,1);
130 
131  // write to double field
132  for(int i = 0; i < numCells; i++) // loop over cells
133  {
134  average(i,0) = 0.0;
135  for(int j = 0; j < numPoints; j++) // loop over IPs
136  average(i,0) += Sacado::ScalarValue<ScalarT>::eval(field(i,j,dim));
137 
138  average(i,0) /= numPoints;
139  }
140 
141  mesh_->setCellFieldData(fieldName+d_mod[dim],blockId,localCellIds,average);
142 
143  }
144  }
145 
146 }
147 
148 } // end panzer_stk
149 
150 #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::MDField< ScalarT, Cell > average
int numPoints
PHX_EVALUATE_FIELDS(ScatterCellAvgQuantity, workset)
Teuchos::RCP< STK_Interface > mesh_
PHX_EVALUATOR_CTOR(ScatterCellAvgQuantity, p)
PHX::MDField< const ScalarT, Cell, IP > field
panzer_stk::STK_Interface::SolutionFieldType VariableField
std::vector< VariableField * > stkFields_