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);
90
91 operator==((phiab/magSf)*nf + Uat);
92}
93
94
95// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96
97Foam::adjointOutletVelocityFvPatchVectorField::
98adjointOutletVelocityFvPatchVectorField
99(
100 const fvPatch& p,
102)
103:
104 fixedValueFvPatchVectorField(p, iF),
106{}
107
108
109Foam::adjointOutletVelocityFvPatchVectorField::
110adjointOutletVelocityFvPatchVectorField
111(
113 const fvPatch& p,
115 const fvPatchFieldMapper& mapper
116)
117:
118 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
119 adjointVectorBoundaryCondition(p, iF, ptf.adjointSolverName_)
120{}
121
122
123Foam::adjointOutletVelocityFvPatchVectorField::
124adjointOutletVelocityFvPatchVectorField
125(
126 const fvPatch& p,
128 const dictionary& dict
129)
130:
131 fixedValueFvPatchVectorField(p, iF),
132 adjointVectorBoundaryCondition(p, iF, dict.get<word>("solverName"))
133{
135 (
136 vectorField("value", dict, p.size())
137 );
138}
139
140
141Foam::adjointOutletVelocityFvPatchVectorField::
142adjointOutletVelocityFvPatchVectorField
143(
146)
147:
148 fixedValueFvPatchVectorField(pivpvf, iF),
150{}
151
152
153// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154
156(
158)
159{
160 assignBoundaryValue();
162}
163
164
166{
168 writeEntry("value", os);
169 os.writeEntry("solverName", adjointSolverName_);
170}
171
172
173// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
174
176(
177 const fvPatchField<vector>& pvf
178)
179{
180 fvPatchField<vector>::operator=(patch().nf()*(patch().nf() & pvf));
181}
182
183
184// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185
186namespace Foam
187{
189 (
192 );
193}
194
195// ************************************************************************* //
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...
void evaluate()
Evaluate boundary conditions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
static const SphericalTensor oneThirdI
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Types of communications.
Definition: UPstream.H:67
Base class for solution control classes.
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Provides the adjoint outlet velocity values (i.e. adjoint velocity in case of a zeroGradient U bounda...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
virtual bool write()
Write the output fields.
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...
Definition: fvPatchField.H:82
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:408
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:68
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
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.
dictionary dict