adjointOutletVelocityFluxFvPatchVectorField.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-2021 PCOpt/NTUA
9  Copyright (C) 2013-2021 FOSS GP
10  Copyright (C) 2019 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 
31 #include "emptyFvPatch.H"
32 #include "fvMatrix.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  fixedValueFvPatchVectorField(p, iF),
46 {}
47 
48 
51 (
53  const fvPatch& p,
55  const fvPatchFieldMapper& mapper
56 )
57 :
58  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
59  adjointVectorBoundaryCondition(p, iF, ptf.adjointSolverName_)
60 {}
61 
62 
65 (
66  const fvPatch& p,
68  const dictionary& dict
69 )
70 :
71  fixedValueFvPatchVectorField(p, iF),
72  adjointVectorBoundaryCondition(p, iF, dict.get<word>("solverName"))
73 {
74  this->readValueEntry(dict, IOobjectOption::MUST_READ);
75 }
76 
77 
80 (
83 )
84 :
85  fixedValueFvPatchVectorField(pivpvf, iF),
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 (
94  fvMatrix<vector>& matrix
95 )
96 {
98  (
100  "adjointOutletVelocityFluxFvPatchVectorField::manipulateMatrix"
101  );
102  vectorField& source = matrix.source();
103  const vectorField& Sf = patch().Sf();
104  const labelList& faceCells = patch().faceCells();
105  const scalarField& magSf = patch().magSf();
106  tmp<vectorField> tvelocitySource(boundaryContrPtr_->velocitySource());
107  const vectorField& velocitySource = tvelocitySource();
108  const fvPatchScalarField& pab = boundaryContrPtr_->pab();
109  const word& fieldName = internalField().name();
110  tmp<tensorField> tgradUab(computePatchGrad<vector>(fieldName));
111  const tensorField& gradUab = tgradUab();
112 
113  // Momentum diffusion coefficient
114  tmp<scalarField> tmomentumDiffusion(boundaryContrPtr_->momentumDiffusion());
115  const scalarField& momentumDiffusion = tmomentumDiffusion();
116 
117  vectorField explDiffusiveFlux
118  (
119  -momentumDiffusion*(gradUab - sphericalTensor::oneThirdI*tr(gradUab))
120  & Sf
121  );
122 
123 // const fvPatchVectorField& Ub = boundaryContrPtr_->Ub();
124 // const fvPatchVectorField& Uab = boundaryContrPtr_->Uab();
125 // vectorField cmFormTerm = (Ub & Uab)*Sf;
126 
127  forAll(faceCells, fI)
128  {
129  const label cI = faceCells[fI];
130  // Contributions from the convection and diffusion term (except from
131  // the transpose part) will be canceled out through the value and
132  // gradient coeffs. The pressure flux will be inserted later through
133  // grad(pa), so it must be canceled out here. Once the typical fluxes
134  // have been canceled out, add the objective flux. velocitySource
135  // includes also fluxes from the adjoint turbulence-dependent terms
136  // found in the adjoint momentum equations.
137  source[cI] +=
138  pab[fI]*Sf[fI]
139 // - cmFormTerm[fI]
140  + explDiffusiveFlux[fI]
141  - velocitySource[fI]*magSf[fI];
142  }
143 }
144 
145 
147 {
148  if (updated())
149  {
150  return;
151  }
152 
153  tmp<vectorField> tnf = patch().nf();
154  const vectorField& nf = tnf();
155 
156  // vectorField Ua = (patchInternalField() & nf) * nf;
157  const fvsPatchScalarField& phia = boundaryContrPtr_->phiab();
158  vectorField Ua((phia/patch().magSf())*nf);
159 
160  operator==(Ua);
162  fixedValueFvPatchVectorField::updateCoeffs();
163 }
164 
165 
168 (
169  const tmp<scalarField>&
170 ) const
171 {
172  return tmp<Field<vector>>::New(this->size(), Zero);
173 }
174 
175 
178 (
179  const tmp<scalarField>&
180 ) const
181 {
182  return tmp<Field<vector>>::New(this->size(), Zero);
183 }
184 
185 
189 {
190  return tmp<Field<vector>>::New(this->size(), Zero);
191 }
192 
193 
197 {
198  return tmp<Field<vector>>::New(this->size(), Zero);
199 }
200 
201 
203 (
204  Ostream& os
205 ) const
206 {
208  os.writeEntry("solverName", adjointSolverName_);
210 }
211 
212 
213 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
214 
215 void Foam::adjointOutletVelocityFluxFvPatchVectorField::operator=
216 (
217  const fvPatchField<vector>& pvf
218 )
219 {
220  fvPatchField<vector>::operator=(patch().nf()*(patch().nf() & pvf));
221 }
222 
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 namespace Foam
227 {
229  (
231  adjointOutletVelocityFluxFvPatchVectorField
232  );
233 }
234 
235 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
fvsPatchField< scalar > fvsPatchScalarField
dictionary dict
fvPatchField< vector > fvPatchVectorField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
adjointOutletVelocityFluxFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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 tmp< Field< vector > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the evaluation of the gradient of this patch...
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
static const SphericalTensor oneThirdI
Macros for easy insertion into run-time selection tables.
virtual tmp< Field< vector > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
#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.
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Engine to manage contributions of the objective functions to the adjoint boundary conditions...
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
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void manipulateMatrix(fvMatrix< vector > &matrix)
add source term in the first cells off the wall due to adjoint WF
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual tmp< Field< vector > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the evaluation of the gradient of this patchFi...
OBJstream os(runTime.globalPath()/outputName)
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Base class for solution control classes.
An outlet boundary condition for patches in which the primal flow exhibits recirculation. Adds the contribution of the objective as an adjoint momentum flux directly to the PDEs, without the need to first compute an adjoint outlet velocity, circumventing thus the division with (almost) zero that manifests in case of primal flow recirculation.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
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.
virtual tmp< Field< vector > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
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)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127