Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherSolution_Tpetra_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef PANZER_GATHER_SOLUTION_TPETRA_IMPL_HPP
12 #define PANZER_GATHER_SOLUTION_TPETRA_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 #include "Phalanx_DataLayout.hpp"
16 
17 #include "Panzer_GlobalIndexer.hpp"
18 #include "Panzer_PureBasis.hpp"
24 #include "Panzer_DOFManager.hpp"
25 
26 #include "Teuchos_FancyOStream.hpp"
27 
28 #include "Tpetra_Vector.hpp"
29 #include "Tpetra_Map.hpp"
30 
31 // **********************************************************************
32 // Specialization: Residual
33 // **********************************************************************
34 
35 template<typename TRAITS,typename LO,typename GO,typename NodeT>
39  const Teuchos::ParameterList& p)
40  : globalIndexer_(indexer)
41  , has_tangent_fields_(false)
42 {
43  typedef std::vector< std::vector<std::string> > vvstring;
44 
46  input.setParameterList(p);
47 
48  const std::vector<std::string> & names = input.getDofNames();
50  const vvstring & tangent_field_names = input.getTangentNames();
51 
52  indexerNames_ = input.getIndexerNames();
53  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
54  globalDataKey_ = input.getGlobalDataKey();
55 
56  // allocate fields
57  gatherFields_.resize(names.size());
58  for (std::size_t fd = 0; fd < names.size(); ++fd) {
59  gatherFields_[fd] =
61  this->addEvaluatedField(gatherFields_[fd]);
62  }
63 
64  // Setup dependent tangent fields if requested
65  if (tangent_field_names.size()>0) {
66  TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
67 
68  has_tangent_fields_ = true;
69  tangentFields_.resize(gatherFields_.size());
70  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
71  tangentFields_[fd].resize(tangent_field_names[fd].size());
72  for (std::size_t i=0; i<tangent_field_names[fd].size(); ++i) {
73  tangentFields_[fd][i] =
74  PHX::MDField<const ScalarT,Cell,NODE>(tangent_field_names[fd][i],basis->functional);
75  this->addDependentField(tangentFields_[fd][i]);
76  }
77  }
78  }
79 
80  // figure out what the first active name is
81  std::string firstName = "<none>";
82  if(names.size()>0)
83  firstName = names[0];
84 
85  std::string n = "GatherSolution (Tpetra): "+firstName+" (Residual)";
86  this->setName(n);
87 }
88 
89 // **********************************************************************
90 template<typename TRAITS,typename LO,typename GO,typename NodeT>
92 postRegistrationSetup(typename TRAITS::SetupData d,
93  PHX::FieldManager<TRAITS>& /* fm */)
94 {
95  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
96 
97  fieldIds_.resize(gatherFields_.size());
98 
99  const Workset & workset_0 = (*d.worksets_)[0];
100  std::string blockId = this->wda(workset_0).block_id;
101  scratch_offsets_.resize(gatherFields_.size());
102 
103  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
104  const std::string& fieldName = indexerNames_[fd];
105  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
106 
107  int fieldNum = fieldIds_[fd];
108  const std::vector<int> & offsets = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
109  scratch_offsets_[fd] = PHX::View<int*>("offsets",offsets.size());
110  Kokkos::deep_copy(scratch_offsets_[fd], Kokkos::View<const int*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(offsets.data(), offsets.size()));
111  }
112 
113  scratch_lids_ = PHX::View<LO**>("lids",gatherFields_[0].extent(0),
114  globalIndexer_->getElementBlockGIDCount(blockId));
115 
116  indexerNames_.clear(); // Don't need this anymore
117 }
118 
119 // **********************************************************************
120 template<typename TRAITS,typename LO,typename GO,typename NodeT>
122 preEvaluate(typename TRAITS::PreEvalData d)
123 {
125 
126  // extract linear object container
127  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject(globalDataKey_));
128 
129  if(tpetraContainer_==Teuchos::null) {
130  // extract linear object container
131  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
132  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(loc);
133  }
134 }
135 
136 // **********************************************************************
137 template<typename TRAITS,typename LO,typename GO,typename NodeT>
139 evaluateFields(typename TRAITS::EvalData workset)
140 {
142 
143  // for convenience pull out some objects from workset
144  std::string blockId = this->wda(workset).block_id;
145  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
146 
148  if (useTimeDerivativeSolutionVector_)
149  x = tpetraContainer_->get_dxdt();
150  else
151  x = tpetraContainer_->get_x();
152 
153  auto x_data = x->getLocalViewDevice(Tpetra::Access::ReadOnly);
154 
155  globalIndexer_->getElementLIDs(this->wda(workset).cell_local_ids_k,scratch_lids_);
156 
157  // NOTE: A reordering of these loops will likely improve performance
158  // The "getGIDFieldOffsets may be expensive. However the
159  // "getElementGIDs" can be cheaper. However the lookup for LIDs
160  // may be more expensive!
161 
162  // gather operation for each cell in workset
163 
164  auto lids = scratch_lids_;
165  for (std::size_t fieldIndex=0; fieldIndex<gatherFields_.size();fieldIndex++) {
166  auto offsets = scratch_offsets_[fieldIndex];
167  auto gather_field = gatherFields_[fieldIndex].get_static_view();
168 
169  Kokkos::parallel_for(localCellIds.size(), KOKKOS_LAMBDA (std::size_t worksetCellIndex) {
170  // loop over basis functions and fill the fields
171  for(std::size_t basis=0;basis<offsets.extent(0);basis++) {
172  int offset = offsets(basis);
173  LO lid = lids(worksetCellIndex,offset);
174 
175  // set the value and seed the FAD object
176  gather_field(worksetCellIndex,basis) = x_data(lid,0);
177  }
178  });
179  }
180 }
181 
182 // **********************************************************************
183 // Specialization: Tangent
184 // **********************************************************************
185 
186 template<typename TRAITS,typename LO,typename GO,typename NodeT>
190  const Teuchos::ParameterList& p)
191  : globalIndexer_(indexer)
192  , has_tangent_fields_(false)
193 {
194  typedef std::vector< std::vector<std::string> > vvstring;
195 
196  GatherSolution_Input input;
197  input.setParameterList(p);
198 
199  const std::vector<std::string> & names = input.getDofNames();
201  const vvstring & tangent_field_names = input.getTangentNames();
202 
203  indexerNames_ = input.getIndexerNames();
204  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
205  globalDataKey_ = input.getGlobalDataKey();
206 
207  // allocate fields
208  gatherFields_.resize(names.size());
209  for (std::size_t fd = 0; fd < names.size(); ++fd) {
210  gatherFields_[fd] =
212  this->addEvaluatedField(gatherFields_[fd]);
213  // Don't allow for sharing so that we can avoid zeroing out the
214  // off-diagonal values of the FAD derivative array.
215  this->addUnsharedField(gatherFields_[fd].fieldTag().clone());
216  }
217 
218  // Setup dependent tangent fields if requested
219  if (tangent_field_names.size()>0) {
220  TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
221 
222  has_tangent_fields_ = true;
223  tangentFields_.resize(gatherFields_.size());
224  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
225  tangentFields_[fd].resize(tangent_field_names[fd].size());
226  for (std::size_t i=0; i<tangent_field_names[fd].size(); ++i) {
227  tangentFields_[fd][i] =
228  PHX::MDField<const RealT,Cell,NODE>(tangent_field_names[fd][i],basis->functional);
229  this->addDependentField(tangentFields_[fd][i]);
230  }
231  }
232  }
233 
234  // figure out what the first active name is
235  std::string firstName = "<none>";
236  if(names.size()>0)
237  firstName = names[0];
238 
239  std::string n = "GatherSolution (Tpetra): "+firstName+" (Tangent)";
240  this->setName(n);
241 }
242 
243 // **********************************************************************
244 template<typename TRAITS,typename LO,typename GO,typename NodeT>
246 postRegistrationSetup(typename TRAITS::SetupData /* d */,
247  PHX::FieldManager<TRAITS>& /* fm */)
248 {
249  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
250 
251  fieldIds_.resize(gatherFields_.size());
252 
253  // Original implementation of tangentFields used vector of
254  // vectors. The inner vectors could have different sizes for each
255  // [fd]. With UVM removal, we need to use a rank 2 view of views. So
256  // we need an extra vector to carry around the inner vector sizes.
257  tangentInnerVectorSizes_ = PHX::View<size_t*>("tangentInnerVectorSizes_",gatherFields_.size());
258  auto tangentInnerVectorSizes_host = Kokkos::create_mirror_view(tangentInnerVectorSizes_);
259  size_t inner_vector_max_size = 0;
260  for (std::size_t fd = 0; fd < tangentFields_.size(); ++fd) {
261  inner_vector_max_size = std::max(inner_vector_max_size,tangentFields_[fd].size());
262  tangentInnerVectorSizes_host(fd) = tangentFields_[fd].size();
263  }
264  Kokkos::deep_copy(tangentInnerVectorSizes_,tangentInnerVectorSizes_host);
265 
266  gatherFieldsVoV_.initialize("GatherSolution_Tpetra<Tangent>::gatherFieldsVoV_",gatherFields_.size());
267  tangentFieldsVoV_.initialize("GatherSolution_Tpetra<Tangent>::tangentFieldsVoV_",gatherFields_.size(),inner_vector_max_size);
268 
269  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
270  const std::string& fieldName = indexerNames_[fd];
271  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
272  gatherFieldsVoV_.addView(gatherFields_[fd].get_static_view(),fd);
273 
274  if (has_tangent_fields_) {
275  for (std::size_t i=0; i<tangentFields_[fd].size(); ++i) {
276  tangentFieldsVoV_.addView(tangentFields_[fd][i].get_static_view(),fd,i);
277  }
278  }
279  }
280 
281  gatherFieldsVoV_.syncHostToDevice();
282  tangentFieldsVoV_.syncHostToDevice();
283 
284  indexerNames_.clear(); // Don't need this anymore
285 }
286 
287 // **********************************************************************
288 template<typename TRAITS,typename LO,typename GO,typename NodeT>
290 preEvaluate(typename TRAITS::PreEvalData d)
291 {
293 
294  // extract linear object container
295  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject(globalDataKey_));
296 
297  if(tpetraContainer_==Teuchos::null) {
298  // extract linear object container
299  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
300  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(loc);
301  }
302 }
303 
304 // **********************************************************************
305 template<typename TRAITS,typename LO,typename GO,typename NodeT>
307 evaluateFields(typename TRAITS::EvalData workset)
308 {
310 
311  // for convenience pull out some objects from workset
312  std::string blockId = this->wda(workset).block_id;
313 
315  if (useTimeDerivativeSolutionVector_)
316  x = tpetraContainer_->get_dxdt();
317  else
318  x = tpetraContainer_->get_x();
319 
320  typedef typename PHX::MDField<ScalarT,Cell,NODE>::array_type::reference_type reference_type;
321  auto cellLocalIdsKokkos = this->wda(workset).getLocalCellIDs();
322  auto lids = globalIndexer_->getLIDs();
323  auto gidFieldOffsetsVoV = Teuchos::rcp_dynamic_cast<const panzer::DOFManager>(globalIndexer_,true)->getGIDFieldOffsetsKokkos(blockId,fieldIds_);
324  auto gidFieldOffsets = gidFieldOffsetsVoV.getViewDevice();
325  auto gatherFieldsDevice = gatherFieldsVoV_.getViewDevice();
326  auto x_view = x->getLocalViewDevice(Tpetra::Access::ReadOnly);
327  auto tangentInnerVectorSizes = this->tangentInnerVectorSizes_;
328 
329  if (has_tangent_fields_) {
330  auto tangentFieldsDevice = tangentFieldsVoV_.getViewDevice();
331  Kokkos::parallel_for("GatherSolutionTpetra<Tangent>",cellLocalIdsKokkos.extent(0),KOKKOS_LAMBDA(const int worksetCellIndex) {
332  for (size_t fieldIndex = 0; fieldIndex < gidFieldOffsets.extent(0); ++fieldIndex) {
333  for(size_t basis=0;basis<gidFieldOffsets(fieldIndex).extent(0);basis++) {
334  int offset = gidFieldOffsets(fieldIndex)(basis);
335  LO lid = lids(cellLocalIdsKokkos(worksetCellIndex),offset);
336  auto gf_ref = (gatherFieldsDevice[fieldIndex])(worksetCellIndex,basis);
337  gf_ref.val() = x_view(lid,0);
338  for (std::size_t i=0; i<tangentInnerVectorSizes(fieldIndex); ++i) {
339  gf_ref.fastAccessDx(i) = tangentFieldsDevice(fieldIndex,i)(worksetCellIndex,basis);
340  }
341  }
342  }
343  });
344  }
345  else {
346  Kokkos::parallel_for("GatherSolutionTpetra<Tangent>",cellLocalIdsKokkos.extent(0),KOKKOS_LAMBDA(const int worksetCellIndex) {
347  for (size_t fieldIndex = 0; fieldIndex < gidFieldOffsets.extent(0); ++fieldIndex) {
348  for(size_t basis=0;basis<gidFieldOffsets(fieldIndex).extent(0);basis++) {
349  int offset = gidFieldOffsets(fieldIndex)(basis);
350  LO lid = lids(cellLocalIdsKokkos(worksetCellIndex),offset);
351  reference_type gf_ref = (gatherFieldsDevice[fieldIndex])(worksetCellIndex,basis);
352  gf_ref.val() = x_view(lid,0);
353  }
354  }
355  });
356  }
357 }
358 
359 // **********************************************************************
360 // Specialization: Jacobian
361 // **********************************************************************
362 
363 template<typename TRAITS,typename LO,typename GO,typename NodeT>
367  const Teuchos::ParameterList& p)
368  : globalIndexer_(indexer)
369 {
370  // typedef std::vector< std::vector<std::string> > vvstring;
371 
372  GatherSolution_Input input;
373  input.setParameterList(p);
374 
375  const std::vector<std::string> & names = input.getDofNames();
377  //const vvstring & tangent_field_names = input.getTangentNames();
378 
379  indexerNames_ = input.getIndexerNames();
380  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
381  globalDataKey_ = input.getGlobalDataKey();
382 
383  gatherSeedIndex_ = input.getGatherSeedIndex();
384  sensitivitiesName_ = input.getSensitivitiesName();
385  disableSensitivities_ = !input.firstSensitivitiesAvailable();
386 
387  gatherFields_.resize(names.size());
388  scratch_offsets_.resize(names.size());
389  for (std::size_t fd = 0; fd < names.size(); ++fd) {
390  PHX::MDField<ScalarT,Cell,NODE> f(names[fd],basis->functional);
391  gatherFields_[fd] = f;
392  this->addEvaluatedField(gatherFields_[fd]);
393  // Don't allow for sharing so that we can avoid zeroing out the
394  // off-diagonal values of the FAD derivative array.
395  this->addUnsharedField(gatherFields_[fd].fieldTag().clone());
396  }
397 
398  // figure out what the first active name is
399  std::string firstName = "<none>";
400  if(names.size()>0)
401  firstName = names[0];
402 
403  // print out convenience
404  if(disableSensitivities_) {
405  std::string n = "GatherSolution (Tpetra, No Sensitivities): "+firstName+" (Jacobian)";
406  this->setName(n);
407  }
408  else {
409  std::string n = "GatherSolution (Tpetra): "+firstName+" (Jacobian) ";
410  this->setName(n);
411  }
412 }
413 
414 // **********************************************************************
415 template<typename TRAITS,typename LO,typename GO,typename NodeT>
417 postRegistrationSetup(typename TRAITS::SetupData d,
418  PHX::FieldManager<TRAITS>& /* fm */)
419 {
420  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
421 
422  fieldIds_.resize(gatherFields_.size());
423 
424  const Workset & workset_0 = (*d.worksets_)[0];
425  std::string blockId = this->wda(workset_0).block_id;
426 
427  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
428  // get field ID from DOF manager
429  const std::string& fieldName = indexerNames_[fd];
430  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
431 
432  int fieldNum = fieldIds_[fd];
433  const std::vector<int> & offsets = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
434  scratch_offsets_[fd] = PHX::View<int*>("offsets",offsets.size());
435  Kokkos::deep_copy(scratch_offsets_[fd], Kokkos::View<const int*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(offsets.data(), offsets.size()));
436  }
437 
438  scratch_lids_ = PHX::View<LO**>("lids",gatherFields_[0].extent(0),
439  globalIndexer_->getElementBlockGIDCount(blockId));
440 
441  indexerNames_.clear(); // Don't need this anymore
442 }
443 
444 // **********************************************************************
445 template<typename TRAITS,typename LO,typename GO,typename NodeT>
447 preEvaluate(typename TRAITS::PreEvalData d)
448 {
449  using Teuchos::RCP;
450  using Teuchos::rcp;
451  using Teuchos::rcp_dynamic_cast;
452 
455 
456  // manage sensitivities
458  if(!disableSensitivities_) {
459  if(d.first_sensitivities_name==sensitivitiesName_)
460  applySensitivities_ = true;
461  else
462  applySensitivities_ = false;
463  }
464  else
465  applySensitivities_ = false;
466 
468 
470 
471  // first try refactored ReadOnly container
472  std::string post = useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X";
473  if(d.gedc->containsDataObject(globalDataKey_+post)) {
474  ged = d.gedc->getDataObject(globalDataKey_+post);
475 
476  RCP<RO_GED> ro_ged = rcp_dynamic_cast<RO_GED>(ged,true);
477 
478  x_vector = ro_ged->getGhostedVector_Tpetra();
479 
480  return;
481  }
482 
483  ged = d.gedc->getDataObject(globalDataKey_);
484 
485  // try to extract linear object container
486  {
487  RCP<LOC> tpetraContainer = rcp_dynamic_cast<LOC>(ged);
488  RCP<LOCPair_GlobalEvaluationData> loc_pair = rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(ged);
489 
490  if(loc_pair!=Teuchos::null) {
491  Teuchos::RCP<LinearObjContainer> loc = loc_pair->getGhostedLOC();
492  // extract linear object container
493  tpetraContainer = rcp_dynamic_cast<LOC>(loc);
494  }
495 
496  if(tpetraContainer!=Teuchos::null) {
497  if (useTimeDerivativeSolutionVector_)
498  x_vector = tpetraContainer->get_dxdt();
499  else
500  x_vector = tpetraContainer->get_x();
501 
502  return; // epetraContainer was found
503  }
504  }
505 
506  // try to extract an EpetraVector_ReadOnly object (this is the last resort!, it throws if not found)
507  {
508  RCP<RO_GED> ro_ged = rcp_dynamic_cast<RO_GED>(ged,true);
509 
510  x_vector = ro_ged->getGhostedVector_Tpetra();
511  }
512 }
513 
514 // **********************************************************************
515 template<typename TRAITS,typename LO,typename GO,typename NodeT>
517 evaluateFields(typename TRAITS::EvalData workset)
518 {
519  // for convenience pull out some objects from workset
520  std::string blockId = this->wda(workset).block_id;
521 
522  double seed_value = 0.0;
523  if (useTimeDerivativeSolutionVector_) {
524  seed_value = workset.alpha;
525  }
526  else if (gatherSeedIndex_<0) {
527  seed_value = workset.beta;
528  }
529  else if(!useTimeDerivativeSolutionVector_) {
530  seed_value = workset.gather_seeds[gatherSeedIndex_];
531  }
532  else {
533  TEUCHOS_ASSERT(false);
534  }
535 
536  // turn off sensitivies: this may be faster if we don't expand the term
537  // but I suspect not because anywhere it is used the full complement of
538  // sensitivies will be needed anyway.
539  if(!applySensitivities_)
540  seed_value = 0.0;
541 
542  // Interface worksets handle DOFs from two element blocks. The
543  // derivative offset for the other element block must be shifted by
544  // the derivative side of my element block.
545  functor_data.dos = 0;
546  if (this->wda.getDetailsIndex() == 1)
547  {
548  // Get the DOF count for my element block.
549  functor_data.dos = globalIndexer_->getElementBlockGIDCount(workset.details(0).block_id);
550  }
551 
552  // switch to a faster assembly
553  bool use_seed = true;
554  if(seed_value==0.0)
555  use_seed = false;
556 
557  globalIndexer_->getElementLIDs(this->wda(workset).cell_local_ids_k,scratch_lids_);
558 
559  // now setup the fuctor_data, and run the parallel_for loop
561 
562  functor_data.x_data = x_vector->getLocalViewDevice(Tpetra::Access::ReadOnly);
563  functor_data.seed_value = seed_value;
564  functor_data.lids = scratch_lids_;
565 
566  // loop over the fields to be gathered
567  for(std::size_t fieldIndex=0;
568  fieldIndex<gatherFields_.size();fieldIndex++) {
569 
570  // setup functor data
571  functor_data.offsets = scratch_offsets_[fieldIndex];
572  functor_data.field = gatherFields_[fieldIndex];
573 
574  if(use_seed)
575  Kokkos::parallel_for(workset.num_cells,*this);
576  else
577  Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device,NoSeed>(0,workset.num_cells),*this);
578  }
580 }
581 
582 // **********************************************************************
583 template<typename TRAITS,typename LO,typename GO,typename NodeT>
584 KOKKOS_INLINE_FUNCTION
586 operator()(const int worksetCellIndex) const
587 {
588  // loop over basis functions and fill the fields
589  for(std::size_t basis=0;basis<functor_data.offsets.extent(0);basis++) {
590  int offset = functor_data.offsets(basis);
591  LO lid = functor_data.lids(worksetCellIndex,offset);
592 
593  // set the value and seed the FAD object
594  if (functor_data.dos == 0)
595  functor_data.field(worksetCellIndex,basis).val() = functor_data.x_data(lid,0);
596  else // Interface conditions need to zero out derivative array
597  functor_data.field(worksetCellIndex,basis) = ScalarT(functor_data.x_data(lid,0));
598 
599  functor_data.field(worksetCellIndex,basis).fastAccessDx(functor_data.dos + offset) = functor_data.seed_value;
600  }
601 }
602 
603 // **********************************************************************
604 template<typename TRAITS,typename LO,typename GO,typename NodeT>
605 KOKKOS_INLINE_FUNCTION
607 operator()(const NoSeed,const int worksetCellIndex) const
608 {
609  // loop over basis functions and fill the fields
610  for(std::size_t basis=0;basis<functor_data.offsets.extent(0);basis++) {
611  int offset = functor_data.offsets(basis);
612  LO lid = functor_data.lids(worksetCellIndex,offset);
613 
614  // set the value and seed the FAD object
615  functor_data.field(worksetCellIndex,basis).val() = functor_data.x_data(lid,0);
616  }
617 }
618 
619 // **********************************************************************
620 
621 #endif
const std::vector< std::vector< std::string > > & getTangentNames() const
Get the name of the tangent fields (tangent only)
const std::vector< std::string > & getDofNames() const
The names of the DOFs to be gathered (all types)
void setParameterList(const Teuchos::ParameterList &pl)
std::string getGlobalDataKey() const
Name of the global evaluation data container to use for the source vector (all types) ...
Teuchos::RCP< const PureBasis > getBasis() const
Basis definiting layout of dof names (all types)
PHX::View< const int * > offsets
Kokkos::View< const LO **, Kokkos::LayoutRight, PHX::Device > lids
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int getGatherSeedIndex() const
What index to use for initializing the seed (Jacobian and Hessian)
Gathers solution values from the Newton solution vector into the nodal fields of the field manager...
bool firstSensitivitiesAvailable()
Are first derivative sensitivities enabled or disabled? (Jacobian and Hessian)
bool useTimeDerivativeSolutionVector() const
Gather a time derivative vector? (all types)
std::string block_id
DEPRECATED - use: getElementBlock()
const std::vector< std::string > & getIndexerNames() const
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< PHX::DataLayout > functional
&lt;Cell,Basis&gt; or &lt;Cell,Basis&gt;
std::string getSensitivitiesName() const
The name of the sensitivities. Enables sensitivities at &quot;preEvaluate&quot; time (Jacobian and Hessian) ...