Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ResponseScatterEvaluator_Functional.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_RESPONSE_SCATTER_EVALUATOR_FUNCTIONAL_HPP
44 #define PANZER_RESPONSE_SCATTER_EVALUATOR_FUNCTIONAL_HPP
45 
46 #include <iostream>
47 #include <string>
48 
49 #include "PanzerDiscFE_config.hpp"
50 #include "Panzer_Dimension.hpp"
51 #include "Panzer_CellData.hpp"
53 #include "Panzer_GlobalIndexer.hpp"
55 
56 #include "Phalanx_Evaluator_Macros.hpp"
57 #include "Phalanx_MDField.hpp"
58 
60 
61 namespace panzer {
62 
64 public:
66 
67  virtual void scatterDerivative(const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
68  panzer::Traits::EvalData workset,
70  const std::vector<Teuchos::ArrayRCP<double> > & dgdx) const = 0;
71 
72 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
73  virtual void scatterHessian(const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
74  panzer::Traits::EvalData workset,
76  const std::vector<Teuchos::ArrayRCP<double> > & d2gdx2) const = 0;
77 #endif
78 };
79 
80 template <typename LO,typename GO>
82 public:
84  {
85  if(globalIndexer!=Teuchos::null)
86  ugis_.push_back(globalIndexer);
87  }
88 
90  : ugis_(ugis) {}
91 
92  void scatterDerivative(const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
93  panzer::Traits::EvalData workset,
95  const std::vector<Teuchos::ArrayRCP<double> > & dgdx) const;
96 
97 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
98  void scatterHessian(const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
99  panzer::Traits::EvalData workset,
101  const std::vector<Teuchos::ArrayRCP<double> > & d2gdx2) const;
102 #endif
103 
104 private:
105 
106  std::vector<Teuchos::RCP<const panzer::GlobalIndexer> > ugis_;
107 };
108 
112 template<typename EvalT, typename Traits>
114  public PHX::EvaluatorDerived<EvalT, Traits> {
115 public:
116 
118  ResponseScatterEvaluator_Functional(const std::string & name,const CellData & cd,
119  const Teuchos::RCP<FunctionalScatterBase> & functionalScatter);
120  ResponseScatterEvaluator_Functional(const std::string & integrandName,const std::string & responseName,const CellData & cd,
121  const Teuchos::RCP<FunctionalScatterBase> & functionalScatter);
122 
123  void evaluateFields(typename Traits::EvalData d);
124 
125  void preEvaluate(typename Traits::PreEvalData d);
126 
127 private:
128  typedef typename EvalT::ScalarT ScalarT;
129 
130  std::string responseName_;
132 
134  PHX::MDField<const ScalarT,panzer::Cell> cellIntegral_; // holds cell integrals
136 };
137 
138 template <typename LO,typename GO>
139 void FunctionalScatter<LO,GO>::scatterDerivative(const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
140  panzer::Traits::EvalData workset,
142  const std::vector<Teuchos::ArrayRCP<double> > & dgdx) const
143 {
144  Kokkos::View<const LO*, PHX::Device> LIDs;
145 
146  // for convenience pull out some objects from workset
147  std::string blockId = wda(workset).block_id;
148 
149  std::vector<int> blockOffsets;
150  computeBlockOffsets(blockId,ugis_,blockOffsets);
151 
152  TEUCHOS_ASSERT(dgdx.size()==ugis_.size());
153 
154  // scatter operation for each cell in workset
155  const std::vector<std::size_t> & localCellIds = wda(workset).cell_local_ids;
156  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
157  std::size_t cellLocalId = localCellIds[worksetCellIndex];
158 
159  for(std::size_t b=0;b<ugis_.size();b++) {
160  int start = blockOffsets[b];
161 
162  LIDs = ugis_[b]->getElementLIDs(cellLocalId);
163 
164  Teuchos::ArrayRCP<double> dgdx_b = dgdx[b];
165 
166  // loop over basis functions
167  for(std::size_t i=0;i<LIDs.size();i++) {
168  dgdx_b[LIDs[i]] += cellIntegral(worksetCellIndex).dx(start+i); // its possible functional is independent of solution value!
169  }
170  }
171  }
172 }
173 
174 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
175 template <typename LO,typename GO>
176 void FunctionalScatter<LO,GO>::scatterHessian(const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
177  panzer::Traits::EvalData workset,
179  const std::vector<Teuchos::ArrayRCP<double> > & d2gdx2) const
180 {
181  Kokkos::View<const LO*, PHX::Device> LIDs;
182 
183  // for convenience pull out some objects from workset
184  std::string blockId = wda(workset).block_id;
185 
186  std::vector<int> blockOffsets;
187  computeBlockOffsets(blockId,ugis_,blockOffsets);
188 
189  TEUCHOS_ASSERT(d2gdx2.size()==ugis_.size());
190 
191  // scatter operation for each cell in workset
192  const std::vector<std::size_t> & localCellIds = wda(workset).cell_local_ids;
193  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
194  std::size_t cellLocalId = localCellIds[worksetCellIndex];
195 
196  for(std::size_t b=0;b<ugis_.size();b++) {
197  int start = blockOffsets[b];
198 
199  LIDs = ugis_[b]->getElementLIDs(cellLocalId);
200 
201  Teuchos::ArrayRCP<double> d2gdx2_b = d2gdx2[b];
202 
203  // loop over basis functions
204  for(std::size_t i=0;i<LIDs.size();i++) {
205  d2gdx2_b[LIDs[i]] += cellIntegral(worksetCellIndex).dx(start+i).dx(0); // its possible functional is independent of solution value!
206  }
207  }
208  }
209 }
210 #endif
211 
212 }
213 
214 #endif
FunctionalScatter(const std::vector< Teuchos::RCP< const panzer::GlobalIndexer > > &ugis)
FunctionalScatter(const Teuchos::RCP< const panzer::GlobalIndexer > &globalIndexer)
ResponseScatterEvaluator_Functional(const std::string &name, const CellData &cd, const Teuchos::RCP< FunctionalScatterBase > &functionalScatter)
A constructor with concrete arguments instead of a parameter list.
void scatterDerivative(const PHX::MDField< const panzer::Traits::Jacobian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &dgdx) const
virtual void scatterDerivative(const PHX::MDField< const panzer::Traits::Jacobian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &dgdx) const =0
void computeBlockOffsets(const std::string &blockId, const std::vector< Teuchos::RCP< GlobalIndexer >> &ugis, std::vector< int > &blockOffsets)
Data for determining cell topology and dimensionality.
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
void scatterHessian(const PHX::MDField< const panzer::Traits::Hessian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &d2gdx2) const
std::vector< Teuchos::RCP< const panzer::GlobalIndexer > > ugis_
#define TEUCHOS_ASSERT(assertion_test)
virtual void scatterHessian(const PHX::MDField< const panzer::Traits::Hessian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &d2gdx2) const =0