Advanced Multi-Physics (AMP)
On-Line Documentation
SourceNonlinearElement.h
Go to the documentation of this file.
1#ifndef included_AMP_SourceNonlinearElement
2#define included_AMP_SourceNonlinearElement
3
4#include <vector>
5
6#include <memory>
7
8// AMP files
9#include "AMP/operators/ElementOperation.h"
10#include "AMP/operators/libmesh/SourcePhysicsModel.h"
11
12// Libmesh headers
14#include "libmesh/libmesh_config.h"
15#undef LIBMESH_ENABLE_REFERENCE_COUNTING
16#include "libmesh/elem.h"
17#include "libmesh/fe_base.h"
18#include "libmesh/fe_type.h"
19#include "libmesh/quadrature_gauss.h"
21
22
23namespace AMP::Operator {
24
30{
31public:
44 explicit SourceNonlinearElement( std::shared_ptr<const ElementOperationParameters> params );
45
48
57 void initializeForCurrentElement( const libMesh::Elem *elem,
58 std::shared_ptr<SourcePhysicsModel> sourceTransportModel );
59
60 void setElementInputVector( const std::vector<std::vector<double>> &elementInputVector )
61 {
62 d_elementInputVector = elementInputVector;
63 }
64
72 void setElementVectors( const std::vector<std::vector<double>> &elementInputVector,
73 const std::vector<std::vector<double>> &elementAuxVector,
74 std::vector<double> &elementOutputVector )
75 {
76 d_elementInputVector = elementInputVector;
77 d_elementAuxVector = elementAuxVector;
78 d_elementOutputVector = &( elementOutputVector );
79 }
80
85 void setElementFlags( const std::string &inputVariableType )
86 {
87 d_isInputType = inputVariableType;
88 }
89
93 void apply() override;
94
95 std::shared_ptr<libMesh::FEBase> getFEBase() { return d_fe; }
96
97 unsigned int getNumberOfGaussPoints() { return ( d_qrule->n_points() ); }
98
99
100protected:
101 std::vector<std::vector<double>> d_elementInputVector;
102
103 std::vector<std::vector<double>> d_elementAuxVector;
104
105 std::vector<double> *d_elementOutputVector;
106
107 std::vector<std::vector<double>> d_elementOtherVectors;
108
109 std::shared_ptr<libMesh::FEType> d_feType;
110
111 std::shared_ptr<libMesh::FEBase> d_fe;
112
113 std::shared_ptr<libMesh::QBase> d_qrule;
114
115 const std::vector<libMesh::Real> *d_JxW;
116
117 const std::vector<std::vector<libMesh::Real>> *d_phi;
118
119 const std::vector<std::vector<libMesh::RealGradient>> *d_dphi;
120
121 std::string d_isInputType;
122
123 const libMesh::Elem *d_elem;
124
126
127 std::shared_ptr<SourcePhysicsModel> d_sourcePhysicsModel;
128
129private:
130};
131} // namespace AMP::Operator
132
133#endif
std::shared_ptr< libMesh::FEBase > d_fe
std::shared_ptr< libMesh::FEBase > getFEBase()
std::vector< std::vector< double > > d_elementInputVector
std::shared_ptr< libMesh::QBase > d_qrule
const std::vector< libMesh::Real > * d_JxW
void setElementInputVector(const std::vector< std::vector< double > > &elementInputVector)
std::vector< std::vector< double > > d_elementAuxVector
const std::vector< std::vector< libMesh::Real > > * d_phi
std::vector< std::vector< double > > d_elementOtherVectors
std::shared_ptr< libMesh::FEType > d_feType
std::shared_ptr< SourcePhysicsModel > d_sourcePhysicsModel
void initializeForCurrentElement(const libMesh::Elem *elem, std::shared_ptr< SourcePhysicsModel > sourceTransportModel)
void setElementVectors(const std::vector< std::vector< double > > &elementInputVector, const std::vector< std::vector< double > > &elementAuxVector, std::vector< double > &elementOutputVector)
void setElementFlags(const std::string &inputVariableType)
const std::vector< std::vector< libMesh::RealGradient > > * d_dphi
SourceNonlinearElement(std::shared_ptr< const ElementOperationParameters > params)
#define DISABLE_WARNINGS
Re-enable warnings.
#define ENABLE_WARNINGS
Suppress all warnings.



Advanced Multi-Physics (AMP)
Oak Ridge National Laboratory
Idaho National Laboratory
Los Alamos National Laboratory
This page automatically produced from the
source code by doxygen
Last updated: Tue Mar 10 2026 13:06:41.
Comments on this page