11 #ifndef __Panzer_GatherTangent_Epetra_impl_hpp__
12 #define __Panzer_GatherTangent_Epetra_impl_hpp__
21 #include "Epetra_Vector.h"
31 #include "Teuchos_Assert.hpp"
34 #include "Thyra_SpmdVectorBase.hpp"
41 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
47 globalIndexer_(indexer),
48 useTimeDerivativeSolutionVector_(false),
49 globalDataKey_(
"Tangent Gather Container")
67 if (p.
isType<
bool>(
"Use Time Derivative Solution Vector"))
69 p.
get<
bool>(
"Use Time Derivative Solution Vector");
70 if (p.
isType<
string>(
"Global Data Key"))
79 MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
89 string firstName(
"<none>");
92 string n(
"GatherTangent (Epetra): " + firstName +
" (" +
93 print<EvalT>() +
")");
102 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
106 typename TRAITS::SetupData ,
109 using std::logic_error;
119 const string& fieldName((*indexerNames_)[fd]);
120 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
124 "GatherTangent_Epetra<Residual>: Could not find field \"" + fieldName +
125 "\" in the global indexer. ");
127 indexerNames_ = null;
135 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
139 typename TRAITS::PreEvalData d)
142 using Teuchos::rcp_dynamic_cast;
145 if (d.gedc->containsDataObject(globalDataKey_))
147 RCP<GED> ged = d.gedc->getDataObject(globalDataKey_);
148 dxdpEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
157 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
161 typename TRAITS::EvalData workset)
168 using Teuchos::ptrFromRef;
170 using Teuchos::rcp_dynamic_cast;
171 using Thyra::SpmdVectorBase;
175 if (dxdpEvRoGed_.is_null())
179 string blockId(this->wda(workset).block_id);
180 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
181 int numCells(localCellIds.size()),
numFields(gatherFields_.size());
188 auto LIDs = globalIndexer_->getLIDs();
189 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
190 Kokkos::deep_copy(LIDs_h, LIDs);
192 for (
int fieldIndex(0); fieldIndex <
numFields; ++fieldIndex)
194 MDField<ScalarT, Cell, NODE>&
field = gatherFields_[fieldIndex];
195 auto field_h = Kokkos::create_mirror_view(field.get_static_view());
196 for (
int cell(0); cell < numCells; ++cell)
198 size_t cellLocalId(localCellIds[cell]);
199 int fieldNum(fieldIds_[fieldIndex]);
200 const vector<int>& elmtOffset =
201 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
202 int numBases(elmtOffset.size());
205 for (
int basis(0); basis < numBases; ++basis)
207 int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
208 field_h(cell, basis) = (*dxdpEvRoGed_)[lid];
211 Kokkos::deep_copy(field.get_static_view(), field_h);
215 #endif // __Panzer_GatherTangent_Epetra_impl_hpp__
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
The fields to be gathered.
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
Create a copy.
T & get(const std::string &name, T def_value)
std::string globalDataKey_
The key identifying the GlobalEvaluationData.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< std::vector< std::string > > indexerNames_
A list of the names of the fields to be gathered.
GatherTangent_Epetra()
Default Constructor (disabled).
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
Post-Registration Setup.
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...
bool isType(const std::string &name) const
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)
void preEvaluate(typename TRAITS::PreEvalData d)
Pre-Evaluate: Sets the tangent vector.
void evaluateFields(typename TRAITS::EvalData d)
Evaluate Fields: Gather operation.
bool useTimeDerivativeSolutionVector_
A flag indicating whether we should be working with or .