Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ScatterResidual_Epetra_Hessian_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_ScatterResidual_Epetra_Hessian_impl_hpp__
12 #define __Panzer_ScatterResidual_Epetra_Hessian_impl_hpp__
13 
14 // only do this if required by the user
15 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
16 
17 #include "Teuchos_RCP.hpp"
18 #include "Teuchos_Assert.hpp"
19 
20 #include "Phalanx_DataLayout.hpp"
21 
22 #include "Epetra_Map.h"
23 #include "Epetra_Vector.h"
24 #include "Epetra_CrsMatrix.h"
25 
26 #include "Panzer_GlobalIndexer.hpp"
28 #include "Panzer_PureBasis.hpp"
33 
34 #include "Phalanx_DataLayout_MDALayout.hpp"
35 
36 #include "Teuchos_FancyOStream.hpp"
37 
38 // **********************************************************************
39 // Specialization: Hessian
40 // **********************************************************************
41 
42 template<typename TRAITS,typename LO,typename GO>
46  const Teuchos::ParameterList& p,
47  bool useDiscreteAdjoint)
48  : globalIndexer_(indexer)
49  , colGlobalIndexer_(cIndexer)
50  , globalDataKey_("Residual Scatter Container")
51  , useDiscreteAdjoint_(useDiscreteAdjoint)
52 {
53  std::string scatterName = p.get<std::string>("Scatter Name");
54  scatterHolder_ =
56 
57  // get names to be evaluated
58  const std::vector<std::string>& names =
59  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Dependent Names"));
60 
61  // grab map from evaluated names to field names
62  fieldMap_ = p.get< Teuchos::RCP< std::map<std::string,std::string> > >("Dependent Map");
63 
65  p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis")->functional;
66 
67  // build the vector of fields that this is dependent on
68  scatterFields_.resize(names.size());
69  for (std::size_t eq = 0; eq < names.size(); ++eq) {
70  scatterFields_[eq] = PHX::MDField<const ScalarT,Cell,NODE>(names[eq],dl);
71 
72  // tell the field manager that we depend on this field
73  this->addDependentField(scatterFields_[eq]);
74  }
75 
76  // this is what this evaluator provides
77  this->addEvaluatedField(*scatterHolder_);
78 
79  if (p.isType<std::string>("Global Data Key"))
80  globalDataKey_ = p.get<std::string>("Global Data Key");
81  if (p.isType<bool>("Use Discrete Adjoint"))
82  useDiscreteAdjoint = p.get<bool>("Use Discrete Adjoint");
83 
84  // discrete adjoint does not work with non-square matrices
85  if(useDiscreteAdjoint)
86  { TEUCHOS_ASSERT(colGlobalIndexer_==globalIndexer_); }
87 
88  this->setName(scatterName+" Scatter Residual Epetra (Jacobian)");
89 }
90 
91 // **********************************************************************
92 template<typename TRAITS,typename LO,typename GO>
94 postRegistrationSetup(typename TRAITS::SetupData /* d */,
95  PHX::FieldManager<TRAITS>& /* fm */)
96 {
97  fieldIds_.resize(scatterFields_.size());
98  // load required field numbers for fast use
99  for(std::size_t fd=0;fd<scatterFields_.size();++fd) {
100  // get field ID from DOF manager
101  std::string fieldName = fieldMap_->find(scatterFields_[fd].fieldTag().name())->second;
102  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
103  }
104 }
105 
106 // **********************************************************************
107 template<typename TRAITS,typename LO,typename GO>
109 preEvaluate(typename TRAITS::PreEvalData d)
110 {
111  // extract linear object container
112  epetraContainer_ = Teuchos::rcp_dynamic_cast<EpetraLinearObjContainer>(d.gedc->getDataObject(globalDataKey_));
113 
114  if(epetraContainer_==Teuchos::null) {
115  // extract linear object container
116  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
117  epetraContainer_ = Teuchos::rcp_dynamic_cast<EpetraLinearObjContainer>(loc);
118  }
119 }
120 
121 // **********************************************************************
122 template<typename TRAITS,typename LO,typename GO>
124 evaluateFields(typename TRAITS::EvalData workset)
125 {
126  using PHX::View;
127  using PHX::Device;
128  using std::vector;
129 
130  std::vector<double> jacRow;
131 
132  bool useColumnIndexer = colGlobalIndexer_!=Teuchos::null;
133 
134  // for convenience pull out some objects from workset
135  std::string blockId = this->wda(workset).block_id;
136  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
137 
138  Teuchos::RCP<Epetra_Vector> r = epetraContainer_->get_f();
139  Teuchos::RCP<Epetra_CrsMatrix> Jac = epetraContainer_->get_A();
140 
142  colGlobalIndexer = useColumnIndexer ? colGlobalIndexer_ : globalIndexer_;
143 
144  // NOTE: A reordering of these loops will likely improve performance
145  // The "getGIDFieldOffsets" may be expensive. However the
146  // "getElementGIDs" can be cheaper. However the lookup for LIDs
147  // may be more expensive!
148 
149  // scatter operation for each cell in workset
150  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
151  std::size_t cellLocalId = localCellIds[worksetCellIndex];
152 
153  auto rLIDs = globalIndexer_->getElementLIDs(cellLocalId);
154  auto initial_cLIDs = colGlobalIndexer->getElementLIDs(cellLocalId);
155  vector<int> cLIDs;
156  for (int i(0); i < static_cast<int>(initial_cLIDs.extent(0)); ++i)
157  cLIDs.push_back(initial_cLIDs(i));
158  if (Teuchos::nonnull(workset.other)) {
159  const std::size_t other_cellLocalId = workset.other->cell_local_ids[worksetCellIndex];
160  auto other_cLIDs = colGlobalIndexer->getElementLIDs(other_cellLocalId);
161  for (int i(0); i < static_cast<int>(other_cLIDs.extent(0)); ++i)
162  cLIDs.push_back(other_cLIDs(i));
163  }
164 
165  // loop over each field to be scattered
166  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++) {
167  int fieldNum = fieldIds_[fieldIndex];
168  const std::vector<int> & elmtOffset = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
169 
170  // loop over the basis functions (currently they are nodes)
171  for(std::size_t rowBasisNum = 0; rowBasisNum < elmtOffset.size(); rowBasisNum++) {
172  const ScalarT scatterField = (scatterFields_[fieldIndex])(worksetCellIndex,rowBasisNum);
173  int rowOffset = elmtOffset[rowBasisNum];
174  int row = rLIDs[rowOffset];
175 
176  // loop over the sensitivity indices: all DOFs on a cell
177  jacRow.resize(scatterField.size());
178 
179  for(int sensIndex=0;sensIndex<scatterField.size();++sensIndex)
180  jacRow[sensIndex] = scatterField.fastAccessDx(sensIndex).fastAccessDx(0);
181 
182  {
183  int err = Jac->SumIntoMyValues(
184  row,
185  std::min(cLIDs.size(), static_cast<size_t>(scatterField.size())),
186  panzer::ptrFromStlVector(jacRow),
187  panzer::ptrFromStlVector(cLIDs));
189  }
190  } // end rowBasisNum
191  } // end fieldIndex
192  }
193 }
194 
195 // **********************************************************************
196 
197 #endif
198 
199 #endif
const Kokkos::View< const panzer::LocalOrdinal *, Kokkos::LayoutRight, PHX::Device > getElementLIDs(panzer::LocalOrdinal localElmtId) const
T & get(const std::string &name, T def_value)
bool nonnull(const std::shared_ptr< T > &p)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
T * ptrFromStlVector(std::vector< T > &v)
Pushes residual values into the residual vector for a Newton-based solve.
int SumIntoMyValues(int MyRow, int NumEntries, const double *Values, const int *Indices)
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)