43 #ifndef __Panzer_GatherTangent_Epetra_impl_hpp__
44 #define __Panzer_GatherTangent_Epetra_impl_hpp__
53 #include "Epetra_Vector.h"
63 #include "Teuchos_Assert.hpp"
66 #include "Thyra_SpmdVectorBase.hpp"
73 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
79 globalIndexer_(indexer),
80 useTimeDerivativeSolutionVector_(false),
81 globalDataKey_(
"Tangent Gather Container")
99 if (p.
isType<
bool>(
"Use Time Derivative Solution Vector"))
101 p.
get<
bool>(
"Use Time Derivative Solution Vector");
102 if (p.
isType<
string>(
"Global Data Key"))
111 MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
121 string firstName(
"<none>");
123 firstName = names[0];
124 string n(
"GatherTangent (Epetra): " + firstName +
" (" +
125 print<EvalT>() +
")");
134 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
138 typename TRAITS::SetupData ,
141 using std::logic_error;
151 const string& fieldName((*indexerNames_)[fd]);
152 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
156 "GatherTangent_Epetra<Residual>: Could not find field \"" + fieldName +
157 "\" in the global indexer. ");
159 indexerNames_ = null;
167 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
171 typename TRAITS::PreEvalData d)
174 using Teuchos::rcp_dynamic_cast;
177 if (d.gedc->containsDataObject(globalDataKey_))
179 RCP<GED> ged = d.gedc->getDataObject(globalDataKey_);
180 dxdpEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged,
true);
189 template<
typename EvalT,
typename TRAITS,
typename LO,
typename GO>
193 typename TRAITS::EvalData workset)
200 using Teuchos::ptrFromRef;
202 using Teuchos::rcp_dynamic_cast;
203 using Thyra::SpmdVectorBase;
207 if (dxdpEvRoGed_.is_null())
211 string blockId(this->wda(workset).block_id);
212 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
213 int numCells(localCellIds.size()),
numFields(gatherFields_.size());
220 for (
int cell(0); cell < numCells; ++cell)
222 size_t cellLocalId(localCellIds[cell]);
223 auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
226 for (
int fieldIndex(0); fieldIndex <
numFields; ++fieldIndex)
228 MDField<ScalarT, Cell, NODE>&
field = gatherFields_[fieldIndex];
229 int fieldNum(fieldIds_[fieldIndex]);
230 const vector<int>& elmtOffset =
231 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
232 int numBases(elmtOffset.size());
235 for (
int basis(0); basis < numBases; ++basis)
237 int offset(elmtOffset[basis]), lid(LIDs[offset]);
238 field(cell, basis) = (*dxdpEvRoGed_)[lid];
244 #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 .