Loading...
Searching...
No Matches
adjointFarFieldTMVar2FvPatchScalarField.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) 2014-2022 PCOpt/NTUA
9 Copyright (C) 2014-2022 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#include "surfaceFields.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
44(
45 const fvPatch& p,
48:
49 fixedValueFvPatchScalarField(p, iF),
51{}
52
53
56(
58 const fvPatch& p,
60 const fvPatchFieldMapper& mapper
62:
63 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
65{}
66
67
70(
71 const fvPatch& p,
73 const dictionary& dict
74)
75:
76 fixedValueFvPatchScalarField(p, iF),
77 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
78{
79 this->readValueEntry(dict, IOobjectOption::MUST_READ);
80}
81
82
85(
88)
89:
90 fixedValueFvPatchScalarField(tppsf, iF),
92{}
93
94
95// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96
98{
99 if (updated())
100 {
101 return;
102 }
103 vectorField nf(patch().nf());
104
105 tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable2Diffusion());
106 const scalarField& nuEff = tnuEff();
107 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
108 const scalarField& magSf = patch().magSf();
109
110 // Patch-adjacent values
111 tmp<scalarField> intf(patchInternalField());
112
113 // Patch deltas
114 const scalarField& delta = patch().deltaCoeffs();
115
116 operator==
117 (
118 pos(phip)
119 *(
120 (nuEff*delta*intf)
121 /(phip/magSf + nuEff*delta)
123 );
124
125 fixedValueFvPatchScalarField::updateCoeffs();
126}
127
128
131(
132 const tmp<scalarField>&
133) const
134{
136
137 // For fixedValue patches
139}
140
141
144(
145 const tmp<scalarField>&
146) const
147{
149
150 // For zeroGradient patches
151 return tmp<Field<scalar>>::New(pos(phip)*(*this));
152}
153
154
156{
158 os.writeEntry("solverName", adjointSolverName_);
160}
161
162
163// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164
166(
168 adjointFarFieldTMVar2FvPatchScalarField
169);
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173} // End namespace Foam
174
175// ************************************************************************* //
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 tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
adjointFarFieldTMVar2FvPatchScalarField(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.
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
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
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.
dimensionedScalar pos(const dimensionedScalar &ds)
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.
adjointBoundaryCondition< scalar > adjointScalarBoundaryCondition
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
dictionary dict
Foam::surfaceFields.