11 #ifndef __Panzer_GatherSolution_Epetra_impl_hpp__
12 #define __Panzer_GatherSolution_Epetra_impl_hpp__
21 #include "Epetra_Vector.h"
33 #include "Teuchos_Assert.hpp"
36 #include "Thyra_SpmdVectorBase.hpp"
43 template<
typename TRAITS,
typename LO,
typename GO>
49 globalIndexer_(indexer),
50 hasTangentFields_(false)
59 using vvstring = std::vector<std::vector<std::string>>;
75 MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
76 this->addEvaluatedField(gatherFields_[fd]);
80 if (tangentFieldNames.size() > 0)
83 hasTangentFields_ =
true;
84 tangentFields_.resize(numFields);
87 int numTangentFields(tangentFieldNames[fd].size());
88 tangentFields_[fd].resize(numTangentFields);
89 for (
int i(0); i < numTangentFields; ++i)
91 tangentFields_[fd][i] =
92 MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
94 this->addDependentField(tangentFields_[fd][i]);
100 string firstName(
"<none>");
102 firstName = names[0];
103 string n(
"GatherSolution (Epetra): " + firstName +
" (Residual)");
112 template<
typename TRAITS,
typename LO,
typename GO>
116 typename TRAITS::SetupData ,
119 using std::logic_error;
128 const string& fieldName(indexerNames_[fd]);
129 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
133 "GatherSolution_Epetra<Residual>: Could not find field \"" +
134 fieldName +
"\" in the global indexer. ")
136 indexerNames_.clear();
144 template<
typename TRAITS,
typename LO,
typename GO>
148 typename TRAITS::PreEvalData d)
152 using Teuchos::rcp_dynamic_cast;
161 string post(useTimeDerivativeSolutionVector_ ?
" - Xdot" :
" - X");
162 if (d.gedc->containsDataObject(globalDataKey_ + post))
164 ged = d.gedc->getDataObject(globalDataKey_ + post);
165 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
170 ged = d.gedc->getDataObject(globalDataKey_);
173 auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
174 auto locPair = rcp_dynamic_cast<LPGED>(ged);
175 if (not locPair.is_null())
177 RCP<LOC> loc = locPair->getGhostedLOC();
178 epetraContainer = rcp_dynamic_cast<ELOC>(loc);
180 if (not epetraContainer.is_null())
182 if (useTimeDerivativeSolutionVector_)
183 x_ = epetraContainer->get_dxdt();
185 x_ = epetraContainer->get_x();
192 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
200 template<
typename TRAITS,
typename LO,
typename GO>
204 typename TRAITS::EvalData workset)
211 using Teuchos::ptrFromRef;
213 using Teuchos::rcp_dynamic_cast;
214 using Thyra::SpmdVectorBase;
217 string blockId(this->wda(workset).block_id);
218 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
219 int numCells(localCellIds.size()),
numFields(gatherFields_.size());
225 auto LIDs = globalIndexer_->getLIDs();
226 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
227 Kokkos::deep_copy(LIDs_h, LIDs);
229 for (
int fieldInd(0); fieldInd <
numFields; ++fieldInd)
231 MDField<ScalarT, Cell, NODE>&
field = gatherFields_[fieldInd];
232 auto field_h = Kokkos::create_mirror_view(field.get_static_view());
233 int fieldNum(fieldIds_[fieldInd]);
234 const vector<int>& elmtOffset =
235 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
236 int numBases(elmtOffset.size());
238 for (
int cell(0); cell < numCells; ++cell)
240 size_t cellLocalId(localCellIds[cell]);
242 for (
int basis(0); basis < numBases; ++basis)
244 int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
246 field_h(cell, basis) = (*xEvRoGed_)[lid];
248 field_h(cell, basis) = (*x_)[lid];
251 Kokkos::deep_copy(field.get_static_view(), field_h);
261 template<
typename TRAITS,
typename LO,
typename GO>
267 globalIndexer_(indexer),
268 hasTangentFields_(false)
277 using vvstring = std::vector<std::vector<std::string>>;
293 MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
294 this->addEvaluatedField(gatherFields_[fd]);
298 if (tangentFieldNames.size() > 0)
301 hasTangentFields_ =
true;
302 tangentFields_.resize(numFields);
305 int numTangentFields(tangentFieldNames[fd].size());
306 tangentFields_[fd].resize(numTangentFields);
307 for (
int i(0); i < numTangentFields; ++i)
309 tangentFields_[fd][i] =
310 MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
312 this->addDependentField(tangentFields_[fd][i]);
318 string firstName(
"<none>");
320 firstName = names[0];
321 string n(
"GatherSolution (Epetra): " + firstName +
" (" +
322 print<EvalT>() +
")");
331 template<
typename TRAITS,
typename LO,
typename GO>
335 typename TRAITS::SetupData ,
338 using std::logic_error;
347 const string& fieldName(indexerNames_[fd]);
348 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
352 "GatherSolution_Epetra<Tangent>: Could not find field \"" + fieldName
353 +
"\" in the global indexer. ")
355 indexerNames_.clear();
363 template<
typename TRAITS,
typename LO,
typename GO>
367 typename TRAITS::PreEvalData d)
371 using Teuchos::rcp_dynamic_cast;
380 string post(useTimeDerivativeSolutionVector_ ?
" - Xdot" :
" - X");
381 if (d.gedc->containsDataObject(globalDataKey_ + post))
383 ged = d.gedc->getDataObject(globalDataKey_ + post);
384 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
389 ged = d.gedc->getDataObject(globalDataKey_);
392 auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
393 auto locPair = rcp_dynamic_cast<LPGED>(ged);
394 if (not locPair.is_null())
396 RCP<LOC> loc = locPair->getGhostedLOC();
397 epetraContainer = rcp_dynamic_cast<ELOC>(loc);
399 if (not epetraContainer.is_null())
401 if (useTimeDerivativeSolutionVector_)
402 x_ = epetraContainer->get_dxdt();
404 x_ = epetraContainer->get_x();
411 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
419 template<
typename TRAITS,
typename LO,
typename GO>
423 typename TRAITS::EvalData workset)
430 using Teuchos::ptrFromRef;
432 using Teuchos::rcp_dynamic_cast;
433 using Thyra::SpmdVectorBase;
436 string blockId(this->wda(workset).block_id);
437 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
438 int numCells(localCellIds.size()),
numFields(gatherFields_.size());
443 auto LIDs = globalIndexer_->getLIDs();
444 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
445 Kokkos::deep_copy(LIDs_h, LIDs);
447 for (
int fieldInd(0); fieldInd <
numFields; ++fieldInd)
449 MDField<ScalarT, Cell, NODE>&
field = gatherFields_[fieldInd];
450 auto field_h = Kokkos::create_mirror_view(field.get_static_view());
451 int fieldNum(fieldIds_[fieldInd]);
452 const vector<int>& elmtOffset =
453 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
454 int numBases(elmtOffset.size());
456 for (
int cell(0); cell < numCells; ++cell)
458 size_t cellLocalId(localCellIds[cell]);
460 for (
int basis(0); basis < numBases; ++basis)
462 int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
464 field_h(cell, basis) = (*xEvRoGed_)[lid];
466 field_h(cell, basis) = (*x_)[lid];
471 if (hasTangentFields_) {
473 int numTangentFields(tangentFields_[fieldInd].size());
474 for (
int i(0); i < numTangentFields; ++i){
475 auto tf = Kokkos::create_mirror_view(tangentFields_[fieldInd][i].get_static_view());
476 Kokkos::deep_copy(tf, tangentFields_[fieldInd][i].get_static_view());
478 for (
int cell(0); cell < numCells; ++cell) {
480 int fieldNum(fieldIds_[fieldInd]);
481 const vector<int>& elmtOffset =
482 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
483 int numBases(elmtOffset.size());
486 for (
int basis(0); basis < numBases; ++basis){
487 field_h(cell, basis).fastAccessDx(i) =
488 tf(cell, basis).val();
493 Kokkos::deep_copy(field.get_static_view(), field_h);
502 template<
typename TRAITS,
typename LO,
typename GO>
508 globalIndexer_(indexer)
533 MDField<ScalarT, Cell, NODE> f(names[fd], basis->functional);
534 gatherFields_[fd] = f;
535 this->addEvaluatedField(gatherFields_[fd]);
539 string firstName(
"<none>"), n(
"GatherSolution (Epetra");
541 firstName = names[0];
542 if (disableSensitivities_)
543 n +=
", No Sensitivities";
544 n +=
"): " + firstName +
" (Jacobian)";
553 template<
typename TRAITS,
typename LO,
typename GO>
557 typename TRAITS::SetupData ,
560 using std::logic_error;
569 const string& fieldName(indexerNames_[fd]);
570 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
574 "GatherSolution_Epetra<Jacobian>: Could not find field \"" +
575 fieldName +
"\" in the global indexer. ")
577 indexerNames_.clear();
585 template<
typename TRAITS,
typename LO,
typename GO>
589 typename TRAITS::PreEvalData d)
593 using Teuchos::rcp_dynamic_cast;
601 applySensitivities_ =
false;
602 if ((not disableSensitivities_ ) and
603 (d.first_sensitivities_name == sensitivitiesName_))
604 applySensitivities_ =
true;
608 string post(useTimeDerivativeSolutionVector_ ?
" - Xdot" :
" - X");
609 if (d.gedc->containsDataObject(globalDataKey_ + post))
611 ged = d.gedc->getDataObject(globalDataKey_ + post);
612 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
617 ged = d.gedc->getDataObject(globalDataKey_);
620 auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
621 auto locPair = rcp_dynamic_cast<LPGED>(ged);
622 if (not locPair.is_null())
624 RCP<LOC> loc = locPair->getGhostedLOC();
625 epetraContainer = rcp_dynamic_cast<ELOC>(loc);
627 if (not epetraContainer.is_null())
629 if (useTimeDerivativeSolutionVector_)
630 x_ = epetraContainer->get_dxdt();
632 x_ = epetraContainer->get_x();
639 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
647 template<
typename TRAITS,
typename LO,
typename GO>
651 typename TRAITS::EvalData workset)
658 using Teuchos::ptrFromRef;
660 using Teuchos::rcp_dynamic_cast;
661 using Thyra::SpmdVectorBase;
664 string blockId(this->wda(workset).block_id);
665 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
666 int numFields(gatherFields_.size()), numCells(localCellIds.size());
670 if (applySensitivities_)
672 if ((useTimeDerivativeSolutionVector_) and (gatherSeedIndex_ < 0))
673 seedValue = workset.alpha;
674 else if (gatherSeedIndex_ < 0)
675 seedValue = workset.beta;
676 else if (not useTimeDerivativeSolutionVector_)
677 seedValue = workset.gather_seeds[gatherSeedIndex_];
685 if (not applySensitivities_)
692 if (this->wda.getDetailsIndex() == 1)
695 dos = globalIndexer_->getElementBlockGIDCount(workset.details(0).block_id);
701 auto LIDs = globalIndexer_->getLIDs();
702 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
703 Kokkos::deep_copy(LIDs_h, LIDs);
705 for (
int fieldInd(0); fieldInd <
numFields; ++fieldInd)
707 MDField<ScalarT, Cell, NODE>&
field = gatherFields_[fieldInd];
708 auto field_h = Kokkos::create_mirror_view(field.get_static_view());
709 int fieldNum(fieldIds_[fieldInd]);
710 const vector<int>& elmtOffset =
711 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
712 int numBases(elmtOffset.size());
714 for (
int cell(0); cell < numCells; ++cell)
716 size_t cellLocalId(localCellIds[cell]);
718 for (
int basis(0); basis < numBases; ++basis)
720 int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
722 field_h(cell, basis) = (*xEvRoGed_)[lid];
724 field_h(cell, basis) = (*x_)[lid];
729 if (applySensitivities_) {
730 int fieldNum(fieldIds_[fieldInd]);
731 const vector<int>& elmtOffset =
732 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
733 int numBases(elmtOffset.size());
736 for (
int cell(0); cell < numCells; ++cell)
739 for (
int basis(0); basis < numBases; ++basis)
742 int offset(elmtOffset[basis]);
744 field_h(cell, basis).fastAccessDx(dos + offset) = seedValue;
748 Kokkos::deep_copy(field.get_static_view(), field_h);
753 #endif // __Panzer_GatherSolution_Epetra_impl_hpp__
Gathers solution values from the Newton solution vector into the nodal fields of the field manager...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
This class provides a boundary exchange communication mechanism for vectors.
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral...
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)