Panzer Version of the Day
Loading...
Searching...
No Matches
Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
panzer::Integrator_GradBasisTimesScalar< EvalT, Traits > Class Template Reference

#include <Panzer_Integrator_GradBasisTimesScalar_decl.hpp>

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

Classes

struct  FieldMultTag
 This empty struct allows us to optimize operator()() depending on the number of field multipliers. More...
 

Public Member Functions

 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.
 
 Integrator_GradBasisTimesScalar (const Teuchos::ParameterList &p)
 ParameterList Constructor.
 
void postRegistrationSetup (typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
 Post-Registration Setup.
 
void evaluateFields (typename Traits::EvalData d)
 Evaluate Fields.
 
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void operator() (const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
 Perform the integration.
 
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void operator() (const FieldMultTag< NUM_FIELD_MULT > &, const size_t &cell) const
 
- Public Member Functions inherited from panzer::EvaluatorWithBaseImpl< Traits >
void setDetailsIndex (const int di)
 An evaluator builder sets the details index.
 
- Public Member Functions inherited from panzer::DomainEvaluator
 DomainEvaluator (DomainType domain=ALL)
 Constructor.
 
virtual ~DomainEvaluator ()=default
 Default destructor.
 
void setDomain (const DomainType domain)
 Set the domain for the evaluator.
 
DomainType getDomain ()
 Get the domain for the evaluator.
 
virtual int cellStartIndex (const panzer::Workset &workset) const
 Returns the starting cell for the specified domain for a given workset.
 
virtual int cellEndIndex (const panzer::Workset &workset) const
 Returns the non-inclusive end cell for the specified domain for a given workset.
 

Private Types

using ScalarT = typename EvalT::ScalarT
 The scalar type.
 
using InnerView = PHX::UnmanagedView< ScalarT ** >
 
using OuterView = PHX::View< InnerView * >
 

Private Member Functions

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

Private Attributes

const panzer::EvaluatorStyle evalStyle_
 An enum determining the behavior of this Evaluator.
 
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > fields_host_
 The fields to which we'll contribute, or in which we'll store, the result of computing this integral.
 
OuterView fields_
 
PHX::MDField< const ScalarT, Cell, IPscalar_
 A field representing the scalar function we're integrating ( $ s $).
 
ScalarT multiplier_
 The scalar multiplier out in front of the integral ( $ M
       $).
 
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
 The scalar multiplier out in front of the integral ( $ M
       $).
 
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.).
 
int numDims_
 The number of dimensions associated with the gradient.
 
std::string basisName_
 The name of the basis we're using.
 
std::size_t basisIndex_
 The index in the Workset bases for our particular BasisIRLayout name.
 
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dimbasis_
 The gradient vector basis information necessary for integration.
 

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_GradBasisTimesScalar< EvalT, Traits >

Definition at line 80 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

Member Typedef Documentation

◆ ScalarT

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

The scalar type.

Definition at line 257 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ InnerView

template<typename EvalT , typename Traits >
using panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::InnerView = PHX::UnmanagedView<ScalarT**>
private

◆ OuterView

template<typename EvalT , typename Traits >
using panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::OuterView = PHX::View<InnerView*>
private

Constructor & Destructor Documentation

◆ Integrator_GradBasisTimesScalar() [1/2]

template<typename EvalT , typename Traits >
panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::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.

Creates an Evaluator to evaluate the integral

\[
  Ma(x)b(x)\cdots\int s(x)\nabla\phi(x)\,dx,
\]

where $ M $ is some constant, $ a(x) $, $ b(x) $, etc., are some fields that depend on position, $ s $ is some scalar function, and $ \phi $ is some vector 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]resNamesThe names of either the contributed or evaluated fields, depending on evalStyle.
[in]scalarThe name of the scalar function being integrated ( $ 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.
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 doesn't match the size of resNames.

Definition at line 71 of file Panzer_Integrator_GradBasisTimesScalar_impl.hpp.

◆ Integrator_GradBasisTimesScalar() [2/2]

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

ParameterList Constructor.

Creates an Evaluator to evaluate the integral

\[
  Ma(x)b(x)\cdots\int s(x)\nabla\phi(x)\,dx,
\]

where $ M $ is some constant, $ a(x) $, $ b(x) $, etc., are some fields that depend on position, $ s $ is some scalar function, and $ \phi $ is some vector 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 Names" type = "std::vector<std::string>" value = (required) />
<Parameter name = "Scalar 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)/>
</ParameterList>
where
  • "Residual Names" are the names for the terms this Evaluator is evaluating,
  • "Scalar Name" is the name of the scalar function being integrated ( $ 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.).

Definition at line 167 of file Panzer_Integrator_GradBasisTimesScalar_impl.hpp.

Member Function Documentation

◆ postRegistrationSetup()

template<typename EvalT , typename Traits >
void panzer::Integrator_GradBasisTimesScalar< 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 multipliers, 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 198 of file Panzer_Integrator_GradBasisTimesScalar_impl.hpp.

◆ evaluateFields()

template<typename EvalT , typename Traits >
void panzer::Integrator_GradBasisTimesScalar< 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 296 of file Panzer_Integrator_GradBasisTimesScalar_impl.hpp.

◆ operator()() [1/2]

template<typename EvalT , typename Traits >
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::operator() ( const FieldMultTag< NUM_FIELD_MULT > &  tag,
const std::size_t &  cell 
) const

Perform the integration.

Generally speaking, for a given cell in the Workset, this routine loops over quadrature points, bases, and dimensions 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.

◆ getValidParameters()

template<typename EvalT , typename TRAITS >
Teuchos::RCP< Teuchos::ParameterList > panzer::Integrator_GradBasisTimesScalar< EvalT, TRAITS >::getValidParameters
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 324 of file Panzer_Integrator_GradBasisTimesScalar_impl.hpp.

◆ operator()() [2/2]

template<typename EvalT , typename Traits >
template<int NUM_FIELD_MULT>
KOKKOS_INLINE_FUNCTION void panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::operator() ( const FieldMultTag< NUM_FIELD_MULT > &  ,
const size_t &  cell 
) const

Member Data Documentation

◆ evalStyle_

template<typename EvalT , typename Traits >
const panzer::EvaluatorStyle panzer::Integrator_GradBasisTimesScalar< 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 267 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ fields_host_

template<typename EvalT , typename Traits >
std::vector<PHX::MDField<ScalarT,Cell,BASIS> > panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::fields_host_
private

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

Definition at line 273 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ fields_

template<typename EvalT , typename Traits >
OuterView panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::fields_
private

◆ scalar_

template<typename EvalT , typename Traits >
PHX::MDField<const ScalarT, Cell, IP> panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::scalar_
private

A field representing the scalar function we're integrating ( $ s $).

Definition at line 282 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ multiplier_

template<typename EvalT , typename Traits >
ScalarT panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::multiplier_
private

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

Definition at line 288 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ fieldMults_

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

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

Definition at line 294 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ kokkosFieldMults_

template<typename EvalT , typename Traits >
PHX::View<PHX::UnmanagedView<const ScalarT**>*> panzer::Integrator_GradBasisTimesScalar< 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 301 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ numDims_

template<typename EvalT , typename Traits >
int panzer::Integrator_GradBasisTimesScalar< EvalT, Traits >::numDims_
private

The number of dimensions associated with the gradient.

Definition at line 306 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ basisName_

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

The name of the basis we're using.

Definition at line 311 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ basisIndex_

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

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

Definition at line 317 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.

◆ basis_

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

The gradient vector basis information necessary for integration.

Definition at line 324 of file Panzer_Integrator_GradBasisTimesScalar_decl.hpp.


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