Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
panzer::Integrator_GradBasisDotVector< EvalT, Traits > Class Template Reference

Computes $ Ma(x)b(x)\cdots\int\vec{s}(x)\cdot\nabla\phi(x)\,dx $. More...

#include <Panzer_Integrator_GradBasisDotVector_decl.hpp>

Inheritance diagram for panzer::Integrator_GradBasisDotVector< EvalT, Traits >:
Inheritance graph
[legend]

Classes

struct  FieldMultTag
 This empty struct allows us to optimize operator()() depending on the number of field multipliers. This is the version that does not use shared memory. More...
 
struct  SharedFieldMultTag
 This empty struct allows us to optimize operator()() depending on the number of field multipliers. This is the shared memory version. More...
 

Public Member Functions

 Integrator_GradBasisDotVector (const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &fluxName, 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. More...
 
 Integrator_GradBasisDotVector (const Teuchos::ParameterList &p)
 ParameterList Constructor. More...
 
void postRegistrationSetup (typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
 Post-Registration Setup. More...
 
void evaluateFields (typename Traits::EvalData d)
 Evaluate Fields. More...
 
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void operator() (const FieldMultTag< NUM_FIELD_MULT > &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
 Perform the integration. More...
 
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void operator() (const SharedFieldMultTag< NUM_FIELD_MULT > &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
 Perform the integration. More...
 
- Public Member Functions inherited from panzer::EvaluatorWithBaseImpl< Traits >
void setDetailsIndex (const int di)
 An evaluator builder sets the details index. More...
 
- Public Member Functions inherited from PHX::EvaluatorWithBaseImpl< Traits >
virtual void evaluateFields (typename Traits::EvalData d) override=0
 
- Public Member Functions inherited from PHX::Evaluator< Traits >
 Evaluator ()
 
virtual ~Evaluator ()
 
virtual void postRegistrationSetup (typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)=0
 
virtual const std::vector
< Teuchos::RCP< FieldTag > > & 
evaluatedFields () const =0
 
virtual const std::vector
< Teuchos::RCP< FieldTag > > & 
contributedFields () const =0
 
virtual const std::vector
< Teuchos::RCP< FieldTag > > & 
dependentFields () const =0
 
virtual const std::vector
< Teuchos::RCP< FieldTag > > & 
unsharedFields () const =0
 
virtual void preEvaluate (typename Traits::PreEvalData d)=0
 
virtual void postEvaluate (typename Traits::PostEvalData d)=0
 
virtual const std::string & getName () const =0
 
virtual void bindField (const PHX::FieldTag &ft, const PHX::any &f)=0
 
virtual PHX::DeviceEvaluator
< Traits > * 
createDeviceEvaluator () const =0
 
virtual void rebuildDeviceEvaluator (PHX::DeviceEvaluator< Traits > *e) const =0
 
virtual void deleteDeviceEvaluator (PHX::DeviceEvaluator< Traits > *e) const =0
 
virtual void printFieldValues (std::ostream &os) const =0
 
- Public Member Functions inherited from panzer::DomainEvaluator
 DomainEvaluator (DomainType domain=ALL)
 Constructor. More...
 
virtual ~DomainEvaluator ()=default
 Default destructor. More...
 
void setDomain (const DomainType domain)
 Set the domain for the evaluator. More...
 
DomainType getDomain ()
 Get the domain for the evaluator. More...
 
virtual int cellStartIndex (const panzer::Workset &workset) const
 Returns the starting cell for the specified domain for a given workset. More...
 
virtual int cellEndIndex (const panzer::Workset &workset) const
 Returns the non-inclusive end cell for the specified domain for a given workset. More...
 

Private Types

using ScalarT = typename EvalT::ScalarT
 The scalar type. More...
 
using scratch_view = Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged >
 Type for shared memory. More...
 

Private Member Functions

Teuchos::RCP
< Teuchos::ParameterList
getValidParameters () const
 Get Valid Parameters. More...
 

Private Attributes

const panzer::EvaluatorStyle evalStyle_
 An enum determining the behavior of this Evaluator. More...
 
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. More...
 
PHX::MDField< const ScalarT,
panzer::Cell, panzer::IP,
panzer::Dim
vector_
 A field representing the vector-valued function we're integrating ( $ \vec{s} $). More...
 
double multiplier_
 The scalar multiplier out in front of the integral ( $ M $). More...
 
std::vector< PHX::MDField
< const ScalarT, panzer::Cell,
panzer::IP > > 
fieldMults_
 The (possibly empty) list of fields that are multipliers out in front of the integral ( $ a(x) $, $ b(x) $, etc.). More...
 
PHX::View< PHX::UnmanagedView
< const ScalarT ** > * > 
kokkosFieldMults_
 The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front of the integral ( $ a(x) $, $ b(x) $, etc.). More...
 
std::string basisName_
 The name of the basis we're using. More...
 
std::size_t basisIndex_
 The index in the Workset bases for our particular BasisIRLayout name. More...
 
PHX::MDField< double,
panzer::Cell, panzer::BASIS,
panzer::IP, panzer::Dim
basis_
 The gradient vector basis information necessary for integration. More...
 
PHX::View< ScalarT * > tmp_
 Temporary used when shared memory is disabled. More...
 

Additional Inherited Members

- Public Types inherited from panzer::DomainEvaluator
enum  DomainType : int {
  OWNED =0, GHOST =1, REAL =2, VIRTUAL =3,
  EXTERNAL =4, ALL =5
}
 Domain types supported by worksets. More...
 
- Protected Attributes inherited from panzer::EvaluatorWithBaseImpl< Traits >
WorksetDetailsAccessor wda
 

Detailed Description

template<typename EvalT, typename Traits>
class panzer::Integrator_GradBasisDotVector< EvalT, Traits >

Computes $ Ma(x)b(x)\cdots\int\vec{s}(x)\cdot\nabla\phi(x)\,dx $.

Evaluates the integral

\[ Ma(x)b(x)\cdots\int\vec{s}(x)\cdot\nabla\phi(x)\,dx, \]

where $ M $ is some constant, $ a(x) $, $ b(x) $, etc., are some fields that depend on position, $ \vec{s} $ is some vector-valued function, and $ \phi $ is some basis.

Definition at line 78 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

Member Typedef Documentation

template<typename EvalT , typename Traits >
using panzer::Integrator_GradBasisDotVector< EvalT, Traits >::ScalarT = typename EvalT::ScalarT
private

The scalar type.

Definition at line 294 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
using panzer::Integrator_GradBasisDotVector< EvalT, Traits >::scratch_view = Kokkos::View<ScalarT* ,typename PHX::DevLayout<ScalarT>::type,typename PHX::exec_space::scratch_memory_space,Kokkos::MemoryUnmanaged>
private

Type for shared memory.

Definition at line 297 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

Constructor & Destructor Documentation

template<typename EvalT , typename Traits >
panzer::Integrator_GradBasisDotVector< EvalT, Traits >::Integrator_GradBasisDotVector ( const panzer::EvaluatorStyle evalStyle,
const std::string &  resName,
const std::string &  fluxName,
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.

Creates an Evaluator to evaluate the integral

\[ Ma(x)b(x)\cdots\int\vec{s}(x)\cdot\nabla\phi(x)\,dx, \]

where $ M $ is some constant, $ a(x) $, $ b(x) $, etc., are some fields that depend on position, $ \vec{s} $ is some vector-valued function, and $ \phi $ is some basis.

Parameters
[in]evalStyleAn enum declaring the behavior of this Evaluator, which is to either:
  • compute and contribute (CONTRIBUTES), or
  • compute and store (EVALUATES).
[in]resNameThe name of either the contributed or evaluated field, depending on evalStyle.
[in]fluxNameThe name of the vector-valued function being integrated ( $ \vec{s} $).
[in]basisThe basis that you'd like to use ( $ \phi $).
[in]irThe integration rule that you'd like to use.
[in]multiplierThe scalar multiplier out in front of the integral you're computing ( $ M $). If not specified, this defaults to 1.
[in]fmNamesA list of names of fields that are multipliers out in front of the integral you're computing ( $ a(x) $, $ b(x) $, etc.). If not specified, this defaults to an empty vector.
[in]vecDLThe vector data layout that you'd like to use. If not specified, this defaults to Teuchos::null and the vector data layout from the given ir is used.
Exceptions
std::invalid_argumentIf any of the inputs are invalid.
std::logic_errorIf the basis supplied does not support the gradient operator, or if the dimension of the space exceeds the dimension of the vector data layout.

Definition at line 67 of file Panzer_Integrator_GradBasisDotVector_impl.hpp.

template<typename EvalT , typename Traits >
panzer::Integrator_GradBasisDotVector< EvalT, Traits >::Integrator_GradBasisDotVector ( const Teuchos::ParameterList p)

ParameterList Constructor.

Creates an Evaluator to evaluate the integral

\[ Ma(x)b(x)\cdots\int\vec{s}(x)\cdot\nabla\phi(x)\,dx, \]

where $ M $ is some constant, $ a(x) $, $ b(x) $, etc., are some fields that depend on position, $ \vec{s} $ is some vector-valued function, and $ \phi $ is some basis.

Note
This constructor exists to preserve the older way of creating an Evaluator with a ParameterList; however, it is strongly advised that you not use this ParameterList Constructor, but rather that you favor the main constructor with its compile-time argument checking instead.
Parameters
[in]pA ParameterList of the form
<ParameterList>
<Parameter name = "Residual Name" type = "std::string" value = (required) />
<Parameter name = "Flux Name" type = "std::string" value = (required) />
<Parameter name = "Basis" type = "RCP<panzer::BasisIRLayout>" value = (required) />
<Parameter name = "IR" type = "RCP<panzer::IntegrationRule>" value = (required) />
<Parameter name = "Multiplier" type = "double" value = (required) />
<Parameter name = "Field Multipliers" type = "RCP<const std::vector<std::string>>" value = null (default)/>
<Parameter name = "Vector Data Layout" type = "RCP<PHX::DataLayout>" value = null (default)/>
</ParameterList>
where
  • "Residual Name" is the name for the term this Evaluator is evaluating,
  • "Flux Name" is the name of the vector-valued function being integrated ( $ \vec{s} $),
  • "Basis" is the basis that you'd like to use ( $ \phi $),
  • "IR" is the integration rule that you'd like to use,
  • "Multiplier" is the scalar multiplier out in front of the integral you're computing ( $ M $), and
  • "Field Multipliers" is an optional list of names of fields that are multipliers out in front of the integral you're computing ( $ a(x) $, $ b(x) $, etc.).
  • "Vector Data Layout" is the vector data layout that you'd like to use.

Definition at line 153 of file Panzer_Integrator_GradBasisDotVector_impl.hpp.

Member Function Documentation

template<typename EvalT , typename Traits >
void panzer::Integrator_GradBasisDotVector< EvalT, Traits >::postRegistrationSetup ( typename Traits::SetupData  d,
PHX::FieldManager< Traits > &  fm 
)

Post-Registration Setup.

Sets the basis index, and gets the PHX::View versions of the field multiplier, if there are any.

Parameters
[in]sdEssentially a list of Worksets, which are collections of cells (elements) that all live on a single process.
[in]fmThis is an unused part of the Evaluator interface.

Definition at line 187 of file Panzer_Integrator_GradBasisDotVector_impl.hpp.

template<typename EvalT , typename Traits >
void panzer::Integrator_GradBasisDotVector< EvalT, Traits >::evaluateFields ( typename Traits::EvalData  d)

Evaluate Fields.

This actually performs the integration by calling operator()() in a Kokkos::parallel_for over the cells in the Workset.

Parameters
[in]worksetThe Workset on which you're going to do the integration.

Definition at line 411 of file Panzer_Integrator_GradBasisDotVector_impl.hpp.

template<typename EvalT , typename Traits >
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void panzer::Integrator_GradBasisDotVector< EvalT, Traits >::operator() ( const FieldMultTag< NUM_FIELD_MULT > &  tag,
const Kokkos::TeamPolicy< PHX::exec_space >::member_type &  team 
) const

Perform the integration.

Generally speaking, for a given cell in the Workset, this routine loops over quadrature points, vector dimensions, and bases to perform the integration, scaling the vector field to be integrated by the multiplier ( $ M $) and any field multipliers ( $ a(x) $, $ b(x) $, etc.).

Note
Optimizations are made for the cases in which we have no field multipliers or only a single one.
Parameters
[in]tagAn indication of the number of field multipliers we have; either 0, 1, or something else.
[in]cellThe cell in the Workset over which to integrate.

Definition at line 225 of file Panzer_Integrator_GradBasisDotVector_impl.hpp.

template<typename EvalT , typename Traits >
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void panzer::Integrator_GradBasisDotVector< EvalT, Traits >::operator() ( const SharedFieldMultTag< NUM_FIELD_MULT > &  tag,
const Kokkos::TeamPolicy< PHX::exec_space >::member_type &  team 
) const

Perform the integration.

Generally speaking, for a given cell in the Workset, this routine loops over quadrature points, vector dimensions, and bases to perform the integration, scaling the vector field to be integrated by the multiplier ( $ M $) and any field multipliers ( $ a(x) $, $ b(x) $, etc.). Uses Shared memory.

Note
Optimizations are made for the cases in which we have no field multipliers or only a single one.
Parameters
[in]tagAn indication of the number of field multipliers we have; either 0, 1, or something else.
[in]cellThe cell in the Workset over which to integrate.

Definition at line 306 of file Panzer_Integrator_GradBasisDotVector_impl.hpp.

template<typename EvalT , typename TRAITS >
Teuchos::RCP< Teuchos::ParameterList > panzer::Integrator_GradBasisDotVector< EvalT, TRAITS >::getValidParameters ( ) const
private

Get Valid Parameters.

Get all the parameters that we support such that the ParameterList Constructor can do some validation of the input ParameterList.

Returns
A ParameterList with all the valid parameters (keys) in it. The values tied to those keys are meaningless default values.

Definition at line 469 of file Panzer_Integrator_GradBasisDotVector_impl.hpp.

Member Data Documentation

template<typename EvalT , typename Traits >
const panzer::EvaluatorStyle panzer::Integrator_GradBasisDotVector< EvalT, Traits >::evalStyle_
private

An enum determining the behavior of this Evaluator.

This Evaluator will compute the result of its integration and then:

  • CONTRIBUTES: contribute it to a specified residual, not saving anything; or
  • EVALUATES: save it under a specified name for future use.

Definition at line 308 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS> panzer::Integrator_GradBasisDotVector< EvalT, Traits >::field_
private

A field to which we'll contribute, or in which we'll store, the result of computing this integral.

Definition at line 314 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
PHX::MDField<const ScalarT, panzer::Cell, panzer::IP, panzer::Dim> panzer::Integrator_GradBasisDotVector< EvalT, Traits >::vector_
private

A field representing the vector-valued function we're integrating ( $ \vec{s} $).

Definition at line 321 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
double panzer::Integrator_GradBasisDotVector< EvalT, Traits >::multiplier_
private

The scalar multiplier out in front of the integral ( $ M $).

Definition at line 327 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
std::vector<PHX::MDField<const ScalarT, panzer::Cell, panzer::IP> > panzer::Integrator_GradBasisDotVector< EvalT, Traits >::fieldMults_
private

The (possibly empty) list of fields that are multipliers out in front of the integral ( $ a(x) $, $ b(x) $, etc.).

Definition at line 334 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
PHX::View<PHX::UnmanagedView<const ScalarT**>* > panzer::Integrator_GradBasisDotVector< EvalT, Traits >::kokkosFieldMults_
private

The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front of the integral ( $ a(x) $, $ b(x) $, etc.).

Definition at line 341 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
std::string panzer::Integrator_GradBasisDotVector< EvalT, Traits >::basisName_
private

The name of the basis we're using.

Definition at line 346 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
std::size_t panzer::Integrator_GradBasisDotVector< EvalT, Traits >::basisIndex_
private

The index in the Workset bases for our particular BasisIRLayout name.

Definition at line 352 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
PHX::MDField<double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim> panzer::Integrator_GradBasisDotVector< EvalT, Traits >::basis_
private

The gradient vector basis information necessary for integration.

Definition at line 359 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.

template<typename EvalT , typename Traits >
PHX::View<ScalarT*> panzer::Integrator_GradBasisDotVector< EvalT, Traits >::tmp_
private

Temporary used when shared memory is disabled.

Definition at line 362 of file Panzer_Integrator_GradBasisDotVector_decl.hpp.


The documentation for this class was generated from the following files: