7 #if !defined(KRATOS_COMPUTE_MATERIAL_DERIVATIVE_SIMPLEX_ELEMENT_H_INCLUDED )
8 #define KRATOS_COMPUTE_MATERIAL_DERIVATIVE_SIMPLEX_ELEMENT_H_INCLUDED
25 #include "utilities/geometry_utilities.h"
59 template<
unsigned int TDim,
60 unsigned int TNumNodes = TDim + 1 >
143 Element(NewId, pGeometry, pProperties)
169 PropertiesType::Pointer pProperties)
const override
175 void CalculateLocalSystem(
MatrixType& rLeftHandSideMatrix,
186 virtual void EquationIdVector(EquationIdVectorType& rResult,
const ProcessInfo& rCurrentProcessInfo)
const override;
193 virtual void GetDofList(DofsVectorType& rElementalDofList,
const ProcessInfo& rCurrentProcessInfo)
const override;
212 int Check(
const ProcessInfo& rCurrentProcessInfo)
const override;
225 virtual std::string
Info()
const override
227 std::stringstream buffer;
228 buffer <<
"ComputeMaterialDerivativeSimplex #" << Id();
233 virtual void PrintInfo(std::ostream& rOStream)
const override
235 rOStream <<
"ComputeMaterialDerivativeSimplex" << TDim <<
"D";
272 virtual void CalculateMassMatrix(
MatrixType& rMassMatrix,
const ProcessInfo& rCurrentProcessInfo)
override;
276 virtual void CalculateLumpedMassMatrix(
MatrixType& rLHSMatrix,
const double Mass);
280 void CalculateWeights(ShapeFunctionDerivativesArrayType& rDN_DX,
Matrix& rNContainer,
Vector& rGaussWeights);
312 void save(
Serializer& rSerializer)
const override
354 ComputeMaterialDerivativeSimplex &
operator=(ComputeMaterialDerivativeSimplex
const& rOther);
357 ComputeMaterialDerivativeSimplex(ComputeMaterialDerivativeSimplex
const& rOther);
375 template<
unsigned int TDim >
383 template<
unsigned int TDim >
388 rOStream << std::endl;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
A post-processing element to recover the Laplacian from the velocity solution.
Definition: calculate_mat_deriv_simplex_element.h:62
ComputeMaterialDerivativeSimplex(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constuctor using geometry and properties.
Definition: calculate_mat_deriv_simplex_element.h:142
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: calculate_mat_deriv_simplex_element.h:93
Node NodeType
Node type (default is: Node)
Definition: calculate_mat_deriv_simplex_element.h:71
ComputeMaterialDerivativeSimplex(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: calculate_mat_deriv_simplex_element.h:123
ComputeMaterialDerivativeSimplex(IndexType NewId=0)
Default constuctor.
Definition: calculate_mat_deriv_simplex_element.h:114
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new element of this type.
Definition: calculate_mat_deriv_simplex_element.h:168
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: calculate_mat_deriv_simplex_element.h:96
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ComputeMaterialDerivativeSimplex)
Pointer definition of ComputeMaterialDerivativeSimplex.
virtual ~ComputeMaterialDerivativeSimplex()
Destructor.
Definition: calculate_mat_deriv_simplex_element.h:147
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: calculate_mat_deriv_simplex_element.h:83
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: calculate_mat_deriv_simplex_element.h:102
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition of nodes container type, redefined from GeometryType.
Definition: calculate_mat_deriv_simplex_element.h:77
Vector VectorType
Vector type for local contributions to the linear system.
Definition: calculate_mat_deriv_simplex_element.h:80
virtual std::string Info() const override
Turn back information as a string.
Definition: calculate_mat_deriv_simplex_element.h:225
virtual void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: calculate_mat_deriv_simplex_element.h:233
ComputeMaterialDerivativeSimplex(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using a geometry object.
Definition: calculate_mat_deriv_simplex_element.h:132
std::vector< std::size_t > EquationIdVectorType
Definition: calculate_mat_deriv_simplex_element.h:89
Geometry< NodeType > GeometryType
Geometry type (using with given NodeType)
Definition: calculate_mat_deriv_simplex_element.h:74
std::size_t SizeType
Definition: calculate_mat_deriv_simplex_element.h:87
Kratos::Matrix ShapeFunctionDerivativesType
Type for a matrix containing the shape function gradients.
Definition: calculate_mat_deriv_simplex_element.h:99
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: calculate_mat_deriv_simplex_element.h:91
std::size_t IndexType
Definition: calculate_mat_deriv_simplex_element.h:85
Base class for all Elements.
Definition: element.h:60
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: element.h:1135
std::size_t IndexType
Definition: flags.h:74
Geometry base class.
Definition: geometry.h:71
This object defines an indexed object.
Definition: indexed_object.h:54
This class defines the node.
Definition: node.h:65
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
F
Definition: hinsberg_optimization.py:144
def load(f)
Definition: ode_solve.py:307