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
35Foam::adjointFarFieldVelocityFvPatchVectorField::
36adjointFarFieldVelocityFvPatchVectorField
37(
38 const fvPatch& p,
40)
41:
42 fixedValueFvPatchVectorField(p, iF),
44{}
45
46
47Foam::adjointFarFieldVelocityFvPatchVectorField::
48adjointFarFieldVelocityFvPatchVectorField
49(
51 const fvPatch& p,
53 const fvPatchFieldMapper& mapper
54)
55:
56 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
57 adjointVectorBoundaryCondition(p, iF, ptf.adjointSolverName_)
58{}
59
60
61Foam::adjointFarFieldVelocityFvPatchVectorField::
62adjointFarFieldVelocityFvPatchVectorField
63(
64 const fvPatch& p,
66 const dictionary& dict
67)
68:
69 fixedValueFvPatchVectorField(p, iF),
70 adjointVectorBoundaryCondition(p, iF, dict.get<word>("solverName"))
71{
73 (
74 vectorField("value", dict, p.size())
75 );
76}
77
78
79Foam::adjointFarFieldVelocityFvPatchVectorField::
80adjointFarFieldVelocityFvPatchVectorField
81(
84)
85:
86 fixedValueFvPatchVectorField(pivpvf, iF),
88{}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
94{
95 if (updated())
96 {
97 return;
98 }
99
100 const scalarField& faceMag = patch().magSf();
101 const vectorField nf(patch().nf());
102
103 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
104
105 scalarField phiOverSurf(phip/faceMag);
106
107 // Ua patch adjacent
108 vectorField Uac(this->patchInternalField());
109
110 // Tangent component of internalField
111 vectorField Uac_t(Uac - nf*(Uac & nf));
112
113 // Inverse distance
114 const scalarField& delta = patch().deltaCoeffs();
115
116 // Objective function and other explicit contributions for
117 // zeroGradient boundaries
118 tmp<vectorField> tsourceVelocity =
119 boundaryContrPtr_->tangentVelocitySource();
120 vectorField& sourceVelocity = tsourceVelocity.ref();
121
122 // Objective function contribution for fixedValue boundaries
123 tmp<vectorField> tsourcePressure =
124 boundaryContrPtr_->normalVelocitySource();
125 vectorField& sourcePressure = tsourcePressure.ref();
126
127 // Momentum diffusion coefficient
128 tmp<scalarField> tmomentumDiffusion =
129 boundaryContrPtr_->momentumDiffusion();
130 scalarField& momentumDiffusion = tmomentumDiffusion.ref();
131
132 scalarField denom(phiOverSurf + momentumDiffusion*delta);
133
134 operator==
135 (
136 // Inlet
137 - neg(phip)*sourcePressure
138
139 // Outlet
140 + pos(phip)
141 *((Uac&nf)*nf + (Uac_t*(momentumDiffusion*delta) - sourceVelocity)/denom)
142 );
143
144 fixedValueFvPatchVectorField::updateCoeffs();
145}
146
147
150(
151 const tmp<scalarField>&
152) const
153{
154 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
155
156 // For fixedValue U patches
157 return tmp<Field<vector>>
158 (
159 new Field<vector>
160 (
162 )
163 );
164}
165
166
169(
170 const tmp<scalarField>&
171) const
172{
173 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
174
175 // For zeroGradient U patches
176 return tmp<Field<vector>>
177 (
178 new Field<vector>
179 (
180 pos(phip)*(*this)
181 )
182 );
183}
184
185
187{
189 writeEntry("value", os);
190 os.writeEntry("solverName", adjointSolverName_);
191}
192
193
194// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195
196namespace Foam
197{
199 (
202 );
203}
204
205// ************************************************************************* //
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...
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
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Base class for solution control classes.
virtual tmp< Field< vector > > valueBoundaryCoeffs(const tmp< scalarField > &) const
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< Field< vector > > valueInternalCoeffs(const tmp< scalarField > &) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
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
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...
Definition: fvsPatchField.H:79
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
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
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
dictionary dict