KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
adjoint_solid_element.h
Go to the documentation of this file.
1 // KRATOS ___| | | |
2 // \___ \ __| __| | | __| __| | | __| _` | |
3 // | | | | | ( | | | | ( | |
4 // _____/ \__|_| \__,_|\___|\__|\__,_|_| \__,_|_| MECHANICS
5 //
6 // License: BSD License
7 // license: StructuralMechanicsApplication/license.txt
8 //
9 // Main authors:
10 //
11 
12 #pragma once
13 
14 // System includes
15 
16 // External include
17 
18 // Project includes
19 #include "includes/define.h"
20 #include "includes/element.h"
23 
24 namespace Kratos
25 {
28 
34 template <class TPrimalElement>
36  : public Element
37 {
38  class ThisExtensions : public AdjointExtensions
39  {
40  Element* mpElement;
41 
42  public:
43  explicit ThisExtensions(Element* pElement);
44 
45  void GetFirstDerivativesVector(std::size_t NodeId,
46  std::vector<IndirectScalar<double>>& rVector,
47  std::size_t Step) override;
48 
49  void GetSecondDerivativesVector(std::size_t NodeId,
50  std::vector<IndirectScalar<double>>& rVector,
51  std::size_t Step) override;
52 
53  void GetAuxiliaryVector(std::size_t NodeId,
54  std::vector<IndirectScalar<double>>& rVector,
55  std::size_t Step) override;
56 
57  void GetFirstDerivativesVariables(std::vector<VariableData const*>& rVariables) const override;
58 
59  void GetSecondDerivativesVariables(std::vector<VariableData const*>& rVariables) const override;
60 
61  void GetAuxiliaryVariables(std::vector<VariableData const*>& rVariables) const override;
62  };
63 
64 public:
67 
69 
73 
74  AdjointSolidElement(IndexType NewId = 0);
75 
76  AdjointSolidElement(IndexType NewId, GeometryType::Pointer pGeometry);
77 
79  GeometryType::Pointer pGeometry,
80  PropertiesType::Pointer pProperties);
81 
85 
86  Element::Pointer Create(IndexType NewId,
87  NodesArrayType const& ThisNodes,
88  PropertiesType::Pointer pProperties) const override;
89 
90  Element::Pointer Create(IndexType NewId,
91  GeometryType::Pointer pGeom,
92  PropertiesType::Pointer pProperties) const override;
93 
94  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
95 
96  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
97 
98  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
99 
100  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
101 
102  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
103 
104  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,
105  const ProcessInfo& rCurrentProcessInfo) override;
106 
107  void CalculateFirstDerivativesLHS(MatrixType& rLeftHandSideMatrix,
108  const ProcessInfo& rCurrentProcessInfo) override;
109 
110  void CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix,
111  const ProcessInfo& rCurrentProcessInfo) override;
112 
113  void GetValuesVector(Vector& rValues, int Step = 0) const override;
114 
116  const ProcessInfo& rCurrentProcessInfo) const override;
117 
118  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
119 
120  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
121 
122  void CalculateSensitivityMatrix(const Variable<array_1d<double, 3>>& rDesignVariable,
123  Matrix& rOutput,
124  const ProcessInfo& rCurrentProcessInfo) override;
125 
127 
128 protected:
129 
130 private:
133 
134  TPrimalElement mPrimalElement;
135 
139 
143  friend class Serializer;
144 
145  void save(Serializer& rSerializer) const override
146  {
148  rSerializer.save("mPrimalElement", mPrimalElement);
149  }
150 
151  void load(Serializer& rSerializer) override
152  {
154  rSerializer.load("mPrimalElement", mPrimalElement);
155  }
157 };
158 
160 
161 } // namespace Kratos.
Interface extensions for adjoint elements and conditions.
Definition: adjoint_extensions.h:37
A template class for creating adjoint elements for solids.
Definition: adjoint_solid_element.h:37
AdjointSolidElement(IndexType NewId=0)
Definition: adjoint_solid_element.cpp:113
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:151
void CalculateSecondDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:212
void GetValuesVector(Vector &rValues, int Step=0) const override
Definition: adjoint_solid_element.cpp:222
void CalculateSensitivityMatrix(const Variable< array_1d< double, 3 >> &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:326
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(AdjointSolidElement)
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Definition: adjoint_solid_element.cpp:281
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: adjoint_solid_element.cpp:312
void FinalizeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:176
void FinalizeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:184
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: adjoint_solid_element.cpp:242
void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:160
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
It creates a new element pointer.
Definition: adjoint_solid_element.cpp:134
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:192
void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:168
void CalculateFirstDerivativesLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_solid_element.cpp:202
Base class for all Elements.
Definition: element.h:60
virtual void GetFirstDerivativesVector(Vector &values, int Step=0) const
Definition: element.h:310
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
virtual void GetSecondDerivativesVector(Vector &values, int Step=0) const
Definition: element.h:320
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
Wrapper for a function which behaves like an arithmetic type.
Definition: indirect_scalar.h:45
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
void load(std::string const &rTag, TDataType &rObject)
Definition: serializer.h:207
void save(std::string const &rTag, std::array< TDataType, TDataSize > const &rObject)
Definition: serializer.h:545
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307