Loading...
Searching...
No Matches
adjointOutletNuaTildaFvPatchScalarField.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#include "fvPatchFieldMapper.H"
33#include "volFields.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
43(
44 const fvPatch& p,
58 const fvPatchFieldMapper& mapper
60:
61 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
63{}
64
65
67(
68 const fvPatch& p,
70 const dictionary& dict
71)
72:
73 fixedValueFvPatchScalarField(p, iF),
74 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
75{
76 this->readValueEntry(dict, IOobjectOption::MUST_READ);
77}
78
79
81(
84)
85:
86 fixedValueFvPatchScalarField(tppsf, iF),
88{}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
94{
95 if (updated())
96 {
97 return;
98 }
99 vectorField nf(patch().nf());
100
101 const fvPatchField<vector>& Ub = boundaryContrPtr_->Ub();
102 tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable1Diffusion());
103 const scalarField& nuEff = tnuEff();
104
105 // Patch-adjacent nuaTilda nuaTildaNei
106 tmp<scalarField> tnuaTildaNei(patchInternalField());
107 const scalarField& nuaTildaNei = tnuaTildaNei();
108
109 const scalarField& delta = patch().deltaCoeffs();
110
111 operator==
112 (
113 (nuEff*delta*nuaTildaNei)
114 /((Ub & nf) + nuEff*delta)
115 );
116
117 fixedValueFvPatchScalarField::updateCoeffs();
118}
119
120
122{
124 os.writeEntry("solverName", adjointSolverName_);
126}
127
128
129// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130
131makePatchTypeField(fvPatchScalarField, adjointOutletNuaTildaFvPatchScalarField);
132
133// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134
135} // End namespace Foam
136
137// ************************************************************************* //
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.
adjointOutletNuaTildaFvPatchScalarField(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