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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 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),
45  adjointVectorBoundaryCondition(p, iF, word::null),
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),
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 {
81  (
82  vectorField("value", dict, p.size())
83  );
84 }
85 
86 
89 (
92 )
93 :
94  fixedValueFvPatchVectorField(pivpvf, iF),
96  kappa_(pivpvf.kappa_),
97  E_(pivpvf.E_)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 (
105  fvMatrix<vector>& matrix
106 )
107 {
108  // Grab ref to the diagonal matrix
109  vectorField& source = matrix.source();
110 
111  // Define boundary condition type
113  SAwallFunctionPatchField;
114 
115  if
116  (
117  isA<SAwallFunctionPatchField>(boundaryContrPtr_->turbulentDiffusivity())
118  && patch().size() != 0
119  )
120  {
121  const tmp<vectorField> tnf = patch().nf();
122  const vectorField& nf = tnf();
123  const scalarField& magSf = patch().magSf();
124 
125  const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
126  const fvPatchField<vector>& Uap = *this;
127 
128  const vectorField Uc(Up.patchInternalField());
129  const vectorField Uc_t(Uc - (Uc & nf)*nf);
130 
131  // By convention, tf has the direction of the tangent PRIMAL velocity
132  // at the first cell off the wall
133  const vectorField tf(Uc_t/mag(Uc_t));
134 
135  tmp<scalarField> tnuw = boundaryContrPtr_->momentumDiffusion();
136  const scalarField& nuw = tnuw();
137  tmp<scalarField> tnu = boundaryContrPtr_->laminarDiffusivity();
138  const scalarField& nu = tnu();
139  tmp<scalarField> tyC = boundaryContrPtr_->wallDistance();
140  const scalarField& yC = tyC();
141 
142  const scalarField magGradU(mag(Up.snGrad()));
143  const scalarField vtau(sqrt(nuw*magGradU));
144  const scalarField uPlus(mag(Uc)/vtau);
145  const scalarField yPlus(yC*vtau/nu);
146  const scalarField kUu(min(kappa_*uPlus, scalar(50)));
147  const scalarField auxA((kappa_/E_)*(exp(kUu)-1 - kUu - 0.5*kUu*kUu));
148  const scalarField auxB(-(1 + auxA)/(yPlus + uPlus*(1 + auxA)));
149 
150  // Tangent components are according to tf
151  tmp<vectorField> tsource = boundaryContrPtr_->normalVelocitySource();
152  const scalarField rt(tsource() & tf);
153  const scalarField Uap_t(Uap & tf);
154 
155  forAll(Up, faceI)
156  {
157  label cellI = patch().faceCells()[faceI];
158  source[cellI] +=
159  2*auxB[faceI]*vtau[faceI]*((rt[faceI] + Uap_t[faceI]))
160  *(Uc[faceI]/mag(Uc[faceI]))*magSf[faceI];
161  }
162  }
163 }
164 
165 
167 {
168  if (updated())
169  {
170  return;
171  }
172 
173  const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
174 
175  // Patch geometry
176  tmp<vectorField> tnf = patch().nf();
177  const vectorField& nf = tnf();
178 
179  // Internal fields
180  vectorField Uac(this->patchInternalField());
182 
183  // Tangent vector based on the direction of Vc
184  vectorField Uc_t(Uc - (Uc & nf)*nf);
185  vectorField tf1(Uc_t/mag(Uc_t));
186 
187  // Tangent vector as the cross product of tf1 x nf
188  vectorField tf2((tf1 ^ nf)/mag(tf1 ^ nf));
189 
190  // Normal adjoint component comes from the objective function
191  tmp<vectorField> tsource = boundaryContrPtr_->normalVelocitySource();
192  vectorField Uan(-(tsource() & nf)*nf);
193 
194  // Tangential adjoint velocity in the t1 direction depends on the primal
195  // wall function used
196  vectorField Uap_t1(patch().size(), Zero);
198  SAwallFunctionPatchField;
199 
200  const fvPatchScalarField& nutb = boundaryContrPtr_->turbulentDiffusivity();
201  if (isA<SAwallFunctionPatchField>(nutb))
202  {
203  Uap_t1 = (Uac & tf1)*tf1;
204  // leaving out second term for now
205  //- (1./delta)*((gradUaC & nf) & tf1)*tf1;
206  }
207  else
208  {
209  Uap_t1 = - (tsource() & tf1)*tf1;
210  }
211 
212  // Adjoint velocity in the t2 direction
213  vectorField Uap_t2(-(tsource() & tf2)*tf2);
214 
215  operator==(Uan + Uap_t1 + Uap_t2);
216 
217  fixedValueFvPatchVectorField::updateCoeffs();
218 }
219 
220 
222 {
224  writeEntry("value", os);
225  os.writeEntry("kappa", kappa_);
226  os.writeEntry("E", E_);
227  os.writeEntry("solverName", adjointSolverName_);
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 namespace Foam
234 {
236  (
239  );
240 }
241 
242 // ************************************************************************* //
Foam::adjointWallVelocityFvPatchVectorField::adjointWallVelocityFvPatchVectorField
adjointWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: adjointWallVelocityFvPatchVectorField.C:39
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::adjointWallVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: adjointWallVelocityFvPatchVectorField.C:221
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
uPlus
scalar uPlus
Definition: evaluateNearWall.H:18
fvMatrix.H
Foam::adjointBoundaryCondition::boundaryContrPtr_
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Definition: adjointBoundaryCondition.H:73
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::adjointWallVelocityFvPatchVectorField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< vector > &matrix)
Definition: adjointWallVelocityFvPatchVectorField.C:104
Foam::adjointBoundaryCondition::adjointSolverName_
word adjointSolverName_
adjointSolver name corresponding to field
Definition: adjointBoundaryCondition.H:65
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::adjointWallVelocityFvPatchVectorField
Adjoint wall velocity boundary condition. If nutUSpaldingWallFunction is employed in the flow solutio...
Definition: adjointWallVelocityFvPatchVectorField.H:68
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::nutUSpaldingWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutUSpaldingWallFunctionFvPatchScalarField.H:145
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::adjointVectorBoundaryCondition
adjointBoundaryCondition< vector > adjointVectorBoundaryCondition
Definition: adjointBoundaryConditionsFwd.H:44
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
adjointWallVelocityFvPatchVectorField.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:445
nutUSpaldingWallFunctionFvPatchScalarField.H
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::adjointWallVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: adjointWallVelocityFvPatchVectorField.C:166
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54