All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
qhdFluxFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on 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 Class
29  Foam::qhdFluxFvPatchScalarField
30 
31 Group
32  grpQGDBCs grpWallBoundaryConditions
33 
34 Description
35  This boundary condition sets the pressure gradient to the provided value
36  such that the flux on the boundary is that specified by the velocity
37  boundary condition.
38 
39 Usage
40  Example of the boundary condition specification:
41  \verbatim
42  <patchName>
43  {
44  type qhdFlux;
45  }
46  \endverbatim
47 
48 SourceFiles
49  qhdFluxFvPatchScalarField.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef qhdFluxFvPatchScalarField_H
54 #define qhdFluxFvPatchScalarField_H
55 
56 #include "fvPatchFields.H"
57 #include "fixedGradientFvPatchFields.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
63 
64 /*---------------------------------------------------------------------------*\
65  Class qhdFluxFvPatchScalarField Declaration
66 \*---------------------------------------------------------------------------*/
67 
69 :
70  public fixedGradientFvPatchScalarField
71 {
72 
73 
74 public:
75 
76  //- Runtime type information
77  TypeName("qhdFlux");
78 
79 
80  // Constructors
81 
82  //- Construct from patch and internal field
84  (
85  const fvPatch&,
86  const DimensionedField<scalar, volMesh>&
87  );
88 
89  //- Construct from patch, internal field and dictionary
91  (
92  const fvPatch&,
93  const DimensionedField<scalar, volMesh>&,
94  const dictionary&
95  );
96 
97  //- Construct by mapping given qhdFluxFvPatchScalarField onto
98  // a new patch
100  (
102  const fvPatch&,
103  const DimensionedField<scalar, volMesh>&,
104  const fvPatchFieldMapper&
105  );
106 
107  //- Construct as copy
109  (
111  );
112 
113  //- Construct and return a clone
114  virtual tmp<fvPatchScalarField> clone() const
115  {
116  return tmp<fvPatchScalarField>
117  (
118  new qhdFluxFvPatchScalarField(*this)
119  );
120  }
121 
122  //- Construct as copy setting internal field reference
124  (
126  const DimensionedField<scalar, volMesh>&
127  );
128 
129  //- Construct and return a clone setting internal field reference
130  virtual tmp<fvPatchScalarField> clone
131  (
132  const DimensionedField<scalar, volMesh>& iF
133  ) const
134  {
135  return tmp<fvPatchScalarField>
136  (
137  new qhdFluxFvPatchScalarField(*this, iF)
138  );
139  }
140 
141 
142  // Member functions
145  //virtual void updateSnGrad(const scalarField& snGradp);
146 
147  //- Update the patch pressure gradient field
148  virtual void updateCoeffs();
149 
150  //- Write
151  virtual void write(Ostream&) const;
152 };
153 
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 } // End namespace Foam
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 //#include "volFields.H"
162 //
163 //namespace Foam
164 //{
165 // template<class GradBC>
166 // inline void setSnGrad
167 // (
168 // volScalarField::Boundary& bf,
169 // const FieldField<fvsPatchField, scalar>& snGrad
170 // )
171 // {
172 // forAll(bf, patchi)
173 // {
174 // if (isA<GradBC>(bf[patchi]))
175 // {
176 // refCast<GradBC>(bf[patchi]).updateSnGrad(snGrad[patchi]);
177 // }
178 // }
179 // }
180 //
181 // template<class GradBC>
182 // inline void setSnGrad
183 // (
184 // volScalarField::Boundary& bf,
185 // const tmp<FieldField<fvsPatchField, scalar>>& tsnGrad
186 // )
187 // {
188 // setSnGrad<GradBC>(bf, tsnGrad());
189 // }
190 //}
191 //
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 #endif
196 
197 // ************************************************************************* //
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
virtual void write(Ostream &) const
Write.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
qhdFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the patch pressure gradient field.
TypeName("qhdFlux")
Runtime type information.