Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Integrator_GradBasisCrossVector_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_GRADBASISCROSSVECTOR_DECL_HPP
12 #define PANZER_EVALUATOR_GRADBASISCROSSVECTOR_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 
34 namespace panzer
35 {
48  template<typename EvalT, typename Traits>
50  :
51  public panzer::EvaluatorWithBaseImpl<Traits>,
52  public PHX::EvaluatorDerived<EvalT, Traits>
53  {
54  public:
55 
103  const std::vector<std::string>& resNames,
104  const std::string& vecName,
105  const panzer::BasisIRLayout& basis,
106  const panzer::IntegrationRule& ir,
107  const double& multiplier = 1,
108  const std::vector<std::string>& fmNames =
109  std::vector<std::string>(),
110  const Teuchos::RCP<PHX::DataLayout>& vecDL = Teuchos::null);
111 
159  const Teuchos::ParameterList& p);
160 
171  void
173  typename Traits::SetupData d,
175 
185  void
187  typename Traits::EvalData d);
188 
193  template<int NUM_FIELD_MULT>
195  {
196  }; // end of struct FieldMultTag
197 
213  template<int NUM_FIELD_MULT>
214  KOKKOS_INLINE_FUNCTION
215  void
216  operator()(
217  const FieldMultTag<NUM_FIELD_MULT>& tag,
218  const std::size_t& cell) const;
219 
220  private:
221 
233  getValidParameters() const;
234 
238  using ScalarT = typename EvalT::ScalarT;
239 
249 
254  std::vector<PHX::MDField<ScalarT,Cell,BASIS>> fields_host_;
255  using InnerView = PHX::UnmanagedView<ScalarT**>;
256  using OuterView = PHX::View<InnerView*>;
258 
264 
269  double multiplier_;
270 
275  std::vector<PHX::MDField<const ScalarT, Cell, IP>> fieldMults_;
276 
282  PHX::View<PHX::UnmanagedView<const ScalarT**>*> kokkosFieldMults_;
283 
287  int numDims_;
288 
293 
297  std::string basisName_;
298 
303  std::size_t basisIndex_;
304 
310 
311  }; // end of class Integrator_GradBasisCrossVector
312 
313 } // end of namespace panzer
314 
315 #endif // PANZER_EVALUATOR_GRADBASISCROSSVECTOR_DECL_HPP
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ...
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
Integrator_GradBasisCrossVector(const panzer::EvaluatorStyle &evalStyle, const std::vector< std::string > &resNames, const std::string &vecName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >(), const Teuchos::RCP< PHX::DataLayout > &vecDL=Teuchos::null)
Main Constructor.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
double multiplier_
The scalar multiplier out in front of the integral ( ).
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
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...
PHX::MDField< const ScalarT, Cell, IP, Dim > vector_
A field representing the vector-valued function we&#39;re integrating ( ).
double multiplier
The scalar multiplier out in front of the integral ( ).
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > basis_
The gradient vector basis information necessary for integration.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
int numDims_
The number of dimensions associated with the vector.
int numGradDims_
The number of dimensions associated with the gradient.
std::string basisName_
The name of the basis we&#39;re using.