43 #ifndef   __Panzer_GatherSolution_Epetra_impl_hpp__ 
   44 #define   __Panzer_GatherSolution_Epetra_impl_hpp__ 
   53 #include "Epetra_Vector.h" 
   65 #include "Teuchos_Assert.hpp" 
   68 #include "Thyra_SpmdVectorBase.hpp" 
   75 template<
typename TRAITS, 
typename LO, 
typename GO>
 
   81   globalIndexer_(indexer),
 
   82   hasTangentFields_(false)
 
   91   using vvstring = std::vector<std::vector<std::string>>;
 
  107       MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
 
  108     this->addEvaluatedField(gatherFields_[fd]);
 
  112   if (tangentFieldNames.size() > 0)
 
  115     hasTangentFields_ = 
true;
 
  116     tangentFields_.resize(numFields);
 
  119       int numTangentFields(tangentFieldNames[fd].size());
 
  120       tangentFields_[fd].resize(numTangentFields);
 
  121       for (
int i(0); i < numTangentFields; ++i)
 
  123         tangentFields_[fd][i] =
 
  124           MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
 
  126         this->addDependentField(tangentFields_[fd][i]);
 
  132   string firstName(
"<none>");
 
  134     firstName = names[0];
 
  135   string n(
"GatherSolution (Epetra): " + firstName + 
" (" +
 
  136     print<EvalT>() + 
")");
 
  145 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  149   typename TRAITS::SetupData ,
 
  152   using std::logic_error;
 
  161     const string& fieldName(indexerNames_[fd]);
 
  162     fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
 
  166       "GatherSolution_Epetra<Residual>: Could not find field \"" +
 
  167       fieldName + 
"\" in the global indexer. ")
 
  169   indexerNames_.clear();
 
  177 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  181   typename TRAITS::PreEvalData d)
 
  185   using Teuchos::rcp_dynamic_cast;
 
  194   string post(useTimeDerivativeSolutionVector_ ? 
" - Xdot" : 
" - X");
 
  195   if (d.gedc->containsDataObject(globalDataKey_ + post))
 
  197     ged       = d.gedc->getDataObject(globalDataKey_ + post);
 
  198     xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, 
true);
 
  203   ged = d.gedc->getDataObject(globalDataKey_);
 
  206     auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
 
  207     auto locPair         = rcp_dynamic_cast<LPGED>(ged);
 
  208     if (not locPair.is_null())
 
  210       RCP<LOC> loc = locPair->getGhostedLOC();
 
  211       epetraContainer = rcp_dynamic_cast<ELOC>(loc);
 
  213     if (not epetraContainer.is_null())
 
  215       if (useTimeDerivativeSolutionVector_)
 
  216         x_ = epetraContainer->get_dxdt();
 
  218         x_ = epetraContainer->get_x();
 
  225   xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, 
true);
 
  233 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  237   typename TRAITS::EvalData workset)
 
  244   using Teuchos::ptrFromRef;
 
  246   using Teuchos::rcp_dynamic_cast;
 
  247   using Thyra::SpmdVectorBase;
 
  250   string blockId(this->wda(workset).block_id);
 
  251   const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
 
  252   int numCells(localCellIds.size()), 
numFields(gatherFields_.size());
 
  260     for (
int cell(0); cell < numCells; ++cell)
 
  262       size_t cellLocalId(localCellIds[cell]);
 
  263       auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
 
  266       for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  268         MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  269         int fieldNum(fieldIds_[fieldInd]);
 
  270         const vector<int>& elmtOffset =
 
  271           globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  272         int numBases(elmtOffset.size());
 
  275         for (
int basis(0); basis < numBases; ++basis)
 
  277           int offset(elmtOffset[basis]), lid(LIDs[offset]);
 
  278           field(cell, basis) = (*xEvRoGed_)[lid];
 
  286     for (
int cell(0); cell < numCells; ++cell)
 
  288       size_t cellLocalId(localCellIds[cell]);
 
  289       auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
 
  292       for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  294         MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  295         int fieldNum(fieldIds_[fieldInd]);
 
  296         const vector<int>& elmtOffset =
 
  297           globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  298         int numBases(elmtOffset.size());
 
  301         for (
int basis(0); basis < numBases; ++basis)
 
  303           int offset(elmtOffset[basis]), lid(LIDs[offset]);
 
  304           field(cell, basis) = (*x_)[lid];
 
  316 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  322   globalIndexer_(indexer),
 
  323   hasTangentFields_(false)
 
  332   using vvstring = std::vector<std::vector<std::string>>;
 
  348       MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
 
  349     this->addEvaluatedField(gatherFields_[fd]);
 
  353   if (tangentFieldNames.size() > 0)
 
  356     hasTangentFields_ = 
