Panzer  Version of the Day
Panzer_STK_ScatterVectorFields_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_VECTOR_FIELDS_IMPL_HPP
44 #define PANZER_STK_SCATTER_VECTOR_FIELDS_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 
54 #include "Panzer_Traits.hpp"
56 
57 #include "Teuchos_FancyOStream.hpp"
58 
59 namespace panzer_stk {
60 
62  mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
63 {
64  TEUCHOS_ASSERT(false);
65 }
66 
67 template <typename EvalT,typename TraitsT>
69 ScatterVectorFields(const std::string & scatterName,
70  const Teuchos::RCP<STK_Interface> mesh,
71  const Teuchos::RCP<const panzer::PointRule> & pointRule,
72  const std::vector<std::string> & names,
73  const std::vector<double> & scaling)
74  : mesh_(mesh)
75  , scaling_(scaling)
76 {
77  using panzer::Cell;
78  using panzer::IP;
79  using panzer::Dim;
80 
81  spatialDimension_ = pointRule->spatial_dimension;
82 
83  // this evaluator assumes you are evaluating at the cell centroid only
84  TEUCHOS_ASSERT(pointRule->num_points==1);
85 
86  // build dependent fields
87  names_ = names;
88  scatterFields_.resize(names.size());
89  for (std::size_t fd = 0; fd < names.size(); ++fd) {
90  scatterFields_[fd] =
91  PHX::MDField<const ScalarT,Cell,IP,Dim>(names_[fd]+"_"+pointRule->getName(),pointRule->dl_vector);
92  this->addDependentField(scatterFields_[fd]);
93  }
94 
95  // setup a dummy field to evaluate
96  PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
97  this->addEvaluatedField(scatterHolder);
98 
99  this->setName(scatterName+": STK-Scatter Vector Fields");
100 }
101 
103 {
104  // this->utils.setFieldData(pointField_,fm);
105 
106  for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd) {
107  // setup the field data object
108  this->utils.setFieldData(scatterFields_[fd],fm);
109  }
110 }
111 
113 {
114  TEUCHOS_ASSERT(false);
115 }
116 
117 template < >
118 void ScatterVectorFields<panzer::Traits::Residual,panzer::Traits>::
119 evaluateFields(panzer::Traits::EvalData workset)
120 {
121  panzer::MDFieldArrayFactory af("",true);
122 
123  std::vector<std::string> dimStrings(3);
124  dimStrings[0] = "X";
125  dimStrings[1] = "Y";
126  dimStrings[2] = "Z";
127 
128  // for convenience pull out some objects from workset
129  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
130  std::string blockId = this->wda(workset).block_id;
131 
132  for(int d=0;d<spatialDimension_;d++) {
133  for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
134  PHX::MDField<const ScalarT,panzer::Cell,panzer::IP,panzer::Dim> & field = scatterFields_[fieldIndex];
135  std::string fieldName = names_[fieldIndex]+dimStrings[d];
136 
137  PHX::MDField<double,panzer::Cell,panzer::NODE> cellValue
138  = af.buildStaticArray<double,panzer::Cell,panzer::NODE>("",field.dimension(0),1);
139 
140  // scaline field value only if the scaling parameter is specified, otherwise use 1.0
141  double scaling = (scaling_.size()>0) ? scaling_[fieldIndex] : 1.0;
142 
143  for(unsigned i=0; i<field.dimension(0);i++)
144  cellValue(i,0) = field(i,0,d);
145 
146  // add in vector value at d^th dimension
147  mesh_->setCellFieldData(fieldName,blockId,localCellIds,cellValue,scaling);
148  }
149  }
150 }
151 
152 } // end panzer_stk
153 
154 #endif
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::Point > > scatterFields_
std::vector< std::string > names_
Teuchos::RCP< const panzer::PointRule > pointRule
Teuchos::RCP< STK_Interface > mesh_
PHX_EVALUATOR_CTOR(ScatterCellAvgQuantity, p)
PHX_EVALUATE_FIELDS(ScatterVectorFields, workset)
PHX::MDField< const ScalarT, Cell, IP > field
PHX_POST_REGISTRATION_SETUP(ScatterVectorFields, d, fm)
ScatterVectorFields(const std::string &scatterName, const Teuchos::RCP< STK_Interface > mesh, const Teuchos::RCP< const panzer::PointRule > &pointRule, const std::vector< std::string > &names, const std::vector< double > &scaling=std::vector< double >())
std::vector< double > scaling_