All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
mQhdFluxFvPatchScalarField.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::mQhdFluxFvPatchScalarField
30 
31 Group
32  grpInletBoundaryConditions 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  Example of the boundary condition specification:
40  \verbatim
41  <patchName>
42  {
43  type mQhdFlux;
44  }
45  \endverbatim
46 
47 See also
48  Foam::fixedGradientFvPatchField
49 
50 SourceFiles
51  mQhdFluxFvPatchScalarField.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef mQhdFluxFvPatchScalarField_H
56 #define mQhdFluxFvPatchScalarField_H
57 
58 #include "fvPatchFields.H"
59 #include "fixedGradientFvPatchFields.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class mQhdFluxFvPatchScalarField Declaration
68 \*---------------------------------------------------------------------------*/
69 
71 :
72  public fixedGradientFvPatchScalarField
73 {
74 
75 
76 public:
77 
78  //- Runtime type information
79  TypeName("mQhdFlux");
80 
81 
82  // Constructors
83 
84  //- Construct from patch and internal field
86  (
87  const fvPatch&,
88  const DimensionedField<scalar, volMesh>&
89  );
90 
91  //- Construct from patch, internal field and dictionary
93  (
94  const fvPatch&,
95  const DimensionedField<scalar, volMesh>&,
96  const dictionary&
97  );
98 
99  //- Construct by mapping given mQhdFluxFvPatchScalarField onto
100  // a new patch
102  (
104  const fvPatch&,
105  const DimensionedField<scalar, volMesh>&,
106  const fvPatchFieldMapper&
107  );
108 
109  //- Construct as copy
111  (
113  );
114 
115  //- Construct and return a clone
116  virtual tmp<fvPatchScalarField> clone() const
117  {
118  return tmp<fvPatchScalarField>
119  (
120  new mQhdFluxFvPatchScalarField(*this)
121  );
122  }
123 
124  //- Construct as copy setting internal field reference
126  (
128  const DimensionedField<scalar, volMesh>&
129  );
130 
131  //- Construct and return a clone setting internal field reference
132  virtual tmp<fvPatchScalarField> clone
133  (
134  const DimensionedField<scalar, volMesh>& iF
135  ) const
136  {
137  return tmp<fvPatchScalarField>
138  (
139  new mQhdFluxFvPatchScalarField(*this, iF)
140  );
141  }
142 
143 
144  // Member functions
147  //virtual void updateSnGrad(const scalarField& snGradp);
148 
149  //- Update the patch pressure gradient field
150  virtual void updateCoeffs();
151 
152  //- Write
153  virtual void write(Ostream&) const;
154 };
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace Foam
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 //#include "volFields.H"
164 //
165 //namespace Foam
166 //{
167 // template<class GradBC>
168 // inline void setSnGrad
169 // (
170 // volScalarField::Boundary& bf,
171 // const FieldField<fvsPatchField, scalar>& snGrad
172 // )
173 // {
174 // forAll(bf, patchi)
175 // {
176 // if (isA<GradBC>(bf[patchi]))
177 // {
178 // refCast<GradBC>(bf[patchi]).updateSnGrad(snGrad[patchi]);
179 // }
180 // }
181 // }
182 //
183 // template<class GradBC>
184 // inline void setSnGrad
185 // (
186 // volScalarField::Boundary& bf,
187 // const tmp<FieldField<fvsPatchField, scalar>>& tsnGrad
188 // )
189 // {
190 // setSnGrad<GradBC>(bf, tsnGrad());
191 // }
192 //}
193 //
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #endif
198 
199 // ************************************************************************* //
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the patch pressure gradient field.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
TypeName("mQhdFlux")
Runtime type information.
mQhdFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...