true;
 
  357     tangentFields_.resize(numFields);
 
  360       int numTangentFields(tangentFieldNames[fd].size());
 
  361       tangentFields_[fd].resize(numTangentFields);
 
  362       for (
int i(0); i < numTangentFields; ++i)
 
  364         tangentFields_[fd][i] =
 
  365           MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
 
  367         this->addDependentField(tangentFields_[fd][i]);
 
  373   string firstName(
"<none>");
 
  375     firstName = names[0];
 
  376   string n(
"GatherSolution (Epetra): " + firstName + 
" (" +
 
  377     print<EvalT>() + 
")");
 
  386 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  390   typename TRAITS::SetupData ,
 
  393   using std::logic_error;
 
  402     const string& fieldName(indexerNames_[fd]);
 
  403     fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
 
  407       "GatherSolution_Epetra<Tangent>:  Could not find field \"" + fieldName
 
  408       + 
"\" in the global indexer. ")
 
  410   indexerNames_.clear();
 
  418 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  422   typename TRAITS::PreEvalData d)
 
  426   using Teuchos::rcp_dynamic_cast;
 
  435   string post(useTimeDerivativeSolutionVector_ ? 
" - Xdot" : 
" - X");
 
  436   if (d.gedc->containsDataObject(globalDataKey_ + post))
 
  438     ged       = d.gedc->getDataObject(globalDataKey_ + post);
 
  439     xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, 
true);
 
  444   ged = d.gedc->getDataObject(globalDataKey_);
 
  447     auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
 
  448     auto locPair         = rcp_dynamic_cast<LPGED>(ged);
 
  449     if (not locPair.is_null())
 
  451       RCP<LOC> loc = locPair->getGhostedLOC();
 
  452       epetraContainer = rcp_dynamic_cast<ELOC>(loc);
 
  454     if (not epetraContainer.is_null())
 
  456       if (useTimeDerivativeSolutionVector_)
 
  457         x_ = epetraContainer->get_dxdt();
 
  459         x_ = epetraContainer->get_x();
 
  466   xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, 
true);
 
  474 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  478   typename TRAITS::EvalData workset)
 
  485   using Teuchos::ptrFromRef;
 
  487   using Teuchos::rcp_dynamic_cast;
 
  488   using Thyra::SpmdVectorBase;
 
  491   string blockId(this->wda(workset).block_id);
 
  492   const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
 
  493   int numCells(localCellIds.size()), 
numFields(gatherFields_.size());
 
  501     for (
int cell(0); cell < numCells; ++cell)
 
  503       size_t cellLocalId(localCellIds[cell]);
 
  504       auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
 
  507       for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  509         MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  510         int fieldNum(fieldIds_[fieldInd]);
 
  511         const vector<int>& elmtOffset =
 
  512           globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  513         int numBases(elmtOffset.size());
 
  516         for (
int basis(0); basis < numBases; ++basis)
 
  518           int offset(elmtOffset[basis]), lid(LIDs[offset]);
 
  519           field(cell, basis) = (*xEvRoGed_)[lid];
 
  527     for (
int cell(0); cell < numCells; ++cell)
 
  529       size_t cellLocalId(localCellIds[cell]);
 
  530       auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
 
  533       for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  535         MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  536         int fieldNum(fieldIds_[fieldInd]);
 
  537         const vector<int>& elmtOffset =
 
  538           globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  539         int numBases(elmtOffset.size());
 
  542         for (
int basis(0); basis < numBases; ++basis)
 
  544           int offset(elmtOffset[basis]), lid(LIDs[offset]);
 
  545           field(cell, basis) = (*x_)[lid];
 
  552   if (hasTangentFields_)
 
  555     for (
int cell(0); cell < numCells; ++cell)
 
  558       for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  560         MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  561         int fieldNum(fieldIds_[fieldInd]);
 
  562         const vector<int>& elmtOffset =
 
  563           globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  564         int numBases(elmtOffset.size());
 
  567         for (
int basis(0); basis < numBases; ++basis)
 
  570           int numTangentFields(tangentFields_[fieldInd].size());
 
  571           for (
int i(0); i < numTangentFields; ++i)
 
  572             field(cell, basis).fastAccessDx(i) =
 
  573               tangentFields_[fieldInd][i](cell, basis).val();
 
  585 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  591   globalIndexer_(indexer)
 
  616     MDField<ScalarT, Cell, NODE> f(names[fd], basis->functional);
 
  617     gatherFields_[fd] = f;
 
  618     this->addEvaluatedField(gatherFields_[fd]);
 
  622   string firstName(
"<none>"), n(
"GatherSolution (Epetra");
 
  624     firstName = names[0];
 
  625   if (disableSensitivities_)
 
  626     n += 
", No Sensitivities";
 
  627   n += 
