Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherSolution_BlockedEpetra_Hessian_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef __Panzer_GatherSolution_BlockedEpetra_Hessian_impl_hpp__
12 #define __Panzer_GatherSolution_BlockedEpetra_Hessian_impl_hpp__
13 
14 // Only do this if required by the user.
15 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
16 
18 //
19 // Include Files
20 //
22 
23 // Panzer
29 #include "Panzer_PureBasis.hpp"
30 #include "Panzer_GlobalIndexer.hpp"
31 
32 // Phalanx
33 #include "Phalanx_DataLayout.hpp"
34 
35 // Teuchos
36 #include "Teuchos_Assert.hpp"
37 #include "Teuchos_FancyOStream.hpp"
38 
39 // Thyra
40 #include "Thyra_ProductVectorBase.hpp"
41 #include "Thyra_SpmdVectorBase.hpp"
42 
44 //
45 // Initializing Constructor (Hessian Specialization)
46 //
48 template<typename TRAITS, typename LO, typename GO>
51  const std::vector<Teuchos::RCP<const GlobalIndexer<LO, int>>>&
52  indexers,
53  const Teuchos::ParameterList& p)
54  :
55  indexers_(indexers)
56 {
57  using panzer::PureBasis;
58  using PHX::MDField;
59  using PHX::print;
60  using std::size_t;
61  using std::string;
62  using std::vector;
63  using Teuchos::RCP;
65  input.setParameterList(p);
66  const vector<string>& names = input.getDofNames();
67  RCP<const PureBasis> basis = input.getBasis();
68  indexerNames_ = input.getIndexerNames();
69  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
70  globalDataKey_ = input.getGlobalDataKey();
71  gatherSeedIndex_ = input.getGatherSeedIndex();
72  sensitivitiesName_ = input.getSensitivitiesName();
73  firstSensitivitiesAvailable_ = input.firstSensitivitiesAvailable();
74  secondSensitivitiesAvailable_ = input.secondSensitivitiesAvailable();
75  sensitivities2ndPrefix_ = input.getSecondSensitivityDataKeyPrefix();
76 
77  // Allocate the fields.
78  int numFields(names.size());
79  gatherFields_.resize(numFields);
80  for (int fd(0); fd < numFields; ++fd)
81  {
82  gatherFields_[fd] =
83  MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
84  this->addEvaluatedField(gatherFields_[fd]);
85  } // end loop over names
86 
87  // Figure out what the first active name is.
88  string firstName("<none>"), n("GatherSolution (BlockedEpetra");
89  if (numFields > 0)
90  firstName = names[0];
91  if (not firstSensitivitiesAvailable_)
92  n += ", No First Sensitivities";
93  n += "): " + firstName + " (" + print<EvalT>() + ")";
94  this->setName(n);
95 } // end of Initializing Constructor (Hessian Specialization)
96 
98 //
99 // postRegistrationSetup() (Hessian Specialization)
100 //
102 template<typename TRAITS, typename LO, typename GO>
103 void
106  typename TRAITS::SetupData /* d */,
107  PHX::FieldManager<TRAITS>& /* fm */)
108 {
109  using std::size_t;
110  using std::string;
111  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
112  int numFields(gatherFields_.size());
113  indexerIds_.resize(numFields);
114  subFieldIds_.resize(numFields);
115  for (int fd(0); fd < numFields; ++fd)
116  {
117  // Get the field ID from the DOF manager.
118  const string& fieldName(indexerNames_[fd]);
119  indexerIds_[fd] = getFieldBlock(fieldName, indexers_);
120  subFieldIds_[fd] = indexers_[indexerIds_[fd]]->getFieldNum(fieldName);
121  TEUCHOS_ASSERT(indexerIds_[fd] >= 0);
122  } // end loop over gatherFields_
123  indexerNames_.clear();
124 } // end of postRegistrationSetup() (Hessian Specialization)
125 
127 //
128 // preEvaluate() (Hessian Specialization)
129 //
131 template<typename TRAITS, typename LO, typename GO>
132 void
135  typename TRAITS::PreEvalData d)
136 {
137  using std::endl;
138  using std::logic_error;
139  using std::string;
140  using Teuchos::RCP;
141  using Teuchos::rcp_dynamic_cast;
142  using Teuchos::typeName;
144  using BELOC = BlockedEpetraLinearObjContainer;
146 
147  // Manage sensitivities.
148  firstApplySensitivities_ = false;
149  if ((firstSensitivitiesAvailable_ ) and
150  (d.first_sensitivities_name == sensitivitiesName_))
151  firstApplySensitivities_ = true;
152  secondApplySensitivities_ = false;
153  if ((secondSensitivitiesAvailable_ ) and
154  (d.second_sensitivities_name == sensitivitiesName_))
155  secondApplySensitivities_ = true;
156 
157  // First try the refactored ReadOnly container.
159  string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
160  if (d.gedc->containsDataObject(globalDataKey_ + post))
161  {
162  ged = d.gedc->getDataObject(globalDataKey_ + post);
163  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
164  }
165 
166  // Otherwise, try the old path.
167  else if (d.gedc->containsDataObject(globalDataKey_))
168  {
169  ged = d.gedc->getDataObject(globalDataKey_);
170 
171  // Try to extract the linear object container.
172  auto roGed = rcp_dynamic_cast<const BVROGED>(ged);
173  auto blockedContainer = rcp_dynamic_cast<const BELOC>(ged);
174  if (not roGed.is_null())
175  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
176  else if (not blockedContainer.is_null())
177  {
178  if (useTimeDerivativeSolutionVector_)
179  x_ = rcp_dynamic_cast<ProductVectorBase<double>>
180  (blockedContainer->get_dxdt());
181  else // if (not useTimeDerivativeSolutionVector_)
182  x_ = rcp_dynamic_cast<ProductVectorBase<double>>
183  (blockedContainer->get_x());
184  } // end if roGed or blockedContainer is non-null
185  } // end if we're doing things the new or old way
186 
187  // Ensure that we actually have something.
188  TEUCHOS_ASSERT((not x_.is_null()) or (not xBvRoGed_.is_null()))
189 
190  // Don't try to extract dx if it's not required.
191  if (not secondApplySensitivities_)
192  return;
193 
194  // Now parse the second derivative direction.
195  if (d.gedc->containsDataObject(sensitivities2ndPrefix_ + globalDataKey_))
196  {
197  ged = d.gedc->getDataObject(sensitivities2ndPrefix_ + globalDataKey_);
198  dxBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
199  } // end if (d.gedc->containsDataObject(...))
200 
201  // Ensure that we actually have something.
202  TEUCHOS_TEST_FOR_EXCEPTION(dxBvRoGed_.is_null(), logic_error,
203  "Cannot find sensitivity vector associated with \"" +
204  sensitivities2ndPrefix_ + globalDataKey_ + "\" and \"" + post + "\".");
205 } // end of preEvaluate() (Hessian Specialization)
206 
208 //
209 // evaluateFields() (Hessian Specialization)
210 //
212 template<typename TRAITS, typename LO, typename GO>
213 void
216  typename TRAITS::EvalData workset)
217 {
218  using PHX::MDField;
219  using std::size_t;
220  using std::string;
221  using std::vector;
222  using Teuchos::ArrayRCP;
223  using Teuchos::ptrFromRef;
224  using Teuchos::RCP;
225  using Teuchos::rcp_dynamic_cast;
227  using Thyra::SpmdVectorBase;
228  using Thyra::VectorBase;
229 
230  // For convenience, pull out some objects from the workset.
231  string blockId(this->wda(workset).block_id);
232  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
233  int numFields(gatherFields_.size()), numCells(localCellIds.size());
234  double seedValue(0);
235  if (firstApplySensitivities_)
236  {
237  if ((useTimeDerivativeSolutionVector_) and (gatherSeedIndex_ < 0))
238  seedValue = workset.alpha;
239  else if (gatherSeedIndex_ < 0)
240  seedValue = workset.beta;
241  else if (not useTimeDerivativeSolutionVector_)
242  seedValue = workset.gather_seeds[gatherSeedIndex_];
243  else
244  TEUCHOS_ASSERT(false);
245  } // end if (firstApplySensitivities_)
246 
247  // Turn off sensitivies. This may be faster if we don't expand the term, but
248  // I suspect not, because anywhere it is used the full complement of
249  // sensitivies will be needed anyway.
250  if (not firstApplySensitivities_)
251  seedValue = 0.0;
252 
253  vector<int> blockOffsets;
254  computeBlockOffsets(blockId, indexers_, blockOffsets);
255  if (x_.is_null())
256  {
257  // Loop over the fields to be gathered.
258  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
259  {
260  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
261  int indexerId(indexerIds_[fieldInd]),
262  subFieldNum(subFieldIds_[fieldInd]);
263 
264  // Grab the local data for inputing.
265  auto xEvRoGed = xBvRoGed_->getGEDBlock(indexerId);
266  auto subRowIndexer = indexers_[indexerId];
267  const vector<int>& elmtOffset =
268  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
269  int numBases(elmtOffset.size());
270 
271  // Gather operation for each cell in the workset.
272  for (int cell(0); cell < numCells; ++cell)
273  {
274  size_t cellLocalId(localCellIds[cell]);
275  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
276 
277  // Loop over the basis functions and fill the fields.
278  for (int basis(0); basis < numBases; ++basis)
279  {
280  int offset(elmtOffset[basis]), lid(LIDs[offset]);
281  field(cell, basis) = (*xEvRoGed)[lid];
282  } // end loop over the basis functions
283  } // end loop over localCellIds
284  } // end loop over the fields to be gathered
285  }
286  else // if (not x_.is_null())
287  {
288  // Loop over the fields to be gathered.
289  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
290  {
291  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
292  int indexerId(indexerIds_[fieldInd]),
293  subFieldNum(subFieldIds_[fieldInd]);
294 
295  // Grab the local data for inputing.
297  rcp_dynamic_cast<SpmdVectorBase<double>>(x_->
298  getNonconstVectorBlock(indexerId))->getLocalData(ptrFromRef(x));
299  auto subRowIndexer = indexers_[indexerId];
300  const vector<int>& elmtOffset =
301  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
302  int numBases(elmtOffset.size());
303 
304  // Gather operation for each cell in the workset.
305  for (int cell(0); cell < numCells; ++cell)
306  {
307  size_t cellLocalId(localCellIds[cell]);
308  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
309 
310  // Loop over the basis functions and fill the fields.
311  for (int basis(0); basis < numBases; ++basis)
312  {
313  int offset(elmtOffset[basis]), lid(LIDs[offset]);
314  field(cell, basis) = x[lid];
315  } // end loop over the basis functions
316  } // end loop over localCellIds
317  } // end loop over the fields to be gathered
318  } // end if (x_.is_null()) or not
319 
320  // Deal with the first sensitivities.
321  if (firstApplySensitivities_)
322  {
323  // Loop over the fields to be gathered.
324  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
325  {
326  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
327  int indexerId(indexerIds_[fieldInd]),
328  subFieldNum(subFieldIds_[fieldInd]);
329  auto subRowIndexer = indexers_[indexerId];
330  const vector<int>& elmtOffset =
331  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
332  int startBlkOffset(blockOffsets[indexerId]),
333  numBases(elmtOffset.size());
334 
335  // Gather operation for each cell in the workset.
336  for (int cell(0); cell < numCells; ++cell)
337  {
338  // Loop over the basis functions and fill the fields.
339  for (int basis(0); basis < numBases; ++basis)
340  {
341  int offset(elmtOffset[basis]);
342  field(cell, basis).fastAccessDx(startBlkOffset + offset) = seedValue;
343  } // end loop over the basis functions
344  } // end loop over localCellIds
345  } // end loop over the fields to be gathered
346  } // end if (firstApplySensitivities_)
347 
348  // Deal with the second sensitivities.
349  if (secondApplySensitivities_)
350  {
351  // Loop over the fields to be gathered.
352  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
353  {
354  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
355  int indexerId(indexerIds_[fieldInd]),
356  subFieldNum(subFieldIds_[fieldInd]);
357 
358  // Grab the local data for inputing.
359  auto dxEvRoGed = dxBvRoGed_->getGEDBlock(indexerId);
360  auto subRowIndexer = indexers_[indexerId];
361  const vector<int>& elmtOffset =
362  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
363  int numBases(elmtOffset.size());
364 
365  // Gather operation for each cell in the workset.
366  for (int cell(0); cell < numCells; ++cell)
367  {
368  size_t cellLocalId(localCellIds[cell]);
369  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
370 
371  // Loop over the basis functions and fill the fields.
372  for (int basis(0); basis < numBases; ++basis)
373  {
374  int offset(elmtOffset[basis]), lid(LIDs[offset]);
375  field(cell, basis).val().fastAccessDx(0) = (*dxEvRoGed)[lid];
376  } // end loop over the basis functions
377  } // end loop over localCellIds
378  } // end loop over the fields to be gathered
379  } // end if (secondApplySensitivities_)
380 } // end of evaluateFields() (Hessian Specialization)
381 
382 #endif // Panzer_BUILD_HESSIAN_SUPPORT
383 
384 #endif // __Panzer_GatherSolution_BlockedEpetra_Hessian_impl_hpp__
const std::vector< std::string > & getDofNames() const
The names of the DOFs to be gathered (all types)
void setParameterList(const Teuchos::ParameterList &pl)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::string getGlobalDataKey() const
Name of the global evaluation data container to use for the source vector (all types) ...
void computeBlockOffsets(const std::string &blockId, const std::vector< Teuchos::RCP< GlobalIndexer >> &ugis, std::vector< int > &blockOffsets)
int numFields
Teuchos::RCP< const PureBasis > getBasis() const
Basis definiting layout of dof names (all types)
This class encapsulates the needs of a gather operation to do a halo exchange for blocked vectors...
int getGatherSeedIndex() const
What index to use for initializing the seed (Jacobian and Hessian)
int getFieldBlock(const std::string &fieldName, const std::vector< Teuchos::RCP< const GlobalIndexer >> &ugis)
bool firstSensitivitiesAvailable()
Are first derivative sensitivities enabled or disabled? (Jacobian and Hessian)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
Post-Registration Setup.
bool useTimeDerivativeSolutionVector() const
Gather a time derivative vector? (all types)
GatherSolution_BlockedEpetra(const Teuchos::ParameterList &p)
Constructor.
virtual void preEvaluate(typename Traits::PreEvalData d)=0
bool secondSensitivitiesAvailable()
Are second derivative sensitivies enabled or disabled (Hessian only)
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
const std::vector< std::string > & getIndexerNames() const
void evaluateFields(typename TRAITS::EvalData d)
Evaluate Fields.
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)
std::string getSecondSensitivityDataKeyPrefix()
What prefix to use for the GEDC for second derivative sensitivity direction (Hessian only) ...
std::string getSensitivitiesName() const
The name of the sensitivities. Enables sensitivities at &quot;preEvaluate&quot; time (Jacobian and Hessian) ...
std::string typeName(const T &t)