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 //
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_BlockedEpetra_impl_hpp__
44 #define __Panzer_GatherSolution_BlockedEpetra_impl_hpp__
45 
47 //
48 // Include Files
49 //
51 
52 // Panzer
56 #include "Panzer_PureBasis.hpp"
60 
61 // Phalanx
62 #include "Phalanx_DataLayout.hpp"
63 
64 // Teuchos
65 #include "Teuchos_Assert.hpp"
66 #include "Teuchos_FancyOStream.hpp"
67 
68 // Thyra
69 #include "Thyra_ProductVectorBase.hpp"
70 #include "Thyra_SpmdVectorBase.hpp"
71 
73 //
74 // Initializing Constructor (Residual Specialization)
75 //
77 template<typename TRAITS, typename LO, typename GO>
81  const std::vector<Teuchos::RCP<const UniqueGlobalIndexer<LO, int>>>&
82  indexers,
83  const Teuchos::ParameterList& p)
84  :
85  indexers_(indexers),
86  hasTangentFields_(false)
87 {
88  using panzer::PureBasis;
89  using PHX::MDField;
90  using PHX::typeAsString;
91  using std::size_t;
92  using std::string;
93  using std::vector;
94  using Teuchos::RCP;
95  using vvstring = std::vector<std::vector<std::string>>;
97  input.setParameterList(p);
98  const vector<string>& names = input.getDofNames();
99  RCP<const PureBasis> basis = input.getBasis();
100  const vvstring& tangentFieldNames = input.getTangentNames();
101  indexerNames_ = input.getIndexerNames();
102  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
103  globalDataKey_ = input.getGlobalDataKey();
104 
105  // Allocate the fields.
106  int numFields(names.size());
107  gatherFields_.resize(numFields);
108  for (int fd(0); fd < numFields; ++fd)
109  {
110  gatherFields_[fd] =
111  MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
112  this->addEvaluatedField(gatherFields_[fd]);
113  } // end loop over names
114 
115  // Setup the dependent tangent fields, if requested.
116  if (tangentFieldNames.size() > 0)
117  {
118  TEUCHOS_ASSERT(gatherFields_.size() == tangentFieldNames.size());
119  hasTangentFields_ = true;
120  tangentFields_.resize(numFields);
121  for (int fd(0); fd < numFields; ++fd)
122  {
123  int numTangentFields(tangentFieldNames[fd].size());
124  tangentFields_[fd].resize(numTangentFields);
125  for (int i(0); i < numTangentFields; ++i)
126  {
127  tangentFields_[fd][i] =
128  MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
129  basis->functional);
130  this->addDependentField(tangentFields_[fd][i]);
131  } // end loop over tangentFieldNames[fd]
132  } // end loop over gatherFields_
133  } // end if we have tangent fields
134 
135  // Figure out what the first active name is.
136  string firstName("<none>");
137  if (numFields > 0)
138  firstName = names[0];
139  string n("GatherSolution (BlockedEpetra): " + firstName + " (" +
140  typeAsString<EvalT>() + ")");
141  this->setName(n);
142 } // end of Initializing Constructor (Residual Specialization)
143 
145 //
146 // postRegistrationSetup() (Residual Specialization)
147 //
149 template<typename TRAITS, typename LO, typename GO>
150 void
151 panzer::
154  typename TRAITS::SetupData /* d */,
155  PHX::FieldManager<TRAITS>& /* fm */)
156 {
157  using std::size_t;
158  using std::string;
159  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
160  int numFields(gatherFields_.size());
161  indexerIds_.resize(numFields);
162  subFieldIds_.resize(numFields);
163  for (int fd(0); fd < numFields; ++fd)
164  {
165  // Get the field ID from the DOF manager.
166  const string& fieldName(indexerNames_[fd]);
167  indexerIds_[fd] = getFieldBlock(fieldName, indexers_);
168  subFieldIds_[fd] = indexers_[indexerIds_[fd]]->getFieldNum(fieldName);
169  TEUCHOS_ASSERT(indexerIds_[fd] >= 0);
170  } // end loop over gatherFields_
171  indexerNames_.clear();
172 } // end of postRegistrationSetup() (Residual Specialization)
173 
175 //
176 // preEvaluate() (Residual Specialization)
177 //
179 template<typename TRAITS, typename LO, typename GO>
180 void
181 panzer::
184  typename TRAITS::PreEvalData d)
185 {
186  using std::logic_error;
187  using std::string;
188  using Teuchos::RCP;
189  using Teuchos::rcp_dynamic_cast;
190  using Teuchos::typeName;
194  using GED = panzer::GlobalEvaluationData;
195 
196  // First try the refactored ReadOnly container.
197  RCP<GED> ged;
198  string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
199  if (d.gedc->containsDataObject(globalDataKey_ + post))
200  {
201  ged = d.gedc->getDataObject(globalDataKey_ + post);
202  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
203  return;
204  } // end of the refactored ReadOnly way
205 
206  // Now try the old path.
207  {
208  ged = d.gedc->getDataObject(globalDataKey_);
209 
210  // Extract the linear object container.
211  auto roGed = rcp_dynamic_cast<const BVROGED>(ged);
212  auto beLoc = rcp_dynamic_cast<const BELOC>(ged);
213  if (not roGed.is_null())
214  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
215  else if (not beLoc.is_null())
216  {
217  if (useTimeDerivativeSolutionVector_)
218  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_dxdt());
219  else // if (not useTimeDerivativeSolutionVector_)
220  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_x());
221  TEUCHOS_TEST_FOR_EXCEPTION(x_.is_null(), logic_error, "Gather " \
222  "Residual: Can't find the x_ vector in GEDC \"" << globalDataKey_ <<
223  "\" (" << post << "). A cast failed for " << ged << ". Type is " <<
224  typeName(*ged));
225  } // end if we have a roGed or beLoc
226  } // end of the old path
227 } // end of preEvaluate() (Residual Specialization)
228 
230 //
231 // evaluateFields() (Residual Specialization)
232 //
234 template<typename TRAITS, typename LO, typename GO>
235 void
236 panzer::
239  typename TRAITS::EvalData workset)
240 {
241  using PHX::MDField;
242  using std::size_t;
243  using std::string;
244  using std::vector;
245  using Teuchos::ArrayRCP;
246  using Teuchos::ptrFromRef;
247  using Teuchos::RCP;
248  using Teuchos::rcp_dynamic_cast;
250  using Thyra::SpmdVectorBase;
251  using Thyra::VectorBase;
252 
253  // For convenience, pull out some objects from the workset.
254  string blockId(this->wda(workset).block_id);
255  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
256  int numFields(gatherFields_.size()), numCells(localCellIds.size());
257 
258  if (x_.is_null())
259  {
260  // Loop over the fields to be gathered.
261  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
262  {
263  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
264  int indexerId(indexerIds_[fieldInd]),
265  subFieldNum(subFieldIds_[fieldInd]);
266 
267  // Grab the local data for inputing.
268  auto xEvRoGed = xBvRoGed_->getGEDBlock(indexerId);
269  auto subRowIndexer = indexers_[indexerId];
270  const vector<int>& elmtOffset =
271  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
272  int numBases(elmtOffset.size());
273 
274  // Gather operation for each cell in the workset.
275  for (int cell(0); cell < numCells; ++cell)
276  {
277  LO cellLocalId = localCellIds[cell];
278  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
279 
280  // Loop over the basis functions and fill the fields.
281  for (int basis(0); basis < numBases; ++basis)
282  {
283  int offset(elmtOffset[basis]), lid(LIDs[offset]);
284  field(cell, basis) = (*xEvRoGed)[lid];
285  } // end loop over the basis functions
286  } // end loop over localCellIds
287  } // end loop over the fields to be gathered
288  }
289  else // if (not x_.is_null())
290  {
291  // Loop over the fields to be gathered.
292  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
293  {
294  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
295  int indexerId(indexerIds_[fieldInd]),
296  subFieldNum(subFieldIds_[fieldInd]);
297 
298  // Grab the local data for inputing.
300  rcp_dynamic_cast<SpmdVectorBase<double>>(x_->
301  getNonconstVectorBlock(indexerId))->getLocalData(ptrFromRef(x));
302  auto subRowIndexer = indexers_[indexerId];
303  const vector<int>& elmtOffset =
304  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
305  int numBases(elmtOffset.size());
306 
307  // Gather operation for each cell in the workset.
308  for (int cell(0); cell < numCells; ++cell)
309  {
310  LO cellLocalId = localCellIds[cell];
311  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
312 
313  // Loop over the basis functions and fill the fields.
314  for (int basis(0); basis < numBases; ++basis)
315  {
316  int offset(elmtOffset[basis]), lid(LIDs[offset]);
317  field(cell, basis) = x[lid];
318  } // end loop over the basis functions
319  } // end loop over localCellIds
320  } // end loop over the fields to be gathered
321  } // end if (x_.is_null()) or not
322 } // end of evaluateFields() (Residual Specialization)
323 
325 //
326 // Initializing Constructor (Tangent Specialization)
327 //
329 template<typename TRAITS, typename LO, typename GO>
332  const std::vector<Teuchos::RCP<const UniqueGlobalIndexer<LO, int>>>&
333  indexers,
334  const Teuchos::ParameterList& p)
335  :
336  indexers_(indexers),
337  hasTangentFields_(false)
338 {
339  using panzer::PureBasis;
340  using PHX::MDField;
341  using PHX::typeAsString;
342  using std::size_t;
343  using std::string;
344  using std::vector;
345  using Teuchos::RCP;
346  using vvstring = std::vector<std::vector<std::string>>;
347  GatherSolution_Input input;
348  input.setParameterList(p);
349  const vector<string>& names = input.getDofNames();
350  RCP<const PureBasis> basis = input.getBasis();
351  const vvstring& tangentFieldNames = input.getTangentNames();
352  indexerNames_ = input.getIndexerNames();
353  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
354  globalDataKey_ = input.getGlobalDataKey();
355 
356  // Allocate the fields.
357  int numFields(names.size());
358  gatherFields_.resize(numFields);
359  for (int fd(0); fd < numFields; ++fd)
360  {
361  gatherFields_[fd] =
362  MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
363  this->addEvaluatedField(gatherFields_[fd]);
364  } // end loop over names
365 
366  // Set up the dependent tangent fields, if requested.
367  if (tangentFieldNames.size() > 0)
368  {
369  TEUCHOS_ASSERT(gatherFields_.size() == tangentFieldNames.size());
370  hasTangentFields_ = true;
371  tangentFields_.resize(numFields);
372  for (int fd(0); fd < numFields; ++fd)
373  {
374  int numTangentFields(tangentFieldNames[fd].size());
375  tangentFields_[fd].resize(numTangentFields);
376  for (int i(0); i < numTangentFields; ++i)
377  {
378  tangentFields_[fd][i] =
379  MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
380  basis->functional);
381  this->addDependentField(tangentFields_[fd][i]);
382  } // end loop over tangentFieldNames
383  } // end loop over gatherFields_
384  } // end if we have tangent fields
385 
386  // Figure out what the first active name is.
387  string firstName("<none>");
388  if (numFields > 0)
389  firstName = names[0];
390  string n("GatherSolution Tangent (BlockedEpetra): " + firstName + " (" +
391  typeAsString<EvalT>() + ")");
392  this->setName(n);
393 } // end of Initializing Constructor (Tangent Specialization)
394 
396 //
397 // postRegistrationSetup() (Tangent Specialization)
398 //
400 template<typename TRAITS, typename LO, typename GO>
401 void
404  typename TRAITS::SetupData /* d */,
405  PHX::FieldManager<TRAITS>& /* fm */)
406 {
407  using std::size_t;
408  using std::string;
409  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
410  int numFields(gatherFields_.size());
411  indexerIds_.resize(numFields);
412  subFieldIds_.resize(numFields);
413  for (int fd(0); fd < numFields; ++fd)
414  {
415  // Get the field ID from the DOF manager.
416  const string& fieldName(indexerNames_[fd]);
417  indexerIds_[fd] = getFieldBlock(fieldName, indexers_);
418  subFieldIds_[fd] = indexers_[indexerIds_[fd]]->getFieldNum(fieldName);
419  TEUCHOS_ASSERT(indexerIds_[fd] >= 0);
420  } // end loop over gatherFields_
421  indexerNames_.clear();
422 } // end of postRegistrationSetup() (Tangent Specialization)
423 
425 //
426 // preEvaluate() (Tangent Specialization)
427 //
429 template<typename TRAITS, typename LO, typename GO>
430 void
433  typename TRAITS::PreEvalData d)
434 {
435  using std::logic_error;
436  using std::string;
437  using Teuchos::RCP;
438  using Teuchos::rcp_dynamic_cast;
439  using Teuchos::typeName;
443  using GED = panzer::GlobalEvaluationData;
444 
445  // First try the refactored ReadOnly container.
446  RCP<GED> ged;
447  string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
448  if (d.gedc->containsDataObject(globalDataKey_ + post))
449  {
450  ged = d.gedc->getDataObject(globalDataKey_ + post);
451  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
452  return;
453  } // end of the refactored ReadOnly way
454 
455  // Now try the old path.
456  {
457  ged = d.gedc->getDataObject(globalDataKey_);
458 
459  // Extract the linear object container.
460  auto roGed = rcp_dynamic_cast<const BVROGED>(ged);
461  auto beLoc = rcp_dynamic_cast<const BELOC>(ged);
462  if (not roGed.is_null())
463  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
464  else if (not beLoc.is_null())
465  {
466  if (useTimeDerivativeSolutionVector_)
467  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_dxdt());
468  else // if (not useTimeDerivativeSolutionVector_)
469  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_x());
470  TEUCHOS_TEST_FOR_EXCEPTION(x_.is_null(), logic_error, "Gather " \
471  "Tangent: Can't find the x_ vector in GEDC \"" << globalDataKey_ <<
472  "\" (" << post << "). A cast failed for " << ged << ". Type is " <<
473  typeName(*ged));
474  } // end if we have a roGed or beLoc
475  } // end of the old path
476 } // end of preEvaluate() (Tangent Specialization)
477 
479 //
480 // evaluateFields() (Tangent Specialization)
481 //
483 template<typename TRAITS, typename LO, typename GO>
484 void
487  typename TRAITS::EvalData workset)
488 {
489  using PHX::MDField;
490  using std::pair;
491  using std::size_t;
492  using std::string;
493  using std::vector;
494  using Teuchos::ArrayRCP;
495  using Teuchos::ptrFromRef;
496  using Teuchos::RCP;
497  using Teuchos::rcp_dynamic_cast;
499  using Thyra::SpmdVectorBase;
500  using Thyra::VectorBase;
501 
502  // For convenience, pull out some objects from the workset.
503  vector<pair<int, GO>> GIDs;
504  string blockId(this->wda(workset).block_id);
505  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
506  int numFields(gatherFields_.size()), numCells(localCellIds.size());
507 
508  if (x_.is_null())
509  {
510  // Loop over the fields to be gathered.
511  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
512  {
513  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
514  int indexerId(indexerIds_[fieldInd]),
515  subFieldNum(subFieldIds_[fieldInd]);
516 
517  // Grab the local data for inputing.
518  auto xEvRoGed = xBvRoGed_->getGEDBlock(indexerId);
519  auto subRowIndexer = indexers_[indexerId];
520  const vector<int>& elmtOffset =
521  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
522  int numBases(elmtOffset.size());
523 
524  // Gather operation for each cell in the workset.
525  for (int cell(0); cell < numCells; ++cell)
526  {
527  LO cellLocalId = localCellIds[cell];
528  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
529 
530  // Loop over the basis functions and fill the fields.
531  for (int basis(0); basis < numBases; ++basis)
532  {
533  int offset(elmtOffset[basis]), lid(LIDs[offset]);
534  field(cell, basis) = (*xEvRoGed)[lid];
535  } // end loop over the basis functions
536  } // end loop over localCellIds
537  } // end loop over the fields to be gathered
538  }
539  else // if (not x_.is_null())
540  {
541  // Loop over the fields to be gathered.
542  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
543  {
544  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
545  int indexerId(indexerIds_[fieldInd]),
546  subFieldNum(subFieldIds_[fieldInd]);
547 
548  // Grab the local data for inputing.
550  rcp_dynamic_cast<SpmdVectorBase<double>>(x_->
551  getNonconstVectorBlock(indexerId))->getLocalData(ptrFromRef(x));
552  auto subRowIndexer = indexers_[indexerId];
553  const vector<int>& elmtOffset =
554  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
555  int numBases(elmtOffset.size());
556 
557  // Gather operation for each cell in the workset.
558  for (int cell(0); cell < numCells; ++cell)
559  {
560  LO cellLocalId = localCellIds[cell];
561  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
562 
563  // Loop over the basis functions and fill the fields.
564  for (int basis(0); basis < numBases; ++basis)
565  {
566  int offset(elmtOffset[basis]), lid(LIDs[offset]);
567  field(cell, basis) = x[lid];
568  } // end loop over the basis functions
569  } // end loop over localCellIds
570  } // end loop over the fields to be gathered
571  } // end if (x_.is_null()) or not
572 
573  // Deal with the tangent fields.
574  if (hasTangentFields_)
575  {
576  // Loop over the fields to be gathered.
577  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
578  {
579  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
580  int indexerId(indexerIds_[fieldInd]),
581  subFieldNum(subFieldIds_[fieldInd]);
582  auto subRowIndexer = indexers_[indexerId];
583  const vector<int>& elmtOffset =
584  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
585  int numBases(elmtOffset.size());
586 
587  // Gather operation for each cell in the workset.
588  for (int cell(0); cell < numCells; ++cell)
589  {
590  // Loop over the basis functions and fill the fields.
591  for (int basis(0); basis < numBases; ++basis)
592  {
593  int numTangentFields(tangentFields_[fieldInd].size());
594  for (int i(0); i < numTangentFields; ++i)
595  field(cell, basis).fastAccessDx(i) =
596  tangentFields_[fieldInd][i](cell, basis).val();
597  } // end loop over the basis functions
598  } // end loop over localCellIds
599  } // end loop over the fields to be gathered
600  } // end if (hasTangentFields_)
601 } // end of evaluateFields() (Tangent Specialization)
602 
604 //
605 // Initializing Constructor (Jacobian Specialization)
606 //
608 template<typename TRAITS, typename LO, typename GO>
609 panzer::
612  const std::vector<Teuchos::RCP<const UniqueGlobalIndexer<LO, int>>>&
613  indexers,
614  const Teuchos::ParameterList& p)
615  :
616  indexers_(indexers)
617 {
618  using panzer::PureBasis;
619  using PHX::MDField;
620  using PHX::typeAsString;
621  using std::size_t;
622  using std::string;
623  using std::vector;
624  using Teuchos::RCP;
625  GatherSolution_Input input;
626  input.setParameterList(p);
627  const vector<string>& names = input.getDofNames();
628  RCP<const PureBasis> basis = input.getBasis();
629  indexerNames_ = input.getIndexerNames();
630  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
631  globalDataKey_ = input.getGlobalDataKey();
632  gatherSeedIndex_ = input.getGatherSeedIndex();
633  sensitivitiesName_ = input.getSensitivitiesName();
634  disableSensitivities_ = not input.firstSensitivitiesAvailable();
635 
636  // Allocate the fields.
637  int numFields(names.size());
638  gatherFields_.resize(numFields);
639  for (int fd(0); fd < numFields; ++fd)
640  {
641  MDField<ScalarT, Cell, NODE> f(names[fd], basis->functional);
642  gatherFields_[fd] = f;
643  this->addEvaluatedField(gatherFields_[fd]);
644  } // end loop over names
645 
646  // Figure out what the first active name is.
647  string firstName("<none>"), n("GatherSolution (BlockedEpetra");
648  if (numFields > 0)
649  firstName = names[0];
650  if (disableSensitivities_)
651  n += ", No Sensitivities";
652  n += "): " + firstName + " (" + typeAsString<EvalT>() + ")";
653  this->setName(n);
654 } // end of Initializing Constructor (Jacobian Specialization)
655 
657 //
658 // postRegistrationSetup() (Jacobian Specialization)
659 //
661 template<typename TRAITS, typename LO, typename GO>
662 void
663 panzer::
666  typename TRAITS::SetupData /* d */,
667  PHX::FieldManager<TRAITS>& /* fm */)
668 {
669  using std::size_t;
670  using std::string;
671  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
672  int numFields(gatherFields_.size());
673  indexerIds_.resize(numFields);
674  subFieldIds_.resize(numFields);
675  for (int fd(0); fd < numFields; ++fd)
676  {
677  // Get the field ID from the DOF manager.
678  const string& fieldName(indexerNames_[fd]);
679  indexerIds_[fd] = getFieldBlock(fieldName, indexers_);
680  subFieldIds_[fd] = indexers_[indexerIds_[fd]]->getFieldNum(fieldName);
681  TEUCHOS_ASSERT(indexerIds_[fd] >= 0);
682  } // end loop over gatherFields_
683  indexerNames_.clear();
684 } // end of postRegistrationSetup() (Jacobian Specialization)
685 
687 //
688 // preEvaluate() (Jacobian Specialization)
689 //
691 template<typename TRAITS, typename LO, typename GO>
692 void
693 panzer::
696  typename TRAITS::PreEvalData d)
697 {
698  using std::endl;
699  using std::logic_error;
700  using std::string;
701  using Teuchos::RCP;
702  using Teuchos::rcp_dynamic_cast;
703  using Teuchos::typeName;
707  using GED = panzer::GlobalEvaluationData;
708 
709  // Manage sensitivities.
710  applySensitivities_ = false;
711  if ((not disableSensitivities_ ) and
712  (d.first_sensitivities_name == sensitivitiesName_))
713  applySensitivities_ = true;
714 
715  // First try the refactored ReadOnly container.
716  RCP<GED> ged;
717  string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
718  if (d.gedc->containsDataObject(globalDataKey_ + post))
719  {
720  ged = d.gedc->getDataObject(globalDataKey_ + post);
721  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
722  return;
723  } // end of the refactored ReadOnly way
724 
725  // Now try the old path.
726  {
727  ged = d.gedc->getDataObject(globalDataKey_);
728 
729  // Extract the linear object container.
730  auto roGed = rcp_dynamic_cast<const BVROGED>(ged);
731  auto beLoc = rcp_dynamic_cast<const BELOC>(ged);
732  if (not roGed.is_null())
733  xBvRoGed_ = rcp_dynamic_cast<BVROGED>(ged, true);
734  else if (not beLoc.is_null())
735  {
736  if (useTimeDerivativeSolutionVector_)
737  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_dxdt());
738  else // if (not useTimeDerivativeSolutionVector_)
739  x_ = rcp_dynamic_cast<ProductVectorBase<double>>(beLoc->get_x());
740  TEUCHOS_TEST_FOR_EXCEPTION(x_.is_null(), logic_error, "Gather " \
741  "Jacobian: Can't find x vector in GEDC \"" << globalDataKey_ <<
742  "\" (" << post << "). A cast failed for " << ged << ". Type is " <<
743  typeName(*ged));
744  } // end if we have a roGed or beLoc
745  } // end of the old path
746 } // end of preEvaluate() (Jacobian Specialization)
747 
749 //
750 // evaluateFields() (Jacobian Specialization)
751 //
753 template<typename TRAITS, typename LO, typename GO>
754 void
755 panzer::
758  typename TRAITS::EvalData workset)
759 {
760  using PHX::MDField;
761  using std::size_t;
762  using std::string;
763  using std::vector;
764  using Teuchos::ArrayRCP;
765  using Teuchos::ptrFromRef;
766  using Teuchos::RCP;
767  using Teuchos::rcp_dynamic_cast;
769  using Thyra::SpmdVectorBase;
770  using Thyra::VectorBase;
771 
772  // For convenience, pull out some objects from the workset.
773  string blockId(this->wda(workset).block_id);
774  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
775  int numFields(gatherFields_.size()), numCells(localCellIds.size());
776  double seedValue(0);
777  if (applySensitivities_)
778  {
779  if ((useTimeDerivativeSolutionVector_) and (gatherSeedIndex_ < 0))
780  seedValue = workset.alpha;
781  else if (gatherSeedIndex_ < 0)
782  seedValue = workset.beta;
783  else if (not useTimeDerivativeSolutionVector_)
784  seedValue = workset.gather_seeds[gatherSeedIndex_];
785  else
786  TEUCHOS_ASSERT(false);
787  } // end if (applySensitivities_)
788 
789  // Turn off sensitivies: This may be faster if we don't expand the term, but
790  // I suspect not, because anywhere it is used the full complement of
791  // sensitivies will be needed anyway.
792  if (not applySensitivities_)
793  seedValue = 0.0;
794 
795  vector<int> blockOffsets;
796  computeBlockOffsets(blockId, indexers_, blockOffsets);
797  int numDerivs(blockOffsets[blockOffsets.size() - 1]);
798 
799  // NOTE: A reordering of these loops will likely improve performance. The
800  // "getGIDFieldOffsets may be expensive. However the "getElementGIDs"
801  // can be cheaper. However the lookup for LIDs may be more expensive!
802  if (x_.is_null())
803  {
804  // Loop over the fields to be gathered.
805  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
806  {
807  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
808  int indexerId(indexerIds_[fieldInd]),
809  subFieldNum(subFieldIds_[fieldInd]);
810 
811  // Grab the local data for inputing.
812  auto xEvRoGed = xBvRoGed_->getGEDBlock(indexerId);
813  auto subRowIndexer = indexers_[indexerId];
814  const vector<int>& elmtOffset =
815  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
816  int startBlkOffset(blockOffsets[indexerId]), numBases(elmtOffset.size());
817 
818  // Gather operation for each cell in the workset.
819  for (int cell(0); cell < numCells; ++cell)
820  {
821  LO cellLocalId = localCellIds[cell];
822  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
823 
824  // Loop over the basis functions and fill the fields.
825  for (int basis(0); basis < numBases; ++basis)
826  {
827  // Set the value and seed the FAD object.
828  int offset(elmtOffset[basis]), lid(LIDs[offset]);
829  field(cell, basis) = ScalarT(numDerivs, (*xEvRoGed)[lid]);
830  field(cell, basis).fastAccessDx(startBlkOffset + offset) = seedValue;
831  } // end loop over the basis functions
832  } // end loop over localCellIds
833  } // end loop over the fields to be gathered
834  }
835  else // if (not x_.is_null())
836  {
837  // Loop over the fields to be gathered.
838  for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
839  {
840  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
841  int indexerId(indexerIds_[fieldInd]),
842  subFieldNum(subFieldIds_[fieldInd]);
843 
844  // Grab the local data for inputing.
846  rcp_dynamic_cast<SpmdVectorBase<double>>(x_->
847  getNonconstVectorBlock(indexerId))->getLocalData(ptrFromRef(x));
848  auto subRowIndexer = indexers_[indexerId];
849  const vector<int>& elmtOffset =
850  subRowIndexer->getGIDFieldOffsets(blockId, subFieldNum);
851  int startBlkOffset(blockOffsets[indexerId]), numBases(elmtOffset.size());
852 
853  // Gather operation for each cell in the workset.
854  for (int cell(0); cell < numCells; ++cell)
855  {
856  LO cellLocalId = localCellIds[cell];
857  auto LIDs = subRowIndexer->getElementLIDs(cellLocalId);
858 
859  // Loop over the basis functions and fill the fields.
860  for (int basis(0); basis < numBases; ++basis)
861  {
862  // Set the value and seed the FAD object.
863  int offset(elmtOffset[basis]), lid(LIDs[offset]);
864  field(cell, basis) = ScalarT(numDerivs, x[lid]);
865  field(cell, basis).fastAccessDx(startBlkOffset + offset) = seedValue;
866  } // end loop over the basis functions
867  } // end loop over localCellIds
868  } // end loop over the fields to be gathered
869  } // end if (x_.is_null()) or not
870 } // end of evaluateFields() (Jacobian Specialization)
871 
872 #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)
int getFieldBlock(const std::string &fieldName, const std::vector< Teuchos::RCP< UniqueGlobalIndexer< LocalOrdinalT, GlobalOrdinalT > > > &ugis)
void setParameterList(const Teuchos::ParameterList &pl)
void computeBlockOffsets(const std::string &blockId, const std::vector< Teuchos::RCP< UniqueGlobalIndexer< LocalOrdinalT, GlobalOrdinalT > > > &ugis, std::vector< int > &blockOffsets)
#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) ...
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)
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.
Gathers solution values from the Newton solution vector into the nodal fields of the field manager...
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)