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 {
75  (
76  vectorField("value", dict, p.size())
77  );
78 }
79 
80 
83 (
86 )
87 :
88  fixedValueFvPatchVectorField(pivpvf, iF),
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 (
97  fvMatrix<vector>& matrix
98 )
99 {
101  (
103  "adjointOutletVelocityFluxFvPatchVectorField::manipulateMatrix"
104  );
105  vectorField& source = matrix.source();
106  const vectorField& Sf = patch().Sf();
107  const labelList& faceCells = patch().faceCells();
108  const scalarField& magSf = patch().magSf();
109  tmp<vectorField> tvelocitySource(boundaryContrPtr_->velocitySource());
110  const vectorField& velocitySource = tvelocitySource();
111  const fvPatchScalarField& pab = boundaryContrPtr_->pab();
112  const word& fieldName = internalField().name();
113  tmp<tensorField> tgradUab(computePatchGrad<vector>(fieldName));
114  const tensorField& gradUab = tgradUab();
115 
116  // Momentum diffusion coefficient
117  tmp<scalarField> tmomentumDiffusion(boundaryContrPtr_->momentumDiffusion());
118  const scalarField& momentumDiffusion = tmomentumDiffusion();
119 
120  vectorField explDiffusiveFlux
121  (
122  -momentumDiffusion*(gradUab - sphericalTensor::oneThirdI*tr(gradUab))
123  & Sf
124  );
125 
126 // const fvPatchVectorField& Ub = boundaryContrPtr_->Ub();
127 // const fvPatchVectorField& Uab = boundaryContrPtr_->Uab();
128 // vectorField cmFormTerm = (Ub & Uab)*Sf;
129 
130  forAll(faceCells, fI)
131  {
132  const label cI = faceCells[fI];
133  // Contributions from the convection and diffusion term (except from
134  // the transpose part) will be canceled out through the value and
135  // gradient coeffs. The pressure flux will be inserted later through
136  // grad(pa), so it must be canceled out here. Once the typical fluxes
137  // have been canceled out, add the objective flux. velocitySource
138  // includes also fluxes from the adjoint turbulence-dependent terms
139  // found in the adjoint momentum equations.
140  source[cI] +=
141  pab[fI]*Sf[fI]
142 // - cmFormTerm[fI]
143  + explDiffusiveFlux[fI]
144  - velocitySource[fI]*magSf[fI];
145  }
146 }
147 
148 
150 {
151  if (updated())
152  {
153  return;
154  }
155 
156  tmp<vectorField> tnf = patch().nf();
157  const vectorField& nf = tnf();
158 
159  // vectorField Ua = (patchInternalField() & nf) * nf;
160  const fvsPatchScalarField& phia = boundaryContrPtr_->phiab();
161  vectorField Ua((phia/patch().magSf())*nf);
162 
163  operator==(Ua);
165  fixedValueFvPatchVectorField::updateCoeffs();
166 }
167 
168 
171 (
172  const tmp<scalarField>&
173 ) const
174 {
175  return tmp<Field<vector>>::New(this->size(), Zero);
176 }
177 
178 
181 (
182  const tmp<scalarField>&
183 ) const
184 {
185  return tmp<Field<vector>>::New(this->size(), Zero);
186 }
187 
188 
192 {
193  return tmp<Field<vector>>::New(this->size(), Zero);
194 }
195 
196 
200 {
201  return tmp<Field<vector>>::New(this->size(), Zero);
202 }
203 
204 
206 (
207  Ostream& os
208 ) const
209 {
211  writeEntry("value", os);
212  os.writeEntry("solverName", adjointSolverName_);
213 }
214 
215 
216 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
217 
218 void Foam::adjointOutletVelocityFluxFvPatchVectorField::operator=
219 (
220  const fvPatchField<vector>& pvf
221 )
222 {
223  fvPatchField<vector>::operator=(patch().nf()*(patch().nf() & pvf));
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 namespace Foam
230 {
232  (
234  adjointOutletVelocityFluxFvPatchVectorField
235  );
236 }
237 
238 // ************************************************************************* //
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
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
fvPatchField< vector > fvPatchVectorField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
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:68
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
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
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
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:413
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:55
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:534
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:352
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:157