fixedWallHeatFluxFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  OpenFOAM: The Open Source CFD Toolbox
7 
8  Copyright (C): 2011-2016 OpenFOAM Foundation
9 
10  OpenFOAM License:
11 
12  OpenFOAM is free software: you can redistribute it and/or modify it
13  under the terms of the GNU General Public License as published by
14  the Free Software Foundation, either version 3 of the License, or
15  (at your option) any later version.
16 
17  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
18  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
19  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
20  for more details.
21 
22  You should have received a copy of the GNU General Public License
23  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 
25 \*---------------------------------------------------------------------------*/
26 
28 #include "addToRunTimeSelectionTable.H"
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40  const fvPatch& p,
41  const DimensionedField<scalar, volMesh>& iF)
42  : fixedGradientFvPatchScalarField(p, iF)
43 {
44 }
45 
49  const fvPatch& p,
50  const DimensionedField<scalar, volMesh>& iF,
51  const fvPatchFieldMapper& mapper)
52  : fixedGradientFvPatchScalarField(tdpvf, p, iF, mapper)
53 {
54 }
55 
58  const fvPatch& p,
59  const DimensionedField<scalar, volMesh>& iF,
60  const dictionary& dict)
61  : fixedGradientFvPatchScalarField(p, iF),
62  heatFlux_(dict.getScalar("heatFlux")),
63  nu_(dict.getScalar("nu")),
64  Pr_(dict.getScalar("Pr")),
65  Prt_(dict.getScalar("Prt")),
66  Cp_(dict.getScalar("Cp"))
67 {
68  fvPatchScalarField::operator=(patchInternalField());
69  gradient() = 0.0;
70 }
71 
75  : fixedGradientFvPatchScalarField(tdpvf)
76 {
77 }
78 
82  const DimensionedField<scalar, volMesh>& iF)
83  : fixedGradientFvPatchScalarField(tdpvf, iF)
84 {
85 }
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 {
90  if (updated())
91  {
92  return;
93  }
94 
95  const label patchi = patch().index();
96 
97  // Retrieve turbulence properties from model
98  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>(
99  IOobject::groupName(
100  turbulenceModel::propertiesName,
101  internalField().group()));
102 
103  const tmp<scalarField> nutw = turbModel.nut(patchi);
104 
105  gradient() = heatFlux_ / (nutw / Prt_ + nu_ / Pr_) / Cp_;
106 
107  fixedGradientFvPatchScalarField::updateCoeffs();
108 }
109 
111 {
112  fvPatchScalarField::write(os);
113  writeEntry("value", os);
114  os.writeEntry("Pr", Pr_);
115  os.writeEntry("Prt", Prt_);
116  os.writeEntry("nu", nu_);
117  os.writeEntry("Cp", Cp_);
118  os.writeEntry("heatFlux", heatFlux_);
119 }
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
124  fvPatchScalarField,
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 } // End namespace Foam
130 
131 // ************************************************************************* //
Foam::fixedWallHeatFluxFvPatchScalarField::fixedWallHeatFluxFvPatchScalarField
fixedWallHeatFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Definition: fixedWallHeatFluxFvPatchScalarField.C:39
fixedWallHeatFluxFvPatchScalarField.H
Foam::fixedWallHeatFluxFvPatchScalarField::write
virtual void write(Ostream &) const
Definition: fixedWallHeatFluxFvPatchScalarField.C:110
p
volScalarField & p
Definition: createRefsPimpleDyM.H:6
Foam
Definition: checkGeometry.C:32
Foam::fixedWallHeatFluxFvPatchScalarField
Definition: fixedWallHeatFluxFvPatchScalarField.H:47
Foam::fixedWallHeatFluxFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Definition: fixedWallHeatFluxFvPatchScalarField.C:88
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, fixedWallHeatFluxFvPatchScalarField)