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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  ajointShapeOptimizationFoam
26 
27 Group
28  grpIncompressibleSolvers
29 
30 Description
31  Steady-state solver for incompressible, turbulent flow of non-Newtonian
32  fluids with optimisation of duct shape by applying "blockage" in regions
33  causing pressure loss as estimated using an adjoint formulation.
34 
35  References:
36  \verbatim
37  "Implementation of a continuous adjoint for topology optimization of
38  ducted flows"
39  C. Othmer,
40  E. de Villiers,
41  H.G. Weller
42  AIAA-2007-3947
43  http://pdf.aiaa.org/preview/CDReadyMCFD07_1379/PV2007_3947.pdf
44  \endverbatim
45 
46  Note that this solver optimises for total pressure loss whereas the
47  above paper describes the method for optimising power-loss.
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include "fvCFD.H"
54 #include "simpleControl.H"
55 #include "fvOptions.H"
56 
57 template<class Type>
58 void zeroCells
59 (
60  GeometricField<Type, fvPatchField, volMesh>& vf,
61  const labelList& cells
62 )
63 {
64  forAll(cells, i)
65  {
66  vf[cells[i]] = pTraits<Type>::zero;
67  }
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73 int main(int argc, char *argv[])
74 {
75  #include "setRootCase.H"
76 
77  #include "createTime.H"
78  #include "createMesh.H"
79 
80  simpleControl simple(mesh);
81 
82  #include "createFields.H"
83  #include "createFvOptions.H"
84  #include "initContinuityErrs.H"
86 
87  turbulence->validate();
88 
89  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90 
91  Info<< "\nStarting time loop\n" << endl;
92 
93  while (simple.loop())
94  {
95  Info<< "Time = " << runTime.timeName() << nl << endl;
96 
97  laminarTransport.lookup("lambda") >> lambda;
98 
99  //alpha +=
100  // mesh.relaxationFactor("alpha")
101  // *(lambda*max(Ua & U, zeroSensitivity) - alpha);
102  alpha +=
103  mesh.fieldRelaxationFactor("alpha")
104  *(min(max(alpha + lambda*(Ua & U), zeroAlpha), alphaMax) - alpha);
105 
107  //zeroCells(alpha, outletCells);
108 
109  // Pressure-velocity SIMPLE corrector
110  {
111  // Momentum predictor
112 
113  tmp<fvVectorMatrix> UEqn
114  (
115  fvm::div(phi, U)
116  + turbulence->divDevReff(U)
117  + fvm::Sp(alpha, U)
118  ==
119  fvOptions(U)
120  );
121 
122  UEqn().relax();
123 
124  fvOptions.constrain(UEqn());
125 
126  solve(UEqn() == -fvc::grad(p));
127 
128  fvOptions.correct(U);
129 
130  volScalarField rAU(1.0/UEqn().A());
131  volVectorField HbyA("HbyA", U);
132  HbyA = rAU*UEqn().H();
133  UEqn.clear();
135  (
136  "phiHbyA",
137  fvc::interpolate(HbyA) & mesh.Sf()
138  );
139  adjustPhi(phiHbyA, U, p);
140 
141  // Non-orthogonal pressure corrector loop
142  while (simple.correctNonOrthogonal())
143  {
144  fvScalarMatrix pEqn
145  (
147  );
148 
149  pEqn.setReference(pRefCell, pRefValue);
150  pEqn.solve();
151 
152  if (simple.finalNonOrthogonalIter())
153  {
154  phi = phiHbyA - pEqn.flux();
155  }
156  }
157 
158  #include "continuityErrs.H"
159 
160  // Explicitly relax pressure for momentum corrector
161  p.relax();
162 
163  // Momentum corrector
164  U = HbyA - rAU*fvc::grad(p);
165  U.correctBoundaryConditions();
166  fvOptions.correct(U);
167  }
168 
169  // Adjoint Pressure-velocity SIMPLE corrector
170  {
171  // Adjoint Momentum predictor
172 
173  volVectorField adjointTransposeConvection((fvc::grad(Ua) & U));
174  //volVectorField adjointTransposeConvection
175  //(
176  // fvc::reconstruct
177  // (
178  // mesh.magSf()*(fvc::snGrad(Ua) & fvc::interpolate(U))
179  // )
180  //);
181 
182  zeroCells(adjointTransposeConvection, inletCells);
183 
184  tmp<fvVectorMatrix> UaEqn
185  (
186  fvm::div(-phi, Ua)
187  - adjointTransposeConvection
188  + turbulence->divDevReff(Ua)
189  + fvm::Sp(alpha, Ua)
190  ==
191  fvOptions(Ua)
192  );
193 
194  UaEqn().relax();
195 
196  fvOptions.constrain(UaEqn());
197 
198  solve(UaEqn() == -fvc::grad(pa));
199 
200  fvOptions.correct(Ua);
201 
202  volScalarField rAUa(1.0/UaEqn().A());
203  volVectorField HbyAa("HbyAa", Ua);
204  HbyAa = rAUa*UaEqn().H();
205  UaEqn.clear();
206  surfaceScalarField phiHbyAa
207  (
208  "phiHbyAa",
209  fvc::interpolate(HbyAa) & mesh.Sf()
210  );
211  adjustPhi(phiHbyAa, Ua, pa);
212 
213  // Non-orthogonal pressure corrector loop
214  while (simple.correctNonOrthogonal())
215  {
216  fvScalarMatrix paEqn
217  (
218  fvm::laplacian(rAUa, pa) == fvc::div(phiHbyAa)
219  );
220 
221  paEqn.setReference(paRefCell, paRefValue);
222  paEqn.solve();
223 
224  if (simple.finalNonOrthogonalIter())
225  {
226  phia = phiHbyAa - paEqn.flux();
227  }
228  }
229 
230  #include "adjointContinuityErrs.H"
231 
232  // Explicitly relax pressure for adjoint momentum corrector
233  pa.relax();
234 
235  // Adjoint momentum corrector
236  Ua = HbyAa - rAUa*fvc::grad(pa);
237  Ua.correctBoundaryConditions();
238  fvOptions.correct(Ua);
239  }
240 
241  laminarTransport.correct();
242  turbulence->correct();
243 
244  runTime.write();
245 
246  Info<< "ExecutionTime = "
247  << runTime.elapsedCpuTime()
248  << " s\n\n" << endl;
249  }
250 
251  Info<< "End\n" << endl;
252 
253  return 0;
254 }
255 
256 
257 // ************************************************************************* //
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:27
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
p
p
Definition: pEqn.H:62
fvOptions.H
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
singlePhaseTransportModel.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
turbulentTransportModel.H
inletCells
const labelList & inletCells
Definition: createFields.H:95
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A
simpleMatrix< scalar > A(Nc)
simple
Simple relative velocity model.
U
U
Definition: pEqn.H:46
createFvOptions.H
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
zeroCells
zeroCells(alpha, inletCells)
solve
rhoEqn solve()
initAdjointContinuityErrs.H
Declare and initialise the cumulative ddjoint continuity error.
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
simpleControl.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
zeroAlpha
dimensionedScalar zeroAlpha("0", dimless/dimTime, 0.0)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
alphaMax
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
adjointContinuityErrs.H
Calculates and prints the continuity errors.
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
setRootCase.H
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
HbyA
HbyA
Definition: pEqn.H:4
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
createMesh.H
createTime.H
fvCFD.H
phiHbyA
phiHbyA
Definition: pEqn.H:21
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:28
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
lambda
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16