Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Integrator_BasisTimesTensorTimesVector.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_Integerator_BasisTimesTensorTimesVector_decl_hpp__
12#define __Panzer_Integerator_BasisTimesTensorTimesVector_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
34namespace panzer
35{
48 template<typename EvalT, typename Traits>
50 :
51 public panzer::EvaluatorWithBaseImpl<Traits>,
52 public PHX::EvaluatorDerived<EvalT, Traits>
53 {
54 public:
55
87 const std::string& resName,
88 const std::string& valName,
89 const panzer::BasisIRLayout& basis,
91 const std::string& tensorName );
92
131 const Teuchos::ParameterList& p);
132
165 const PHX::FieldTag& resTag,
166 const PHX::FieldTag& valTag,
167 const BasisDescriptor& bd,
168 const IntegrationDescriptor& id,
169 const PHX::FieldTag& tensorTag);
170
182 void
184 typename Traits::SetupData sd,
186
196 void
198 typename Traits::EvalData workset);
199
212 KOKKOS_INLINE_FUNCTION
213 void operator()(const std::size_t& cell) const;
214
215 private:
216
227 Teuchos::RCP<Teuchos::ParameterList>
228 getValidParameters() const;
229
233 using ScalarT = typename EvalT::ScalarT;
234
244
250
255
260
265 PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS> field_;
266
271 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
273
277 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP, panzer::Dim, panzer::Dim>
279
284 PHX::View<const ScalarT****> kokkosTensor_;
285
286
291
296
300 std::string basisName_;
301
306 std::size_t basisIndex_;
307
311 PHX::MDField<const double, panzer::Cell, panzer::BASIS, panzer::IP,
313
315 PHX::View<ScalarT*> tmp_;
316
317 }; // end of class Integrator_BasisTimesTensorTimesVector
318
319} // end of namespace panzer
320
321#endif // __Panzer_Integerator_BasisTimesTensorTimesVector_decl_hpp__
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
IntegrationDescriptor id_
The IntegrationDescriptor for the quadrature to use.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
PHX::View< ScalarT * > tmp_
Scratch space for caching temporary values in the kokkos kernel.
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
KOKKOS_INLINE_FUNCTION void operator()(const std::size_t &cell) const
Perform the integration.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field_
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
PHX::View< const ScalarT **** > kokkosTensor_
The PHX::View representation of the tensor fields that are multiplied out in front of the integral ( ...
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
bool useDescriptors_
A flag indicating whether or not to use the descriptor interface.
PHX::MDField< const double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > basis_
The vector basis information necessary for integration.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we're integrating ( ).
BasisDescriptor bd_
The BasisDescriptor for the basis to use.
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.