Loading...
Searching...
No Matches
adjointOutletVelocityFvPatchVectorField.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 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
32#include "emptyFvPatch.H"
33
34// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
35
36void Foam::adjointOutletVelocityFvPatchVectorField::assignBoundaryValue()
37{
38 const scalarField& magSf = patch().magSf();
39 tmp<vectorField> tnf(patch().nf());
40 const vectorField& nf = tnf();
41
42 // Primal normal velocity
43 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
44 const scalarField phiOverSurf(phip/magSf);
45
46 // Ua patch adjacent
47 vectorField Uac(this->patchInternalField());
48
49 // Tangent component of internalField
50 vectorField Uac_t(Uac - nf*(Uac & nf));
51
52 // Adjoint normal velocity
53 const fvsPatchField<scalar>& phiab = boundaryContrPtr_->phiab();
54
55 // Inverse distance
56 const scalarField& delta = patch().deltaCoeffs();
57
58 // Objective function and other explicit contributions
59 tmp<vectorField> tsource(boundaryContrPtr_->tangentVelocitySource());
60 const vectorField& source = tsource();
61
62 // Momentum diffusion coefficient
63 tmp<scalarField> tmomentumDiffusion
64 (
65 boundaryContrPtr_->momentumDiffusion()
66 );
67 const scalarField& momentumDiffusion = tmomentumDiffusion();
68
69 // Part of the diffusive flux related to div(nuEff*dev(grad(Ua).T()))
70 const word& fieldName = internalField().name();
71 tmp<tensorField> tgradUaf(computePatchGrad<vector>(fieldName));
72 const tensorField& gradUaf = tgradUaf();
73 const vectorField explDiffusiveFlux
74 (
75 momentumDiffusion
76 *(gradUaf - sphericalTensor::oneThirdI*tr(gradUaf)) & nf
77 );
78 const vectorField explDiffusiveFlux_t
79 (
80 explDiffusiveFlux - (explDiffusiveFlux & nf)*nf
81 );
82
83 // Auxiliary quantities
84 scalarField nd(momentumDiffusion*delta);
85 // Denominator. Susceptible to zero values in case of back flow
86 // Should use adjointOutletVelocityFlux in such cases.
87 scalarField denom(phiOverSurf + nd);
88
89 vectorField Uat((nd*Uac_t - explDiffusiveFlux_t - source)/denom);
91 operator==((phiab/magSf)*nf + Uat);
92}
93
94
95// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96
99(
100 const fvPatch& p,
103:
104 fixedValueFvPatchVectorField(p, iF),
106{}
107
108
111(
113 const fvPatch& p,
115 const fvPatchFieldMapper& mapper
117:
118 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
120{}
121
122
125(
126 const fvPatch& p,
128 const dictionary& dict
129)
130:
131 fixedValueFvPatchVectorField(p, iF),
132 adjointVectorBoundaryCondition(p, iF, dict.get<word>("solverName"))
133{
134 this->readValueEntry(dict, IOobjectOption::MUST_READ);
135}
136
137
140(
143)
144:
145 fixedValueFvPatchVectorField(pivpvf, iF),
147{}
148
149
150// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151
153(
156{
157 assignBoundaryValue();
159}
160
161
163{
165 os.writeEntry("solverName", adjointSolverName_);
167}
168
169
170// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
171
172void Foam::adjointOutletVelocityFvPatchVectorField::operator=
173(
174 const fvPatchField<vector>& pvf
175)
176{
177 fvPatchField<vector>::operator=(patch().nf()*(patch().nf() & pvf));
179
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183namespace Foam
184{
186 (
188 adjointOutletVelocityFvPatchVectorField
189 );
190}
191
192// ************************************************************************* //
scalar delta
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
commsTypes
Communications types.
Definition UPstream.H:81
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
tmp< Field< typename Foam::outerProduct< Foam::vector, Type2 >::type > > computePatchGrad(word name)
Provides the adjoint outlet velocity values (i.e. adjoint velocity in case of a zeroGradient U bounda...
adjointOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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 > &)
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
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.
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.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
adjointBoundaryCondition< vector > adjointVectorBoundaryCondition
fvPatchField< vector > fvPatchVectorField
dictionary dict