Loading...
Searching...
No Matches
regularisationPDE.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) 2021 PCOpt/NTUA
9 Copyright (C) 2021 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\*---------------------------------------------------------------------------*/
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35 defineTypeNameAndDebug(regularisationPDE, 0);
36 defineRunTimeSelectionTable(regularisationPDE, dictionary);
37}
38
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
43(
44 fvScalarMatrix& bTildaEqn,
45 const bool isTopoField,
46 const scalar minSetValue
47)
48{
49 const fvMesh& mesh = bTildaEqn.psi().internalField().mesh();
51 DynamicList<scalar> values(mesh.nCells()/100);
52 setValues(mesh, cells, values, isTopoField, minSetValue);
53 bTildaEqn.setValues(cells, values);
54}
55
56
58(
59 const fvMesh& mesh,
61 DynamicList<scalar>& values,
62 bool isTopoField,
63 const scalar minSetValue
64)
65{
66 const labelList& IOcells = mesh.cellZones()[zones_.IOzoneID()];
67 cells.append(IOcells);
68 values.append(scalarField(IOcells.size(), minSetValue));
69
70 forAll(zones_.fixedPorousZoneIDs(), zI)
71 {
72 const label cellZoneID = zones_.fixedPorousZoneIDs()[zI];
73 const labelList& zoneCells = mesh.cellZones()[cellZoneID];
74 cells.append(zoneCells);
75 values.append
76 (
78 (
79 zoneCells.size(),
80 isTopoField ?
81 zones_.fixedPorousValues()[zI] : scalar(0)
82 )
83 );
84 }
85
86 for (label cellZoneID : zones_.fixedZeroPorousZoneIDs())
87 {
88 const labelList& zoneCells = mesh.cellZones()[cellZoneID];
89 cells.append(zoneCells);
90 values.append(scalarField(zoneCells.size(), minSetValue));
91 }
92}
93
94
96{
97 scalar averageVol(gAverage(mesh_.V().field()));
98 const Vector<label>& geometricD = mesh_.geometricD();
99 const boundBox& bounds = mesh_.bounds();
100 forAll(geometricD, iDir)
101 {
102 if (geometricD[iDir] == -1)
103 {
104 averageVol /= bounds.span()[iDir];
105 }
106 }
107 scalar radius = pow(averageVol, scalar(1)/scalar(mesh_.nGeometricD()));
108 Info<< "Computed a mean radius of " << radius << endl;
109 return radius;
110}
111
112
113// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114
115Foam::regularisationPDE::regularisationPDE
116(
117 const fvMesh& mesh,
118 const dictionary& dict,
119 const topOZones& zones
120)
121:
122 mesh_(mesh),
123 dict_(dict),
124 zones_(zones),
125 growFromWalls_(dict.getOrDefault<bool>("growFromWalls", false))
126{}
127
128
129// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
130
132(
133 const fvMesh& mesh,
134 const dictionary& dict,
135 const topOZones& zones
136)
137{
138 const word modelType =
139 dict.getOrDefault<word>("regularisationPDE", "Helmholtz");
140
141 auto* ctorPtr = dictionaryConstructorTable(modelType);
142
143 Info<< "regularisationPDE type " << modelType << endl;
144
145 if (!ctorPtr)
146 {
148 (
149 dict,
150 "regularisationPDE",
151 modelType,
152 *dictionaryConstructorTablePtr_
153 ) << exit(FatalIOError);
154 }
155
156 return autoPtr<regularisationPDE>
157 (
158 ctorPtr(mesh, dict, zones)
159 );
160}
161
162
163// ************************************************************************* //
const Mesh & mesh() const noexcept
Return const reference to mesh.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
void setValues(const labelUList &cellLabels, const Type &value)
Set solution in given cells to the specified value and eliminate the corresponding equations from the...
Definition fvMatrix.C:971
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base class for selecting the regulatisation PDE.
bool growFromWalls_
Whether to apply a fixedValue BC or zeroGradient one to alphaTilda, when regularisation is performed.
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.
const topOZones & zones_
Cell zones related to topology optimisation.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
const cellShapeList & cells
Namespace for bounding specifications. At the moment, mostly for tables.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
List< label > labelList
A List of labels.
Definition List.H:62
fvMatrix< scalar > fvScalarMatrix
messageStream Info
Information stream (stdout output on master, null elsewhere).
DynamicID< cellZoneMesh > cellZoneID
Foam::cellZoneID.
Definition ZoneIDs.H:38
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299