Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Response_ExtremeValue_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_Response_ExtremeValue_impl_hpp__
12 #define __Panzer_Response_ExtremeValue_impl_hpp__
13 
14 #include "Teuchos_Comm.hpp"
15 #include "Teuchos_CommHelpers.hpp"
16 #include "Teuchos_dyn_cast.hpp"
17 
18 #include "PanzerDiscFE_config.hpp"
19 #ifdef PANZER_HAVE_EPETRA_STACK
20 #include "Epetra_LocalMap.h"
21 #endif
22 
23 #include "Sacado_Traits.hpp"
24 
25 namespace panzer {
26 
27 template <typename EvalT>
30 {
31  double locValue = Sacado::scalarValue(value);
32  double glbValue = 0.0;
33 
34  // do global summation
35  if(useMax_)
36  Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_MAX, static_cast<Thyra::Ordinal>(1), &locValue,&glbValue);
37  else
38  Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_MIN, static_cast<Thyra::Ordinal>(1), &locValue,&glbValue);
39 
40  value = glbValue;
41 
42 #ifdef PANZER_HAVE_EPETRA_STACK
43  // built data in vectors
44  if(this->useEpetra()) {
45  // use epetra
46  this->getEpetraVector()[0] = glbValue;
47  }
48  else
49 #endif
50  {
51  // use thyra
52  TEUCHOS_ASSERT(this->useThyra());
53 
54  this->getThyraVector()[0] = glbValue;
55  }
56 }
57 
58 template < >
61 {
62  using Teuchos::rcp_dynamic_cast;
63 
64  Teuchos::RCP<Thyra::MultiVectorBase<double> > dgdx_unique = getDerivative();
65 
66  uniqueContainer_ = linObjFactory_->buildLinearObjContainer();
67  Teuchos::rcp_dynamic_cast<ThyraObjContainer<double> >(uniqueContainer_)->set_x_th(dgdx_unique->col(0));
68 
69  linObjFactory_->ghostToGlobalContainer(*ghostedContainer_,*uniqueContainer_,LinearObjContainer::X);
70 
71  uniqueContainer_ = Teuchos::null;
72 }
73 
74 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
75 template < >
78 {
79  using Teuchos::rcp_dynamic_cast;
80 
81  Teuchos::RCP<Thyra::MultiVectorBase<double> > dgdx_unique = getDerivative();
82 
83  uniqueContainer_ = linObjFactory_->buildLinearObjContainer();
84  Teuchos::rcp_dynamic_cast<ThyraObjContainer<double> >(uniqueContainer_)->set_x_th(dgdx_unique->col(0));
85 
86  linObjFactory_->ghostToGlobalContainer(*ghostedContainer_,*uniqueContainer_,LinearObjContainer::X);
87 
88  uniqueContainer_ = Teuchos::null;
89 }
90 #endif
91 
92 template < >
95 {
96  const int n = value.size();
97  const int num_deriv = this->numDeriv();
98  TEUCHOS_ASSERT(n == 0 || n == num_deriv);
99  ScalarT glbValue = ScalarT(num_deriv, 0.0);
100 
101  // do global min/max on value
102  if(useMax_)
103  Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_MAX, Thyra::Ordinal(1), &value.val(), &glbValue.val());
104  else
105  Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_MIN, Thyra::Ordinal(1), &value.val(), &glbValue.val());
106 
107  // find the minimum processor who's local value == the global min/max value
108  if (num_deriv > 0) {
109  int locProc = value.val() == glbValue.val() ? this->getComm()->getRank() : this->getComm()->getSize();
110  int glbProc = 0;
111  Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_MIN, Thyra::Ordinal(1), &locProc, &glbProc);
112 
113  // now broadcast the derivatives from proc glbProc
114  Teuchos::broadcast(*this->getComm(), glbProc, Thyra::Ordinal(n), &glbValue.fastAccessDx(0));
115  }
116 
117  value = glbValue;
118 
119  // copy data in vectors
120 #ifdef PANZER_HAVE_EPETRA_STACK
121  if(this->useEpetra()) {
122  // use epetra
123  Epetra_MultiVector& deriv = this->getEpetraMultiVector();
124  for (int i=0; i<num_deriv; ++i)
125  deriv[i][0] = glbValue.dx(i);
126  }
127  else
128 #endif
129  {
130  // use thyra
131  TEUCHOS_ASSERT(this->useThyra());
132  Thyra::ArrayRCP< Thyra::ArrayRCP<double> > deriv = this->getThyraMultiVector();
133  for (int i=0; i<num_deriv; ++i)
134  deriv[i][0] = glbValue.dx(i);
135  }
136 }
137 
138 // Do nothing unless derivatives are actually required
139 template <typename EvalT>
142 
143 // derivatives are required for
144 template < >
147 {
148  setDerivativeVectorSpace(soln_vs);
149 }
150 
151 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
152 // derivatives are required for
153 template < >
156 {
157  setDerivativeVectorSpace(soln_vs);
158 }
159 #endif
160 
161 // Do nothing unless derivatives are required
162 template <typename EvalT>
164 adjustForDirichletConditions(const GlobalEvaluationData & /* localBCRows */, const GlobalEvaluationData & /* globalBCRows */) { }
165 
166 // Do nothing unless derivatives are required
167 template < >
170 {
171  linObjFactory_->adjustForDirichletConditions(Teuchos::dyn_cast<const LinearObjContainer>(localBCRows),
172  Teuchos::dyn_cast<const LinearObjContainer>(globalBCRows),
173  *ghostedContainer_,true,true);
174 }
175 
176 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
177 // Do nothing unless derivatives are required
178 template < >
181 {
182  linObjFactory_->adjustForDirichletConditions(Teuchos::dyn_cast<const LinearObjContainer>(localBCRows),
183  Teuchos::dyn_cast<const LinearObjContainer>(globalBCRows),
184  *ghostedContainer_,true,true);
185 }
186 #endif
187 
188 }
189 
190 #endif
virtual void scatterResponse()
This simply does global summation, then shoves the result into a vector.
void setSolnVectorSpace(const Teuchos::RCP< const Thyra::VectorSpaceBase< double > > &soln_vs)
Set solution vector space.
void adjustForDirichletConditions(const GlobalEvaluationData &localBCRows, const GlobalEvaluationData &globalBCRows)
#define TEUCHOS_ASSERT(assertion_test)