noFluxFvPatchScalarField.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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8License
9 This file is part of OpenFOAM.
10
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
15
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19 for more details.
20
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23
24\*---------------------------------------------------------------------------*/
25
27#include "addToRunTimeSelectionTable.H"
28#include "linear.H"
29
30// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31
33(
34 const fvPatch& p,
35 const DimensionedField<scalar, volMesh>& iF
36)
37 :
38 fixedGradientFvPatchScalarField(p, iF),
39 rhorAUfName_("rhorAUf"),
40 phigName_("phig")
41{}
42
44(
45 const noFlux& ptf,
46 const fvPatch& p,
47 const DimensionedField<scalar, volMesh>& iF,
48 const fvPatchFieldMapper& mapper
49)
50 :
51 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
52 rhorAUfName_(ptf.rhorAUfName_),
53 phigName_(ptf.phigName_)
54{
55 patchType() = ptf.patchType();
56
57 // Map gradient. Set unmapped values and overwrite with mapped ptf
58 gradient() = 0.0;
59 mapper(gradient(), ptf.gradient());
60
61 // Evaluate the value field from the gradient if the internal field is valid
62 if (notNull(iF) && iF.size())
63 {
64 scalarField::operator=
65 (
66 // patchInternalField() + gradient()/patch().deltaCoeffs()
67 // ***HGW Hack to avoid the construction of mesh.deltaCoeffs
68 // which fails for AMI patches for some mapping operations
69 patchInternalField() + gradient()*(patch().nf() & patch().delta())
70 );
71 }
72 else
73 {
74 // Enforce mapping of values so we have a valid starting value. This
75 // constructor is used when reconstructing fields
76 mapper(*this, ptf);
77 }
78}
79
81(
82 const fvPatch& p,
83 const DimensionedField<scalar, volMesh>& iF,
84 const dictionary& dict
85)
86 :
87 fixedGradientFvPatchScalarField(p, iF),
88 rhorAUfName_(dict.lookupOrDefault<word>("rhorAUf", "rhorAUf")),
89 phigName_(dict.lookupOrDefault<word>("phig", "phig"))
90{
91 if (dict.found("value") && dict.found("gradient"))
92 {
93 fvPatchField<scalar>::operator=
94 (
95 scalarField("value", dict, p.size())
96 );
97 gradient() = scalarField("gradient", dict, p.size());
98 }
99 else
100 {
101 fvPatchField<scalar>::operator=(patchInternalField());
102 gradient() = 0.0;
103 }
104}
105
107(
108 const noFlux& ptf
109)
110 :
111 fixedGradientFvPatchScalarField(ptf),
112 rhorAUfName_(ptf.rhorAUfName_),
113 phigName_(ptf.phigName_)
114{}
115
117(
118 const noFlux& ptf,
119 const DimensionedField<scalar, volMesh>& iF
120)
121 :
122 fixedGradientFvPatchScalarField(ptf, iF),
123 rhorAUfName_(ptf.rhorAUfName_),
124 phigName_(ptf.phigName_)
125{}
126
127// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128
130{
131 if (updated())
132 {
133 return;
134 }
135
136 const fvsPatchField<scalar>& rhorAUf=
137 patch().lookupPatchField<surfaceScalarField, scalar>(rhorAUfName_);
138
139 const fvsPatchField<scalar>& phig=
140 patch().lookupPatchField<surfaceScalarField, scalar>(phigName_);
141
142 gradient()=phig/rhorAUf/(patch().magSf());
143
144 fixedGradientFvPatchScalarField::updateCoeffs();
145}
146
147void Foam::noFlux::write(Ostream& os) const
148{
149 fixedGradientFvPatchScalarField::write(os);
150 writeEntry(os, "value", *this);
151}
152
153// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154
155namespace Foam
156{
158(
159 fvPatchScalarField,
160 noFlux
161);
162}
163
164// ************************************************************************* //
volScalarField & p
Definition: createFields.H:10
surfaceScalarField phig("phig",(fvc::interpolate(rho) *rhorAUf *g) &mesh.Sf())
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
virtual void write(Ostream &) const
noFlux(const fvPatch &, const DimensionedField< scalar, volMesh > &)
virtual void updateCoeffs()
makePatchTypeField(fvPatchScalarField, HydrothermalHeatFlux)