"): " + firstName + 
" (" + print<EvalT>() + 
")";
 
  636 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  640   typename TRAITS::SetupData ,
 
  643   using std::logic_error;
 
  652     const string& fieldName(indexerNames_[fd]);
 
  653     fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
 
  657       "GatherSolution_Epetra<Jacobian>: Could not find field \"" +
 
  658       fieldName + 
"\" in the global indexer. ")
 
  660   indexerNames_.clear();
 
  668 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  672   typename TRAITS::PreEvalData d)
 
  676   using Teuchos::rcp_dynamic_cast;
 
  684   applySensitivities_ = 
false;
 
  685   if ((not disableSensitivities_                       )  and
 
  686       (d.first_sensitivities_name == sensitivitiesName_))
 
  687     applySensitivities_ = 
true;
 
  691   string post(useTimeDerivativeSolutionVector_ ? 
" - Xdot" : 
" - X");
 
  692   if (d.gedc->containsDataObject(globalDataKey_ + post))
 
  694     ged       = d.gedc->getDataObject(globalDataKey_ + post);
 
  695     xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, 
true);
 
  700   ged = d.gedc->getDataObject(globalDataKey_);
 
  703     auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
 
  704     auto locPair         = rcp_dynamic_cast<LPGED>(ged);
 
  705     if (not locPair.is_null())
 
  707       RCP<LOC> loc = locPair->getGhostedLOC();
 
  708       epetraContainer = rcp_dynamic_cast<ELOC>(loc);
 
  710     if (not epetraContainer.is_null())
 
  712       if (useTimeDerivativeSolutionVector_)
 
  713         x_ = epetraContainer->get_dxdt();
 
  715         x_ = epetraContainer->get_x();
 
  722   xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, 
true);
 
  730 template<
typename TRAITS, 
typename LO, 
typename GO>
 
  734   typename TRAITS::EvalData workset)
 
  741   using Teuchos::ptrFromRef;
 
  743   using Teuchos::rcp_dynamic_cast;
 
  744   using Thyra::SpmdVectorBase;
 
  747   string blockId(this->wda(workset).block_id);
 
  748   const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
 
  749   int numFields(gatherFields_.size()), numCells(localCellIds.size());
 
  753   if (applySensitivities_)
 
  755     if ((useTimeDerivativeSolutionVector_) and (gatherSeedIndex_ < 0))
 
  756       seedValue = workset.alpha;
 
  757     else if (gatherSeedIndex_ < 0)
 
  758       seedValue = workset.beta;
 
  759     else if (not useTimeDerivativeSolutionVector_)
 
  760       seedValue = workset.gather_seeds[gatherSeedIndex_];
 
  768   if (not applySensitivities_)
 
  775   if (this->wda.getDetailsIndex() == 1)
 
  778     dos = globalIndexer_->getElementBlockGIDCount(workset.details(0).block_id);
 
  787     for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  789       MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  790       int fieldNum(fieldIds_[fieldInd]);
 
  791       const vector<int>& elmtOffset =
 
  792         globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  793       int numBases(elmtOffset.size());
 
  796       for (
int cell(0); cell < numCells; ++cell)
 
  798         size_t cellLocalId(localCellIds[cell]);
 
  799   auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
 
  802         for (
int basis(0); basis < numBases; ++basis)
 
  804           int offset(elmtOffset[basis]), lid(LIDs[offset]);
 
  805           field(cell, basis) = (*xEvRoGed_)[lid];
 
  813     for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  815       MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  816       int fieldNum(fieldIds_[fieldInd]);
 
  817       const vector<int>& elmtOffset =
 
  818         globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  819       int numBases(elmtOffset.size());
 
  822       for (
int cell(0); cell < numCells; ++cell)
 
  824         size_t cellLocalId(localCellIds[cell]);
 
  825   auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
 
  828         for (
int basis(0); basis < numBases; ++basis)
 
  830           int offset(elmtOffset[basis]), lid(LIDs[offset]);
 
  831           field(cell, basis) = (*x_)[lid];
 
  838   if (applySensitivities_)
 
  841     for (
int fieldInd(0); fieldInd < 
numFields; ++fieldInd)
 
  843       MDField<ScalarT, Cell, NODE>& 
field = gatherFields_[fieldInd];
 
  844       int fieldNum(fieldIds_[fieldInd]);
 
  845       const vector<int>& elmtOffset =
 
  846         globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
 
  847       int numBases(elmtOffset.size());
 
  850       for (
int cell(0); cell < numCells; ++cell)
 
  853         for (
int basis(0); basis < numBases; ++basis)
 
  856           int offset(elmtOffset[basis]);
 
  857           field(cell, basis).fastAccessDx(dos + offset) = seedValue;
 
  864 #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)