12 #if !defined(KRATOS_FLUX_CONDITION_H_INCLUDED )
13 #define KRATOS_FLUX_CONDITION_H_INCLUDED
38 namespace FluxConditionInternals
43 template<
unsigned int TNodeNumber >
53 mNodalFluxes(TNodeNumber,0.0)
66 for (
unsigned int g = 0; g <
NumGauss; g++)
68 mIntegrationWeights[g] = DetJ[g]*IntegrationPoints[g].Weight();
71 for (
unsigned int i = 0;
i < TNodeNumber;
i++)
73 mNodalFluxes[
i] = rGeometry[
i].FastGetSolutionStepValue(rFluxVar);
82 double N(
unsigned int i)
const
84 return mShapeFunctionValues(mGaussPoint,
i);
89 return mNodalFluxes[
i];
94 double flux = mNodalFluxes[0]*mShapeFunctionValues(mGaussPoint,0);
95 for (
unsigned int i = 1;
i < TNodeNumber;
i++)
97 flux += mNodalFluxes[
i]*mShapeFunctionValues(mGaussPoint,
i);
104 return mIntegrationWeights[mGaussPoint];
111 unsigned int mGaussPoint;
115 Matrix mShapeFunctionValues;
117 Vector mIntegrationWeights;
133 template<
unsigned int TNodeNumber >
157 Properties::Pointer pProperties);
166 Condition::Pointer
Create(
169 Properties::Pointer pProperties)
const override;
171 Condition::Pointer
Create(
173 GeometryType::Pointer pGeom,
174 Properties::Pointer pProperties)
const override;
187 const ProcessInfo& rCurrentProcessInfo)
const override;
191 const ProcessInfo& CurrentProcessInfo)
const override;
202 std::vector<double>& rValues,
212 std::vector<Vector>& rValues,
217 std::vector<Matrix>& rValues,
225 std::string
Info()
const override;
228 void PrintInfo(std::ostream& rOStream)
const override;
231 void PrintData(std::ostream& rOStream)
const override;
261 void save(
Serializer& rSerializer)
const override;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
std::size_t IndexType
Definition: flags.h:74
A basic Neumann condition for convection-diffusion problems.
Definition: flux_condition.h:135
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, Properties::Pointer pProperties) const override
Definition: flux_condition.cpp:45
Condition::VectorType VectorType
Definition: flux_condition.h:144
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(FluxCondition)
Pointer definition of FluxCondition.
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: flux_condition.cpp:82
Condition::MatrixType MatrixType
Definition: flux_condition.h:143
GeometryData::IntegrationMethod GetIntegrationMethod() const override
Definition: flux_condition.cpp:166
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 > > &rVariable, std::vector< array_1d< double, 3 > > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Definition: flux_condition.cpp:172
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: flux_condition.cpp:282
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: flux_condition.cpp:276
FluxCondition()
Definition: flux_condition.cpp:343
void CalculateNormal(array_1d< double, 3 > &An)
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: flux_condition.cpp:63
~FluxCondition() override
Destructor.
Definition: flux_condition.cpp:38
void AddIntegrationPointRHSContribution(VectorType &rRightHandSideVector, const FluxConditionInternals::IntegrationData< TNodeNumber > &rData)
Definition: flux_condition.cpp:292
std::string Info() const override
Turn back information as a string.
Definition: flux_condition.cpp:268
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: flux_condition.cpp:112
void GetDofList(DofsVectorType &ConditionalDofList, const ProcessInfo &CurrentProcessInfo) const override
Definition: flux_condition.cpp:139
Definition: flux_condition.h:45
double GaussPointFlux() const
Definition: flux_condition.h:92
void SetCurrentGaussPoint(unsigned int g)
Definition: flux_condition.h:77
IntegrationData(Geometry< Node > &rGeometry, const Variable< double > &rFluxVar)
Definition: flux_condition.h:48
double IntegrationWeight() const
Definition: flux_condition.h:102
double NodalFlux(unsigned int i) const
Definition: flux_condition.h:87
unsigned int NumGauss
Definition: flux_condition.h:107
double N(unsigned int i) const
Definition: flux_condition.h:82
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
const Matrix & ShapeFunctionsValues() const
Definition: geometry.h:3393
Vector & DeterminantOfJacobian(Vector &rResult) const
Definition: geometry.h:3154
const IntegrationPointsArrayType & IntegrationPoints() const
Definition: geometry.h:2284
SizeType IntegrationPointsNumber() const
Definition: geometry.h:2257
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
def load(f)
Definition: ode_solve.py:307
integer i
Definition: TensorModule.f:17