Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_DOF.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_DOF_DECL_HPP
12#define PANZER_EVALUATOR_DOF_DECL_HPP
13
14#include <string>
15#include "Phalanx_Evaluator_Macros.hpp"
16#include "Phalanx_MDField.hpp"
17#include "PanzerDiscFE_config.hpp"
19
20namespace panzer {
21
23template<typename EvalT, typename TRAITS>
24class DOF : public panzer::EvaluatorWithBaseImpl<TRAITS>,
25 public PHX::EvaluatorDerived<EvalT, TRAITS> {
26public:
27
28 DOF(const Teuchos::ParameterList& p);
29
37 DOF(const PHX::FieldTag & input,
38 const PHX::FieldTag & output,
39 const panzer::BasisDescriptor & bd,
41
42 void postRegistrationSetup(typename TRAITS::SetupData d,
44
45 void evaluateFields(typename TRAITS::EvalData d);
46
47private:
48
49 typedef typename EvalT::ScalarT ScalarT;
50
54
55 PHX::MDField<const ScalarT,Cell,Point> dof_basis;
56
57 PHX::MDField<ScalarT,Cell,Point> dof_ip_scalar;
58 PHX::MDField<ScalarT,Cell,Point,Dim> dof_ip_vector;
59
60 std::string basis_name;
61 std::size_t basis_index;
62
64};
65
69template<typename TRAITS>
70class DOF<typename TRAITS::Jacobian,TRAITS> :
71 public panzer::EvaluatorWithBaseImpl<TRAITS>,
72 public PHX::EvaluatorDerived<typename TRAITS::Jacobian, TRAITS> {
73public:
74
75 DOF(const Teuchos::ParameterList& p);
76
77 DOF(const PHX::FieldTag & input,
78 const PHX::FieldTag & output,
79 const panzer::BasisDescriptor & bd,
81
82 void postRegistrationSetup(typename TRAITS::SetupData d,
84
85 void preEvaluate(typename TRAITS::PreEvalData d);
86
87 void evaluateFields(typename TRAITS::EvalData d);
88
89private:
90
92
96
97 PHX::MDField<const ScalarT,Cell,Point> dof_basis;
98
99 PHX::MDField<ScalarT,Cell,Point> dof_ip_scalar;
100 PHX::MDField<ScalarT,Cell,Point,Dim> dof_ip_vector;
101
102 std::string basis_name;
103 std::size_t basis_index;
104
107 PHX::View<int*> offsets_array;
108 std::string sensitivities_name; // This sets which gather operations have sensitivities
109 // and thus which DOF operations can use accelerated jacobians
110
112};
113
114}
115
116#endif
DOF(const PHX::FieldTag &input, const PHX::FieldTag &output, const panzer::BasisDescriptor &bd, const panzer::IntegrationDescriptor &id)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
PHX::MDField< const ScalarT, Cell, Point > dof_basis
void evaluateFields(typename TRAITS::EvalData d)
PHX::MDField< ScalarT, Cell, Point > dof_ip_scalar
PHX::MDField< ScalarT, Cell, Point, Dim > dof_ip_vector
DOF(const Teuchos::ParameterList &p)
void preEvaluate(typename TRAITS::PreEvalData d)
panzer::Traits::Jacobian::ScalarT ScalarT
Interpolates basis DOF values to IP DOF values.
PHX::MDField< const ScalarT, Cell, Point > dof_basis
panzer::BasisDescriptor bd_
EvalT::ScalarT ScalarT
PHX::MDField< ScalarT, Cell, Point, Dim > dof_ip_vector
bool is_vector_basis
std::size_t basis_index
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
panzer::IntegrationDescriptor id_
void evaluateFields(typename TRAITS::EvalData d)
bool use_descriptors_
std::string basis_name
PHX::MDField< ScalarT, Cell, Point > dof_ip_scalar
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.