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_base_potential_flow_element.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 //
11 // Main authors: Marc Nunez, based on A. Geiser, M. Fusseder, I. Lopez and R. Rossi work
12 //
13 
14 #if !defined(KRATOS_ADJOINT_BASE_POTENTIAL_FLOW_ELEMENT_H_INCLUDED )
15 #define KRATOS_ADJOINT_BASE_POTENTIAL_FLOW_ELEMENT_H_INCLUDED
16 
17 
18 // Project includes
19 #include "includes/element.h"
20 
21 namespace Kratos
22 {
23 
24 template <class TPrimalElement>
26 {
27 public:
28 
31 
32  typedef Element BaseType;
33  static constexpr int NumNodes = TPrimalElement::TNumNodes;
34  static constexpr int Dim = TPrimalElement::TDim;
35 
40 
44 
49  : Element(NewId),
50  mpPrimalElement(Kratos::make_intrusive<TPrimalElement>(NewId))
51  {};
52 
54  GeometryType::Pointer pGeometry)
55  : Element(NewId, pGeometry),
56  mpPrimalElement(Kratos::make_intrusive<TPrimalElement>(NewId, pGeometry))
57  {
58  }
59 
61  GeometryType::Pointer pGeometry,
62  PropertiesType::Pointer pProperties)
63  : Element(NewId, pGeometry, pProperties),
64  mpPrimalElement(Kratos::make_intrusive<TPrimalElement>(NewId, pGeometry, pProperties))
65  {
66  }
71 
76 
80 
83  {
84  BaseType::operator=(rOther);
85  Flags::operator =(rOther);
86  return *this;
87  }
88 
92 
93  IntegrationMethod GetIntegrationMethod() const override;
94 
95  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
96 
97  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
98 
99  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
100 
101  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
102  VectorType& rRightHandSideVector,
103  const ProcessInfo& rCurrentProcessInfo) override;
104 
105  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,
106  const ProcessInfo& rCurrentProcessInfo) override;
107 
108  void CalculateRightHandSide(VectorType& rRightHandSideVector,
109  const ProcessInfo& rCurrentProcessInfo) override;
110 
111  void CalculateOnIntegrationPoints(const Variable<int>& rVariable,
112  std::vector<int>& rValues,
113  const ProcessInfo& rCurrentProcessInfo) override;
114 
115  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
116  std::vector<double>& rValues,
117  const ProcessInfo& rCurrentProcessInfo) override;
118 
120  std::vector< array_1d<double,3> >& rValues,
121  const ProcessInfo& rCurrentProcessInfo) override;
122 
123  void GetValuesVector(Vector& rValues, int Step=0) const override;
124 
125  void CalculateSensitivityMatrix(const Variable<double>& rDesignVariable,
126  Matrix& rOutput,
127  const ProcessInfo& rCurrentProcessInfo) override;
128 
129  void CalculateSensitivityMatrix(const Variable<array_1d<double,3> >& rDesignVariable,
130  Matrix& rOutput,
131  const ProcessInfo& rCurrentProcessInfo) override;
132 
133  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& CurrentProcessInfo) const override;
134 
135  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& CurrentProcessInfo) const override;
136 
137  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
138 
139  std::string Info() const override;
140 
141  void PrintInfo(std::ostream& rOStream) const override;
142 
143  void PrintData(std::ostream& rOStream) const override;
144 
145  Element::Pointer pGetPrimalElement();
146 
147 protected:
148 
149  Element::Pointer mpPrimalElement;
150 
151  void GetValuesOnSplitElement(Vector& split_element_values, const array_1d<double,NumNodes>& distances) const;
152 
153 private:
154 
155  friend class Serializer;
156 
157  void save(Serializer& rSerializer) const override;
158 
159  void load(Serializer& rSerializer) override;
160 
161 }; // Class AdjointBasePotentialFlowElement
162 
163 
164 } // namespace Kratos.
165 
166 #endif
Definition: adjoint_base_potential_flow_element.h:26
IntegrationMethod GetIntegrationMethod() const override
Definition: adjoint_base_potential_flow_element.cpp:25
Element BaseType
Definition: adjoint_base_potential_flow_element.h:32
Element::Pointer mpPrimalElement
Definition: adjoint_base_potential_flow_element.h:149
~AdjointBasePotentialFlowElement() override
Definition: adjoint_base_potential_flow_element.h:75
void CalculateSensitivityMatrix(const Variable< double > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:151
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: adjoint_base_potential_flow_element.cpp:318
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &CurrentProcessInfo) const override
Definition: adjoint_base_potential_flow_element.cpp:163
std::string Info() const override
Turn back information as a string.
Definition: adjoint_base_potential_flow_element.cpp:304
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: adjoint_base_potential_flow_element.cpp:271
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(AdjointBasePotentialFlowElement)
void CalculateOnIntegrationPoints(const Variable< int > &rVariable, std::vector< int > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:79
Element::Pointer pGetPrimalElement()
Definition: adjoint_base_potential_flow_element.cpp:324
AdjointBasePotentialFlowElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: adjoint_base_potential_flow_element.h:60
AdjointBasePotentialFlowElement(AdjointBasePotentialFlowElement const &rOther)
Definition: adjoint_base_potential_flow_element.h:70
AdjointBasePotentialFlowElement(IndexType NewId=0)
Definition: adjoint_base_potential_flow_element.h:48
AdjointBasePotentialFlowElement & operator=(AdjointBasePotentialFlowElement const &rOther)
Assignment operator.
Definition: adjoint_base_potential_flow_element.h:82
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:62
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:31
void GetValuesVector(Vector &rValues, int Step=0) const override
Definition: adjoint_base_potential_flow_element.cpp:103
void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:37
void GetValuesOnSplitElement(Vector &split_element_values, const array_1d< double, NumNodes > &distances) const
Definition: adjoint_base_potential_flow_element.cpp:332
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: adjoint_base_potential_flow_element.cpp:312
void FinalizeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:46
static constexpr int NumNodes
Definition: adjoint_base_potential_flow_element.h:33
static constexpr int Dim
Definition: adjoint_base_potential_flow_element.h:34
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:71
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: adjoint_base_potential_flow_element.cpp:52
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &CurrentProcessInfo) const override
Definition: adjoint_base_potential_flow_element.cpp:218
AdjointBasePotentialFlowElement(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: adjoint_base_potential_flow_element.h:53
Base class for all Elements.
Definition: element.h:60
Element & operator=(Element const &rOther)
Assignment operator.
Definition: element.h:179
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Matrix MatrixType
Definition: element.h:90
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
Flags & operator=(Flags const &rOther)
Assignment operator.
Definition: flags.h:151
IntegrationMethod
Definition: geometry_data.h:76
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
intrusive_ptr< C > make_intrusive(Args &&...args)
Definition: smart_pointers.h:36
def load(f)
Definition: ode_solve.py:307