43#ifndef PANZER_GATHER_SOLUTION_TPETRA_IMPL_HPP
44#define PANZER_GATHER_SOLUTION_TPETRA_IMPL_HPP
46#include "Teuchos_Assert.hpp"
47#include "Phalanx_DataLayout.hpp"
58#include "Teuchos_FancyOStream.hpp"
60#include "Tpetra_Vector.hpp"
61#include "Tpetra_Map.hpp"
67template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
70 const Teuchos::RCP<const panzer::GlobalIndexer> & indexer,
71 const Teuchos::ParameterList& p)
75 typedef std::vector< std::vector<std::string> > vvstring;
80 const std::vector<std::string> & names = input.
getDofNames();
81 Teuchos::RCP<const panzer::PureBasis> basis = input.
getBasis();
90 for (std::size_t fd = 0; fd < names.size(); ++fd) {
92 PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
97 if (tangent_field_names.size()>0) {
98 TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
100 has_tangent_fields_ = true;
101 tangentFields_.resize(gatherFields_.size());
102 for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
103 tangentFields_[fd].resize(tangent_field_names[fd].size());
104 for (std::size_t i=0; i<tangent_field_names[fd].size(); ++i) {
105 tangentFields_[fd][i] =
106 PHX::MDField<const ScalarT,Cell,NODE>(tangent_field_names[fd][i],basis->functional);
107 this->addDependentField(tangentFields_[fd][i]);
113 std::string firstName =
"<none>";
115 firstName = names[0];
117 std::string n =
"GatherSolution (Tpetra): "+firstName+
" (Residual)";
122template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
131 const Workset & workset_0 = (*d.worksets_)[0];
132 std::string blockId = this->
wda(workset_0).block_id;
140 const std::vector<int> & offsets =
globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
142 Kokkos::deep_copy(
scratch_offsets_[fd], Kokkos::View<const int*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(offsets.data(), offsets.size()));
152template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
163 Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(
globalDataKey_),
true)->getGhostedLOC();
169template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
176 std::string blockId = this->
wda(workset).block_id;
177 const std::vector<std::size_t> & localCellIds = this->
wda(workset).cell_local_ids;
179 Teuchos::RCP<typename LOC::VectorType> x;
185 auto x_data = x->getLocalViewDevice(Tpetra::Access::ReadOnly);
197 for (std::size_t fieldIndex=0; fieldIndex<
gatherFields_.size();fieldIndex++) {
199 auto gather_field =
gatherFields_[fieldIndex].get_static_view();
201 Kokkos::parallel_for(localCellIds.size(), KOKKOS_LAMBDA (std::size_t worksetCellIndex) {
203 for(std::size_t basis=0;basis<offsets.extent(0);basis++) {
204 int offset = offsets(basis);
205 LO lid = lids(worksetCellIndex,offset);
208 gather_field(worksetCellIndex,basis) = x_data(lid,0);
218template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
221 const Teuchos::RCP<const panzer::GlobalIndexer> & indexer,
222 const Teuchos::ParameterList& p)
226 typedef std::vector< std::vector<std::string> > vvstring;
231 const std::vector<std::string> & names = input.
getDofNames();
232 Teuchos::RCP<const panzer::PureBasis> basis = input.
getBasis();
241 for (std::size_t fd = 0; fd < names.size(); ++fd) {
243 PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
251 if (tangent_field_names.size()>0) {
252 TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
254 has_tangent_fields_ = true;
255 tangentFields_.resize(gatherFields_.size());
256 for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
257 tangentFields_[fd].resize(tangent_field_names[fd].size());
258 for (std::size_t i=0; i<tangent_field_names[fd].size(); ++i) {
259 tangentFields_[fd][i] =
260 PHX::MDField<const RealT,Cell,NODE>(tangent_field_names[fd][i],basis->functional);
261 this->addDependentField(tangentFields_[fd][i]);
267 std::string firstName =
"<none>";
269 firstName = names[0];
271 std::string n =
"GatherSolution (Tpetra): "+firstName+
" (Tangent)";
276template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
291 size_t inner_vector_max_size = 0;
293 inner_vector_max_size = std::max(inner_vector_max_size,
tangentFields_[fd].size());
320template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
331 Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(
globalDataKey_),
true)->getGhostedLOC();
337template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
344 std::string blockId = this->
wda(workset).block_id;
346 Teuchos::RCP<typename LOC::VectorType> x;
352 typedef typename PHX::MDField<ScalarT,Cell,NODE>::array_type::reference_type reference_type;
353 auto cellLocalIdsKokkos = this->
wda(workset).getLocalCellIDs();
355 auto gidFieldOffsetsVoV = Teuchos::rcp_dynamic_cast<const panzer::DOFManager>(
globalIndexer_,
true)->getGIDFieldOffsetsKokkos(blockId,
fieldIds_);
356 auto gidFieldOffsets = gidFieldOffsetsVoV.getViewDevice();
358 auto x_view = x->getLocalViewDevice(Tpetra::Access::ReadOnly);
363 Kokkos::parallel_for(
"GatherSolutionTpetra<Tangent>",cellLocalIdsKokkos.extent(0),KOKKOS_LAMBDA(
const int worksetCellIndex) {
364 for (
size_t fieldIndex = 0; fieldIndex < gidFieldOffsets.extent(0); ++fieldIndex) {
365 for(
size_t basis=0;basis<gidFieldOffsets(fieldIndex).extent(0);basis++) {
366 int offset = gidFieldOffsets(fieldIndex)(basis);
367 LO lid = lids(cellLocalIdsKokkos(worksetCellIndex),offset);
368 auto gf_ref = (gatherFieldsDevice[fieldIndex])(worksetCellIndex,basis);
369 gf_ref.val() = x_view(lid,0);
370 for (std::size_t i=0; i<tangentInnerVectorSizes(fieldIndex); ++i) {
371 gf_ref.fastAccessDx(i) = tangentFieldsDevice(fieldIndex,i)(worksetCellIndex,basis);
378 Kokkos::parallel_for(
"GatherSolutionTpetra<Tangent>",cellLocalIdsKokkos.extent(0),KOKKOS_LAMBDA(
const int worksetCellIndex) {
379 for (
size_t fieldIndex = 0; fieldIndex < gidFieldOffsets.extent(0); ++fieldIndex) {
380 for(
size_t basis=0;basis<gidFieldOffsets(fieldIndex).extent(0);basis++) {
381 int offset = gidFieldOffsets(fieldIndex)(basis);
382 LO lid = lids(cellLocalIdsKokkos(worksetCellIndex),offset);
383 reference_type gf_ref = (gatherFieldsDevice[fieldIndex])(worksetCellIndex,basis);
384 gf_ref.val() = x_view(lid,0);
395template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
398 const Teuchos::RCP<const panzer::GlobalIndexer> & indexer,
399 const Teuchos::ParameterList& p)
407 const std::vector<std::string> & names = input.
getDofNames();
408 Teuchos::RCP<const panzer::PureBasis> basis = input.
getBasis();
421 for (std::size_t fd = 0; fd < names.size(); ++fd) {
422 PHX::MDField<ScalarT,Cell,NODE> f(names[fd],basis->functional);
431 std::string firstName =
"<none>";
433 firstName = names[0];
437 std::string n =
"GatherSolution (Tpetra, No Sensitivities): "+firstName+
" (Jacobian)";
441 std::string n =
"GatherSolution (Tpetra): "+firstName+
" (Jacobian) ";
447template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
456 const Workset & workset_0 = (*d.worksets_)[0];
457 std::string blockId = this->
wda(workset_0).block_id;
477template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
483 using Teuchos::rcp_dynamic_cast;
501 RCP<GlobalEvaluationData> ged;
508 RCP<RO_GED> ro_ged = rcp_dynamic_cast<RO_GED>(ged,
true);
510 x_vector = ro_ged->getGhostedVector_Tpetra();
519 RCP<LOC> tpetraContainer = rcp_dynamic_cast<LOC>(ged);
520 RCP<LOCPair_GlobalEvaluationData> loc_pair = rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(ged);
522 if(loc_pair!=Teuchos::null) {
523 Teuchos::RCP<LinearObjContainer> loc = loc_pair->getGhostedLOC();
525 tpetraContainer = rcp_dynamic_cast<LOC>(loc);
528 if(tpetraContainer!=Teuchos::null) {
530 x_vector = tpetraContainer->get_dxdt();
532 x_vector = tpetraContainer->get_x();
540 RCP<RO_GED> ro_ged = rcp_dynamic_cast<RO_GED>(ged,
true);
542 x_vector = ro_ged->getGhostedVector_Tpetra();
547template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
552 std::string blockId = this->
wda(workset).block_id;
565 TEUCHOS_ASSERT(
false);
578 if (this->
wda.getDetailsIndex() == 1)
585 bool use_seed =
true;
599 for(std::size_t fieldIndex=0;
607 Kokkos::parallel_for(workset.num_cells,*
this);
609 Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device,NoSeed>(0,workset.num_cells),*
this);
611 functor_data.x_data = Kokkos::View<const double**, Kokkos::LayoutLeft,PHX::Device>();
615template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
616KOKKOS_INLINE_FUNCTION
621 for(std::size_t basis=0;basis<
functor_data.offsets.extent(0);basis++) {
636template<
typename TRAITS,
typename LO,
typename GO,
typename NodeT>
637KOKKOS_INLINE_FUNCTION
642 for(std::size_t basis=0;basis<
functor_data.offsets.extent(0);basis++) {
WorksetDetailsAccessor wda
Teuchos::RCP< const panzer::GlobalIndexer > globalIndexer_
bool useTimeDerivativeSolutionVector_
std::string sensitivitiesName_
std::vector< int > fieldIds_
struct panzer::GatherSolution_Tpetra< panzer::Traits::Jacobian, TRAITS, LO, GO, NodeT >::@206057163070241233110272241152374302267306142275 functor_data
std::string globalDataKey_
panzer::Traits::Jacobian::ScalarT ScalarT
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
PHX::View< int ** > scratch_lids_
std::vector< std::string > indexerNames_
PHX::View< const int * > offsets
std::vector< PHX::View< int * > > scratch_offsets_
bool disableSensitivities_
Teuchos::RCP< typename TpetraLinearObjContainer< double, LO, GO, NodeT >::VectorType > x_vector
bool useTimeDerivativeSolutionVector_
std::vector< PHX::View< int * > > scratch_offsets_
Teuchos::RCP< const panzer::GlobalIndexer > globalIndexer_
std::vector< int > fieldIds_
PHX::View< int ** > scratch_lids_
std::vector< std::string > indexerNames_
std::string globalDataKey_
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
Teuchos::RCP< const TpetraLinearObjContainer< double, LO, GO, NodeT > > tpetraContainer_
PHX::ViewOfViews3< 2, PHX::View< const RealT ** > > tangentFieldsVoV_
std::vector< std::string > indexerNames_
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
std::vector< int > fieldIds_
Teuchos::RCP< const TpetraLinearObjContainer< double, LO, GO, NodeT > > tpetraContainer_
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
PHX::ViewOfViews3< 1, PHX::View< ScalarT ** > > gatherFieldsVoV_
bool useTimeDerivativeSolutionVector_
PHX::View< size_t * > tangentInnerVectorSizes_
std::string globalDataKey_
std::vector< std::vector< PHX::MDField< const RealT, Cell, NODE > > > tangentFields_
Teuchos::RCP< const panzer::GlobalIndexer > globalIndexer_
Gathers solution values from the Newton solution vector into the nodal fields of the field manager.