Loading...
Searching...
No Matches
adjointInletVelocityFvPatchVectorField.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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
39(
40 const fvPatch& p,
43:
44 fixedValueFvPatchVectorField(p, iF),
46{}
47
48
51(
53 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 if (updated())
95 {
96 return;
97 }
98
99 // Objective function contribution
100 tmp<vectorField> tsource = boundaryContrPtr_->normalVelocitySource();
101 const vectorField& source = tsource();
102
103 operator==(-source);
104
105 fixedValueFvPatchVectorField::updateCoeffs();
106}
107
108
111(
112 const tmp<scalarField>&
113) const
114{
115 return tmp<Field<vector>>::New(this->size(), pTraits<vector>::one);
116}
117
118
121(
123) const
124{
125 return tmp<Field<vector>>::New(this->size(), Zero);
126}
127
128
130{
132 os.writeEntry("solverName", adjointSolverName_);
135
136
137// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138
139namespace Foam
140{
142 (
144 adjointInletVelocityFvPatchVectorField
145 );
146}
147
148// ************************************************************************* //
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
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
adjointInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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()
Add explicit sink to zero adjoint velocity tangential motion at the cells next to the inlet.
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
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
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.
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.
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
adjointBoundaryCondition< vector > adjointVectorBoundaryCondition
fvPatchField< vector > fvPatchVectorField
dictionary dict