Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_GatherSolution_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_GatherSolution_Epetra_Hessian_impl_hpp__
12#define __Panzer_GatherSolution_Epetra_Hessian_impl_hpp__
13
14// Only do this if required by the user.
15#ifdef Panzer_BUILD_HESSIAN_SUPPORT
16
18//
19// Include Files
20//
22
23// Epetra
24#include "Epetra_Map.h"
25#include "Epetra_Vector.h"
26
27// Panzer
32#include "Panzer_PureBasis.hpp"
34
35// Phalanx
36#include "Phalanx_DataLayout.hpp"
37
38// Teuchos
39#include "Teuchos_Assert.hpp"
40#include "Teuchos_FancyOStream.hpp"
41
42// Thyra
43#include "Thyra_SpmdVectorBase.hpp"
44
46//
47// Initializing Constructor (Hessian Specialization)
48//
50template<typename TRAITS, typename LO, typename GO>
53 const Teuchos::RCP<const panzer::GlobalIndexer>& indexer,
54 const Teuchos::ParameterList& p)
55 :
56 globalIndexer_(indexer)
57{
59 using PHX::MDField;
60 using PHX::print;
61 using std::size_t;
62 using std::string;
63 using std::vector;
64 using Teuchos::RCP;
66 input.setParameterList(p);
67 const vector<string>& names = input.getDofNames();
68 RCP<const PureBasis> basis = input.getBasis();
69 indexerNames_ = input.getIndexerNames();
70 useTimeDerivativeSolutionVector_ = input.useTimeDerivativeSolutionVector();
71 globalDataKey_ = input.getGlobalDataKey();
72 gatherSeedIndex_ = input.getGatherSeedIndex();
73 sensitivitiesName_ = input.getSensitivitiesName();
74 firstSensitivitiesAvailable_ = input.firstSensitivitiesAvailable();
75 secondSensitivitiesAvailable_ = input.secondSensitivitiesAvailable();
76 sensitivities2ndPrefix_ = input.getSecondSensitivityDataKeyPrefix();
77
78 // Allocate the fields.
79 int numFields(names.size());
80 gatherFields_.resize(numFields);
81 for (int fd(0); fd < numFields; ++fd)
82 {
83 gatherFields_[fd] =
84 MDField<ScalarT, Cell, NODE>(names[fd], basis->functional);
85 this->addEvaluatedField(gatherFields_[fd]);
86 } // end loop over names
87
88 // Figure out what the first active name is.
89 string firstName("<none>"), n("Gather Solution (Epetra");
90 if (numFields > 0)
91 firstName = names[0];
92 if (not firstSensitivitiesAvailable_)
93 n += ", No First Sensitivities";
94 n += "): " + firstName + " (" + print<EvalT>() + ") ";
95 this->setName(n);
96} // end of Initializing Constructor (Hessian Specialization)
97
99//
100// postRegistrationSetup() (Hessian Specialization)
101//
103template<typename TRAITS, typename LO, typename GO>
104void
107 typename TRAITS::SetupData /* d */,
109{
110 using std::logic_error;
111 using std::size_t;
112 using std::string;
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 "GatherSolution_Epetra<Hessian>: Could not find field \"" + fieldName
125 + "\" in the global indexer. ")
126 } // end loop over gatherFields_
127 indexerNames_.clear();
128} // end of postRegistrationSetup() (Hessian Specialization)
129
131//
132// preEvaluate() (Hessian Specialization)
133//
135template<typename TRAITS, typename LO, typename GO>
136void
139 typename TRAITS::PreEvalData d)
140{
141 using std::logic_error;
142 using std::string;
143 using Teuchos::RCP;
144 using Teuchos::rcp_dynamic_cast;
148 using LOC = panzer::LinearObjContainer;
150
151 // Manage sensitivities.
152 firstApplySensitivities_ = false;
153 if ((firstSensitivitiesAvailable_ ) and
154 (d.first_sensitivities_name == sensitivitiesName_))
155 firstApplySensitivities_ = true;
156 secondApplySensitivities_ = false;
157 if ((secondSensitivitiesAvailable_ ) and
158 (d.second_sensitivities_name == sensitivitiesName_))
159 secondApplySensitivities_ = true;
160
161 // First try refactored ReadOnly container.
162 RCP<GED> ged;
163 string post(useTimeDerivativeSolutionVector_ ? " - Xdot" : " - X");
164 if (d.gedc->containsDataObject(globalDataKey_ + post))
165 {
166 ged = d.gedc->getDataObject(globalDataKey_ + post);
167 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
168 }
169
170 // Otherwise, try the old path.
171 else if (d.gedc->containsDataObject(globalDataKey_))
172 {
173 ged = d.gedc->getDataObject(globalDataKey_);
174
175 // Try to extract the linear object container.
176 xEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged);
177 auto epetraContainer = rcp_dynamic_cast<ELOC>(ged);
178
179 // Handle the LOCPair case.
180 {
181 auto locPair = rcp_dynamic_cast<LPGED>(ged);
182 if (not locPair.is_null())
183 {
184 RCP<LOC> loc = locPair->getGhostedLOC();
185 epetraContainer = rcp_dynamic_cast<ELOC>(loc);
186 } // end if (not locPair.is_null())
187 } // end of the LOCPair case
188 if ((xEvRoGed_.is_null()) and (not epetraContainer.is_null()))
189 {
190 if (useTimeDerivativeSolutionVector_)
191 x_ = epetraContainer->get_dxdt();
192 else // if (not useTimeDerivativeSolutionVector_)
193 x_ = epetraContainer->get_x();
194 } // end if ((xEvRoGed_.is_null()) and (not epetraContainer.is_null()))
195 } // end if we're doing things the new or old way
196
197 // Ensure that we actually have something.
198 TEUCHOS_TEST_FOR_EXCEPTION((x_.is_null()) and (xEvRoGed_.is_null()),
199 logic_error, "GatherSolution_Epetra_Hessian::preEvaluate(): Unable to " \
200 "find solution vector.")
201
202 // Don't try to extract dx if it's not required.
203 if (not secondApplySensitivities_)
204 return;
205
206 // Now parse the second derivative direction.
207 if (d.gedc->containsDataObject(sensitivities2ndPrefix_ + globalDataKey_))
208 {
209 ged = d.gedc->getDataObject(sensitivities2ndPrefix_ + globalDataKey_);
210 dxEvRoGed_ = rcp_dynamic_cast<EVROGED>(ged, true);
211
212 // Ensure that we actually have something.
213 TEUCHOS_TEST_FOR_EXCEPTION(dxEvRoGed_.is_null(), logic_error, "Cannot " \
214 "find sensitivity vector associated with \"" + sensitivities2ndPrefix_ +
215 globalDataKey_ + "\" and \"" + post + "\".")
216 } // end of parsing the second derivative direction
217} // end of preEvaluate() (Hessian Specialization)
218
220//
221// evaluateFields() (Hessian Specialization)
222//
224template<typename TRAITS, typename LO, typename GO>
225void
228 typename TRAITS::EvalData workset)
229{
230 using PHX::MDField;
231 using std::size_t;
232 using std::string;
233 using std::vector;
234 using Teuchos::ArrayRCP;
235 using Teuchos::ptrFromRef;
236 using Teuchos::RCP;
237 using Teuchos::rcp_dynamic_cast;
238 using Thyra::SpmdVectorBase;
239
240 // For convenience, pull out some objects from the workset.
241 string blockId(this->wda(workset).block_id);
242 const vector<size_t>& localCellIds = this->wda(workset).cell_local_ids;
243 int numFields(gatherFields_.size()), numCells(localCellIds.size());
244
245 // Set a sensitivity seed value.
246 double seedValue(0);
247 if (firstApplySensitivities_)
248 {
249 if ((useTimeDerivativeSolutionVector_) and (gatherSeedIndex_ < 0))
250 seedValue = workset.alpha;
251 else if (gatherSeedIndex_ < 0)
252 seedValue = workset.beta;
253 else if (not useTimeDerivativeSolutionVector_)
254 seedValue = workset.gather_seeds[gatherSeedIndex_];
255 else
256 TEUCHOS_ASSERT(false);
257 } // end if (firstApplySensitivities_)
258
259 // Turn off sensitivies. This may be faster if we don't expand the term, but
260 // I suspect not, because anywhere it is used the full complement of
261 // sensitivies will be needed anyway.
262 if (not firstApplySensitivities_)
263 seedValue = 0.0;
264
265 // Interface worksets handle DOFs from two element blocks. The derivative
266 // offset for the other element block must be shifted by the derivative side
267 // of my element block.
268 int dos(0);
269 if (this->wda.getDetailsIndex() == 1)
270 {
271 // Get the DOF count for my element block.
272 dos = globalIndexer_->getElementBlockGIDCount(workset.details(0).block_id);
273 } // end if (this->wda.getDetailsIndex() == 1)
274
275 if (x_.is_null())
276 {
277 // Loop over the fields to be gathered.
278 for (int fieldIndex(0); fieldIndex < numFields; ++fieldIndex)
279 {
280 MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldIndex];
281 int fieldNum(fieldIds_[fieldIndex]);
282 const vector<int>& elmtOffset =
283 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
284 int numBases(elmtOffset.size());
285
286 // Gather operation for each cell in the workset.
287 for (int cell(0); cell < numCells; ++cell)
288 {
289 size_t cellLocalId(localCellIds[cell]);
290 auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
291
292 // Loop over the basis functions and fill the fields.
293 for (int basis(0); basis < numBases; ++basis)
294 {
295 int offset(elmtOffset[basis]), lid(LIDs[offset]);
296 field(cell, basis) = (*xEvRoGed_)[lid];
297 } // end loop over the basis functions
298 } // end loop over localCellIds
299 } // end loop over the fields to be gathered
300 }
301 else // if (not x_.is_null())
302 {
303 // Loop over the fields to be gathered.
304 for (int fieldIndex(0); fieldIndex < numFields; ++fieldIndex)
305 {
306 MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldIndex];
307 int fieldNum(fieldIds_[fieldIndex]);
308 const vector<int>& elmtOffset =
309 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
310 int numBases(elmtOffset.size());
311
312 // Gather operation for each cell in the workset.
313 for (int cell(0); cell < numCells; ++cell)
314 {
315 size_t cellLocalId(localCellIds[cell]);
316 auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
317
318 // Loop over the basis functions and fill the fields.
319 for (int basis(0); basis < numBases; ++basis)
320 {
321 int offset(elmtOffset[basis]), lid(LIDs[offset]);
322 field(cell, basis) = (*x_)[lid];
323 } // end loop over the basis functions
324 } // end loop over localCellIds
325 } // end loop over the fields to be gathered
326 } // end if (x_.is_null()) or not
327
328 // Deal with the first sensitivities.
329 if (firstApplySensitivities_)
330 {
331 // Loop over the fields to be gathered.
332 for (int fieldIndex(0); fieldIndex < numFields; ++fieldIndex)
333 {
334 MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldIndex];
335 int fieldNum(fieldIds_[fieldIndex]);
336 const vector<int>& elmtOffset =
337 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
338 int numBases(elmtOffset.size());
339
340 // Gather operation for each cell in the workset.
341 for (int cell(0); cell < numCells; ++cell)
342 {
343 // Loop over the basis functions and fill the fields.
344 for (int basis(0); basis < numBases; ++basis)
345 {
346 int offset(elmtOffset[basis]);
347 field(cell, basis).fastAccessDx(dos + offset) = seedValue;
348 } // end loop over the basis functions
349 } // end loop over localCellIds
350 } // end loop over the fields to be gathered
351 } // end if (firstApplySensitivities_)
352
353 // Deal with the second sensitivities.
354 if (secondApplySensitivities_)
355 {
356 // Loop over the fields to be gathered.
357 for (int fieldIndex(0); fieldIndex < numFields; ++fieldIndex)
358 {
359 MDField<ScalarT, Cell, NODE>& field = gatherFields_[fieldIndex];
360 int fieldNum(fieldIds_[fieldIndex]);
361 const vector<int>& elmtOffset =
362 globalIndexer_->getGIDFieldOffsets(blockId, fieldNum);
363 int numBases(elmtOffset.size());
364
365 // Gather operation for each cell in the workset.
366 for (int cell(0); cell < numCells; ++cell)
367 {
368 size_t cellLocalId(localCellIds[cell]);
369 auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
370
371 // Loop over the basis functions and fill the fields.
372 for (int basis(0); basis < numBases; ++basis)
373 {
374 int offset(elmtOffset[basis]), lid(LIDs[offset]);
375 field(cell, basis).val().fastAccessDx(0) = (*dxEvRoGed_)[lid];
376 } // end loop over the basis functions
377 } // end loop over localCellIds
378 } // end loop over the fields to be gathered
379 } // end if (secondApplySensitivities_)
380} // end of evaluateFields() (Hessian Specialization)
381
382#endif // Panzer_BUILD_HESSIAN_SUPPORT
383
384#endif // __Panzer_GatherSolution_Epetra_Hessian_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)
std::string getSecondSensitivityDataKeyPrefix()
What prefix to use for the GEDC for second derivative sensitivity direction (Hessian 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)
bool secondSensitivitiesAvailable()
Are second derivative sensitivies enabled or disabled (Hessian only)
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.