Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherTangent_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_GatherTangent_Epetra_impl_hpp__
12 #define __Panzer_GatherTangent_Epetra_impl_hpp__
13 
15 //
16 // Include Files
17 //
19 
20 // Epetra
21 #include "Epetra_Vector.h"
22 
23 // Panzer
26 #include "Panzer_PureBasis.hpp"
27 #include "Panzer_GlobalIndexer.hpp"
29 
30 // Teuchos
31 #include "Teuchos_Assert.hpp"
32 
33 // Thyra
34 #include "Thyra_SpmdVectorBase.hpp"
35 
37 //
38 // Initializing Constructor
39 //
41 template<typename EvalT, typename TRAITS, typename LO, typename GO>
45  const Teuchos::ParameterList& p)
46  :
47  globalIndexer_(indexer),
48  useTimeDerivativeSolutionVector_(false),
49  globalDataKey_("Tangent Gather Container")
50 {
51  using panzer::PureBasis;
52  using PHX::MDField;
53  using PHX::print;
54  using std::size_t;
55  using std::string;
56  using std::vector;
57  using Teuchos::RCP;
58 
59  // Get the necessary information from the ParameterList.
60  const vector<string>& names = *(p.get<RCP<vector<string>>>("DOF Names"));
61  indexerNames_ = p.get<RCP<vector<string>>>("Indexer Names");
63  if (p.isType<RCP<PureBasis>>("Basis"))
64  basis = p.get<RCP<PureBasis>>("Basis");
65  else // if (not p.isType<RCP<PureBasis>>("Basis"))
66  basis = p.get<RCP<const PureBasis>>("Basis");
67  if (p.isType<bool>("Use Time Derivative Solution Vector"))
69  p.get<bool>("Use Time Derivative Solution Vector");
70  if (p.isType<string>("Global Data Key"))
71  globalDataKey_ = p.get<string>("Global Data Key");
72 
73  // Allocate fields.
74  int numFields(names.size());
75  gatherFields_.resize(numFields);
76  for (int fd(0); fd < numFields; ++fd)
77  {
78  gatherFields_[fd] =
79  MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
80  this->addEvaluatedField(gatherFields_[fd]);
81 
82  // This fixes the case of dxdpEvRoGed_ being null and no
83  // operations performed during evalaute. Keeps the field with
84  // initial zero state.
85  this->addUnsharedField(gatherFields_[fd].fieldTag().clone());
86  } // end loop over names
87 
88  // Figure out what the first active name is.
89  string firstName("<none>");
90  if (numFields > 0)
91  firstName = names[0];
92  string n("GatherTangent (Epetra): " + firstName + " (" +
93  print<EvalT>() + ")");
94  this->setName(n);
95 } // end of Initializing Constructor
96 
98 //
99 // postRegistrationSetup()
100 //
102 template<typename EvalT, typename TRAITS, typename LO, typename GO>
103 void
106  typename TRAITS::SetupData /* d */,
107  PHX::FieldManager<TRAITS>& /* fm */)
108 {
109  using std::logic_error;
110  using std::size_t;
111  using std::string;
112  using Teuchos::null;
113  TEUCHOS_ASSERT(gatherFields_.size() == indexerNames_->size());
114  int numFields(gatherFields_.size());
115  fieldIds_.resize(numFields);
116  for (int fd(0); fd < numFields; ++fd)
117  {
118  // Get the field ID from the DOF manager.
119  const string& fieldName((*indexerNames_)[fd]);
120  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
121 
122  // This is the error return code; raise the alarm.
123  TEUCHOS_TEST_FOR_EXCEPTION(fieldIds_[fd] == -1, logic_error,
124  "GatherTangent_Epetra<Residual>: Could not find field \"" + fieldName +
125  "\" in the global indexer. ");
126  } // end loop over gatherFields_
127  indexerNames_ = null;
128 } // end of postRegistrationSetup()
129 
131 //
132 // preEvaluate()
133 //
135 template<typename EvalT, typename TRAITS, typename LO, typename GO>
136 void
139  typename TRAITS::PreEvalData d)
140 {
141  using Teuchos::RCP;
142  using Teuchos::rcp_dynamic_cast;
144  using GED = GlobalEvaluationData;
145  if (d.gedc->containsDataObject(globalDataKey_))
146  {
147  RCP<GED> ged = d.gedc->getDataObject(globalDataKey_);
148  dxdpEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
149  }
150 } // end of preEvaluate()
151 
153 //
154 // evaluateFields()
155 //
157 template<typename EvalT, typename TRAITS, typename LO, typename GO>
158 void
161  typename TRAITS::EvalData workset)
162 {
163  using PHX::MDField;
164  using std::size_t;
165  using std::string;
166  using std::vector;
167  using Teuchos::ArrayRCP;
168  using Teuchos::ptrFromRef;
169  using Teuchos::RCP;
170  using Teuchos::rcp_dynamic_cast;
171  using Thyra::SpmdVectorBase;
172 
173  // If no global evaluation data container was set, then this evaluator
174  // becomes a no-op.
175  if (dxdpEvRoGed_.is_null())
176  return;
177 
178  // For convenience, pull out some objects from the workset.
179  string blockId(this->wda(workset).block_id);
180  const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
181  int numCells(localCellIds.size()), numFields(gatherFields_.size());
182 
183  // NOTE: A reordering of these loops will likely improve performance. The
184  // "getGIDFieldOffsets may be expensive. However the "getElementGIDs"
185  // can be cheaper. However the lookup for LIDs may be more expensive!
186 
187  // Gather operation for each cell in the workset.
188  auto LIDs = globalIndexer_->getLIDs();
189  auto LIDs_h = Kokkos::create_mirror_view(LIDs);
190  Kokkos::deep_copy(LIDs_h, LIDs);
191  // Loop over the fields to be gathered.
192  for (int fieldIndex(0); fieldIndex < numFields; ++fieldIndex)
193  {
194  MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldIndex];
195  auto field_h = Kokkos::create_mirror_view(field.get_static_view());
196  for (int cell(0); cell < numCells; ++cell)
197  {
198  size_t cellLocalId(localCellIds[cell]);
199  int fieldNum(fieldIds_[fieldIndex]);
200  const vector<int>& elmtOffset =
201  globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
202  int numBases(elmtOffset.size());
203 
204  // Loop over the basis functions and fill the fields.
205  for (int basis(0); basis < numBases; ++basis)
206  {
207  int offset(elmtOffset[basis]), lid(LIDs_h(cellLocalId, offset));
208  field_h(cell, basis) = (*dxdpEvRoGed_)[lid];
209  } // end loop over the basis functions
210  } // end loop over the cells in the workset
211  Kokkos::deep_copy(field.get_static_view(), field_h);
212  } // end loop over the fields to be gathered
213 } // end of evaluateFields()
214 
215 #endif // __Panzer_GatherTangent_Epetra_impl_hpp__
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
The fields to be gathered.
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
Create a copy.
T & get(const std::string &name, T def_value)
std::string globalDataKey_
The key identifying the GlobalEvaluationData.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< std::vector< std::string > > indexerNames_
A list of the names of the fields to be gathered.
GatherTangent_Epetra()
Default Constructor (disabled).
int numFields
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
Post-Registration Setup.
This class provides a boundary exchange communication mechanism for vectors.
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
bool isType(const std::string &name) const
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)
void preEvaluate(typename TRAITS::PreEvalData d)
Pre-Evaluate: Sets the tangent vector.
void evaluateFields(typename TRAITS::EvalData d)
Evaluate Fields: Gather operation.
bool useTimeDerivativeSolutionVector_
A flag indicating whether we should be working with or .