adjointOutletKaFvPatchScalarField.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-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "fvPatchFieldMapper.H"
32#include "volFields.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
42(
43 const fvPatch& p,
45)
46:
47 fixedValueFvPatchScalarField(p, iF),
49{}
50
51
53(
55 const fvPatch& p,
57 const fvPatchFieldMapper& mapper
58)
59:
60 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
61 adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
62{}
63
64
66(
67 const fvPatch& p,
69 const dictionary& dict
70)
71:
72 fixedValueFvPatchScalarField(p, iF),
73 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
74{
76 (
77 scalarField("value", dict, p.size())
78 );
79}
80
81
83(
86)
87:
88 fixedValueFvPatchScalarField(tppsf, iF),
90{}
91
92
93// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94
96{
97 if (updated())
98 {
99 return;
100 }
101 vectorField nf(patch().nf());
102
103 const fvPatchField<vector>& Ub = boundaryContrPtr_->Ub();
104 tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable1Diffusion());
105 const scalarField& nuEff = tnuEff();
106
107 // Patch-adjacent ka
108 tmp<scalarField> tkaNei(patchInternalField());
109 const scalarField& kaNei = tkaNei();
110
111 const scalarField& delta = patch().deltaCoeffs();
112
113 // Source from the objective
114 tmp<scalarField> source = boundaryContrPtr_->adjointTMVariable1Source();
115
116 operator==
117 (
118 (nuEff*delta*kaNei - source)
119 /((Ub & nf) + nuEff*delta)
120 );
121
122 fixedValueFvPatchScalarField::updateCoeffs();
123}
124
125
127{
129 writeEntry("value", os);
130 os.writeEntry("solverName", adjointSolverName_);
131}
132
133
134// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135
137
138// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140} // End namespace Foam
141
142// ************************************************************************* //
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.
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
word adjointSolverName_
adjointSolver name corresponding to field
virtual void updateCoeffs()
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: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
A class for managing temporary objects.
Definition: tmp.H:65
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.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dictionary dict