Loading...
Searching...
No Matches
fieldRegularisation.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fieldRegularisation.H"
31#include "wallFvPatch.H"
32#include "fvm.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38 defineTypeNameAndDebug(fieldRegularisation, 1);
39}
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43Foam::fieldRegularisation::fieldRegularisation
44(
45 fvMesh& mesh,
46 const scalarField& alpha,
47 const topOZones& zones,
48 const dictionary& dict
49)
50:
51 mesh_(mesh),
52 dict_(dict),
53 zones_(zones),
54 regularise_(dict.getOrDefault<bool>("regularise", false)),
55 project_(dict.getOrDefault<bool>("project", regularise_)),
56 radius_(regularisationRadius::New(mesh, dict, false)),
57 alpha_(alpha),
58 alphaTilda_
59 (
60 regularise_
61 ? new volScalarField
62 (
64 (
65 "alphaTilda",
66 mesh_.time().timeName(),
67 mesh_,
68 IOobject::READ_IF_PRESENT,
69 IOobject::AUTO_WRITE
70 ),
71 mesh_,
73 fixedValueFvPatchScalarField::typeName
74 )
75 : nullptr
76 ),
77 sharpenFunction_
78 (
79 project_ ?
81 nullptr
82 ),
83 regularisationPDE_(regularisationPDE::New(mesh, dict, zones)),
84 betaArg_(regularise_ ? alphaTilda_().primitiveField() : alpha_),
85 growFromWalls_(dict.getOrDefault<bool>("growFromWalls", false)),
86 beta_
87 (
89 (
90 "beta",
91 mesh_.time().timeName(),
92 mesh_,
93 IOobject::READ_IF_PRESENT,
94 IOobject::AUTO_WRITE
95 ),
96 mesh_,
98 fvPatchFieldBase::zeroGradientType()
99 )
100{
102 << "Regularise " << Switch(regularise_) << nl
103 << "Project " << Switch(project_) << endl;
104}
105
106
107// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108
110{
111 if (regularise_)
112 {
113 regularise(alpha_, alphaTilda_(), true);
114 }
115
116 if (project_)
117 {
118 sharpenFunction_->interpolate(betaArg_, beta_.primitiveFieldRef());
119 }
120 else
121 {
122 beta_.primitiveFieldRef() = betaArg_;
123 }
124
125 beta_.correctBoundaryConditions();
126}
127
128
130(
131 const scalarField& source,
132 scalarField& result,
133 const bool isTopoField,
134 const regularisationRadius& radius
136{
138 regularise(alphaTilda_(), source, result, isTopoField, radius);
139}
140
141
143(
144 const scalarField& source,
145 scalarField& result,
146 const bool isTopoField
148{
150 regularise(alphaTilda_(), source, result, isTopoField, radius_());
151}
152
153
155{
156 // Add dBeta/dBetaArg
157 if (project_)
158 {
159 sens *= sharpenFunction_->derivative(betaArg_);
160 }
161 // Add part due to regularisation
162 if (regularise_)
163 {
164 // Solve the adjoint to the regularisation equation
165 regularise(sens, sens, false);
166 }
167
168 // Add volume
169 sens *= mesh_.V();
170}
171
172
173// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool growFromWalls_
Whether to apply a fixedValue BC or zeroGradient one to alphaTilda, when regularisation is performed.
autoPtr< regularisationRadius > radius_
Smoothing radius.
bool regularise_
Perform regulaisation on alpha before inputing it on beta?
virtual void updateBeta()
Update the beta field.
bool project_
Perform the projection (sharpening) step?
volScalarField beta_
Beta is the field used for all interpolations between fluid and solid in topology optimisation.
autoPtr< topOInterpolationFunction > sharpenFunction_
Function used to sharpen the field after regularisation.
void postProcessSens(scalarField &sens)
Update part of fieldRegularisation to the sensitivitiy derivatives.
autoPtr< volScalarField > alphaTilda_
The regularised alpha field, if regulatisation is performed.
void regularise(const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius)
Regularise an externally provided radius.
autoPtr< regularisationPDE > regularisationPDE_
PDE used for the regularisation.
const topOZones & zones_
Cell zones related to topology optimisation.
const scalarField & betaArg_
Argument of the beta field.
const scalarField & alpha_
Alpha field (design variables of topology optimisation).
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Template invariant parts for fvPatchField.
Base class for selecting the regulatisation PDE.
Base class for selecting the regulatisation radius.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict