adjointWallVelocityFvPatchVectorField.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 PCOpt/NTUA
9  Copyright (C) 2013-2023 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 
32 #include "fvMatrix.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  fixedValueFvPatchVectorField(p, iF),
46  kappa_(0.41),
47  E_(9.8)
48 {}
49 
50 
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
61  adjointVectorBoundaryCondition(p, iF, ptf.adjointSolverName_),
62  kappa_(ptf.kappa_),
63  E_(ptf.E_)
64 {}
65 
66 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
75  fixedValueFvPatchVectorField(p, iF),
76  adjointVectorBoundaryCondition(p, iF, dict.get<word>("solverName")),
77  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
78  E_(dict.getOrDefault<scalar>("E", 9.8))
79 {
80  this->readValueEntry(dict, IOobjectOption::MUST_READ);
81 }
82 
83 
86 (
89 )
90 :
91  fixedValueFvPatchVectorField(pivpvf, iF),
93  kappa_(pivpvf.kappa_),
94  E_(pivpvf.E_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 (
102  fvMatrix<vector>& matrix
103 )
104 {
106  (
108  "adjointWallVelocityFvPatchVectorField::manipulateMatrix"
109  );
110  // Grab ref to the diagonal matrix
111  vectorField& source = matrix.source();
112 
113  // Define boundary condition type
115  SAwallFunctionPatchField;
116 
117  tmp<fvPatchScalarField> nutPatch(boundaryContrPtr_->turbulentDiffusivity());
118  if
119  (
120  isA<SAwallFunctionPatchField>(nutPatch())
121  && patch().size() != 0
122  )
123  {
124  const tmp<vectorField> tnf = patch().nf();
125  const vectorField& nf = tnf();
126  const scalarField& magSf = patch().magSf();
127 
128  const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
129  const fvPatchField<vector>& Uap = *this;
130 
131  const vectorField Uc(Up.patchInternalField());
132  const vectorField Uc_t(Uc - (Uc & nf)*nf);
133 
134  // By convention, tf has the direction of the tangent PRIMAL velocity
135  // at the first cell off the wall
136  const vectorField tf(Uc_t/mag(Uc_t));
137 
138  tmp<scalarField> tnuw = boundaryContrPtr_->momentumDiffusion();
139  const scalarField& nuw = tnuw();
140  tmp<scalarField> tnu = boundaryContrPtr_->laminarDiffusivity();
141  const scalarField& nu = tnu();
142  tmp<scalarField> tyC = boundaryContrPtr_->wallDistance();
143  const scalarField& yC = tyC();
144 
145  const scalarField magGradU(mag(Up.snGrad()));
146  const scalarField vtau(sqrt(nuw*magGradU));
147  const scalarField uPlus(mag(Uc)/vtau);
148  const scalarField yPlus(yC*vtau/nu);
149  const scalarField kUu(min(kappa_*uPlus, scalar(50)));
150  const scalarField auxA((kappa_/E_)*(exp(kUu)-1 - kUu - 0.5*kUu*kUu));
151  const scalarField auxB(-(1 + auxA)/(yPlus + uPlus*(1 + auxA)));
152 
153  // Tangent components are according to tf
154  tmp<vectorField> tsource = boundaryContrPtr_->normalVelocitySource();
155  const scalarField rt(tsource() & tf);
156  const scalarField Uap_t(Uap & tf);
157 
158  forAll(Up, faceI)
159  {
160  label cellI = patch().faceCells()[faceI];
161  source[cellI] +=
162  2*auxB[faceI]*vtau[faceI]*((rt[faceI] + Uap_t[faceI]))
163  *(Uc[faceI]/mag(Uc[faceI]))*magSf[faceI];
164  }
165  }
166 }
167 
168 
170 {
171  if (updated())
172  {
173  return;
174  }
175 
176  const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
177 
178  // Patch geometry
179  tmp<vectorField> tnf = patch().nf();
180  const vectorField& nf = tnf();
181 
182  // Internal fields
183  vectorField Uac(this->patchInternalField());
184  vectorField Uc(Up.patchInternalField());
185 
186  // Tangent vector based on the direction of Vc
187  vectorField Uc_t(Uc - (Uc & nf)*nf);
188  vectorField tf1(Uc_t/mag(Uc_t));
189 
190  // Tangent vector as the cross product of tf1 x nf
191  vectorField tf2((tf1 ^ nf)/mag(tf1 ^ nf));
192 
193  // Normal adjoint component comes from the objective function
194  tmp<vectorField> tsource = boundaryContrPtr_->normalVelocitySource();
195  vectorField Uan(-(tsource() & nf)*nf);
196 
197  // Tangential adjoint velocity in the t1 direction depends on the primal
198  // wall function used
199  vectorField Uap_t1(patch().size(), Zero);
201  SAwallFunctionPatchField;
202 
203  tmp<fvPatchScalarField> nutb(boundaryContrPtr_->turbulentDiffusivity());
204  if (isA<SAwallFunctionPatchField>(nutb()))
205  {
206  Uap_t1 = (Uac & tf1)*tf1;
207  // leaving out second term for now
208  //- (1./delta)*((gradUaC & nf) & tf1)*tf1;
209  }
210  else
211  {
212  Uap_t1 = - (tsource() & tf1)*tf1;
213  }
214 
215  // Adjoint velocity in the t2 direction
216  vectorField Uap_t2(-(tsource() & tf2)*tf2);
217 
218  operator==(Uan + Uap_t1 + Uap_t2);
220  fixedValueFvPatchVectorField::updateCoeffs();
221 }
222 
223 
225 {
227  os.writeEntry("kappa", kappa_);
228  os.writeEntry("E", E_);
229  os.writeEntry("solverName", adjointSolverName_);
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 namespace Foam
237 {
239  (
241  adjointWallVelocityFvPatchVectorField
242  );
243 }
244 
245 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Adjoint wall velocity boundary condition. If nutUSpaldingWallFunction is employed in the flow solutio...
scalar uPlus
dictionary dict
fvPatchField< vector > fvPatchVectorField
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:213
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
#define addProfiling(Name, Descr)
Define profiling trigger with specified name and description string.
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on velocity (i.e. U). Using Spalding&#39;s law gives a continuous nut profile to the wall.
dimensionedScalar exp(const dimensionedScalar &ds)
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Engine to manage contributions of the objective functions to the adjoint boundary conditions...
virtual void manipulateMatrix(fvMatrix< vector > &matrix)
In case of High-Re runs based on the nutUSpaldingWallFunction add source terms in the first cell cent...
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:221
OBJstream os(runTime.globalPath()/outputName)
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
adjointWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
scalar yPlus
Base class for solution control classes.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
volScalarField & nu
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127