objectiveMoment.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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019-2020 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "objectiveMoment.H"
31 #include "createZeroField.H"
32 #include "wallFvPatch.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 namespace objectives
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(objectiveMoment, 0);
47 (
48  objectiveIncompressible,
49  objectiveMoment,
50  dictionary
51 );
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const fvMesh& mesh,
59  const dictionary& dict,
60  const word& adjointSolverName,
61  const word& primalSolverName
62 )
63 :
64  objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
65  momentPatches_
66  (
67  mesh_.boundaryMesh().patchSet
68  (
69  dict.get<wordRes>("patches")
70  ).sortedToc()
71  ),
72  momentDirection_(dict.get<vector>("direction")),
73  rotationCentre_(dict.get<vector>("rotationCenter")),
74  Aref_(dict.get<scalar>("Aref")),
75  lRef_(dict.get<scalar>("lRef")),
76  rhoInf_(dict.get<scalar>("rhoInf")),
77  UInf_(dict.get<scalar>("UInf")),
78  invDenom_(2./(rhoInf_*UInf_*UInf_*Aref_*lRef_)),
79  stressXPtr_
80  (
81  Foam::createZeroFieldPtr<vector>
82  (
83  mesh_, "stressX", dimLength/sqr(dimTime)
84  )
85  ),
86  stressYPtr_
87  (
88  Foam::createZeroFieldPtr<vector>
89  (
90  mesh_, "stressY", dimLength/sqr(dimTime)
91  )
92  ),
93  stressZPtr_
94  (
95  Foam::createZeroFieldPtr<vector>
96  (
97  mesh_, "stressZ", dimLength/sqr(dimTime)
98  )
99  ),
100  devReff_(vars_.turbulence()->devReff()())
101 {
102  // Sanity check and print info
103  if (momentPatches_.empty())
104  {
106  << "No valid patch name on which to minimize " << type() << endl
107  << exit(FatalError);
108  }
109  if (debug)
110  {
111  Info<< "Minimizing " << type() << " in patches:" << endl;
112  for (const label patchI : momentPatches_)
113  {
114  Info<< "\t " << mesh_.boundary()[patchI].name() << endl;
115  }
116  }
117 
118  // Allocate boundary field pointers
119  bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
120  bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
121  bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
122  bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 {
130  vector pressureMoment(Zero);
131  vector viscousMoment(Zero);
132  vector cumulativeMoment(Zero);
133 
134  // Update field here and use the same value for all functions
135  const volScalarField& p = vars_.pInst();
136  devReff_ = vars_.turbulence()->devReff()();
137 
138  for (const label patchI : momentPatches_)
139  {
140  const fvPatch& patch = mesh_.boundary()[patchI];
141  const vectorField& Sf = patch.Sf();
142  vectorField dx(patch.Cf() - rotationCentre_);
143  pressureMoment += gSum
144  (
145  rhoInf_*(dx ^ Sf)*p.boundaryField()[patchI]
146  );
147 
148  // Viscous term calculated using the full tensor derivative
149  viscousMoment += gSum
150  (
151  rhoInf_*(dx^(devReff_.boundaryField()[patchI] & Sf))
152  );
153  }
154 
155  cumulativeMoment = pressureMoment + viscousMoment;
156 
157  scalar moment = cumulativeMoment & momentDirection_;
158  scalar Cm = moment*invDenom_;
159  DebugInfo<<
160  "Moment|Coeff " << moment << "|" << Cm << endl;
161  J_ = Cm;
162  return Cm;
163 }
164 
165 
167 {
168  if (computeMeanFields_)
169  {
170  const volVectorField& U = vars_.U();
172  turbVars = vars_.RASModelVariables();
173  const singlePhaseTransportModel& lamTransp = vars_.laminarTransport();
174 
175  devReff_ = turbVars->devReff(lamTransp, U)();
176  }
177 }
178 
179 
181 {
182  for (const label patchI : momentPatches_)
183  {
184  const fvPatch& patch = mesh_.boundary()[patchI];
185  vectorField dx(patch.Cf() - rotationCentre_);
186  bdJdpPtr_()[patchI] = (momentDirection_ ^ dx)*invDenom_*rhoInf_;
187  }
188 }
189 
190 
192 {
193  const volScalarField& p = vars_.p();
194 
195  for (const label patchI : momentPatches_)
196  {
197  const fvPatch& patch = mesh_.boundary()[patchI];
198  const vectorField dx(patch.Cf() - rotationCentre_);
199  bdSdbMultPtr_()[patchI] =
200  (
201  (
202  rhoInf_*
203  (
204  (momentDirection_^dx) &
205  (
206  devReff_.boundaryField()[patchI]
207  )
208  )
209  )
210  + rhoInf_ * (momentDirection_^dx) * p.boundaryField()[patchI]
211  )
212  *invDenom_;
213  }
214 }
215 
216 
218 {
219  const volScalarField& p = vars_.p();
220  const volVectorField& U = vars_.U();
221 
223  turbVars = vars_.RASModelVariables();
224  const singlePhaseTransportModel& lamTransp = vars_.laminarTransport();
225  volScalarField nuEff(lamTransp.nu() + turbVars->nutRef());
226  volTensorField gradU(fvc::grad(U));
227 
228  // Explicitly correct the boundary gradient to get rid of the
229  // tangential component
230  forAll(mesh_.boundary(), patchI)
231  {
232  const fvPatch& patch = mesh_.boundary()[patchI];
233  if (isA<wallFvPatch>(patch))
234  {
235  tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
236  const vectorField& nf = tnf();
237  gradU.boundaryFieldRef()[patchI] =
238  nf * U.boundaryField()[patchI].snGrad();
239  }
240  }
241 
242  volTensorField stress(nuEff*(gradU + T(gradU)));
243 
244  stressXPtr_().replace(0, stress.component(0));
245  stressXPtr_().replace(1, stress.component(1));
246  stressXPtr_().replace(2, stress.component(2));
247 
248  stressYPtr_().replace(0, stress.component(3));
249  stressYPtr_().replace(1, stress.component(4));
250  stressYPtr_().replace(2, stress.component(5));
251 
252  stressZPtr_().replace(0, stress.component(6));
253  stressZPtr_().replace(1, stress.component(7));
254  stressZPtr_().replace(2, stress.component(8));
255 
256  volTensorField gradStressX(fvc::grad(stressXPtr_()));
257  volTensorField gradStressY(fvc::grad(stressYPtr_()));
258  volTensorField gradStressZ(fvc::grad(stressZPtr_()));
259 
260  volVectorField gradp(fvc::grad(p));
261 
262  for (const label patchI : momentPatches_)
263  {
264  const fvPatch& patch = mesh_.boundary()[patchI];
265  tmp<vectorField> tnf = patch.nf();
266  const vectorField& nf = tnf();
267  vectorField dx(patch.Cf() - rotationCentre_);
268  vectorField aux(momentDirection_^dx);
269  bdxdbMultPtr_()[patchI] =
270  (
271  (
272  (
273  -(aux.component(0) * gradStressX.boundaryField()[patchI])
274  -(aux.component(1) * gradStressY.boundaryField()[patchI])
275  -(aux.component(2) * gradStressZ.boundaryField()[patchI])
276  ) & nf
277  )
278  + (momentDirection_ & (dx^nf))*gradp.boundaryField()[patchI]
279  )
280  *invDenom_*rhoInf_;
281  }
282 }
283 
284 
286 {
287  const volScalarField& p = vars_.p();
288 
289  for (const label patchI : momentPatches_)
290  {
291  const fvPatch& patch = mesh_.boundary()[patchI];
292  tmp<vectorField> tnf = patch.nf();
293  const vectorField& nf = tnf();
294  const vectorField dx(patch.Cf() - rotationCentre_);
295  const vectorField force
296  (
297  rhoInf_
298  *(
299  ((p.boundaryField()[patchI]*nf)
300  + (devReff_.boundaryField()[patchI] & nf))
301  )
302  );
303  bdxdbDirectMultPtr_()[patchI] =
304  (force^momentDirection_)*invDenom_*rhoInf_;
305  }
306 }
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 } // End namespace objectives
312 } // End namespace Foam
313 
314 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::objectiveIncompressible::vars_
const incompressibleVars & vars_
Definition: objectiveIncompressible.H:62
Foam::objectives::addToRunTimeSelectionTable
addToRunTimeSelectionTable(objectiveIncompressible, objectiveForce, dictionary)
wallFvPatch.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::objectives::objectiveMoment::update_dxdbMultiplier
void update_dxdbMultiplier()
Update delta(x)/delta b multiplier.
Definition: objectiveMoment.C:217
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::singlePhaseTransportModel
A simple single-phase transport model based on viscosityModel.
Definition: singlePhaseTransportModel.H:57
objectiveMoment.H
Foam::incompressibleVars::p
const volScalarField & p() const
Return const reference to pressure.
Definition: incompressibleVars.C:305
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::incompressibleVars::RASModelVariables
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Definition: incompressibleVars.C:444
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::objectives::defineTypeNameAndDebug
defineTypeNameAndDebug(objectiveForce, 0)
Foam::Field< vector >
Foam::objective::J_
scalar J_
Objective function value and weight.
Definition: objective.H:77
createZeroField.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::incompressibleVars::laminarTransport
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
Definition: incompressibleVars.C:418
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::incompressibleVars::U
const volVectorField & U() const
Return const reference to velocity.
Definition: incompressibleVars.C:331
Foam::objectiveIncompressible
Abstract base class for objective functions in incompressible flows.
Definition: objectiveIncompressible.H:54
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::objectives::objectiveMoment::update_dSdbMultiplier
void update_dSdbMultiplier()
Update delta(n dS)/delta b multiplier.
Definition: objectiveMoment.C:191
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::objectives::objectiveMoment::objectiveMoment
objectiveMoment(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
Definition: objectiveMoment.C:57
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressibleVars::turbulence
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Definition: incompressibleVars.C:431
Foam::objectives::objectiveMoment::update_meanValues
void update_meanValues()
Update mean drag and lift values.
Definition: objectiveMoment.C:166
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
Foam::objectives::objectiveMoment::update_dxdbDirectMultiplier
void update_dxdbDirectMultiplier()
Definition: objectiveMoment.C:285
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::objective::computeMeanFields_
bool computeMeanFields_
Definition: objective.H:72
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::singlePhaseTransportModel::nu
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: singlePhaseTransportModel.C:72
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::objectives::objectiveMoment::J
scalar J()
Return the objective function value.
Definition: objectiveMoment.C:128
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::objectives::objectiveMoment::update_boundarydJdp
void update_boundarydJdp()
Update values to be added to the adjoint wall velocity.
Definition: objectiveMoment.C:180
Foam::objective::mesh_
const fvMesh & mesh_
Definition: objective.H:67
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::objectiveIncompressible::bdJdpPtr_
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
Definition: objectiveIncompressible.H:88
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::incompressibleVars::pInst
const volScalarField & pInst() const
Return const reference to pressure.
Definition: incompressibleVars.C:382