Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherSolution_BlockedEpetra_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_impl_hpp__
12 #define __Panzer_GatherSolution_BlockedEpetra_impl_hpp__
13 
15 //
16 // Include Files
17 //
19 
20 // Panzer
24 #include "Panzer_PureBasis.hpp"
25 #include "Panzer_GlobalIndexer.hpp"
28 
29 // Phalanx
30 #include "Phalanx_DataLayout.hpp"
31 
32 // Teuchos
33 #include "Teuchos_Assert.hpp"
34 #include "Teuchos_FancyOStream.hpp"
35 
36 // Thyra
37 #include "Thyra_ProductVectorBase.hpp"
38 #include "Thyra_SpmdVectorBase.hpp"
39 
41 //
42 // Initializing Constructor (Residual Specialization)
43 //
45 template<typename TRAITS, typename LO, typename GO>
49  const std::vector<Teuchos::RCP<const GlobalIndexer>>&
50  indexers,
51  const Teuchos::ParameterList& p)
52  :
53  indexers_(indexers),
54  hasTangentFields_(false)
55 {
56  using panzer::PureBasis;
57  using PHX::MDField;
58  using PHX::print;
59  using std::size_t;
60  using std::string;
61  using std::vector;
62  using Teuchos::RCP;
63  using vvstring = std::vector<std::vector<std::string>>;
65  input.setParameterList(p);
66  const vector<string>& names = input.getDofNames();
67  RCP<const PureBasis> basis = input.getBasis();
68  const vvstring& tangentFieldNames = input.getTangentNames();
69  indexerNames_ = input.getIndexerNames();
70  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
71  globalDataKey_ = input.getGlobalDataKey();
72 
73  // Allocate the fields.
74  int numFields(names.size());
75  gatherFields_.resize(numFields);
76  for (int fd(0); fd < numFields; ++fd)
77  {
78  gatherFields_[fd] =
79  MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
80  this->addEvaluatedField(gatherFields_[fd]);
81  } // end loop over names
82 
83  // Setup the dependent tangent fields, if requested.
84  if (tangentFieldNames.size() > 0)
85  {
86  TEUCHOS_ASSERT(gatherFields_.size() == tangentFieldNames.size());
87  hasTangentFields_ = true;
88  tangentFields_.resize(numFields);
89  for (int fd(0); fd < numFields; ++fd)
90  {
91  int numTangentFields(tangentFieldNames[fd].size());
92  tangentFields_[fd].resize(numTangentFields);
93  for (int i(0); i < numTangentFields; ++i)
94  {
95  tangentFields_[fd][i] =
96  MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
97  basis->functional);
98  this->addDependentField(tangentFields_[fd][i]);
99  } // end loop over tangentFieldNames[fd]
100  } // end loop over gatherFields_
101  } // end if we have tangent fields
102 
103  // Figure out what the first active name is.
104  string firstName("<none>");
105  if (numFields > 0)
106  firstName = names[0];
107  string n("GatherSolution (BlockedEpetra): " + firstName + " (" +
108  print<EvalT>() + ")");
109  this->setName(n);
110 } // end of Initializing Constructor (Residual Specialization)
111 
113 //
114 // postRegistrationSetup() (Residual Specialization)
115 //
117 template<typename TRAITS, typename LO, typename GO>
118 void
119 panzer::
122  typename TRAITS::SetupData /* d */,
123  PHX::FieldManager<TRAITS>& /* fm */)
124 {
125  using std::size_t;
126  using std::string;
127  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
128  int numFields(gatherFields_.size());
129  indexerIds_.resize(numFields);
130  subFieldIds_.resize(numFields);
131  for (int fd(0); fd < numFields; ++fd)
132  {
133  // Get the field ID from the DOF manager.
134  const string& fieldName(indexerNames_[fd]);
135  indexerIds_[fd] = getFieldBlock(fieldName, indexers_);
136  subFieldIds_[fd] = indexers_[indexerIds_[fd]]->getFieldNum(fieldName);
137  TEUCHOS_ASSERT(indexerIds_[fd] >= 0);
138  } // end loop over gatherFields_
139  indexerNames_.clear();
140 } // end of postRegistrationSetup() (Residual Specialization)
141 
143 //
144 // preEvaluate() (Residual Specialization)
145 //
147 template<typename TRAITS, typename LO, typename GO>
148 void
149 panzer::
152  typename TRAITS::PreEvalData d)
153 {
154  using std::logic_error;
155  using std::string;
156  using Teuchos::RCP;
157  using Teuchos::rcp_dynamic_cast;
158  using Teuchos::typeName;
162  using GED = panzer::GlobalEvaluationData;
163 
164  // First try the refactored ReadOnly container.
165  RCP<GED> ged;
166  string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
167  if (d.gedc->containsDataObject(globalDataKey_ + post))
168  {
169  ged = d.gedc->getDataObject(globalDataKey_ + post);
170  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
171  return;
172  } // end of the refactored ReadOnly way
173 
174  // Now try the old path.
175  {
176  ged = d.gedc->getDataObject(globalDataKey_);
177 
178  // Extract the linear object container.
179  auto roGed = rcp_dynamic_cast<const BVROGED>(ged);
180  auto beLoc = rcp_dynamic_cast<const BELOC>(ged);
181  if (not roGed.is_null())
182  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
183  else if (not beLoc.is_null())
184  {
185  if (useTimeDerivativeSolutionVector_)
186  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_dxdt());
187  else // if (not useTimeDerivativeSolutionVector_)
188  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_x());
189  TEUCHOS_TEST_FOR_EXCEPTION(x_.is_null(), logic_error, "Gather " \
190  "Residual: Can't find the x_ vector in GEDC \"" << globalDataKey_ <<
191  "\" (" << post << "). A cast failed for " << ged << ". Type is " <<
192  typeName(*ged));
193  } // end if we have a roGed or beLoc
194  } // end of the old path
195 } // end of preEvaluate() (Residual Specialization)
196 
198 //
199 // evaluateFields() (Residual Specialization)
200 //
202 template<typename TRAITS, typename LO, typename GO>
203 void
204 panzer::
207  typename TRAITS::EvalData workset)
208 {
209  using PHX::MDField;
210  using std::size_t;
211  using std::string;
212  using std::vector;
213  using Teuchos::ArrayRCP;
214  using Teuchos::ptrFromRef;
215  using Teuchos::RCP;
216  using Teuchos::rcp_dynamic_cast;
218  using Thyra::SpmdVectorBase;
219  using Thyra::VectorBase;
220 
221  // For convenience, pull out some objects from the workset.
222  string blockId(this->wda(workset).block_id);
223  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
224  int numFields(gatherFields_.size()), numCells(localCellIds.size());
225 
226  // Loop over the fields to be gathered.
227  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
228  {
229  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
230  auto field_h = Kokkos::create_mirror_view(field.get_static_view());
231 
232  int indexerId(indexerIds_[fieldInd]),
233  subFieldNum(subFieldIds_[fieldInd]);
234 
235  // Grab the local data for inputing.
238 
239  if(not x_.is_null()) {
240  rcp_dynamic_cast<SpmdVectorBase<double>>(x_->
241  getNonconstVectorBlock(indexerId))->getLocalData(ptrFromRef(x));
242  }
243  else {
244  xEvRoGed = xBvRoGed_->getGEDBlock(indexerId);
245  }
246 
247  auto subRowIndexer = indexers_[indexerId];
248  const vector<int>& elmtOffset =
249  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
250  int numBases(elmtOffset.size());
251 
252  auto LIDs = subRowIndexer->getLIDs();
253  auto LIDs_h = Kokkos::create_mirror_view(LIDs);
254  Kokkos::deep_copy(LIDs_h, LIDs);
255 
256  // Gather operation for each cell in the workset.
257  for (int cell(0); cell < numCells; ++cell)
258  {
259  LO cellLocalId = localCellIds[cell];
260 
261  // Loop over the basis functions and fill the fields.
262  for (int basis(0); basis < numBases; ++basis)
263  {
264  int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
265  if(x_.is_null())
266  field_h(cell, basis) = (*xEvRoGed)[lid];
267  else
268  field_h(cell, basis) = x[lid];
269  } // end loop over the basis functions
270  } // end loop over localCellIds
271  Kokkos::deep_copy(field.get_static_view(), field_h);
272  } // end loop over the fields to be gathered
273 } // end of evaluateFields() (Residual Specialization)
274 
276 //
277 // Initializing Constructor (Tangent Specialization)
278 //
280 template<typename TRAITS, typename LO, typename GO>
283  const std::vector<Teuchos::RCP<const GlobalIndexer>>&
284  indexers,
285  const Teuchos::ParameterList& p)
286  :
287  indexers_(indexers),
288  hasTangentFields_(false)
289 {
290  using panzer::PureBasis;
291  using PHX::MDField;
292  using PHX::print;
293  using std::size_t;
294  using std::string;
295  using std::vector;
296  using Teuchos::RCP;
297  using vvstring = std::vector<std::vector<std::string>>;
298  GatherSolution_Input input;
299  input.setParameterList(p);
300  const vector<string>& names = input.getDofNames();
301  RCP<const PureBasis> basis = input.getBasis();
302  const vvstring& tangentFieldNames = input.getTangentNames();
303  indexerNames_ = input.getIndexerNames();
304  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
305  globalDataKey_ = input.getGlobalDataKey();
306 
307  // Allocate the fields.
308  int numFields(names.size());
309  gatherFields_.resize(numFields);
310  for (int fd(0); fd < numFields; ++fd)
311  {
312  gatherFields_[fd] =
313  MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
314  this->addEvaluatedField(gatherFields_[fd]);
315  } // end loop over names
316 
317  // Set up the dependent tangent fields, if requested.
318  if (tangentFieldNames.size() > 0)
319  {
320  TEUCHOS_ASSERT(gatherFields_.size() == tangentFieldNames.size());
321  hasTangentFields_ = true;
322  tangentFields_.resize(numFields);
323  for (int fd(0); fd < numFields; ++fd)
324  {
325  int numTangentFields(tangentFieldNames[fd].size());
326  tangentFields_[fd].resize(numTangentFields);
327  for (int i(0); i < numTangentFields; ++i)
328  {
329  tangentFields_[fd][i] =
330  MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
331  basis->functional);
332  this->addDependentField(tangentFields_[fd][i]);
333  } // end loop over tangentFieldNames
334  } // end loop over gatherFields_
335  } // end if we have tangent fields
336 
337  // Figure out what the first active name is.
338  string firstName("<none>");
339  if (numFields > 0)
340  firstName = names[0];
341  string n("GatherSolution Tangent (BlockedEpetra): " + firstName + " (" +
342  print<EvalT>() + ")");
343  this->setName(n);
344 } // end of Initializing Constructor (Tangent Specialization)
345 
347 //
348 // postRegistrationSetup() (Tangent Specialization)
349 //
351 template<typename TRAITS, typename LO, typename GO>
352 void
355  typename TRAITS::SetupData /* d */,
356  PHX::FieldManager<TRAITS>& /* fm */)
357 {
358  using std::size_t;
359  using std::string;
360  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
361  int numFields(gatherFields_.size());
362  indexerIds_.resize(numFields);
363  subFieldIds_.resize(numFields);
364  for (int fd(0); fd < numFields; ++fd)
365  {
366  // Get the field ID from the DOF manager.
367  const string& fieldName(indexerNames_[fd]);
368  indexerIds_[fd] = getFieldBlock(fieldName, indexers_);
369  subFieldIds_[fd] = indexers_[indexerIds_[fd]]->getFieldNum(fieldName);
370  TEUCHOS_ASSERT(indexerIds_[fd] >= 0);
371  } // end loop over gatherFields_
372  indexerNames_.clear();
373 } // end of postRegistrationSetup() (Tangent Specialization)
374 
376 //
377 // preEvaluate() (Tangent Specialization)
378 //
380 template<typename TRAITS, typename LO, typename GO>
381 void
384  typename TRAITS::PreEvalData d)
385 {
386  using std::logic_error;
387  using std::string;
388  using Teuchos::RCP;
389  using Teuchos::rcp_dynamic_cast;
390  using Teuchos::typeName;
394  using GED = panzer::GlobalEvaluationData;
395 
396  // First try the refactored ReadOnly container.
397  RCP<GED> ged;
398  string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
399  if (d.gedc->containsDataObject(globalDataKey_ + post))
400  {
401  ged = d.gedc->getDataObject(globalDataKey_ + post);
402  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
403  return;
404  } // end of the refactored ReadOnly way
405 
406  // Now try the old path.
407  {
408  ged = d.gedc->getDataObject(globalDataKey_);
409 
410  // Extract the linear object container.
411  auto roGed = rcp_dynamic_cast<const BVROGED>(ged);
412  auto beLoc = rcp_dynamic_cast<const BELOC>(ged);
413  if (not roGed.is_null())
414  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
415  else if (not beLoc.is_null())
416  {
417  if (useTimeDerivativeSolutionVector_)
418  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_dxdt());
419  else // if (not useTimeDerivativeSolutionVector_)
420  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_x());
421  TEUCHOS_TEST_FOR_EXCEPTION(x_.is_null(), logic_error, "Gather " \
422  "Tangent: Can't find the x_ vector in GEDC \"" << globalDataKey_ <<
423  "\" (" << post << "). A cast failed for " << ged << ". Type is " <<
424  typeName(*ged));
425  } // end if we have a roGed or beLoc
426  } // end of the old path
427 } // end of preEvaluate() (Tangent Specialization)
428 
430 //
431 // evaluateFields() (Tangent Specialization)
432 //
434 template<typename TRAITS, typename LO, typename GO>
435 void
438  typename TRAITS::EvalData workset)
439 {
440  using PHX::MDField;
441  using std::pair;
442  using std::size_t;
443  using std::string;
444  using std::vector;
445  using Teuchos::ArrayRCP;
446  using Teuchos::ptrFromRef;
447  using Teuchos::RCP;
448  using Teuchos::rcp_dynamic_cast;
450  using Thyra::SpmdVectorBase;
451  using Thyra::VectorBase;
452 
453  // For convenience, pull out some objects from the workset.
454  vector<pair<int, GO>> GIDs;
455  string blockId(this->wda(workset).block_id);
456  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
457  int numFields(gatherFields_.size()), numCells(localCellIds.size());
458 
459  if (x_.is_null())
460  {
461  // Loop over the fields to be gathered.
462  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
463  {
464  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
465  int indexerId(indexerIds_[fieldInd]),
466  subFieldNum(subFieldIds_[fieldInd]);
467 
468  // Grab the local data for inputing.
469  auto xEvRoGed = xBvRoGed_->getGEDBlock(indexerId);
470  auto subRowIndexer = indexers_[indexerId];
471  const vector<int>& elmtOffset =
472  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
473  int numBases(elmtOffset.size());
474 
475  // Gather operation for each cell in the workset.
476  for (int cell(0); cell < numCells; ++cell)
477  {
478  LO cellLocalId = localCellIds[cell];
479  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
480 
481  // Loop over the basis functions and fill the fields.
482  for (int basis(0); basis < numBases; ++basis)
483  {
484  int offset(elmtOffset[basis]), lid(LIDs[offset]);
485  field(cell, basis) = (*xEvRoGed)[lid];
486  } // end loop over the basis functions
487  } // end loop over localCellIds
488  } // end loop over the fields to be gathered
489  }
490  else // if (not x_.is_null())
491  {
492  // Loop over the fields to be gathered.
493  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
494  {
495  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
496  int indexerId(indexerIds_[fieldInd]),
497  subFieldNum(subFieldIds_[fieldInd]);
498 
499  // Grab the local data for inputing.
501  rcp_dynamic_cast<SpmdVectorBase<double>>(x_->
502  getNonconstVectorBlock(indexerId))->getLocalData(ptrFromRef(x));
503  auto subRowIndexer = indexers_[indexerId];
504  const vector<int>& elmtOffset =
505  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
506  int numBases(elmtOffset.size());
507 
508  // Gather operation for each cell in the workset.
509  for (int cell(0); cell < numCells; ++cell)
510  {
511  LO cellLocalId = localCellIds[cell];
512  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
513 
514  // Loop over the basis functions and fill the fields.
515  for (int basis(0); basis < numBases; ++basis)
516  {
517  int offset(elmtOffset[basis]), lid(LIDs[offset]);
518  field(cell, basis) = x[lid];
519  } // end loop over the basis functions
520  } // end loop over localCellIds
521  } // end loop over the fields to be gathered
522  } // end if (x_.is_null()) or not
523 
524  // Deal with the tangent fields.
525  if (hasTangentFields_)
526  {
527  // Loop over the fields to be gathered.
528  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
529  {
530  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
531  int indexerId(indexerIds_[fieldInd]),
532  subFieldNum(subFieldIds_[fieldInd]);
533  auto subRowIndexer = indexers_[indexerId];
534  const vector<int>& elmtOffset =
535  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
536  int numBases(elmtOffset.size());
537 
538  // Gather operation for each cell in the workset.
539  for (int cell(0); cell < numCells; ++cell)
540  {
541  // Loop over the basis functions and fill the fields.
542  for (int basis(0); basis < numBases; ++basis)
543  {
544  int numTangentFields(tangentFields_[fieldInd].size());
545  for (int i(0); i < numTangentFields; ++i)
546  field(cell, basis).fastAccessDx(i) =
547  tangentFields_[fieldInd][i](cell, basis).val();
548  } // end loop over the basis functions
549  } // end loop over localCellIds
550  } // end loop over the fields to be gathered
551  } // end if (hasTangentFields_)
552 } // end of evaluateFields() (Tangent Specialization)
553 
555 //
556 // Initializing Constructor (Jacobian Specialization)
557 //
559 template<typename TRAITS, typename LO, typename GO>
560 panzer::
563  const std::vector<Teuchos::RCP<const GlobalIndexer>>&
564  indexers,
565  const Teuchos::ParameterList& p)
566  :
567  indexers_(indexers)
568 {
569  using panzer::PureBasis;
570  using PHX::MDField;
571  using PHX::print;
572  using std::size_t;
573  using std::string;
574  using std::vector;
575  using Teuchos::RCP;
576  GatherSolution_Input input;
577  input.setParameterList(p);
578  const vector<string>& names = input.getDofNames();
579  RCP<const PureBasis> basis = input.getBasis();
580  indexerNames_ = input.getIndexerNames();
581  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
582  globalDataKey_ = input.getGlobalDataKey();
583  gatherSeedIndex_ = input.getGatherSeedIndex();
584  sensitivitiesName_ = input.getSensitivitiesName();
585  disableSensitivities_ = not input.firstSensitivitiesAvailable();
586 
587  // Allocate the fields.
588  int numFields(names.size());
589  gatherFields_.resize(numFields);
590  for (int fd(0); fd < numFields; ++fd)
591  {
592  MDField<ScalarT, Cell, NODE> f(names[fd], basis->functional);
593  gatherFields_[fd] = f;
594  this->addEvaluatedField(gatherFields_[fd]);
595  } // end loop over names
596 
597  // Figure out what the first active name is.
598  string firstName("<none>"), n("GatherSolution (BlockedEpetra");
599  if (numFields > 0)
600  firstName = names[0];
601  if (disableSensitivities_)
602  n += ", No Sensitivities";
603  n += "): " + firstName + " (" + print<EvalT>() + ")";
604  this->setName(n);
605 } // end of Initializing Constructor (Jacobian Specialization)
606 
608 //
609 // postRegistrationSetup() (Jacobian Specialization)
610 //
612 template<typename TRAITS, typename LO, typename GO>
613 void
614 panzer::
617  typename TRAITS::SetupData /* d */,
618  PHX::FieldManager<TRAITS>& /* fm */)
619 {
620  using std::size_t;
621  using std::string;
622  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
623  int numFields(gatherFields_.size());
624  indexerIds_.resize(numFields);
625  subFieldIds_.resize(numFields);
626  for (int fd(0); fd < numFields; ++fd)
627  {
628  // Get the field ID from the DOF manager.
629  const string& fieldName(indexerNames_[fd]);
630  indexerIds_[fd] = getFieldBlock(fieldName, indexers_);
631  subFieldIds_[fd] = indexers_[indexerIds_[fd]]->getFieldNum(fieldName);
632  TEUCHOS_ASSERT(indexerIds_[fd] >= 0);
633  } // end loop over gatherFields_
634  indexerNames_.clear();
635 } // end of postRegistrationSetup() (Jacobian Specialization)
636 
638 //
639 // preEvaluate() (Jacobian Specialization)
640 //
642 template<typename TRAITS, typename LO, typename GO>
643 void
644 panzer::
647  typename TRAITS::PreEvalData d)
648 {
649  using std::endl;
650  using std::logic_error;
651  using std::string;
652  using Teuchos::RCP;
653  using Teuchos::rcp_dynamic_cast;
654  using Teuchos::typeName;
658  using GED = panzer::GlobalEvaluationData;
659 
660  // Manage sensitivities.
661  applySensitivities_ = false;
662  if ((not disableSensitivities_ ) and
663  (d.first_sensitivities_name == sensitivitiesName_))
664  applySensitivities_ = true;
665 
666  // First try the refactored ReadOnly container.
667  RCP<GED> ged;
668  string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
669  if (d.gedc->containsDataObject(globalDataKey_ + post))
670  {
671  ged = d.gedc->getDataObject(globalDataKey_ + post);
672  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
673  return;
674  } // end of the refactored ReadOnly way
675 
676  // Now try the old path.
677  {
678  ged = d.gedc->getDataObject(globalDataKey_);
679 
680  // Extract the linear object container.
681  auto roGed = rcp_dynamic_cast<const BVROGED>(ged);
682  auto beLoc = rcp_dynamic_cast<const BELOC>(ged);
683  if (not roGed.is_null())
684  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
685  else if (not beLoc.is_null())
686  {
687  if (useTimeDerivativeSolutionVector_)
688  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_dxdt());
689  else // if (not useTimeDerivativeSolutionVector_)
690  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_x());
691  TEUCHOS_TEST_FOR_EXCEPTION(x_.is_null(), logic_error, "Gather " \
692  "Jacobian: Can't find x vector in GEDC \"" << globalDataKey_ <<
693  "\" (" << post << "). A cast failed for " << ged << ". Type is " <<
694  typeName(*ged));
695  } // end if we have a roGed or beLoc
696  } // end of the old path
697 } // end of preEvaluate() (Jacobian Specialization)
698 
700 //
701 // evaluateFields() (Jacobian Specialization)
702 //
704 template<typename TRAITS, typename LO, typename GO>
705 void
706 panzer::
709  typename TRAITS::EvalData workset)
710 {
711  using PHX::MDField;
712  using std::size_t;
713  using std::string;
714  using std::vector;
715  using Teuchos::ArrayRCP;
716  using Teuchos::ptrFromRef;
717  using Teuchos::RCP;
718  using Teuchos::rcp_dynamic_cast;
720  using Thyra::SpmdVectorBase;
721  using Thyra::VectorBase;
722 
723  // For convenience, pull out some objects from the workset.
724  string blockId(this->wda(workset).block_id);
725  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
726  int numFields(gatherFields_.size()), numCells(localCellIds.size());
727  double seedValue(0);
728  if (applySensitivities_)
729  {
730  if ((useTimeDerivativeSolutionVector_) and (gatherSeedIndex_ < 0))
731  seedValue = workset.alpha;
732  else if (gatherSeedIndex_ < 0)
733  seedValue = workset.beta;
734  else if (not useTimeDerivativeSolutionVector_)
735  seedValue = workset.gather_seeds[gatherSeedIndex_];
736  else
737  TEUCHOS_ASSERT(false);
738  } // end if (applySensitivities_)
739 
740  // Turn off sensitivies: This may be faster if we don't expand the term, but
741  // I suspect not, because anywhere it is used the full complement of
742  // sensitivies will be needed anyway.
743  if (not applySensitivities_)
744  seedValue = 0.0;
745 
746  vector<int> blockOffsets;
747  computeBlockOffsets(blockId, indexers_, blockOffsets);
748  int numDerivs(blockOffsets[blockOffsets.size() - 1]);
749 
750  // NOTE: A reordering of these loops will likely improve performance. The
751  // "getGIDFieldOffsets may be expensive. However the "getElementGIDs"
752  // can be cheaper. However the lookup for LIDs may be more expensive!
753 
754  // Loop over the fields to be gathered.
755  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
756  {
757  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
758  auto field_h = Kokkos::create_mirror_view(field.get_view());
759 
760  int indexerId(indexerIds_[fieldInd]), subFieldNum(subFieldIds_[fieldInd]);
761 
762  // Grab the local data for inputing.
765  if(not x_.is_null()) {
766  rcp_dynamic_cast<SpmdVectorBase<double>>(x_->getNonconstVectorBlock(indexerId))->getLocalData(ptrFromRef(x));
767  }else {
768  xEvRoGed = xBvRoGed_->getGEDBlock(indexerId);
769  }
770 
771  auto subRowIndexer = indexers_[indexerId];
772  const vector<int>& elmtOffset =
773  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
774  int startBlkOffset(blockOffsets[indexerId]), numBases(elmtOffset.size());
775 
776  auto LIDs = subRowIndexer->getLIDs();
777  auto LIDs_h = Kokkos::create_mirror_view(LIDs);
778  Kokkos::deep_copy(LIDs_h, LIDs);
779 
780  // Gather operation for each cell in the workset.
781  for (int cell(0); cell < numCells; ++cell)
782  {
783  LO cellLocalId = localCellIds[cell];
784 
785  // Loop over the basis functions and fill the fields.
786  for (int basis(0); basis < numBases; ++basis)
787  {
788  // Set the value and seed the FAD object.
789  int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
790  if(x_.is_null())
791  field_h(cell, basis) = ScalarT(numDerivs, (*xEvRoGed)[lid]);
792  else
793  field_h(cell, basis) = ScalarT(numDerivs, x[lid]);
794 
795  field_h(cell, basis).fastAccessDx(startBlkOffset + offset) = seedValue;
796  } // end loop over the basis functions
797  } // end loop over localCellIds
798  Kokkos::deep_copy(field.get_static_view(), field_h);
799  } // end loop over the fields to be gathered
800 } // end of evaluateFields() (Jacobian Specialization)
801 
802 #endif // __Panzer_GatherSolution_BlockedEpetra_impl_hpp__
const std::vector< std::vector< std::string > > & getTangentNames() const
Get the name of the tangent fields (tangent only)
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
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 getSensitivitiesName() const
The name of the sensitivities. Enables sensitivities at &quot;preEvaluate&quot; time (Jacobian and Hessian) ...
std::string typeName(const T &t)