Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Integrator_GradBasisTimesScalar_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_EVALUATOR_GRADBASISTIMESSCALAR_DECL_HPP
12 #define PANZER_EVALUATOR_GRADBASISTIMESSCALAR_DECL_HPP
13 
15 //
16 // Include Files
17 //
19 
20 // C++
21 #include <string>
22 
23 // Kokkos
24 #include "Kokkos_DynRankView.hpp"
25 
26 // Panzer
29 
30 // Phalanx
31 #include "Phalanx_Evaluator_Derived.hpp"
32 #include "Phalanx_MDField.hpp"
33 
45 namespace panzer
46 {
47  template<typename EvalT, typename Traits>
49  :
50  public panzer::EvaluatorWithBaseImpl<Traits>,
51  public PHX::EvaluatorDerived<EvalT, Traits>
52  {
53  public:
54 
93  const std::vector<std::string>& resNames,
94  const std::string& scalar,
95  const panzer::BasisIRLayout& basis,
96  const panzer::IntegrationRule& ir,
97  const double& multiplier = 1,
98  const std::vector<std::string>& fmNames =
99  std::vector<std::string>());
100 
145  const Teuchos::ParameterList& p);
146 
157  void
159  typename Traits::SetupData d,
161 
171  void
173  typename Traits::EvalData d);
174 
179  template<int NUM_FIELD_MULT>
181  {
182  }; // end of struct FieldMultTag
183 
200  template<int NUM_FIELD_MULT>
201  KOKKOS_INLINE_FUNCTION
202  void
203  operator()(
204  const FieldMultTag<NUM_FIELD_MULT>& tag,
205  const std::size_t& cell) const;
206 
207  private:
208 
220  getValidParameters() const;
221 
225  using ScalarT = typename EvalT::ScalarT;
226 
236 
241  std::vector<PHX::MDField<ScalarT,Cell,BASIS>> fields_host_;
242  using InnerView = PHX::UnmanagedView<ScalarT**>;
243  using OuterView = PHX::View<InnerView*>;
245 
251 
257 
262  std::vector<PHX::MDField<const ScalarT, Cell, IP>> fieldMults_;
263 
269  PHX::View<PHX::UnmanagedView<const ScalarT**>*> kokkosFieldMults_;
270 
274  int numDims_;
275 
279  std::string basisName_;
280 
285  std::size_t basisIndex_;
286 
293 
294  }; // end of class Integrator_GradBasisTimesScalar
295 
296 } // end of namespace panzer
297 
298 #endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_DECL_HPP
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we&#39;re integrating ( ).
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > basis_
The gradient vector basis information necessary for integration.
EvaluatorStyle
An indication of how an Evaluator will behave.
std::string basisName_
The name of the basis we&#39;re using.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
double multiplier
The scalar multiplier out in front of the integral ( ).
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
ScalarT multiplier_
The scalar multiplier out in front of the integral ( ).
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The scalar multiplier out in front of the integral ( ).
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > fields_host_
The fields to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
int numDims_
The number of dimensions associated with the gradient.
Integrator_GradBasisTimesScalar(const panzer::EvaluatorStyle &evalStyle, const std::vector< std::string > &resNames, const std::string &scalar, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >())
Main Constructor.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.