43#ifndef PANZER_GATHER_SOLUTION_BLOCKED_EPETRA_IMPL_HPP
44#define PANZER_GATHER_SOLUTION_BLOCKED_EPETRA_IMPL_HPP
46#include "Teuchos_Assert.hpp"
47#include "Phalanx_DataLayout.hpp"
52#include "Panzer_TpetraLinearObjFactory.hpp"
57#include "Teuchos_FancyOStream.hpp"
59#include "Thyra_SpmdVectorBase.hpp"
60#include "Thyra_ProductVectorBase.hpp"
62#include "Tpetra_Map.hpp"
64template <
typename EvalT,
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
67 const Teuchos::RCP<const BlockedDOFManager> & indexer,
68 const Teuchos::ParameterList& p)
70 const std::vector<std::string>& names =
71 *(p.get< Teuchos::RCP< std::vector<std::string> > >(
"DOF Names"));
73 Teuchos::RCP<panzer::PureBasis> basis =
74 p.get< Teuchos::RCP<panzer::PureBasis> >(
"Basis");
76 for (std::size_t fd = 0; fd < names.size(); ++fd) {
77 PHX::MDField<ScalarT,Cell,NODE> field = PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
78 this->addEvaluatedField(field.fieldTag());
81 this->setName(
"Gather Solution");
88template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
91 const Teuchos::RCP<const BlockedDOFManager> & indexer,
92 const Teuchos::ParameterList& p)
96 typedef std::vector< std::vector<std::string> > vvstring;
101 const std::vector<std::string> & names = input.
getDofNames();
102 Teuchos::RCP<const panzer::PureBasis> basis = input.
getBasis();
111 for (std::size_t fd = 0; fd < names.size(); ++fd) {
113 PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
118 if (tangent_field_names.size()>0) {
119 TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
121 has_tangent_fields_ = true;
122 tangentFields_.resize(gatherFields_.size());
123 for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
124 tangentFields_[fd].resize(tangent_field_names[fd].size());
125 for (std::size_t i=0; i<tangent_field_names[fd].size(); ++i) {
126 tangentFields_[fd][i] =
127 PHX::MDField<const ScalarT,Cell,NODE>(tangent_field_names[fd][i],basis->functional);
128 this->addDependentField(tangentFields_[fd][i]);
134 std::string firstName =
"<none>";
136 firstName = names[0];
138 std::string n =
"GatherSolution (BlockedTpetra): "+firstName+
" (Residual)";
143template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
150 const Workset & workset_0 = (*d.worksets_)[0];
151 const std::string blockId = this->
wda(workset_0).block_id;
157 int maxElementBlockGIDCount = -1;
161 const int globalFieldNum =
globalIndexer_->getFieldNum(fieldName);
167 fieldOffsets_[fd] = PHX::View<int*>(
"GatherSolution_BlockedTpetra(Residual):fieldOffsets",offsets.size());
168 auto hostFieldOffsets = Kokkos::create_mirror_view(
fieldOffsets_[fd]);
169 for(std::size_t i=0; i < offsets.size(); ++i)
170 hostFieldOffsets(i) = offsets[i];
173 maxElementBlockGIDCount = std::max(
fieldGlobalIndexers_[fd]->getElementBlockGIDCount(blockId),maxElementBlockGIDCount);
179 worksetLIDs_ = PHX::View<LO**>(
"GatherSolution_BlockedTpetra(Residual):worksetLIDs",
181 maxElementBlockGIDCount);
187template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
196template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
201 using Teuchos::rcp_dynamic_cast;
205 const PHX::View<const int*>& localCellIds = this->
wda(workset).cell_local_ids_k;
207 RCP<ProductVectorBase<ScalarT>> thyraBlockSolution;
209 thyraBlockSolution = rcp_dynamic_cast<ProductVectorBase<ScalarT>>(
blockedContainer_->get_dxdt(),
true);
211 thyraBlockSolution = rcp_dynamic_cast<ProductVectorBase<ScalarT>>(
blockedContainer_->get_x(),
true);
214 int currentWorksetLIDSubBlock = -1;
215 for (std::size_t fieldIndex = 0; fieldIndex <
gatherFields_.size(); fieldIndex++) {
218 const std::string blockId = this->
wda(workset).block_id;
224 const auto& tpetraSolution = *((rcp_dynamic_cast<Thyra::TpetraVector<ScalarT,LO,GO,NodeT>>(thyraBlockSolution->getNonconstVectorBlock(
productVectorBlockIndex_[fieldIndex]),
true))->getTpetraVector());
225 const auto& kokkosSolution = tpetraSolution.getLocalViewDevice(Tpetra::Access::ReadOnly);
232 Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device>(0,workset.num_cells), KOKKOS_LAMBDA (
const int& cell) {
233 for(
int basis=0; basis < static_cast<int>(fieldOffsets.size()); ++basis) {
234 const int lid = worksetLIDs(cell,fieldOffsets(basis));
235 fieldValues(cell,basis) = kokkosSolution(lid,0);
246template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
249 const Teuchos::RCP<const BlockedDOFManager> & indexer,
250 const Teuchos::ParameterList& p)
254 typedef std::vector< std::vector<std::string> > vvstring;
259 const std::vector<std::string> & names = input.
getDofNames();
260 Teuchos::RCP<const panzer::PureBasis> basis = input.
getBasis();
269 for (std::size_t fd = 0; fd < names.size(); ++fd) {
271 PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
276 if (tangent_field_names.size()>0) {
277 TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
279 has_tangent_fields_ = true;
280 tangentFields_.resize(gatherFields_.size());
281 for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
282 tangentFields_[fd].resize(tangent_field_names[fd].size());
283 for (std::size_t i=0; i<tangent_field_names[fd].size(); ++i) {
284 tangentFields_[fd][i] =
285 PHX::MDField<const ScalarT,Cell,NODE>(tangent_field_names[fd][i],basis->functional);
286 this->addDependentField(tangentFields_[fd][i]);
292 std::string firstName =
"<none>";
294 firstName = names[0];
296 std::string n =
"GatherSolution (BlockedTpetra): "+firstName+
" (Tangent)";
301template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
320template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
329template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
334 using Teuchos::ArrayRCP;
335 using Teuchos::ptrFromRef;
336 using Teuchos::rcp_dynamic_cast;
339 using Thyra::SpmdVectorBase;
342 Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
343 out.setShowProcRank(
true);
344 out.setOutputToRootOnly(-1);
346 std::vector<std::pair<int,GO> > GIDs;
347 std::vector<LO> LIDs;
350 std::string blockId = this->
wda(workset).block_id;
351 const std::vector<std::size_t> & localCellIds = this->
wda(workset).cell_local_ids;
353 Teuchos::RCP<ProductVectorBase<double> > x;
360 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
361 LO cellLocalId = localCellIds[worksetCellIndex];
363 gidIndexer_->getElementGIDsPair(cellLocalId,GIDs,blockId);
366 LIDs.resize(GIDs.size());
367 for(std::size_t i=0;i<GIDs.size();i++) {
371 LIDs[i] = x_map->getLocalElement(GIDs[i].second);
375 Teuchos::ArrayRCP<const double> local_x;
376 for (std::size_t fieldIndex=0; fieldIndex<
gatherFields_.size();fieldIndex++) {
378 int indexerId =
gidIndexer_->getFieldBlock(fieldNum);
381 RCP<SpmdVectorBase<double> > block_x = rcp_dynamic_cast<SpmdVectorBase<double> >(x->getNonconstVectorBlock(indexerId));
382 block_x->getLocalData(ptrFromRef(local_x));
384 const std::vector<int> & elmtOffset =
gidIndexer_->getGIDFieldOffsets(blockId,fieldNum);
387 for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
388 int offset = elmtOffset[basis];
389 int lid = LIDs[offset];
392 (
gatherFields_[fieldIndex])(worksetCellIndex,basis) = local_x[lid];
394 (
gatherFields_[fieldIndex])(worksetCellIndex,basis).val() = local_x[lid];
396 (
gatherFields_[fieldIndex])(worksetCellIndex,basis).fastAccessDx(i) =
408template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
411 const Teuchos::RCP<const BlockedDOFManager> & indexer,
412 const Teuchos::ParameterList& p)
418 const std::vector<std::string> & names = input.
getDofNames();
419 Teuchos::RCP<const panzer::PureBasis> basis = input.
getBasis();
428 for (std::size_t fd = 0; fd < names.size(); ++fd) {
429 PHX::MDField<ScalarT,Cell,NODE> f(names[fd],basis->functional);
435 std::string firstName =
"<none>";
437 firstName = names[0];
441 std::string n =
"GatherSolution (BlockedTpetra, No Sensitivities): "+firstName+
" (Jacobian)";
445 std::string n =
"GatherSolution (BlockedTpetra): "+firstName+
" ("+PHX::print<EvalT>()+
") ";
451template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
458 const Workset & workset_0 = (*d.worksets_)[0];
459 const std::string blockId = this->
wda(workset_0).block_id;
464 int maxElementBlockGIDCount = -1;
468 const int globalFieldNum =
globalIndexer_->getFieldNum(fieldName);
471 fieldIds_[fd] = subGlobalIndexer->getFieldNum(fieldName);
473 const std::vector<int>& offsets = subGlobalIndexer->getGIDFieldOffsets(blockId,
fieldIds_[fd]);
474 fieldOffsets_[fd] = PHX::View<int*>(
"GatherSolution_BlockedTpetra(Jacobian):fieldOffsets",offsets.size());
475 auto hostOffsets = Kokkos::create_mirror_view(
fieldOffsets_[fd]);
476 for (std::size_t i=0; i < offsets.size(); ++i)
477 hostOffsets(i) = offsets[i];
479 maxElementBlockGIDCount = std::max(subGlobalIndexer->getElementBlockGIDCount(blockId),maxElementBlockGIDCount);
485 worksetLIDs_ = PHX::View<LO**>(
"ScatterResidual_BlockedTpetra(Residual):worksetLIDs",
487 maxElementBlockGIDCount);
490 const auto& blockGlobalIndexers =
globalIndexer_->getFieldDOFManagers();
491 const int numBlocks =
static_cast<int>(
globalIndexer_->getFieldDOFManagers().size());
492 blockOffsets_ = PHX::View<LO*>(
"GatherSolution_BlockedTpetra(Jacobian):blockOffsets_",
494 const auto hostBlockOffsets = Kokkos::create_mirror_view(
blockOffsets_);
495 for (
int blk=0;blk<numBlocks;++blk) {
497 hostBlockOffsets(blk) = blockOffset;
499 hostBlockOffsets(numBlocks) = hostBlockOffsets(numBlocks-1) + blockGlobalIndexers[blockGlobalIndexers.size()-1]->getElementBlockGIDCount(blockId);
505template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
514template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
519 using Teuchos::rcp_dynamic_cast;
523 const auto& localCellIds = this->
wda(workset).cell_local_ids_k;
526 RCP<ProductVectorBase<double>> blockedSolution;
528 blockedSolution = rcp_dynamic_cast<ProductVectorBase<double> >(
blockedContainer_->get_dxdt());
529 seedValue = workset.alpha;
532 blockedSolution = rcp_dynamic_cast<ProductVectorBase<double> >(
blockedContainer_->get_x());
533 seedValue = workset.beta;
543 int currentWorksetLIDSubBlock = -1;
544 for (std::size_t fieldIndex = 0; fieldIndex <
gatherFields_.size(); fieldIndex++) {
548 const std::string blockId = this->
wda(workset).block_id;
550 blockIndexer->getElementLIDs(localCellIds,
worksetLIDs_,num_dofs);
555 const auto& subblockSolution = *((rcp_dynamic_cast<Thyra::TpetraVector<RealType,LO,GO,NodeT>>(blockedSolution->getNonconstVectorBlock(blockRowIndex),
true))->getTpetraVector());
556 const auto kokkosSolution = subblockSolution.getLocalViewDevice(Tpetra::Access::ReadOnly);
559 const PHX::View<const int*> fieldOffsets =
fieldOffsets_[fieldIndex];
561 const PHX::View<ScalarT**> fieldValues =
gatherFields_[fieldIndex].get_static_view();
563 auto blockOffsets_h = Kokkos::create_mirror_view(blockOffsets);
564 Kokkos::deep_copy(blockOffsets_h, blockOffsets);
565 const int blockStart = blockOffsets_h(blockRowIndex);
567 Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device>(0,workset.num_cells), KOKKOS_LAMBDA (
const int& cell) {
568 for (
int basis=0; basis < static_cast<int>(fieldOffsets.size()); ++basis) {
569 const int rowLID = worksetLIDs(cell,fieldOffsets(basis));
570 fieldValues(cell,basis).zero();
571 fieldValues(cell,basis).val() = kokkosSolution(rowLID,0);
572 fieldValues(cell,basis).fastAccessDx(blockStart+fieldOffsets(basis)) = seedValue;
WorksetDetailsAccessor wda
bool disableSensitivities_
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
std::vector< int > fieldIds_
std::vector< std::string > indexerNames_
std::vector< int > productVectorBlockIndex_
PHX::View< LO ** > worksetLIDs_
Local indices for unknowns.
Teuchos::RCP< const BlockedDOFManager > globalIndexer_
bool useTimeDerivativeSolutionVector_
TRAITS::RealType RealType
Teuchos::RCP< const BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > > blockedContainer_
PHX::View< LO * > blockOffsets_
The offset values of the blocked DOFs per element. Size of number of blocks in the product vector + 1...
std::string globalDataKey_
std::vector< PHX::View< int * > > fieldOffsets_
Offset into the cell lids for each field. Size of number of fields to scatter.
bool useTimeDerivativeSolutionVector_
std::vector< PHX::View< int * > > fieldOffsets_
Offset into the cell lids for each field.
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
Teuchos::RCP< const BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > > blockedContainer_
std::string globalDataKey_
std::vector< int > productVectorBlockIndex_
Teuchos::RCP< const BlockedDOFManager > globalIndexer_
Maps the local (field,element,basis) triplet to a global ID for.
std::vector< std::string > indexerNames_
std::vector< Teuchos::RCP< const panzer::GlobalIndexer > > fieldGlobalIndexers_
Vector of global indexers, one for each field to gather, respectively.
std::vector< int > fieldIds_
Field IDs in the local product vector block (not global field id)
PHX::View< LO ** > worksetLIDs_
Local indices for unknowns.
Teuchos::RCP< const BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > > blockedContainer_
Teuchos::RCP< const BlockedDOFManager > gidIndexer_
std::vector< std::string > indexerNames_
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
bool useTimeDerivativeSolutionVector_
std::string globalDataKey_
std::vector< int > fieldIds_
std::vector< std::vector< PHX::MDField< const ScalarT, Cell, NODE > > > tangentFields_
Gathers solution values from the Newton solution vector into the nodal fields of the field manager.
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &vm)
void evaluateFields(typename TRAITS::EvalData d)
GatherSolution_BlockedTpetra(const Teuchos::RCP< const BlockedDOFManager > &indexer)