43 #ifndef PANZER_RESPONSE_SCATTER_EVALUATOR_FUNCTIONAL_HPP
44 #define PANZER_RESPONSE_SCATTER_EVALUATOR_FUNCTIONAL_HPP
49 #include "PanzerDiscFE_config.hpp"
56 #include "Phalanx_Evaluator_Macros.hpp"
57 #include "Phalanx_MDField.hpp"
67 virtual void scatterDerivative(
const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
72 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
73 virtual void scatterHessian(
const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
80 template <
typename LO,
typename GO>
85 if(globalIndexer!=Teuchos::null)
86 ugis_.push_back(globalIndexer);
92 void scatterDerivative(
const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
97 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
98 void scatterHessian(
const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
106 std::vector<Teuchos::RCP<const panzer::GlobalIndexer> >
ugis_;
112 template<
typename EvalT,
typename Traits>
114 public PHX::EvaluatorDerived<EvalT, Traits> {
138 template <
typename LO,
typename GO>
144 Kokkos::View<const LO*, PHX::Device> LIDs;
147 std::string blockId = wda(workset).block_id;
149 std::vector<int> blockOffsets;
155 const std::vector<std::size_t> & localCellIds = wda(workset).cell_local_ids;
156 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
157 std::size_t cellLocalId = localCellIds[worksetCellIndex];
159 for(std::size_t b=0;b<ugis_.size();b++) {
160 int start = blockOffsets[b];
162 LIDs = ugis_[b]->getElementLIDs(cellLocalId);
167 for(std::size_t i=0;i<LIDs.size();i++) {
168 dgdx_b[LIDs[i]] += cellIntegral(worksetCellIndex).dx(start+i);
174 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
175 template <
typename LO,
typename GO>
181 Kokkos::View<const LO*, PHX::Device> LIDs;
184 std::string blockId = wda(workset).block_id;
186 std::vector<int> blockOffsets;
192 const std::vector<std::size_t> & localCellIds = wda(workset).cell_local_ids;
193 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
194 std::size_t cellLocalId = localCellIds[worksetCellIndex];
196 for(std::size_t b=0;b<ugis_.size();b++) {
197 int start = blockOffsets[b];
199 LIDs = ugis_[b]->getElementLIDs(cellLocalId);
204 for(std::size_t i=0;i<LIDs.size();i++) {
205 d2gdx2_b[LIDs[i]] += cellIntegral(worksetCellIndex).dx(start+i).dx(0);
FunctionalScatter(const std::vector< Teuchos::RCP< const panzer::GlobalIndexer > > &ugis)
FunctionalScatter(const Teuchos::RCP< const panzer::GlobalIndexer > &globalIndexer)
virtual ~FunctionalScatterBase()
ResponseScatterEvaluator_Functional(const std::string &name, const CellData &cd, const Teuchos::RCP< FunctionalScatterBase > &functionalScatter)
A constructor with concrete arguments instead of a parameter list.
void scatterDerivative(const PHX::MDField< const panzer::Traits::Jacobian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &dgdx) const
Teuchos::RCP< FunctionalScatterBase > scatterObj_
virtual void scatterDerivative(const PHX::MDField< const panzer::Traits::Jacobian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &dgdx) const =0
void preEvaluate(typename Traits::PreEvalData d)
void computeBlockOffsets(const std::string &blockId, const std::vector< Teuchos::RCP< GlobalIndexer >> &ugis, std::vector< int > &blockOffsets)
std::string responseName_
PHX::MDField< const ScalarT, panzer::Cell > cellIntegral_
Data for determining cell topology and dimensionality.
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
void evaluateFields(typename Traits::EvalData d)
Teuchos::RCP< PHX::FieldTag > scatterHolder_
void scatterHessian(const PHX::MDField< const panzer::Traits::Hessian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &d2gdx2) const
std::vector< Teuchos::RCP< const panzer::GlobalIndexer > > ugis_
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< Response_Functional< EvalT > > responseObj_
virtual void scatterHessian(const PHX::MDField< const panzer::Traits::Hessian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &d2gdx2) const =0