Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
12License
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,
55 const fvPatchFieldMapper& mapper
57:
58 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
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 labelUList& 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
161
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
204 Ostream& os
205) const
206{
208 os.writeEntry("solverName", adjointSolverName_);
210}
211
212
213// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
214
215void Foam::adjointOutletVelocityFluxFvPatchVectorField::operator=
216(
217 const fvPatchField<vector>& pvf
218)
219{
220 fvPatchField<vector>::operator=(patch().nf()*(patch().nf() & pvf));
222
223
224// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225
226namespace Foam
227{
229 (
231 adjointOutletVelocityFluxFvPatchVectorField
232 );
233}
234
235// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static const SphericalTensor oneThirdI
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
tmp< Field< typename Foam::outerProduct< Foam::vector, Type2 >::type > > computePatchGrad(word name)
An outlet boundary condition for patches in which the primal flow exhibits recirculation....
virtual tmp< Field< vector > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the evaluation of the gradient of this patchFi...
virtual tmp< Field< vector > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
adjointOutletVelocityFluxFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void manipulateMatrix(fvMatrix< vector > &matrix)
add source term in the first cells off the wall due to adjoint WF
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< Field< vector > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the evaluation of the gradient of this patch...
virtual tmp< Field< vector > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Field< Type > & source() noexcept
Definition fvMatrix.H:535
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
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.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
adjointBoundaryCondition< vector > adjointVectorBoundaryCondition
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< vector > fvPatchVectorField
fvsPatchField< scalar > fvsPatchScalarField
fvPatchField< scalar > fvPatchScalarField
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299