Loading...
Searching...
No Matches
adjointFarFieldVelocityFvPatchVectorField.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
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
37(
38 const fvPatch& p,
53 const fvPatchFieldMapper& mapper
55:
56 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
58{}
59
60
63(
64 const fvPatch& p,
66 const dictionary& dict
67)
68:
69 fixedValueFvPatchVectorField(p, iF),
70 adjointVectorBoundaryCondition(p, iF, dict.get<word>("solverName"))
71{
72 this->readValueEntry(dict, IOobjectOption::MUST_READ);
73}
74
75
78(
81)
82:
83 fixedValueFvPatchVectorField(pivpvf, iF),
85{}
86
87
88// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89
91{
92 if (updated())
93 {
94 return;
95 }
96
97 const scalarField& faceMag = patch().magSf();
98 const vectorField nf(patch().nf());
99
100 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
101
102 scalarField phiOverSurf(phip/faceMag);
103
104 // Ua patch adjacent
105 vectorField Uac(this->patchInternalField());
106
107 // Tangent component of internalField
108 vectorField Uac_t(Uac - nf*(Uac & nf));
109
110 // Inverse distance
111 const scalarField& delta = patch().deltaCoeffs();
112
113 // Objective function and other explicit contributions for
114 // zeroGradient boundaries
115 tmp<vectorField> tsourceVelocity =
116 boundaryContrPtr_->tangentVelocitySource();
117 vectorField& sourceVelocity = tsourceVelocity.ref();
118
119 // Objective function contribution for fixedValue boundaries
120 tmp<vectorField> tsourcePressure =
121 boundaryContrPtr_->normalVelocitySource();
122 vectorField& sourcePressure = tsourcePressure.ref();
123
124 // Momentum diffusion coefficient
125 tmp<scalarField> tmomentumDiffusion =
126 boundaryContrPtr_->momentumDiffusion();
127 scalarField& momentumDiffusion = tmomentumDiffusion.ref();
128
129 scalarField denom(phiOverSurf + momentumDiffusion*delta);
130
131 operator==
132 (
133 // Inlet
134 - neg(phip)*sourcePressure
135
136 // Outlet
137 + pos(phip)
138 *((Uac&nf)*nf + (Uac_t*(momentumDiffusion*delta) - sourceVelocity)/denom)
139 );
140
141 fixedValueFvPatchVectorField::updateCoeffs();
142}
143
144
147(
148 const tmp<scalarField>&
149) const
150{
151 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
152
153 // For fixedValue U patches
154 return tmp<Field<vector>>
155 (
156 new Field<vector>
157 (
159 )
160 );
161}
162
163
166(
167 const tmp<scalarField>&
168) const
169{
170 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
171
172 // For zeroGradient U patches
173 return tmp<Field<vector>>
174 (
175 new Field<vector>
177 pos(phip)*(*this)
178 )
179 );
180}
181
182
184{
186 os.writeEntry("solverName", adjointSolverName_);
189
190
191// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192
193namespace Foam
194{
196 (
198 adjointFarFieldVelocityFvPatchVectorField
199 );
200}
201
202// ************************************************************************* //
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...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
virtual tmp< Field< vector > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
adjointFarFieldVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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,...
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
adjointBoundaryCondition< vector > adjointVectorBoundaryCondition
fvPatchField< vector > fvPatchVectorField
dictionary dict