objectiveForce.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-2023, 2022 PCOpt/NTUA
9  Copyright (C) 2013-2023, 2022 FOSS GP
10  Copyright (C) 2019-2023 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 "objectiveForce.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 namespace objectives
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(objectiveForce, 0);
45 (
46  objectiveIncompressible,
47  objectiveForce,
49 );
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const fvMesh& mesh,
57  const dictionary& dict,
58  const word& adjointSolverName,
59  const word& primalSolverName
60 )
61 :
62  objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
63  forcePatches_
64  (
65  mesh_.boundaryMesh().patchSet
66  (
67  dict.get<wordRes>("patches")
68  ).sortedToc()
69  ),
70  forceDirection_(dict.get<vector>("direction")),
71  Aref_(dict.get<scalar>("Aref")),
72  rhoInf_(dict.get<scalar>("rhoInf")),
73  UInf_(dict.get<scalar>("UInf"))
74 {
75  // Sanity check and print info
76  if (forcePatches_.empty())
77  {
79  << "No valid patch name on which to minimize " << type() << endl
80  << exit(FatalError);
81  }
82  if (debug)
83  {
84  Info<< "Minimizing " << type() << " in patches:" << endl;
85  for (const label patchI : forcePatches_)
86  {
87  Info<< "\t " << mesh_.boundary()[patchI].name() << endl;
88  }
89  }
90 
91  // Allocate boundary field pointers
92  bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
93  bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
94  bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
95  bdJdnutPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
96  bdJdGradUPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 scalar objectiveForce::J()
103 {
104  vector pressureForce(Zero);
105  vector viscousForce(Zero);
106  vector cumulativeForce(Zero);
107 
108 
109  const volScalarField& p = vars_.pInst();
112 
113  volSymmTensorField devReff(turbulence->devReff());
114 
115  for (const label patchI : forcePatches_)
116  {
117  const vectorField& Sf = mesh_.Sf().boundaryField()[patchI];
118  pressureForce += gSum(Sf*p.boundaryField()[patchI]);
119  viscousForce += gSum(devReff.boundaryField()[patchI] & Sf);
120  }
121 
122  cumulativeForce = pressureForce + viscousForce;
123 
124  scalar force = cumulativeForce & forceDirection_;
125 
126  // Intentionally not using denom - derived might implement virtual denom()
127  // function differently
128  scalar Cforce = force/(0.5*UInf_*UInf_*Aref_);
129 
130  DebugInfo
131  << "Force|Coeff " << force << "|" << Cforce << endl;
133  J_ = Cforce;
134 
135  return Cforce;
136 }
137 
138 
140 {
141  for (const label patchI : forcePatches_)
142  {
143  bdJdpPtr_()[patchI] = forceDirection_/denom();
144  }
145 }
146 
147 
149 {
150  // Compute contributions with mean fields, if present
151  const volScalarField& p = vars_.p();
152  const volVectorField& U = vars_.U();
155  const singlePhaseTransportModel& lamTransp = vars_.laminarTransport();
156 
157  tmp<volSymmTensorField> tdevReff = turbVars->devReff(lamTransp, U);
158  const volSymmTensorField& devReff = tdevReff();
159 
160  for (const label patchI : forcePatches_)
161  {
162  bdSdbMultPtr_()[patchI] =
163  (
164  (
165  forceDirection_& devReff.boundaryField()[patchI]
166  )
167  + (forceDirection_)*p.boundaryField()[patchI]
168  )
169  /denom();
170  }
171 }
172 
173 
175 {
176  const volScalarField& p = vars_.p();
177  const volVectorField& U = vars_.U();
178 
180  turbVars = vars_.RASModelVariables();
181  const singlePhaseTransportModel& lamTransp = vars_.laminarTransport();
182 
183  // We only need to modify the boundaryField of gradU locally.
184  // If grad(U) is cached then
185  // a. The .ref() call fails since the tmp is initialised from a
186  // const ref
187  // b. we would be changing grad(U) for all other places in the code
188  // that need it
189  // So, always allocate new memory and avoid registering the new field
190  tmp<volTensorField> tgradU =
191  volTensorField::New("gradULocal", fvc::grad(U));
192  volTensorField& gradU = tgradU.ref();
193  volTensorField::Boundary& gradUbf = gradU.boundaryFieldRef();
194 
195  // Explicitly correct the boundary gradient to get rid of
196  // the tangential component
197  forAll(mesh_.boundary(), patchI)
198  {
199  const fvPatch& patch = mesh_.boundary()[patchI];
200  if (isA<wallFvPatch>(patch))
201  {
202  tmp<vectorField> tnf = patch.nf();
203  gradUbf[patchI] = tnf*U.boundaryField()[patchI].snGrad();
204  }
205  }
206 
207  // Term coming from gradp
208  tmp<volVectorField> tgradp(fvc::grad(p));
209  const volVectorField& gradp = tgradp.cref();
210  for (const label patchI : forcePatches_)
211  {
212  bdxdbMultPtr_()[patchI] =
213  (forceDirection_ & mesh_.boundary()[patchI].nf())
214  *gradp.boundaryField()[patchI]/denom();
215  }
216  tgradp.clear();
217 
218  // Term coming from stresses
219  tmp<volScalarField> tnuEff = lamTransp.nu() + turbVars->nut();
220  tmp<volSymmTensorField> tstress = tnuEff*twoSymm(tgradU);
221  const volSymmTensorField& stress = tstress.cref();
222  autoPtr<volVectorField> ptemp
223  (Foam::createZeroFieldPtr<vector>( mesh_, "temp", sqr(dimVelocity)));
224  volVectorField& temp = ptemp.ref();
225 
226  for (label idir = 0; idir < pTraits<vector>::nComponents; ++idir)
227  {
228  unzipRow(stress, idir, temp);
229  volTensorField gradStressDir(fvc::grad(temp));
230  for (const label patchI : forcePatches_)
231  {
232  const fvPatch& patch = mesh_.boundary()[patchI];
233  tmp<vectorField> tnf = patch.nf();
234  bdxdbMultPtr_()[patchI] -=
235  forceDirection_.component(idir)
236  *(gradStressDir.boundaryField()[patchI] & tnf)/denom();
237  }
238  }
239 }
240 
241 
243 {
244  const volVectorField& U = vars_.U();
246 
247  for (const label patchI : forcePatches_)
248  {
249  const fvPatch& patch = mesh_.boundary()[patchI];
250  tmp<vectorField> tnf = patch.nf();
251  bdJdnutPtr_()[patchI] =
252  - ((devGradU.boundaryField()[patchI] & forceDirection_) & tnf)
253  /denom();
254  }
255 }
256 
257 
259 {
262  const singlePhaseTransportModel& lamTransp = vars_.laminarTransport();
263  volScalarField nuEff(lamTransp.nu() + turbVars->nut());
264  for (const label patchI : forcePatches_)
265  {
266  const fvPatch& patch = mesh_.boundary()[patchI];
267  const vectorField& Sf = patch.Sf();
268  bdJdGradUPtr_()[patchI] =
269  - nuEff.boundaryField()[patchI]
271  }
272 }
273 
275 scalar objectiveForce::denom() const
276 {
277  return 0.5*UInf_*UInf_*Aref_;
278 }
279 
280 
282 {
283  return forceDirection_;
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace objectives
290 } // End namespace Foam
291 
292 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
const volScalarField & pInst() const
Return const reference to pressure.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:171
const surfaceVectorField & Sf() const
Return cell face area vectors.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void update_boundarydJdGradU()
Update dJ/dGradU multiplier.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
objectiveForce(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
virtual void update_boundarydJdnut()
Update dJ/dnut multiplier.
void unzipRow(const FieldField< Field, SymmTensor< Cmpt >> &input, const direction idx, FieldField< Field, Vector< Cmpt >> &result)
Extract a symmTensor field field row (x,y,z) == (0,1,2)
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
virtual void update_boundarydJdp()
Update values to be added to the adjoint wall velocity.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
virtual void update_dSdbMultiplier()
Update delta(n dS)/delta b multiplier.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
dynamicFvMesh & mesh
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Type gSum(const FieldField< Field, Type > &f)
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
A class for handling words, derived from Foam::string.
Definition: word.H:63
defineTypeNameAndDebug(objectivePartialVolume, 1)
const volScalarField & p() const
Return const reference to pressure.
virtual scalar denom() const
Return denominator, without density.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const fvMesh & mesh_
Definition: objective.H:63
scalar J_
Objective function value and weight.
Definition: objective.H:76
#define DebugInfo
Report an information message using Foam::Info.
const volVectorField & U() const
Return const reference to velocity.
static tmp< GeometricField< tensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< tensor >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTable.H:337
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
virtual const vector & forceDirection() const
Return force direction.
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
const incompressibleVars & vars_
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void update_dxdbMultiplier()
Update delta(x)/delta b multiplier.
const std::string patch
OpenFOAM patch number as a std::string.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
virtual scalar J()
Return the objective function value.
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:161
A simple single-phase transport model based on viscosityModel.
Abstract base class for objective functions in incompressible flows.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity