Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_GatherSolution_Epetra_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_GatherSolution_Epetra_impl_hpp__
12#define __Panzer_GatherSolution_Epetra_impl_hpp__
13
15//
16// Include Files
17//
19
20// Epetra
21#include "Epetra_Vector.h"
22
23// Panzer
28#include "Panzer_PureBasis.hpp"
31
32// Teuchos
33#include "Teuchos_Assert.hpp"
34
35// Thyra
36#include "Thyra_SpmdVectorBase.hpp"
37
39//
40// Initializing Constructor (Residual Specialization)
41//
43template<typename TRAITS, typename LO, typename GO>
46 const Teuchos::RCP<const panzer::GlobalIndexer>& indexer,
47 const Teuchos::ParameterList& p)
48 :
49 globalIndexer_(indexer),
50 hasTangentFields_(false)
51{
53 using PHX::MDField;
54 using PHX::print;
55 using std::size_t;
56 using std::vector;
57 using std::string;
58 using Teuchos::RCP;
59 using vvstring = std::vector<std::vector<std::string>>;
61 input.setParameterList(p);
62 const vector<string>& names = input.getDofNames();
63 RCP<const PureBasis> basis = input.getBasis();
64 const vvstring& tangentFieldNames = input.getTangentNames();
65 indexerNames_ = input.getIndexerNames();
66 useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
67 globalDataKey_ = input.getGlobalDataKey();
68
69 // Allocate the fields.
70 int numFields(names.size());
71 gatherFields_.resize(numFields);
72 for (int fd(0); fd < numFields; ++fd)
73 {
74 gatherFields_[fd] =
75 MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
76 this->addEvaluatedField(gatherFields_[fd]);
77 } // end loop over names
78
79 // Set up the dependent tangent fields, if requested.
80 if (tangentFieldNames.size() > 0)
81 {
82 TEUCHOS_ASSERT(gatherFields_.size() == tangentFieldNames.size())
83 hasTangentFields_ = true;
84 tangentFields_.resize(numFields);
85 for (int fd(0); fd < numFields; ++fd)
86 {
87 int numTangentFields(tangentFieldNames[fd].size());
88 tangentFields_[fd].resize(numTangentFields);
89 for (int i(0); i < numTangentFields; ++i)
90 {
91 tangentFields_[fd][i] =
92 MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
93 basis->functional);
94 this->addDependentField(tangentFields_[fd][i]);
95 } // end loop over tangentFieldNames[fd]
96 } // end loop over gatherFields
97 } // end if (tangentFieldNames.size() > 0)
98
99 // Figure out what the first active name is.
100 string firstName("<none>");
101 if (numFields > 0)
102 firstName = names[0];
103 string n("GatherSolution (Epetra): " + firstName + " (Residual)");
104 this->setName(n);
105} // end of Initializing Constructor (Residual Specialization)
106
108//
109// postRegistrationSetup() (Residual Specialization)
110//
112template<typename TRAITS, typename LO, typename GO>
113void
116 typename TRAITS::SetupData /* d */,
118{
119 using std::logic_error;
120 using std::size_t;
121 using std::string;
122 TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size())
123 int numFields(gatherFields_.size());
124 fieldIds_.resize(numFields);
125 for (int fd(0); fd < numFields; ++fd)
126 {
127 // Get the field ID from the DOF manager.
128 const string& fieldName(indexerNames_[fd]);
129 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
130
131 // This is the error return code; raise the alarm.
132 TEUCHOS_TEST_FOR_EXCEPTION(fieldIds_[fd] == -1, logic_error,
133 "GatherSolution_Epetra<Residual>: Could not find field \"" +
134 fieldName + "\" in the global indexer. ")
135 } // end loop over gatherFields_
136 indexerNames_.clear();
137} // end of postRegistrationSetup() (Residual Specialization)
138
140//
141// preEvaluate() (Residual Specialization)
142//
144template<typename TRAITS, typename LO, typename GO>
145void
148 typename TRAITS::PreEvalData d)
149{
150 using std::string;
151 using Teuchos::RCP;
152 using Teuchos::rcp_dynamic_cast;
156 using LOC = panzer::LinearObjContainer;
158
159 // First try the refactored ReadOnly container.
160 RCP<GED> ged;
161 string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
162 if (d.gedc->containsDataObject(globalDataKey_ + post))
163 {
164 ged = d.gedc->getDataObject(globalDataKey_ + post);
165 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
166 return;
167 } // end of the refactored ReadOnly way
168
169 // Now try the old path.
170 ged = d.gedc->getDataObject(globalDataKey_);
171 {
172 // Try to extract the linear object container.
173 auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
174 auto locPair = rcp_dynamic_cast<LPGED>(ged);
175 if (not locPair.is_null())
176 {
177 RCP<LOC> loc = locPair->getGhostedLOC();
178 epetraContainer = rcp_dynamic_cast<ELOC>(loc);
179 } // end if (not locPair.is_null())
180 if (not epetraContainer.is_null())
181 {
182 if (useTimeDerivativeSolutionVector_)
183 x_ = epetraContainer->get_dxdt();
184 else // if (not useTimeDerivativeSolutionVector_)
185 x_ = epetraContainer->get_x();
186 return;
187 } // end if (not epetraContainer.is_null())
188 } // end of the old path
189
190 // As a last resort, try to extract an EpetraVector_ReadOnly object. This
191 // will throw an exception if it doesn't work.
192 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
193} // end of preEvaluate() (Residual Specialization)
194
196//
197// evaluateFields() (Residual Specialization)
198//
200template<typename TRAITS, typename LO, typename GO>
201void
204 typename TRAITS::EvalData workset)
205{
206 using PHX::MDField;
207 using std::size_t;
208 using std::string;
209 using std::vector;
210 using Teuchos::ArrayRCP;
211 using Teuchos::ptrFromRef;
212 using Teuchos::RCP;
213 using Teuchos::rcp_dynamic_cast;
214 using Thyra::SpmdVectorBase;
215
216 // For convenience, pull out some objects from the workset.
217 string blockId(this->wda(workset).block_id);
218 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
219 int numCells(localCellIds.size()), numFields(gatherFields_.size());
220
221 // NOTE: A reordering of these loops will likely improve performance. The
222 // "getGIDFieldOffsets may be expensive. However the "getElementGIDs"
223 // can be cheaper. However the lookup for LIDs may be more expensive!
224
225 auto LIDs = globalIndexer_->getLIDs();
226 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
227 Kokkos::deep_copy(LIDs_h, LIDs);
228 // Loop over the fields to be gathered.
229 for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
230 {
231 MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
232 auto field_h = Kokkos::create_mirror_view(field.get_static_view());
233 int fieldNum(fieldIds_[fieldInd]);
234 const vector<int>& elmtOffset =
235 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
236 int numBases(elmtOffset.size());
237 // Gather operation for each cell in the workset.
238 for (int cell(0); cell < numCells; ++cell)
239 {
240 size_t cellLocalId(localCellIds[cell]);
241 // Loop over the basis functions and fill the fields.
242 for (int basis(0); basis < numBases; ++basis)
243 {
244 int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
245 if (x_.is_null())
246 field_h(cell, basis) = (*xEvRoGed_)[lid];
247 else
248 field_h(cell, basis) = (*x_)[lid];
249 } // end loop over the basis functions
250 } // end loop over localCellIds
251 Kokkos::deep_copy(field.get_static_view(), field_h);
252 } // end loop over the fields to be gathered
253
254} // end of evaluateFields() (Residual Specialization)
255
257//
258// Initializing Constructor (Tangent Specialization)
259//
261template<typename TRAITS, typename LO, typename GO>
264 const Teuchos::RCP<const panzer::GlobalIndexer>& indexer,
265 const Teuchos::ParameterList& p)
266 :
267 globalIndexer_(indexer),
268 hasTangentFields_(false)
269{
270 using panzer::PureBasis;
271 using PHX::MDField;
272 using PHX::print;
273 using std::size_t;
274 using std::string;
275 using std::vector;
276 using Teuchos::RCP;
277 using vvstring = std::vector<std::vector<std::string>>;
279 input.setParameterList(p);
280 const vector<string>& names = input.getDofNames();
281 RCP<const PureBasis> basis = input.getBasis();
282 const vvstring& tangentFieldNames = input.getTangentNames();
283 indexerNames_ = input.getIndexerNames();
284 useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
285 globalDataKey_ = input.getGlobalDataKey();
286
287 // Allocate the fields.
288 int numFields(names.size());
289 gatherFields_.resize(numFields);
290 for (int fd(0); fd < numFields; ++fd)
291 {
292 gatherFields_[fd] =
293 MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
294 this->addEvaluatedField(gatherFields_[fd]);
295 } // end loop over names
296
297 // Set up the dependent tangent fields, if requested.
298 if (tangentFieldNames.size() > 0)
299 {
300 TEUCHOS_ASSERT(gatherFields_.size() == tangentFieldNames.size())
301 hasTangentFields_ = true;
302 tangentFields_.resize(numFields);
303 for (int fd(0); fd < numFields; ++fd)
304 {
305 int numTangentFields(tangentFieldNames[fd].size());
306 tangentFields_[fd].resize(numTangentFields);
307 for (int i(0); i < numTangentFields; ++i)
308 {
309 tangentFields_[fd][i] =
310 MDField<const ScalarT, Cell, NODE>(tangentFieldNames[fd][i],
311 basis->functional);
312 this->addDependentField(tangentFields_[fd][i]);
313 } // end loop over tangeng_field_names[fd]
314 } // end loop over gatherFields_
315 } // end if (tangentFieldNames.size() > 0)
316
317 // Figure out what the first active name is.
318 string firstName("<none>");
319 if (numFields > 0)
320 firstName = names[0];
321 string n("GatherSolution (Epetra): " + firstName + " (" +
322 print<EvalT>() + ")");
323 this->setName(n);
324} // end of Initializing Constructor (Tangent Specialization)
325
327//
328// postRegistrationSetup() (Tangent Specialization)
329//
331template<typename TRAITS, typename LO, typename GO>
332void
335 typename TRAITS::SetupData /* d */,
337{
338 using std::logic_error;
339 using std::size_t;
340 using std::string;
341 TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size())
342 int numFields(gatherFields_.size());
343 fieldIds_.resize(numFields);
344 for (int fd(0); fd < numFields; ++fd)
345 {
346 // Get the field ID from the DOF manager.
347 const string& fieldName(indexerNames_[fd]);
348 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
349
350 // This is the error return code; raise the alarm.
351 TEUCHOS_TEST_FOR_EXCEPTION(fieldIds_[fd] == -1, logic_error,
352 "GatherSolution_Epetra<Tangent>: Could not find field \"" + fieldName
353 + "\" in the global indexer. ")
354 } // end loop over gatherFields_
355 indexerNames_.clear();
356} // end of postRegistrationSetup() (Tangent Specialization)
357
359//
360// preEvaluate() (Tangent Specialization)
361//
363template<typename TRAITS, typename LO, typename GO>
364void
367 typename TRAITS::PreEvalData d)
368{
369 using std::string;
370 using Teuchos::RCP;
371 using Teuchos::rcp_dynamic_cast;
375 using LOC = panzer::LinearObjContainer;
377
378 // First try the refactored ReadOnly container.
379 RCP<GED> ged;
380 string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
381 if (d.gedc->containsDataObject(globalDataKey_ + post))
382 {
383 ged = d.gedc->getDataObject(globalDataKey_ + post);
384 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
385 return;
386 } // end of the refactored ReadOnly way
387
388 // Now try the old path.
389 ged = d.gedc->getDataObject(globalDataKey_);
390 {
391 // Try to extract the linear object container.
392 auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
393 auto locPair = rcp_dynamic_cast<LPGED>(ged);
394 if (not locPair.is_null())
395 {
396 RCP<LOC> loc = locPair->getGhostedLOC();
397 epetraContainer = rcp_dynamic_cast<ELOC>(loc);
398 } // end if (not locPair.is_null())
399 if (not epetraContainer.is_null())
400 {
401 if (useTimeDerivativeSolutionVector_)
402 x_ = epetraContainer->get_dxdt();
403 else // if (not useTimeDerivativeSolutionVector_)
404 x_ = epetraContainer->get_x();
405 return;
406 } // end if (not epetraContainer.is_null())
407 } // end of the old path
408
409 // As a last resort, try to extract an EpetraVector_ReadOnly object. This
410 // will throw an exception if it doesn't work.
411 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
412} // end of preEvaluate() (Tangent Specialization)
413
415//
416// evaluateFields() (Tangent Specialization)
417//
419template<typename TRAITS, typename LO, typename GO>
420void
423 typename TRAITS::EvalData workset)
424{
425 using PHX::MDField;
426 using std::size_t;
427 using std::string;
428 using std::vector;
429 using Teuchos::ArrayRCP;
430 using Teuchos::ptrFromRef;
431 using Teuchos::RCP;
432 using Teuchos::rcp_dynamic_cast;
433 using Thyra::SpmdVectorBase;
434
435 // For convenience, pull out some objects from the workset.
436 string blockId(this->wda(workset).block_id);
437 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
438 int numCells(localCellIds.size()), numFields(gatherFields_.size());
439
440 // NOTE: A reordering of these loops will likely improve performance. The
441 // "getGIDFieldOffsets may be expensive. However the "getElementGIDs"
442 // can be cheaper. However the lookup for LIDs may be more expensive!
443 auto LIDs = globalIndexer_->getLIDs();
444 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
445 Kokkos::deep_copy(LIDs_h, LIDs);
446 // Loop over the fields to be gathered.
447 for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
448 {
449 MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
450 auto field_h = Kokkos::create_mirror_view(field.get_static_view());
451 int fieldNum(fieldIds_[fieldInd]);
452 const vector<int>& elmtOffset =
453 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
454 int numBases(elmtOffset.size());
455 // Gather operation for each cell in the workset.
456 for (int cell(0); cell < numCells; ++cell)
457 {
458 size_t cellLocalId(localCellIds[cell]);
459 // Loop over the basis functions and fill the fields.
460 for (int basis(0); basis < numBases; ++basis)
461 {
462 int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
463 if (x_.is_null())
464 field_h(cell, basis) = (*xEvRoGed_)[lid];
465 else
466 field_h(cell, basis) = (*x_)[lid];
467 } // end loop over the basis functions
468 } // end loop over localCellIds
469
470 // Deal with the tangent fields.
471 if (hasTangentFields_) {
472 // Loop over the tangent fields and fill them in.
473 int numTangentFields(tangentFields_[fieldInd].size());
474 for (int i(0); i < numTangentFields; ++i){
475 auto tf = Kokkos::create_mirror_view(tangentFields_[fieldInd][i].get_static_view());
476 Kokkos::deep_copy(tf, tangentFields_[fieldInd][i].get_static_view());
477 // Loop over the cells in the workset.
478 for (int cell(0); cell < numCells; ++cell) {
479 // Loop over the fields to be gathered.
480 int fieldNum2(fieldIds_[fieldInd]);
481 const vector<int>& elmtOffset2 =
482 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum2);
483 int numBases2(elmtOffset2.size());
484
485 // Loop over the basis functions.
486 for (int basis(0); basis < numBases2; ++basis){
487 field_h(cell, basis).fastAccessDx(i) =
488 tf(cell, basis).val();
489 } // end loop over the basis functions
490 } // end loop over the cells in the workset
491 } // end loop over numTangentFields
492 } // end if (hasTangentFields_)
493 Kokkos::deep_copy(field.get_static_view(), field_h);
494 } // end loop over the fields to be gathered
495} // end of evaluateFields() (Tangent Specialization)
496
498//
499// Initializing Constructor (Jacobian Specialization)
500//
502template<typename TRAITS, typename LO, typename GO>
505 const Teuchos::RCP<const panzer::GlobalIndexer>& indexer,
506 const Teuchos::ParameterList& p)
507 :
508 globalIndexer_(indexer)
509{
510 using panzer::PureBasis;
511 using PHX::MDField;
512 using PHX::print;
513 using std::size_t;
514 using std::string;
515 using std::vector;
516 using Teuchos::RCP;
518 input.setParameterList(p);
519 const vector<string>& names = input.getDofNames();
520 RCP<const PureBasis> basis = input.getBasis();
521 indexerNames_ = input.getIndexerNames();
522 useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
523 globalDataKey_ = input.getGlobalDataKey();
524 gatherSeedIndex_ = input.getGatherSeedIndex();
525 sensitivitiesName_ = input.getSensitivitiesName();
526 disableSensitivities_ = not input.firstSensitivitiesAvailable();
527
528 // Allocate the fields.
529 int numFields(names.size());
530 gatherFields_.resize(numFields);
531 for (int fd(0); fd < numFields; ++fd)
532 {
533 MDField<ScalarT, Cell, NODE> f(names[fd], basis->functional);
534 gatherFields_[fd] = f;
535 this->addEvaluatedField(gatherFields_[fd]);
536 } // end loop over names
537
538 // Figure out what the first active name is.
539 string firstName("<none>"), n("GatherSolution (Epetra");
540 if (numFields > 0)
541 firstName = names[0];
542 if (disableSensitivities_)
543 n += ", No Sensitivities";
544 n += "): " + firstName + " (Jacobian)";
545 this->setName(n);
546} // end of Initializing Constructor (Jacobian Specialization)
547
549//
550// postRegistrationSetup() (Jacobian Specialization)
551//
553template<typename TRAITS, typename LO, typename GO>
554void
557 typename TRAITS::SetupData /* d */,
559{
560 using std::logic_error;
561 using std::size_t;
562 using std::string;
563 TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_.size())
564 int numFields(gatherFields_.size());
565 fieldIds_.resize(numFields);
566 for (int fd(0); fd < numFields; ++fd)
567 {
568 // Get the field ID from the DOF manager.
569 const string& fieldName(indexerNames_[fd]);
570 fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
571
572 // This is the error return code; raise the alarm.
573 TEUCHOS_TEST_FOR_EXCEPTION(fieldIds_[fd] == -1, logic_error,
574 "GatherSolution_Epetra<Jacobian>: Could not find field \"" +
575 fieldName + "\" in the global indexer. ")
576 } // end loop over gatherFields_
577 indexerNames_.clear();
578} // end of postRegistrationSetup() (Jacobian Specialization)
579
581//
582// preEvaluate() (Jacobian Specialization)
583//
585template<typename TRAITS, typename LO, typename GO>
586void
589 typename TRAITS::PreEvalData d)
590{
591 using std::string;
592 using Teuchos::RCP;
593 using Teuchos::rcp_dynamic_cast;
597 using LOC = panzer::LinearObjContainer;
599
600 // Manage sensitivities.
601 applySensitivities_ = false;
602 if ((not disableSensitivities_ ) and
603 (d.first_sensitivities_name == sensitivitiesName_))
604 applySensitivities_ = true;
605
606 // First try the refactored ReadOnly container.
607 RCP<GED> ged;
608 string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
609 if (d.gedc->containsDataObject(globalDataKey_ + post))
610 {
611 ged = d.gedc->getDataObject(globalDataKey_ + post);
612 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
613 return;
614 } // end of the refactored ReadOnly way
615
616 // Now try the old path.
617 ged = d.gedc->getDataObject(globalDataKey_);
618 {
619 // Try to extract the linear object container.
620 auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
621 auto locPair = rcp_dynamic_cast<LPGED>(ged);
622 if (not locPair.is_null())
623 {
624 RCP<LOC> loc = locPair->getGhostedLOC();
625 epetraContainer = rcp_dynamic_cast<ELOC>(loc);
626 } // end if (not locPair.is_null())
627 if (not epetraContainer.is_null())
628 {
629 if (useTimeDerivativeSolutionVector_)
630 x_ = epetraContainer->get_dxdt();
631 else // if (not useTimeDerivativeSolutionVector_)
632 x_ = epetraContainer->get_x();
633 return;
634 } // end if (not epetraContainer.is_null())
635 } // end of the old path
636
637 // As a last resort, try to extract an EpetraVector_ReadOnly object. This
638 // will throw an exception if it doesn't work.
639 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
640} // end of preEvaluate() (Jacobian Specialization)
641
643//
644// evaluateFields() (Jacobian Specialization)
645//
647template<typename TRAITS, typename LO, typename GO>
648void
651 typename TRAITS::EvalData workset)
652{
653 using PHX::MDField;
654 using std::size_t;
655 using std::string;
656 using std::vector;
657 using Teuchos::ArrayRCP;
658 using Teuchos::ptrFromRef;
659 using Teuchos::RCP;
660 using Teuchos::rcp_dynamic_cast;
661 using Thyra::SpmdVectorBase;
662
663 // For convenience, pull out some objects from the workset.
664 string blockId(this->wda(workset).block_id);
665 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
666 int numFields(gatherFields_.size()), numCells(localCellIds.size());
667
668 // Set a sensitivity seed value.
669 double seedValue(0);
670 if (applySensitivities_)
671 {
672 if ((useTimeDerivativeSolutionVector_) and (gatherSeedIndex_ < 0))
673 seedValue = workset.alpha;
674 else if (gatherSeedIndex_ < 0)
675 seedValue = workset.beta;
676 else if (not useTimeDerivativeSolutionVector_)
677 seedValue = workset.gather_seeds[gatherSeedIndex_];
678 else
679 TEUCHOS_ASSERT(false)
680 } // end if (applySensitivities_)
681
682 // Turn off sensitivies. This may be faster if we don't expand the term, but
683 // I suspect not, because anywhere it is used the full complement of
684 // sensitivies will be needed anyway.
685 if (not applySensitivities_)
686 seedValue = 0.0;
687
688 // Interface worksets handle DOFs from two element blocks. The derivative
689 // offset for the other element block must be shifted by the derivative side
690 // of my element block.
691 int dos(0);
692 if (this->wda.getDetailsIndex() == 1)
693 {
694 // Get the DOF count for my element block.
695 dos = globalIndexer_->getElementBlockGIDCount(workset.details(0).block_id);
696 } // end if (this->wda.getDetailsIndex() == 1)
697
698 // NOTE: A reordering of these loops will likely improve performance. The
699 // "getGIDFieldOffsets may be expensive. However the "getElementGIDs"
700 // can be cheaper. However the lookup for LIDs may be more expensive!
701 auto LIDs = globalIndexer_->getLIDs();
702 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
703 Kokkos::deep_copy(LIDs_h, LIDs);
704 // Loop over the fields to be gathered.
705 for (int fieldInd(0); fieldInd < numFields; ++fieldInd)
706 {
707 MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldInd];
708 auto field_h = Kokkos::create_mirror_view(field.get_static_view());
709 int fieldNum(fieldIds_[fieldInd]);
710 const vector<int>& elmtOffset =
711 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
712 int numBases(elmtOffset.size());
713 // Gather operation for each cell in the workset.
714 for (int cell(0); cell < numCells; ++cell)
715 {
716 size_t cellLocalId(localCellIds[cell]);
717 // Loop over the basis functions and fill the fields.
718 for (int basis(0); basis < numBases; ++basis)
719 {
720 int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
721 if (x_.is_null())
722 field_h(cell, basis) = (*xEvRoGed_)[lid];
723 else
724 field_h(cell, basis) = (*x_)[lid];
725 } // end loop over the basis functions
726 } // end loop over localCellIds
727
728 // Deal with the sensitivities.
729 if (applySensitivities_) {
730 int fieldNum2(fieldIds_[fieldInd]);
731 const vector<int>& elmtOffset2 =
732 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum2);
733 int numBases2(elmtOffset2.size());
734
735 // Gather operation for each cell in the workset.
736 for (int cell(0); cell < numCells; ++cell)
737 {
738 // Loop over the basis functions.
739 for (int basis(0); basis < numBases2; ++basis)
740 {
741 // Seed the FAD object.
742 int offset(elmtOffset2[basis]);
743
744 field_h(cell, basis).fastAccessDx(dos + offset) = seedValue;
745 } // end loop over the basis functions
746 } // end loop over localCellIds
747 } // end if (applySensitivities_)
748 Kokkos::deep_copy(field.get_static_view(), field_h);
749 } // end loop over the fields to be gathered
750
751} // end of evaluateFields() (Jacobian Specialization)
752
753#endif // __Panzer_GatherSolution_Epetra_impl_hpp__
int numFields
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.
This class provides a boundary exchange communication mechanism for vectors.
Gathers solution values from the Newton solution vector into the nodal fields of the field manager.
std::string getGlobalDataKey() const
Name of the global evaluation data container to use for the source vector (all types)
void setParameterList(const Teuchos::ParameterList &pl)
int getGatherSeedIndex() const
What index to use for initializing the seed (Jacobian and Hessian)
const std::vector< std::vector< std::string > > & getTangentNames() const
Get the name of the tangent fields (tangent only)
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
std::string getSensitivitiesName() const
The name of the sensitivities. Enables sensitivities at "preEvaluate" time (Jacobian and Hessian)
Teuchos::RCP< const PureBasis > getBasis() const
Basis definiting layout of dof names (all types)
const std::vector< std::string > & getDofNames() const
The names of the DOFs to be gathered (all types)
Description and data layouts associated with a particular basis.