Panzer  Version of the Day
Panzer_GatherSolution_Epetra_Hessian_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_GatherSolution_Epetra_Hessian_impl_hpp__
44 #define __Panzer_GatherSolution_Epetra_Hessian_impl_hpp__
45 
46 // only do this if required by the user
47 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
48 
49 #include "Teuchos_Assert.hpp"
50 #include "Phalanx_DataLayout.hpp"
51 
53 #include "Panzer_PureBasis.hpp"
58 
59 #include "Teuchos_FancyOStream.hpp"
60 
61 #include "Epetra_Vector.h"
62 #include "Epetra_Map.h"
63 
64 namespace panzer {
65 
66 // **********************************************************************
67 // Specialization: Hessian
68 // **********************************************************************
69 
70 template<typename TRAITS,typename LO,typename GO>
73  const Teuchos::RCP<const panzer::UniqueGlobalIndexer<LO,GO> > & indexer,
74  const Teuchos::ParameterList& p)
75  : globalIndexer_(indexer)
76 {
77  typedef std::vector< std::vector<std::string> > vvstring;
78 
79  GatherSolution_Input input;
80  input.setParameterList(p);
81 
82  const std::vector<std::string> & names = input.getDofNames();
83  Teuchos::RCP<const panzer::PureBasis> basis = input.getBasis();
84 
85  indexerNames_ = input.getIndexerNames();
86  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
87  globalDataKey_ = input.getGlobalDataKey();
88 
89  gatherSeedIndex_ = input.getGatherSeedIndex();
90  sensitivitiesName_ = input.getSensitivitiesName();
91  firstSensitivitiesAvailable_ = input.firstSensitivitiesAvailable();
92 
93  secondSensitivitiesAvailable_ = input.secondSensitivitiesAvailable();
94  sensitivities2ndPrefix_ = input.getSecondSensitivityDataKeyPrefix();
95 
96  // allocate fields
97  gatherFields_.resize(names.size());
98  for (std::size_t fd = 0; fd < names.size(); ++fd) {
99  PHX::MDField<ScalarT,Cell,NODE> f(names[fd],basis->functional);
100  gatherFields_[fd] = f;
101  this->addEvaluatedField(gatherFields_[fd]);
102  }
103 
104  // figure out what the first active name is
105  std::string firstName = "<none>";
106  if(names.size()>0)
107  firstName = names[0];
108 
109  // print out convenience
110  if(!firstSensitivitiesAvailable_) {
111  std::string n = "GatherSolution (Epetra, No First Sensitivities): "+firstName+" (Hessian)";
112  this->setName(n);
113  }
114  else {
115  std::string n = "GatherSolution (Epetra): "+firstName+" ("+PHX::typeAsString<EvalT>()+") ";
116  this->setName(n);
117  }
118 }
119 
120 // **********************************************************************
121 template<typename TRAITS,typename LO,typename GO>
123 postRegistrationSetup(typename TRAITS::SetupData d,
125 {
126  // globalIndexer_ = d.globalIndexer_;
127  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
128 
129  fieldIds_.resize(gatherFields_.size());
130 
131  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
132  // get field ID from DOF manager
133  const std::string& fieldName = indexerNames_[fd];
134  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
135 
136  // this is the error return code, raise the alarm
137  if(fieldIds_[fd]==-1) {
138  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
139  "GatherSolution_Epetra<Hessian>: Could not find field \"" + fieldName + "\" in the global indexer. ");
140  // wouldn't it be nice to print more information???
141  }
142 
143  // setup the field data object
144  this->utils.setFieldData(gatherFields_[fd],fm);
145  }
146 
147  indexerNames_.clear(); // Don't need this anymore
148 }
149 
150 // **********************************************************************
151 template<typename TRAITS,typename LO,typename GO>
153 preEvaluate(typename TRAITS::PreEvalData d)
154 {
155  using Teuchos::RCP;
156  using Teuchos::rcp_dynamic_cast;
157 
158  // manage sensitivities
160  if(firstSensitivitiesAvailable_) {
161  if(d.first_sensitivities_name==sensitivitiesName_)
162  firstApplySensitivities_ = true;
163  else
164  firstApplySensitivities_ = false;
165  }
166  else
167  firstApplySensitivities_ = false;
168 
169  if(secondSensitivitiesAvailable_) {
170  if(d.second_sensitivities_name==sensitivitiesName_)
171  secondApplySensitivities_ = true;
172  else
173  secondApplySensitivities_ = false;
174  }
175  else
176  secondApplySensitivities_ = false;
177 
179 
180  RCP<GlobalEvaluationData> ged;
181 
182  // first try refactored ReadOnly container
183  std::string post = useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X";
184  if(d.gedc.containsDataObject(globalDataKey_+post)) {
185  ged = d.gedc.getDataObject(globalDataKey_+post);
186 
187  RCP<EpetraVector_ReadOnly_GlobalEvaluationData> ro_ged = rcp_dynamic_cast<EpetraVector_ReadOnly_GlobalEvaluationData>(ged,true);
188 
189  x_ = ro_ged->getGhostedVector_Epetra();
190  }
191  else if(d.gedc.containsDataObject(globalDataKey_)) {
192 
193  ged = d.gedc.getDataObject(globalDataKey_);
194 
195  // try to extract linear object container
196  RCP<EpetraVector_ReadOnly_GlobalEvaluationData> ro_ged = rcp_dynamic_cast<EpetraVector_ReadOnly_GlobalEvaluationData>(ged);
197  RCP<EpetraLinearObjContainer> epetraContainer = rcp_dynamic_cast<EpetraLinearObjContainer>(ged);
198 
199  // handle LOCPair case
200  {
201  RCP<LOCPair_GlobalEvaluationData> loc_pair = rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(ged);
202 
203  if(loc_pair!=Teuchos::null) {
204  Teuchos::RCP<LinearObjContainer> loc = loc_pair->getGhostedLOC();
205  // extract linear object container
206  epetraContainer = rcp_dynamic_cast<EpetraLinearObjContainer>(loc);
207  }
208  }
209 
210  if(ro_ged!=Teuchos::null) {
211  RCP<EpetraVector_ReadOnly_GlobalEvaluationData> ro_ged = rcp_dynamic_cast<EpetraVector_ReadOnly_GlobalEvaluationData>(ged,true);
212 
213  x_ = ro_ged->getGhostedVector_Epetra();
214  }
215  else if(epetraContainer!=Teuchos::null) {
216  if (useTimeDerivativeSolutionVector_)
217  x_ = epetraContainer->get_dxdt();
218  else
219  x_ = epetraContainer->get_x();
220  }
221  } // end "else if" contains(globalDataKey)
222 
223  // post condition
224  TEUCHOS_ASSERT(x_!=Teuchos::null); // someone has to find the x_ vector
225 
226  // don't try to extract the dx, if its not required
227  if(!secondApplySensitivities_) {
228  dx_ = Teuchos::null; // do set it to null though!
229  return;
230  }
231 
232  // get second derivative perturbation
234 
235  // now parse the second derivative direction
236  if(d.gedc.containsDataObject(sensitivities2ndPrefix_+globalDataKey_)) {
237  ged = d.gedc.getDataObject(sensitivities2ndPrefix_+globalDataKey_);
238 
239  RCP<EpetraVector_ReadOnly_GlobalEvaluationData> ro_ged = rcp_dynamic_cast<EpetraVector_ReadOnly_GlobalEvaluationData>(ged,true);
240 
241  dx_ = ro_ged->getGhostedVector_Epetra();
242  }
243 
244  // post conditions
245  TEUCHOS_TEST_FOR_EXCEPTION(dx_==Teuchos::null,std::logic_error,
246  "Cannot find sensitivity vector associated with \""+sensitivities2ndPrefix_+globalDataKey_+"\" and \""+post+"\""); // someone has to find the dx_ vector
247 }
248 
249 // **********************************************************************
250 template<typename TRAITS,typename LO,typename GO>
252 evaluateFields(typename TRAITS::EvalData workset)
253 {
254  // for convenience pull out some objects from workset
255  std::string blockId = this->wda(workset).block_id;
256  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
257 
258  double seed_value = 0.0;
259  if(firstApplySensitivities_) {
260  if (useTimeDerivativeSolutionVector_ && gatherSeedIndex_<0) {
261  seed_value = workset.alpha;
262  }
263  else if (gatherSeedIndex_<0) {
264  seed_value = workset.beta;
265  }
266  else if(!useTimeDerivativeSolutionVector_) {
267  seed_value = workset.gather_seeds[gatherSeedIndex_];
268  }
269  else {
270  TEUCHOS_ASSERT(false);
271  }
272  }
273 
274  Epetra_Vector & x = *x_;
275 
276  // turn off sensitivies: this may be faster if we don't expand the term
277  // but I suspect not because anywhere it is used the full complement of
278  // sensitivies will be needed anyway.
279  if(!firstApplySensitivities_)
280  seed_value = 0.0;
281 
282  // Interface worksets handle DOFs from two element blocks. The derivative
283  // offset for the other element block must be shifted by the derivative side
284  // of my element block.
285  int dos = 0;
286  if (this->wda.getDetailsIndex() == 1) {
287  // Get the DOF count for my element block.
288  dos = globalIndexer_->getElementBlockGIDCount(workset.details(0).block_id);
289  }
290 
291  // loop over the fields to be gathered
292  for(std::size_t fieldIndex=0;
293  fieldIndex<gatherFields_.size();fieldIndex++) {
294  PHX::MDField<ScalarT,Cell,NODE> & field = gatherFields_[fieldIndex];
295  int fieldNum = fieldIds_[fieldIndex];
296  const std::vector<int> & elmtOffset = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
297 
298  // gather operation for each cell in workset
299  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
300  std::size_t cellLocalId = localCellIds[worksetCellIndex];
301 
302  const std::vector<int> & LIDs = globalIndexer_->getElementLIDs(cellLocalId);
303 
304  if(!firstApplySensitivities_) {
305  // loop over basis functions and fill the fields
306  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
307  int offset = elmtOffset[basis];
308  int lid = LIDs[offset];
309 
310  // set the value and seed the FAD object
311  field(worksetCellIndex,basis) = x[lid];
312  }
313  }
314  else {
315  // loop over basis functions and fill the fields
316  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
317  int offset = elmtOffset[basis];
318  int lid = LIDs[offset];
319 
320  // set the value and seed the FAD object
321  field(worksetCellIndex,basis).val() = x[lid];
322  field(worksetCellIndex,basis).fastAccessDx(dos + offset) = seed_value;
323  }
324  }
325 
326  // this is the direction to use for the second derivative
327  if(secondApplySensitivities_) {
328  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
329  int offset = elmtOffset[basis];
330  int lid = LIDs[offset];
331  field(worksetCellIndex,basis).val().fastAccessDx(0) = (*dx_)[lid];
332  }
333  }
334  }
335  }
336 }
337 
338 // **********************************************************************
339 
340 } // end namespace panzer
341 
342 #endif
343 
344 #endif
Gathers solution values from the Newton solution vector into the nodal fields of the field manager...
PHX::MDField< const ScalarT > input
PHX::MDField< const ScalarT, Cell, IP > field
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.