Loading...
Searching...
No Matches
regularisationPDE.H
Go to the documentation of this file.
1
2/*---------------------------------------------------------------------------*\
3 ========= |
4 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / O peration |
6 \\ / A nd | www.openfoam.com
7 \\/ M anipulation |
8-------------------------------------------------------------------------------
9 Copyright (C) 2021 PCOpt/NTUA
10 Copyright (C) 2021 FOSS GP
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28
29Class
30 Foam::regularisationPDE
31
32Description
33 Base class for selecting the regulatisation PDE
34
35
36SourceFiles
37 regularisationPDE.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef regularisationPDE_H
42#define regularisationPDE_H
43
44#include "fvMesh.H"
45#include "topOZones.H"
46#include "fvMatrices.H"
49
50// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52namespace Foam
53{
55/*---------------------------------------------------------------------------*\
56 Class regularisationPDE Declaration
57\*---------------------------------------------------------------------------*/
58
59class regularisationPDE
60{
61
62private:
63
64 // Private Member Functions
65
66 //- No copy construct
67 regularisationPDE(const regularisationPDE&) = delete;
68
69 //- No copy assignment
70 void operator=(const regularisationPDE&) = delete;
71
72
73protected:
74
75 // Protected data
77 const fvMesh& mesh_;
79 const dictionary dict_;
80
81 //- Cell zones related to topology optimisation
82 const topOZones& zones_;
84 //- Whether to apply a fixedValue BC or zeroGradient one to alphaTilda,
85 //- when regularisation is performed
86 bool growFromWalls_;
87
88
89 // Protected Member Functions
90
91 //- Set fixed bTilda values
92 void setValues
93 (
94 fvScalarMatrix& bTildaEqn,
95 const bool isTopoField,
96 const scalar minSetValue = Zero
97 );
98
99 //- Gather the cells with constant values from all constrained zones
100 void setValues
101 (
102 const fvMesh& mesh,
104 DynamicList<scalar>& values,
105 bool isTopoField,
106 const scalar minSetValue = Zero
107 );
108
109 //- Compute smoothing radius, if not directly given
110 scalar computeRadius();
111
112
113public:
114
115 //- Runtime type information
116 TypeName("regularisationPDE");
117
118
119 // Declare run-time constructor selection table
120
122 (
123 autoPtr,
124 regularisationPDE,
126 (
127 const fvMesh& mesh,
128 const dictionary& dict,
129 const topOZones& zones
130 ),
131 (mesh, dict, zones)
132 );
133
134
135 // Constructors
136
137 //- Construct from components
138 regularisationPDE
139 (
140 const fvMesh& mesh,
141 const dictionary& dict,
142 const topOZones& zones
143 );
144
145
146 // Selectors
147
148 //- Construct and return the selected regularisationPDE
150 (
151 const fvMesh& mesh,
152 const dictionary& dict,
153 const topOZones& zones
154 );
155
156
157 //- Destructor
158 virtual ~regularisationPDE() = default;
159
160
161 // Member Functions
162
163 //- Regularise field (or sensitivities)
164 virtual void regularise
165 (
166 const volScalarField& aTilda,
167 const scalarField& source,
168 scalarField& result,
169 const bool isTopoField,
170 const regularisationRadius& radius,
171 const scalar minSetValue = Zero,
172 const bool fixATildaValues = true
173 ) = 0;
174};
176
177// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178
179} // End namespace Foam
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183#endif
184
185// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
bool growFromWalls_
Whether to apply a fixedValue BC or zeroGradient one to alphaTilda, when regularisation is performed.
declareRunTimeSelectionTable(autoPtr, regularisationPDE, dictionary,(const fvMesh &mesh, const dictionary &dict, const topOZones &zones),(mesh, dict, zones))
scalar computeRadius()
Compute smoothing radius, if not directly given.
void setValues(fvScalarMatrix &bTildaEqn, const bool isTopoField, const scalar minSetValue=Zero)
Set fixed bTilda values.
static autoPtr< regularisationPDE > New(const fvMesh &mesh, const dictionary &dict, const topOZones &zones)
Construct and return the selected regularisationPDE.
virtual void regularise(const volScalarField &aTilda, const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius, const scalar minSetValue=Zero, const bool fixATildaValues=true)=0
Regularise field (or sensitivities).
virtual ~regularisationPDE()=default
Destructor.
const topOZones & zones_
Cell zones related to topology optimisation.
TypeName("regularisationPDE")
Runtime type information.
Base class for selecting the regulatisation radius.
dynamicFvMesh & mesh
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const cellShapeList & cells
Namespace for OpenFOAM.
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68