Loading...
Searching...
No Matches
adjointShapeOptimizationFoam.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) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Application
27 ajointShapeOptimizationFoam
28
29Group
30 grpIncompressibleSolvers
31
32Description
33 Steady-state solver for incompressible, turbulent flow of non-Newtonian
34 fluids with optimisation of duct shape by applying "blockage" in regions
35 causing pressure loss as estimated using an adjoint formulation.
36
37 References:
38 \verbatim
39 "Implementation of a continuous adjoint for topology optimization of
40 ducted flows"
41 C. Othmer,
42 E. de Villiers,
43 H.G. Weller
44 AIAA-2007-3947
45 http://pdf.aiaa.org/preview/CDReadyMCFD07_1379/PV2007_3947.pdf
46 \endverbatim
47
48 Note that this solver optimises for total pressure loss whereas the
49 above paper describes the method for optimising power-loss.
50
51\*---------------------------------------------------------------------------*/
52
53#include "fvCFD.H"
56#include "simpleControl.H"
57#include "fvOptions.H"
58
59template<class Type>
60void zeroCells
61(
62 GeometricField<Type, fvPatchField, volMesh>& vf,
63 const labelUList& cells
64)
65{
66 UIndirectList<Type>(vf.primitiveField(), cells) = Zero;
67}
68
69
70// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71
72int main(int argc, char *argv[])
73{
74 argList::addNote
75 (
76 "Steady-state solver for incompressible, turbulent flow"
77 " of non-Newtonian fluids with duct shape optimisation"
78 " by applying 'blockage' in regions causing pressure loss"
79 );
80
81 #include "postProcess.H"
82
83 #include "addCheckCaseOptions.H"
84 #include "setRootCaseLists.H"
85 #include "createTime.H"
86 #include "createMesh.H"
87 #include "createControl.H"
88 #include "createFields.H"
89 #include "initContinuityErrs.H"
91
92 turbulence->validate();
93
94 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95
96 Info<< "\nStarting time loop\n" << endl;
97
98 while (simple.loop())
99 {
100 Info<< "Time = " << runTime.timeName() << nl << endl;
101
102 //alpha +=
103 // mesh.relaxationFactor("alpha")
104 // *(lambda*max(Ua & U, zeroSensitivity) - alpha);
105 alpha +=
106 mesh.fieldRelaxationFactor("alpha")
107 *(min(max(alpha + lambda*(Ua & U), zeroAlpha), alphaMax) - alpha);
108
110 //zeroCells(alpha, outletCells);
111
112 // Pressure-velocity SIMPLE corrector
113 {
114 // Momentum predictor
115
116 tmp<fvVectorMatrix> tUEqn
117 (
118 fvm::div(phi, U)
119 + turbulence->divDevReff(U)
120 + fvm::Sp(alpha, U)
121 ==
122 fvOptions(U)
123 );
124 fvVectorMatrix& UEqn = tUEqn.ref();
125
126 UEqn.relax();
127
128 fvOptions.constrain(UEqn);
129
130 solve(UEqn == -fvc::grad(p));
131
132 fvOptions.correct(U);
133
134 volScalarField rAU(1.0/UEqn.A());
136 tUEqn.clear();
137 surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA));
138 adjustPhi(phiHbyA, U, p);
139
140 // Update the pressure BCs to ensure flux consistency
142
143 // Non-orthogonal pressure corrector loop
144 while (simple.correctNonOrthogonal())
145 {
146 fvScalarMatrix pEqn
147 (
148 fvm::laplacian(rAU, p) == fvc::div(phiHbyA)
149 );
150
151 pEqn.setReference(pRefCell, pRefValue);
152 pEqn.solve();
153
154 if (simple.finalNonOrthogonalIter())
155 {
156 phi = phiHbyA - pEqn.flux();
157 }
158 }
159
160 #include "continuityErrs.H"
161
162 // Explicitly relax pressure for momentum corrector
163 p.relax();
164
165 // Momentum corrector
166 U = HbyA - rAU*fvc::grad(p);
167 U.correctBoundaryConditions();
168 fvOptions.correct(U);
169 }
170
171 // Adjoint Pressure-velocity SIMPLE corrector
172 {
173 // Adjoint Momentum predictor
174
175 volVectorField adjointTransposeConvection((fvc::grad(Ua) & U));
176 //volVectorField adjointTransposeConvection
177 //(
178 // fvc::reconstruct
179 // (
180 // mesh.magSf()*fvc::dotInterpolate(fvc::snGrad(Ua), U)
181 // )
182 //);
183
184 zeroCells(adjointTransposeConvection, inletCells);
185
186 tmp<fvVectorMatrix> tUaEqn
187 (
188 fvm::div(-phi, Ua)
189 - adjointTransposeConvection
190 + turbulence->divDevReff(Ua)
191 + fvm::Sp(alpha, Ua)
192 ==
193 fvOptions(Ua)
194 );
195 fvVectorMatrix& UaEqn = tUaEqn.ref();
196
197 UaEqn.relax();
198
199 fvOptions.constrain(UaEqn);
200
201 solve(UaEqn == -fvc::grad(pa));
202
203 fvOptions.correct(Ua);
204
205 volScalarField rAUa(1.0/UaEqn.A());
206 volVectorField HbyAa("HbyAa", Ua);
207 HbyAa = rAUa*UaEqn.H();
208 tUaEqn.clear();
209 surfaceScalarField phiHbyAa("phiHbyAa", fvc::flux(HbyAa));
210 adjustPhi(phiHbyAa, Ua, pa);
211
212 // Non-orthogonal pressure corrector loop
213 while (simple.correctNonOrthogonal())
214 {
216 (
217 fvm::laplacian(rAUa, pa) == fvc::div(phiHbyAa)
218 );
219
220 paEqn.setReference(paRefCell, paRefValue);
221 paEqn.solve();
222
223 if (simple.finalNonOrthogonalIter())
224 {
225 phia = phiHbyAa - paEqn.flux();
226 }
227 }
228
229 #include "adjointContinuityErrs.H"
230
231 // Explicitly relax pressure for adjoint momentum corrector
232 pa.relax();
233
234 // Adjoint momentum corrector
235 Ua = HbyAa - rAUa*fvc::grad(pa);
236 Ua.correctBoundaryConditions();
237 fvOptions.correct(Ua);
238 }
239
240 laminarTransport.correct();
241 turbulence->correct();
242
243 runTime.write();
244
245 runTime.printExecutionTime(Info);
246 }
247
248 Info<< "End\n" << endl;
249
250 return 0;
251}
252
253
254// ************************************************************************* //
Required Classes.
Calculates and prints the continuity errors.
fv::options & fvOptions
const scalar pRefValue
const label pRefCell
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
U
Definition pEqn.H:72
volScalarField & p
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
fvVectorMatrix & UEqn
Definition UEqn.H:13
phiHbyA
Definition pcEqn.H:73
HbyA
Definition pcEqn.H:74
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
compressible::turbulenceModel & turbulence
const cellShapeList & cells
adjustPhi(phiHbyA, U, p_rgh)
Declare and initialise the cumulative ddjoint continuity error.
tmp< volScalarField > rAU
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
fvMatrix< vector > fvVectorMatrix
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
fvScalarMatrix paEqn(fvm::d2dt2(pa) - sqr(c0) *fvc::laplacian(pa))
Execute application functionObjects to post-process existing results.
const dictionary & simple
volScalarField & alpha
CEqn solve()
singlePhaseTransportModel laminarTransport(U, phi)
dimensionedScalar zeroAlpha(dimless/dimTime, Zero)
const labelUList & inletCells
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
zeroCells(alpha, inletCells)
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)