Loading...
Searching...
No Matches
adjointOutletWaFvPatchScalarField.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,
46:
47 fixedValueFvPatchScalarField(p, iF),
49{}
50
51
53(
55 const fvPatch& p,
57 const fvPatchFieldMapper& mapper
59:
60 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
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{
75 this->readValueEntry(dict, IOobjectOption::MUST_READ);
76}
77
78
80(
83)
84:
85 fixedValueFvPatchScalarField(tppsf, iF),
87{}
88
89
90// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91
93{
94 if (updated())
95 {
96 return;
97 }
98 vectorField nf(patch().nf());
99
100 const fvPatchField<vector>& Ub = boundaryContrPtr_->Ub();
101 tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable2Diffusion());
102 const scalarField& nuEff = tnuEff();
103
104 // Patch-adjacent wa
105 tmp<scalarField> twaNei = patchInternalField();
106 const scalarField& waNei = twaNei();
107
108 const scalarField& delta = patch().deltaCoeffs();
109
110 // Source from the objective
111 tmp<scalarField> source = boundaryContrPtr_->adjointTMVariable2Source();
112
113 operator==
114 (
115 (nuEff*delta*waNei - source)
116 /((Ub & nf) + nuEff*delta)
117 );
118
119 fixedValueFvPatchScalarField::updateCoeffs();
120}
121
122
124{
126 os.writeEntry("solverName", adjointSolverName_);
128}
129
130
131// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132
133makePatchTypeField(fvPatchScalarField, adjointOutletWaFvPatchScalarField);
134
135// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136
137} // End namespace Foam
138
139// ************************************************************************* //
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...
@ 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 void updateCoeffs()
Update the coefficients associated with the patch field.
adjointOutletWaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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 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.
adjointBoundaryCondition< scalar > adjointScalarBoundaryCondition
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
fvPatchField< scalar > fvPatchScalarField
dictionary dict