Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Response_Functional_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_Functional_impl_hpp__
12#define __Panzer_Response_Functional_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
25namespace panzer {
26
27template <typename EvalT>
30{
31 double locValue = Sacado::scalarValue(value);
32 double glbValue = 0.0;
33
34 // do global summation
35 Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_SUM, static_cast<Thyra::Ordinal>(1), &locValue,&glbValue);
36
37 value = glbValue;
38
39 // built data in vectors
40#ifdef PANZER_HAVE_EPETRA_STACK
41 if(this->useEpetra()) {
42 // use epetra
43 this->getEpetraVector()[0] = glbValue;
44 }
45 else
46 #endif
47 {
48 // use thyra
49 TEUCHOS_ASSERT(this->useThyra());
50
51 this->getThyraVector()[0] = glbValue;
52 }
53}
54
55template < >
58{
59 using Teuchos::rcp_dynamic_cast;
60
61 Teuchos::RCP<Thyra::MultiVectorBase<double> > dgdx_unique = getDerivative();
62
63 // if its null, don't do anything
64 if(dgdx_unique==Teuchos::null)
65 return;
66
67 uniqueContainer_ = linObjFactory_->buildLinearObjContainer();
68 Teuchos::rcp_dynamic_cast<ThyraObjContainer<double> >(uniqueContainer_)->set_x_th(dgdx_unique->col(0));
69
70 linObjFactory_->ghostToGlobalContainer(*ghostedContainer_,*uniqueContainer_,LinearObjContainer::X);
71
72 uniqueContainer_ = Teuchos::null;
73}
74
75#ifdef Panzer_BUILD_HESSIAN_SUPPORT
76template < >
79{
80 using Teuchos::rcp_dynamic_cast;
81
82 Teuchos::RCP<Thyra::MultiVectorBase<double> > dgdx_unique = getDerivative();
83
84 // if its null, don't do anything
85 if(dgdx_unique==Teuchos::null)
86 return;
87
88 uniqueContainer_ = linObjFactory_->buildLinearObjContainer();
89 Teuchos::rcp_dynamic_cast<ThyraObjContainer<double> >(uniqueContainer_)->set_x_th(dgdx_unique->col(0));
90
91 linObjFactory_->ghostToGlobalContainer(*ghostedContainer_,*uniqueContainer_,LinearObjContainer::X);
92
93 uniqueContainer_ = Teuchos::null;
94}
95#endif
96
97template < >
100{
101 const int n = value.size();
102 const int num_deriv = this->numDeriv();
103 TEUCHOS_ASSERT(n == 0 || n == num_deriv);
104 ScalarT glbValue = ScalarT(num_deriv, 0.0);
105
106 // do global summation -- it is possible to do the reduceAll() on the Fad's directly, but it is somewhat
107 // complicated for DFad (due to temporaries that might get created). Since this is just a sum, it is
108 // easier to do the reduction for each value and derivative component.
109 Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_SUM, Thyra::Ordinal(1), &value.val(), &glbValue.val());
110 if (num_deriv > 0)
111 Teuchos::reduceAll(*this->getComm(), Teuchos::REDUCE_SUM, Thyra::Ordinal(n), value.dx(), &glbValue.fastAccessDx(0));
112
113 value = glbValue;
114
115 // copy data in vectors
116#ifdef PANZER_HAVE_EPETRA_STACK
117 if(this->useEpetra()) {
118 // use epetra
119 Epetra_MultiVector& deriv = this->getEpetraMultiVector();
120 for (int i=0; i<num_deriv; ++i)
121 deriv[i][0] = glbValue.dx(i);
122 }
123 else
124#endif
125 {
126 // use thyra
127 TEUCHOS_ASSERT(this->useThyra());
128 Thyra::ArrayRCP< Thyra::ArrayRCP<double> > deriv = this->getThyraMultiVector();
129 for (int i=0; i<num_deriv; ++i)
130 deriv[i][0] = glbValue.dx(i);
131 }
132}
133
134// Do nothing unless derivatives are actually required
135template <typename EvalT>
137setSolnVectorSpace(const Teuchos::RCP<const Thyra::VectorSpaceBase<double> > & /* soln_vs */) { }
138
139// derivatives are required for
140template < >
142setSolnVectorSpace(const Teuchos::RCP<const Thyra::VectorSpaceBase<double> > & soln_vs)
143{
144 setDerivativeVectorSpace(soln_vs);
145}
146
147// derivatives are required for
148#ifdef Panzer_BUILD_HESSIAN_SUPPORT
149template < >
151setSolnVectorSpace(const Teuchos::RCP<const Thyra::VectorSpaceBase<double> > & soln_vs)
152{
153 setDerivativeVectorSpace(soln_vs);
154}
155#endif
156
157// Do nothing unless derivatives are required
158template <typename EvalT>
160adjustForDirichletConditions(const GlobalEvaluationData & /* localBCRows */, const GlobalEvaluationData & /* globalBCRows */) { }
161
162// Do nothing unless derivatives are required
163template < >
165adjustForDirichletConditions(const GlobalEvaluationData & localBCRows,const GlobalEvaluationData & globalBCRows)
166{
167 linObjFactory_->adjustForDirichletConditions(Teuchos::dyn_cast<const LinearObjContainer>(localBCRows),
168 Teuchos::dyn_cast<const LinearObjContainer>(globalBCRows),
169 *ghostedContainer_,true,true);
170}
171
172#ifdef Panzer_BUILD_HESSIAN_SUPPORT
173// Do nothing unless derivatives are required
174template < >
176adjustForDirichletConditions(const GlobalEvaluationData & localBCRows,const GlobalEvaluationData & globalBCRows)
177{
178 linObjFactory_->adjustForDirichletConditions(Teuchos::dyn_cast<const LinearObjContainer>(localBCRows),
179 Teuchos::dyn_cast<const LinearObjContainer>(globalBCRows),
180 *ghostedContainer_,true,true);
181}
182#endif
183
184}
185
186#endif
virtual void scatterResponse()
This simply does global summation, then shoves the result into a vector.
void adjustForDirichletConditions(const GlobalEvaluationData &localBCRows, const GlobalEvaluationData &globalBCRows)
void setSolnVectorSpace(const Teuchos::RCP< const Thyra::VectorSpaceBase< double > > &soln_vs)
Set solution vector space.