All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
adjointShapeOptimisationFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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  adjointShapeOptimisationFoam
26 
27 Description
28  Steady-state solver for incompressible, turbulent flow of non-Newtonian
29  fluids with optimisation of duct shape by applying "blockage" in regions
30  causing pressure loss as estimated using an adjoint formulation.
31 
32  References:
33  \verbatim
34  "Implementation of a continuous adjoint for topology optimisation of
35  ducted flows"
36  C. Othmer,
37  E. de Villiers,
38  H.G. Weller
39  AIAA-2007-3947
40  http://pdf.aiaa.org/preview/CDReadyMCFD07_1379/PV2007_3947.pdf
41  \endverbatim
42 
43  Note that this solver optimises for total pressure loss whereas the
44  above paper describes the method for optimising power-loss.
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "fvCFD.H"
51 #include "simpleControl.H"
52 #include "pressureReference.H"
53 #include "fvModels.H"
54 #include "fvConstraints.H"
55 
56 template<class Type>
57 void zeroCells
58 (
59  GeometricField<Type, fvPatchField, volMesh>& vf,
60  const labelList& cells
61 )
62 {
63  forAll(cells, i)
64  {
65  vf[cells[i]] = Zero;
66  }
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 int main(int argc, char *argv[])
73 {
74  #include "postProcess.H"
75 
76  #include "setRootCaseLists.H"
77  #include "createTime.H"
78  #include "createMesh.H"
79  #include "createControl.H"
80  #include "createFields.H"
81  #include "initContinuityErrs.H"
83 
84  turbulence->validate();
85 
86  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
87 
88  Info<< "\nStarting time loop\n" << endl;
89 
90  while (simple.loop(runTime))
91  {
92  Info<< "Time = " << runTime.timeName() << nl << endl;
93 
94  laminarTransport.lookup("lambda") >> lambda;
95 
96  // alpha +=
97  // mesh.relaxationFactor("alpha")
98  // *(lambda*max(Ua & U, zeroSensitivity) - alpha);
99  alpha +=
100  mesh.fieldRelaxationFactor("alpha")
101  *(min(max(alpha + lambda*(Ua & U), zeroAlpha), alphaMax) - alpha);
102 
104  // zeroCells(alpha, outletCells);
105 
106  // Pressure-velocity SIMPLE corrector
107  {
108  fvModels.correct();
109 
110  // Momentum predictor
111 
112  tmp<fvVectorMatrix> tUEqn
113  (
114  fvm::div(phi, U)
115  + turbulence->divDevSigma(U)
116  + fvm::Sp(alpha, U)
117  ==
118  fvModels.source(U)
119  );
120  fvVectorMatrix& UEqn = tUEqn.ref();
121 
122  UEqn.relax();
123 
124  fvConstraints.constrain(UEqn);
125 
126  solve(UEqn == -fvc::grad(p));
127 
129 
130  volScalarField rAU(1.0/UEqn.A());
131  volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
132  tUEqn.clear();
134  adjustPhi(phiHbyA, U, p);
135 
136  // Update the pressure BCs to ensure flux consistency
138 
139  // Non-orthogonal pressure corrector loop
140  while (simple.correctNonOrthogonal())
141  {
142  fvScalarMatrix pEqn
143  (
145  );
146 
147  pEqn.setReference
148  (
149  pressureReference.refCell(),
150  pressureReference.refValue()
151  );
152 
153  pEqn.solve();
154 
155  if (simple.finalNonOrthogonalIter())
156  {
157  phi = phiHbyA - pEqn.flux();
158  }
159  }
160 
161  #include "continuityErrs.H"
162 
163  // Explicitly relax pressure for momentum corrector
164  p.relax();
165 
166  // Momentum corrector
167  U = HbyA - rAU*fvc::grad(p);
168  U.correctBoundaryConditions();
170  }
171 
172  // Adjoint Pressure-velocity SIMPLE corrector
173  {
174  // Adjoint Momentum predictor
175 
176  volVectorField adjointTransposeConvection((fvc::grad(Ua) & U));
177  // volVectorField adjointTransposeConvection
178  //(
179  // fvc::reconstruct
180  // (
181  // mesh.magSf()*fvc::dotInterpolate(fvc::snGrad(Ua), U)
182  // )
183  //);
184 
185  zeroCells(adjointTransposeConvection, inletCells);
186 
187  tmp<fvVectorMatrix> tUaEqn
188  (
189  fvm::div(-phi, Ua)
190  - adjointTransposeConvection
191  + turbulence->divDevSigma(Ua)
192  + fvm::Sp(alpha, Ua)
193  ==
194  fvModels.source(Ua)
195  );
196  fvVectorMatrix& UaEqn = tUaEqn.ref();
197 
198  UaEqn.relax();
199 
200  fvConstraints.constrain(UaEqn);
201 
202  solve(UaEqn == -fvc::grad(pa));
203 
205 
206  volScalarField rAUa(1.0/UaEqn.A());
207  volVectorField HbyAa("HbyAa", Ua);
208  HbyAa = rAUa*UaEqn.H();
209  tUaEqn.clear();
210  surfaceScalarField phiHbyAa("phiHbyAa", fvc::flux(HbyAa));
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();
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 // ************************************************************************* //
pressureReference & pressureReference
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const labelList & inletCells
Definition: createFields.H:94
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
zeroCells(alpha, inletCells)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
dimensionedScalar zeroAlpha("0", dimless/dimTime, 0.0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
const cellShapeList & cells
Foam::fvConstraints & fvConstraints
phi
Definition: correctPhi.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
static const zero Zero
Definition: zero.H:97
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
Declare and initialise the cumulative adjoint continuity error.
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvModels.source(rho, U))
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
phiHbyA
Definition: pEqn.H:32
Foam::fvModels & fvModels
rhoEqn solve()
U
Definition: pEqn.H:72
tmp< volScalarField > rAU
volVectorField & HbyA
Definition: pEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Calculates and prints the continuity errors.
messageStream Info
volScalarField & p
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
singlePhaseTransportModel laminarTransport(U, phi)
simpleControl simple(mesh)