adjointFarFieldNuaTildaFvPatchScalarField.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#include "surfaceFields.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43adjointFarFieldNuaTildaFvPatchScalarField::
44adjointFarFieldNuaTildaFvPatchScalarField
45(
46 const fvPatch& p,
48)
49:
50 fixedValueFvPatchScalarField(p, iF),
52{}
53
54
55adjointFarFieldNuaTildaFvPatchScalarField::
56adjointFarFieldNuaTildaFvPatchScalarField
57(
59 const fvPatch& p,
61 const fvPatchFieldMapper& mapper
62)
63:
64 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
65 adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
66{}
67
68
69adjointFarFieldNuaTildaFvPatchScalarField::
70adjointFarFieldNuaTildaFvPatchScalarField
71(
72 const fvPatch& p,
74 const dictionary& dict
75)
76:
77 fixedValueFvPatchScalarField(p, iF),
78 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
79{
81 (
82 scalarField("value", dict, p.size())
83 );
84}
85
86
87adjointFarFieldNuaTildaFvPatchScalarField::
88adjointFarFieldNuaTildaFvPatchScalarField
89(
92)
93:
94 fixedValueFvPatchScalarField(tppsf, iF),
96{}
97
98
99// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100
102{
103 if (updated())
104 {
105 return;
106 }
107 vectorField nf(patch().nf());
108
109 const fvPatchField<vector>& Ub = boundaryContrPtr_->Ub();
110 tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable1Diffusion());
111 const scalarField& nuEff = tnuEff();
112 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
113
114 // Patch-adjacent nuaTilda nuaTildaNei
115 tmp<scalarField> tnuaTildaNei(patchInternalField());
116 const scalarField& nuaTildaNei = tnuaTildaNei();
117
118 // Patch deltas
119 const scalarField& delta = patch().deltaCoeffs();
120
121 operator==
122 (
123 pos(phip)
124 *(
125 (nuEff*delta*nuaTildaNei)
126 /((Ub & nf) + nuEff*delta)
127 )
128 );
129
130 fixedValueFvPatchScalarField::updateCoeffs();
131}
132
133
136(
137 const tmp<scalarField>&
138) const
139{
140 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
141
142 // For fixedValue nuTilda patches
144}
145
146
149(
150 const tmp<scalarField>&
151) const
152{
153 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
154
155 // For zeroGradient nuTilda patches
156 return tmp<Field<scalar>>::New(pos(phip)*(*this));
157}
158
159
161{
163 writeEntry("value", os);
164 os.writeEntry("solverName", adjointSolverName_);
165}
166
167
168// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169
171(
174);
175
176// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177
178} // End namespace Foam
179
180// ************************************************************************* //
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 tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
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
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
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.
dimensionedScalar pos(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar neg(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
Foam::surfaceFields.