Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Integrator_GradBasisDotTensorTimesVector_decl.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_Integrator_GradBasisDotTensorTimesVector_decl_hpp__
12 #define __Panzer_Integrator_GradBasisDotTensorTimesVector_decl_hpp__
13 
15 //
16 // Include Files
17 //
19 
20 // C++
21 #include <string>
22 
23 // Panzer
26 
27 // Phalanx
28 #include "Phalanx_Evaluator_Derived.hpp"
29 #include "Phalanx_MDField.hpp"
30 
31 namespace panzer
32 {
45  template<typename EvalT, typename Traits>
47  :
48  public panzer::EvaluatorWithBaseImpl<Traits>,
49  public PHX::EvaluatorDerived<EvalT, Traits>
50  {
51  public:
52 
89  const std::string& resName,
90  const std::string& fluxName,
91  const panzer::BasisIRLayout& basis,
92  const panzer::IntegrationRule& ir,
93  const std::string& tensorName,
94  const Teuchos::RCP<PHX::DataLayout>& vecDL = Teuchos::null);
95 
137  const Teuchos::ParameterList& p);
138 
149  void
151  typename Traits::SetupData d,
153 
163  void
165  typename Traits::EvalData d);
166 
173  struct FieldMultTag {};
174 
181 
198  // template<int NUM_FIELD_MULT>
199  KOKKOS_INLINE_FUNCTION
200  void
201  operator()(
202  const FieldMultTag& tag,
203  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const;
204 
221  // template<int NUM_FIELD_MULT>
222  KOKKOS_INLINE_FUNCTION
223  void
224  operator()(
225  const SharedFieldMultTag& tag,
226  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const;
227 
228  private:
229 
241  getValidParameters() const;
242 
246  using ScalarT = typename EvalT::ScalarT;
247 
249  using scratch_view = Kokkos::View<ScalarT* ,typename PHX::DevLayout<ScalarT>::type,typename PHX::exec_space::scratch_memory_space,Kokkos::MemoryUnmanaged>;
250 
261 
267 
274 
280 
285  PHX::View<const ScalarT****> kokkosTensor_;
286 
290  std::string basisName_;
291 
296  std::size_t basisIndex_;
297 
304 
306  PHX::View<ScalarT*> tmp_;
307 
308  }; // end of class Integrator_GradBasisDotTensorTimesVector
309 
310 } // end of namespace panzer
311 
312 #endif // __Panzer_Integrator_GradBasisDotTensorTimesVector_decl_hpp__
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
Perform the integration.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim, panzer::Dim > tensor_
The tensor field( ).
EvaluatorStyle
An indication of how an Evaluator will behave.
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > basis_
The gradient vector basis information necessary for integration.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field_
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
PHX::View< ScalarT * > tmp_
Temporary used when shared memory is disabled.
Integrator_GradBasisDotTensorTimesVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &fluxName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const std::string &tensorName, const Teuchos::RCP< PHX::DataLayout > &vecDL=Teuchos::null)
Main Constructor.
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we&#39;re integrating ( ).
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
PHX::View< const ScalarT **** > kokkosTensor_
The PHX::View representation of the tensor fields that are multiplied out in front of the integral ( ...