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 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_GATHER_SOLUTION_TPETRA_IMPL_HPP
44 #define PANZER_GATHER_SOLUTION_TPETRA_IMPL_HPP
45 
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
48 
49 #include "Panzer_GlobalIndexer.hpp"
50 #include "Panzer_PureBasis.hpp"
56 
57 #include "Teuchos_FancyOStream.hpp"
58 
59 #include "Tpetra_Vector.hpp"
60 #include "Tpetra_Map.hpp"
61 
62 // **********************************************************************
63 // Specialization: Residual
64 // **********************************************************************
65 
66 template<typename TRAITS,typename LO,typename GO,typename NodeT>
70  const Teuchos::ParameterList& p)
71  : globalIndexer_(indexer)
72  , has_tangent_fields_(false)
73 {
74  typedef std::vector< std::vector<std::string> > vvstring;
75 
77  input.setParameterList(p);
78 
79  const std::vector<std::string> & names = input.getDofNames();
81  const vvstring & tangent_field_names = input.getTangentNames();
82 
83  indexerNames_ = input.getIndexerNames();
84  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
85  globalDataKey_ = input.getGlobalDataKey();
86 
87  // allocate fields
88  gatherFields_.resize(names.size());
89  for (std::size_t fd = 0; fd < names.size(); ++fd) {
90  gatherFields_[fd] =
91  PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
92  this->addEvaluatedField(gatherFields_[fd]);
93  }
94 
95  // Setup dependent tangent fields if requested
96  if (tangent_field_names.size()>0) {
97  TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
98 
99  has_tangent_fields_ = true;
100  tangentFields_.resize(gatherFields_.size());
101  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
102  tangentFields_[fd].resize(tangent_field_names[fd].size());
103  for (std::size_t i=0; i<tangent_field_names[fd].size(); ++i) {
104  tangentFields_[fd][i] =
105  PHX::MDField<const ScalarT,Cell,NODE>(tangent_field_names[fd][i],basis->functional);
106  this->addDependentField(tangentFields_[fd][i]);
107  }
108  }
109  }
110 
111  // figure out what the first active name is
112  std::string firstName = "<none>";
113  if(names.size()>0)
114  firstName = names[0];
115 
116  std::string n = "GatherSolution (Tpetra): "+firstName+" (Residual)";
117  this->setName(n);
118 }
119 
120 // **********************************************************************
121 template<typename TRAITS,typename LO,typename GO,typename NodeT>
123 postRegistrationSetup(typename TRAITS::SetupData d,
124  PHX::FieldManager<TRAITS>& /* fm */)
125 {
126  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
127 
128  fieldIds_.resize(gatherFields_.size());
129 
130  const Workset & workset_0 = (*d.worksets_)[0];
131  std::string blockId = this->wda(workset_0).block_id;
132  scratch_offsets_.resize(gatherFields_.size());
133 
134  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
135  const std::string& fieldName = indexerNames_[fd];
136  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
137 
138  int fieldNum = fieldIds_[fd];
139  const std::vector<int> & offsets = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
140  scratch_offsets_[fd] = Kokkos::View<int*,PHX::Device>("offsets",offsets.size());
141  for(std::size_t i=0;i<offsets.size();i++)
142  scratch_offsets_[fd](i) = offsets[i];
143  }
144 
145  scratch_lids_ = Kokkos::View<LO**,PHX::Device>("lids",gatherFields_[0].extent(0),
146  globalIndexer_->getElementBlockGIDCount(blockId));
147 
148  indexerNames_.clear(); // Don't need this anymore
149 }
150 
151 // **********************************************************************
152 template<typename TRAITS,typename LO,typename GO,typename NodeT>
154 preEvaluate(typename TRAITS::PreEvalData d)
155 {
157 
158  // extract linear object container
159  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject(globalDataKey_));
160 
161  if(tpetraContainer_==Teuchos::null) {
162  // extract linear object container
163  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
164  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(loc);
165  }
166 }
167 
168 // **********************************************************************
169 template<typename TRAITS,typename LO,typename GO,typename NodeT>
171 evaluateFields(typename TRAITS::EvalData workset)
172 {
174 
175  // for convenience pull out some objects from workset
176  std::string blockId = this->wda(workset).block_id;
177  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
178 
180  if (useTimeDerivativeSolutionVector_)
181  x = tpetraContainer_->get_dxdt();
182  else
183  x = tpetraContainer_->get_x();
184 
185  auto x_data = x->template getLocalView<PHX::Device>();
186 
187  globalIndexer_->getElementLIDs(this->wda(workset).cell_local_ids_k,scratch_lids_);
188 
189  // NOTE: A reordering of these loops will likely improve performance
190  // The "getGIDFieldOffsets may be expensive. However the
191  // "getElementGIDs" can be cheaper. However the lookup for LIDs
192  // may be more expensive!
193 
194  // gather operation for each cell in workset
195 
196  auto lids = scratch_lids_;
197  for (std::size_t fieldIndex=0; fieldIndex<gatherFields_.size();fieldIndex++) {
198  auto offsets = scratch_offsets_[fieldIndex];
199  auto gather_field = gatherFields_[fieldIndex];
200 
201  Kokkos::parallel_for(localCellIds.size(), KOKKOS_LAMBDA (std::size_t worksetCellIndex) {
202  // loop over basis functions and fill the fields
203  for(std::size_t basis=0;basis<offsets.extent(0);basis++) {
204  int offset = offsets(basis);
205  LO lid = lids(worksetCellIndex,offset);
206 
207  // set the value and seed the FAD object
208  gather_field(worksetCellIndex,basis) = x_data(lid,0);
209  }
210  });
211  }
212 }
213 
214 // **********************************************************************
215 // Specialization: Tangent
216 // **********************************************************************
217 
218 template<typename TRAITS,typename LO,typename GO,typename NodeT>
222  const Teuchos::ParameterList& p)
223  : globalIndexer_(indexer)
224  , has_tangent_fields_(false)
225 {
226  typedef std::vector< std::vector<std::string> > vvstring;
227 
228  GatherSolution_Input input;
229  input.setParameterList(p);
230 
231  const std::vector<std::string> & names = input.getDofNames();
233  const vvstring & tangent_field_names = input.getTangentNames();
234 
235  indexerNames_ = input.getIndexerNames();
236  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
237  globalDataKey_ = input.getGlobalDataKey();
238 
239  // allocate fields
240  gatherFields_.resize(names.size());
241  for (std::size_t fd = 0; fd < names.size(); ++fd) {
242  gatherFields_[fd] =
243  PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
244  this->addEvaluatedField(gatherFields_[fd]);
245  // Don't allow for sharing so that we can avoid zeroing out the
246  // off-diagonal values of the FAD derivative array.
247  this->addUnsharedField(gatherFields_[fd].fieldTag().clone());
248  }
249 
250  // Setup dependent tangent fields if requested
251  if (tangent_field_names.size()>0) {
252  TEUCHOS_ASSERT(gatherFields_.size() == tangent_field_names.size());
253 
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]);
262  }
263  }
264  }
265 
266  // figure out what the first active name is
267  std::string firstName = "<none>";
268  if(names.size()>0)
269  firstName = names[0];
270 
271  std::string n = "GatherSolution (Tpetra): "+firstName+" (Tangent)";
272  this->setName(n);
273 }
274 
275 // **********************************************************************
276 template<typename TRAITS,typename LO,typename GO,typename NodeT>
278 postRegistrationSetup(typename TRAITS::SetupData /* d */,
279  PHX::FieldManager<TRAITS>& /* fm */)
280 {
281  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
282 
283  fieldIds_.resize(gatherFields_.size());
284 
285  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
286  const std::string& fieldName = indexerNames_[fd];
287  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
288  }
289 
290  indexerNames_.clear(); // Don't need this anymore
291 }
292 
293 // **********************************************************************
294 template<typename TRAITS,typename LO,typename GO,typename NodeT>
296 preEvaluate(typename TRAITS::PreEvalData d)
297 {
299 
300  // extract linear object container
301  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject(globalDataKey_));
302 
303  if(tpetraContainer_==Teuchos::null) {
304  // extract linear object container
305  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
306  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(loc);
307  }
308 }
309 
310 // **********************************************************************
311 template<typename TRAITS,typename LO,typename GO,typename NodeT>
313 evaluateFields(typename TRAITS::EvalData workset)
314 {
316 
317  // for convenience pull out some objects from workset
318  std::string blockId = this->wda(workset).block_id;
319  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
320 
322  if (useTimeDerivativeSolutionVector_)
323  x = tpetraContainer_->get_dxdt();
324  else
325  x = tpetraContainer_->get_x();
326 
327  Teuchos::ArrayRCP<const double> x_array = x->get1dView();
328 
329  // NOTE: A reordering of these loops will likely improve performance
330  // The "getGIDFieldOffsets may be expensive. However the
331  // "getElementGIDs" can be cheaper. However the lookup for LIDs
332  // may be more expensive!
333 
334  typedef typename PHX::MDField<ScalarT,Cell,NODE>::array_type::reference_type reference_type;
335 
336  if (has_tangent_fields_) {
337  // gather operation for each cell in workset
338  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
339  std::size_t cellLocalId = localCellIds[worksetCellIndex];
340 
341  auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
342 
343  // loop over the fields to be gathered
344  for (std::size_t fieldIndex=0; fieldIndex<gatherFields_.size();fieldIndex++) {
345  int fieldNum = fieldIds_[fieldIndex];
346  const std::vector<int> & elmtOffset = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
347 
348  const std::vector< PHX::MDField<const RealT,Cell,NODE> >& tf_ref =
349  tangentFields_[fieldIndex];
350  const std::size_t num_tf = tf_ref.size();
351 
352  // loop over basis functions and fill the fields
353  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
354  int offset = elmtOffset[basis];
355  LO lid = LIDs[offset];
356  reference_type gf_ref = (gatherFields_[fieldIndex])(worksetCellIndex,basis);
357  gf_ref.val() = x_array[lid];
358  for (std::size_t i=0; i<num_tf; ++i)
359  gf_ref.fastAccessDx(i) = tf_ref[i](worksetCellIndex,basis);
360  }
361  }
362  }
363  }
364  else {
365  // gather operation for each cell in workset
366  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
367  std::size_t cellLocalId = localCellIds[worksetCellIndex];
368 
369  auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
370 
371  // loop over the fields to be gathered
372  for (std::size_t fieldIndex=0; fieldIndex<gatherFields_.size();fieldIndex++) {
373  int fieldNum = fieldIds_[fieldIndex];
374  const std::vector<int> & elmtOffset = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
375 
376  // loop over basis functions and fill the fields
377  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
378  int offset = elmtOffset[basis];
379  LO lid = LIDs[offset];
380  reference_type gf_ref = (gatherFields_[fieldIndex])(worksetCellIndex,basis);
381  gf_ref.val() = x_array[lid];
382  }
383  }
384  }
385  }
386 }
387 
388 // **********************************************************************
389 // Specialization: Jacobian
390 // **********************************************************************
391 
392 template<typename TRAITS,typename LO,typename GO,typename NodeT>
396  const Teuchos::ParameterList& p)
397  : globalIndexer_(indexer)
398 {
399  // typedef std::vector< std::vector<std::string> > vvstring;
400 
401  GatherSolution_Input input;
402  input.setParameterList(p);
403 
404  const std::vector<std::string> & names = input.getDofNames();
406  //const vvstring & tangent_field_names = input.getTangentNames();
407 
408  indexerNames_ = input.getIndexerNames();
409  useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
410  globalDataKey_ = input.getGlobalDataKey();
411 
412  gatherSeedIndex_ = input.getGatherSeedIndex();
413  sensitivitiesName_ = input.getSensitivitiesName();
414  disableSensitivities_ = !input.firstSensitivitiesAvailable();
415 
416  gatherFields_.resize(names.size());
417  scratch_offsets_.resize(names.size());
418  for (std::size_t fd = 0; fd < names.size(); ++fd) {
419  PHX::MDField<ScalarT,Cell,NODE> f(names[fd],basis->functional);
420  gatherFields_[fd] = f;
421  this->addEvaluatedField(gatherFields_[fd]);
422  // Don't allow for sharing so that we can avoid zeroing out the
423  // off-diagonal values of the FAD derivative array.
424  this->addUnsharedField(gatherFields_[fd].fieldTag().clone());
425  }
426 
427  // figure out what the first active name is
428  std::string firstName = "<none>";
429  if(names.size()>0)
430  firstName = names[0];
431 
432  // print out convenience
433  if(disableSensitivities_) {
434  std::string n = "GatherSolution (Tpetra, No Sensitivities): "+firstName+" (Jacobian)";
435  this->setName(n);
436  }
437  else {
438  std::string n = "GatherSolution (Tpetra): "+firstName+" ("+PHX::print<EvalT>()+") ";
439  this->setName(n);
440  }
441 }
442 
443 // **********************************************************************
444 template<typename TRAITS,typename LO,typename GO,typename NodeT>
446 postRegistrationSetup(typename TRAITS::SetupData d,
447  PHX::FieldManager<TRAITS>& /* fm */)
448 {
449  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size());
450 
451  fieldIds_.resize(gatherFields_.size());
452 
453  const Workset & workset_0 = (*d.worksets_)[0];
454  std::string blockId = this->wda(workset_0).block_id;
455 
456  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
457  // get field ID from DOF manager
458  const std::string& fieldName = indexerNames_[fd];
459  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
460 
461  int fieldNum = fieldIds_[fd];
462  const std::vector<int> & offsets = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
463  scratch_offsets_[fd] = Kokkos::View<int*,PHX::Device>("offsets",offsets.size());
464  for(std::size_t i=0;i<offsets.size();i++)
465  scratch_offsets_[fd](i) = offsets[i];
466  }
467 
468  scratch_lids_ = Kokkos::View<LO**,PHX::Device>("lids",gatherFields_[0].extent(0),
469  globalIndexer_->getElementBlockGIDCount(blockId));
470 
471  indexerNames_.clear(); // Don't need this anymore
472 }
473 
474 // **********************************************************************
475 template<typename TRAITS,typename LO,typename GO,typename NodeT>
477 preEvaluate(typename TRAITS::PreEvalData d)
478 {
479  using Teuchos::RCP;
480  using Teuchos::rcp;
481  using Teuchos::rcp_dynamic_cast;
482 
485 
486  // manage sensitivities
488  if(!disableSensitivities_) {
489  if(d.first_sensitivities_name==sensitivitiesName_)
490  applySensitivities_ = true;
491  else
492  applySensitivities_ = false;
493  }
494  else
495  applySensitivities_ = false;
496 
498 
500 
501  // first try refactored ReadOnly container
502  std::string post = useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X";
503  if(d.gedc->containsDataObject(globalDataKey_+post)) {
504  ged = d.gedc->getDataObject(globalDataKey_+post);
505 
506  RCP<RO_GED> ro_ged = rcp_dynamic_cast<RO_GED>(ged,true);
507 
508  x_vector = ro_ged->getGhostedVector_Tpetra();
509 
510  x_vector->template sync<PHX::Device>();
511 
512  return;
513  }
514 
515  ged = d.gedc->getDataObject(globalDataKey_);
516 
517  // try to extract linear object container
518  {
519  RCP<LOC> tpetraContainer = rcp_dynamic_cast<LOC>(ged);
520  RCP<LOCPair_GlobalEvaluationData> loc_pair = rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(ged);
521 
522  if(loc_pair!=Teuchos::null) {
523  Teuchos::RCP<LinearObjContainer> loc = loc_pair->getGhostedLOC();
524  // extract linear object container
525  tpetraContainer = rcp_dynamic_cast<LOC>(loc);
526  }
527 
528  if(tpetraContainer!=Teuchos::null) {
529  if (useTimeDerivativeSolutionVector_)
530  x_vector = tpetraContainer->get_dxdt();
531  else
532  x_vector = tpetraContainer->get_x();
533 
534  x_vector->template sync<PHX::Device>();
535 
536  return; // epetraContainer was found
537  }
538  }
539 
540  // try to extract an EpetraVector_ReadOnly object (this is the last resort!, it throws if not found)
541  {
542  RCP<RO_GED> ro_ged = rcp_dynamic_cast<RO_GED>(ged,true);
543 
544  x_vector = ro_ged->getGhostedVector_Tpetra();
545  }
546 
547  x_vector->template sync<PHX::Device>();
548 }
549 
550 // **********************************************************************
551 template<typename TRAITS,typename LO,typename GO,typename NodeT>
553 evaluateFields(typename TRAITS::EvalData workset)
554 {
555  // for convenience pull out some objects from workset
556  std::string blockId = this->wda(workset).block_id;
557 
558  double seed_value = 0.0;
559  if (useTimeDerivativeSolutionVector_) {
560  seed_value = workset.alpha;
561  }
562  else if (gatherSeedIndex_<0) {
563  seed_value = workset.beta;
564  }
565  else if(!useTimeDerivativeSolutionVector_) {
566  seed_value = workset.gather_seeds[gatherSeedIndex_];
567  }
568  else {
569  TEUCHOS_ASSERT(false);
570  }
571 
572  // turn off sensitivies: this may be faster if we don't expand the term
573  // but I suspect not because anywhere it is used the full complement of
574  // sensitivies will be needed anyway.
575  if(!applySensitivities_)
576  seed_value = 0.0;
577 
578  // switch to a faster assembly
579  bool use_seed = true;
580  if(seed_value==0.0)
581  use_seed = false;
582 
583  globalIndexer_->getElementLIDs(this->wda(workset).cell_local_ids_k,scratch_lids_);
584 
585  // now setup the fuctor_data, and run the parallel_for loop
587 
588  functor_data.x_data = x_vector->template getLocalView<PHX::Device>();
589  functor_data.seed_value = seed_value;
590  functor_data.lids = scratch_lids_;
591 
592  // loop over the fields to be gathered
593  for(std::size_t fieldIndex=0;
594  fieldIndex<gatherFields_.size();fieldIndex++) {
595 
596  // setup functor data
597  functor_data.offsets = scratch_offsets_[fieldIndex];
598  functor_data.field = gatherFields_[fieldIndex];
599 
600  if(use_seed)
601  Kokkos::parallel_for(workset.num_cells,*this);
602  else
603  Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device,NoSeed>(0,workset.num_cells),*this);
604  }
605 }
606 
607 // **********************************************************************
608 template<typename TRAITS,typename LO,typename GO,typename NodeT>
609 KOKKOS_INLINE_FUNCTION
611 operator()(const int worksetCellIndex) const
612 {
613  // loop over basis functions and fill the fields
614  for(std::size_t basis=0;basis<functor_data.offsets.extent(0);basis++) {
615  int offset = functor_data.offsets(basis);
616  LO lid = functor_data.lids(worksetCellIndex,offset);
617 
618  // set the value and seed the FAD object
619  functor_data.field(worksetCellIndex,basis).val() = functor_data.x_data(lid,0);
620  functor_data.field(worksetCellIndex,basis).fastAccessDx(offset) = functor_data.seed_value;
621  }
622 }
623 
624 // **********************************************************************
625 template<typename TRAITS,typename LO,typename GO,typename NodeT>
626 KOKKOS_INLINE_FUNCTION
628 operator()(const NoSeed,const int worksetCellIndex) const
629 {
630  // loop over basis functions and fill the fields
631  for(std::size_t basis=0;basis<functor_data.offsets.extent(0);basis++) {
632  int offset = functor_data.offsets(basis);
633  LO lid = functor_data.lids(worksetCellIndex,offset);
634 
635  // set the value and seed the FAD object
636  functor_data.field(worksetCellIndex,basis).val() = functor_data.x_data(lid,0);
637  }
638 }
639 
640 // **********************************************************************
641 
642 #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)
Kokkos::View< const LO **, 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)
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) ...
Kokkos::View< const int *, PHX::Device > offsets