HydrothermalHeatFluxFvPatchScalarField.H
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
24Class
25 Foam::HydrothermalHeatFlux
26
27Description
28 Fixed heat flux boundary condition to specify temperature gradient. Input heat flux [W/m^2].
29 \heading Patch usage
30 \table
31 Property | Description | Required | Default value
32 q | heat flux value | yes |
33 kr | variable name for thermal conductivity in transport dict | no | kr
34 shape | special shape | no | fixed
35 \endtable
36 Example of the boundary condition specification:
37 \verbatim
38 myPatch
39 {
40 type HydrothermalHeatFlux;
41 q uniform 5; //W/m^2
42 kr kr;
43 value uniform 0; // place holder
44 shape gaussian2d; //gaussian2d, gaussian3d
45 x0 0; //for gaussian shape: peak position
46 z0 0; //for gaussian3d shape: peak position
47 qmax 5; //for gaussian shape
48 qmin 0.05; //for gaussian shape
49 c 500; //for gaussian shape: heaf width
50 }
51 \endverbatim
52
53SourceFiles
54 HydrothermalHeatFlux.C
55
56\*---------------------------------------------------------------------------*/
57
58
59#ifndef HydrothermalHeatFlux_H
60#define HydrothermalHeatFlux_H
61
62#include "fvPatchFields.H"
63#include "fixedGradientFvPatchFields.H"
64
65// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67namespace Foam
68{
69
70/*---------------------------------------------------------------------------*\
71 Class HydrothermalHeatFluxFvPatch Declaration
72 \*---------------------------------------------------------------------------*/
73
75:
76 public fixedGradientFvPatchScalarField
77{
78 // Private data
79
80 //- Name of the velocity and permeability fields used to calculate the wall BC in darcy's law
81 word krName_; //Thermal conductivity: default is "kr"
82 scalarField q_; //heat flux [W/m^2]
83 // --------- Special heat flux boundary condition --------------
84 word shape_; //Gaussian
85 // 1. Gaussian shape
93 scalar x1_, x2_, q1_, q2_; //variables for linear heat flux boundary
94public:
95
96 //- Runtime type information
97 TypeName("hydrothermalHeatFlux");
98
99
100 // Constructors
101
102 //- Construct from patch and internal field
104 (
105 const fvPatch&,
106 const DimensionedField<scalar, volMesh>&
107 );
108
109 //- Construct from patch, internal field and dictionary
111 (
112 const fvPatch&,
113 const DimensionedField<scalar, volMesh>&,
114 const dictionary&
115 );
116
117 //- Construct by mapping given
118 // HydrothermalHeatFlux onto a new patch
120 (
122 const fvPatch&,
123 const DimensionedField<scalar, volMesh>&,
124 const fvPatchFieldMapper&
125 );
126
127 //- Construct as copy
129 (
131 );
132
133 //- Construct and return a clone
134 virtual tmp<fvPatchScalarField> clone() const
135 {
136 return tmp<fvPatchScalarField>
137 (
138 new HydrothermalHeatFlux(*this)
139 );
140 }
141
142 //- Construct as copy setting internal field reference
144 (
146 const DimensionedField<scalar, volMesh>&
147 );
148
149 //- Construct and return a clone setting internal field reference
150 virtual tmp<fvPatchScalarField> clone
151 (
152 const DimensionedField<scalar, volMesh>& iF
153 ) const
154 {
155 return tmp<fvPatchScalarField>
156 (
157 new HydrothermalHeatFlux(*this, iF)
158 );
159 }
160
161 // Member functions
162
163 //- Update the coefficients associated with the patch field
164 virtual void updateCoeffs();
165
166 //- Write
167 virtual void write(Ostream&) const;
168};
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172} // End namespace Foam
173
174// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175
176#endif
177
178// ************************************************************************* //
virtual tmp< fvPatchScalarField > clone(const DimensionedField< scalar, volMesh > &iF) const
virtual tmp< fvPatchScalarField > clone() const
HydrothermalHeatFlux(const fvPatch &, const DimensionedField< scalar, volMesh > &)
TypeName("hydrothermalHeatFlux")