MRFZoneDF.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  Description:
7  A modified version of MRF that allows changing the rotation speed at
8  the runtime to enable derivatives wrt the rotation speed
9 
10  This class is modified from OpenFOAM's source code
11  src/finiteVolume/cfdTools/general/MRF
12 
13  OpenFOAM: The Open Source CFD Toolbox
14 
15  Copyright (C): 2011-2016 OpenFOAM Foundation
16 
17  OpenFOAM License:
18 
19  OpenFOAM is free software: you can redistribute it and/or modify it
20  under the terms of the GNU General Public License as published by
21  the Free Software Foundation, either version 3 of the License, or
22  (at your option) any later version.
23 
24  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
25  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
26  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
27  for more details.
28 
29  You should have received a copy of the GNU General Public License
30  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef MRFZoneDF_H
35 #define MRFZoneDF_H
36 
37 #include "dictionary.H"
38 #include "wordList.H"
39 #include "labelList.H"
40 #include "dimensionedScalar.H"
41 #include "dimensionedVector.H"
42 #include "volFieldsFwd.H"
43 #include "surfaceFields.H"
44 #include "fvMatricesFwd.H"
45 #include "mapPolyMesh.H"
46 #include "Function1.H"
47 #include "autoPtr.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 // Forward declaration of classes
55 class fvMesh;
56 
57 /*---------------------------------------------------------------------------*\
58  Class MRFZoneDF Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 class MRFZoneDF
62 {
63  // Private data
64 
65  //- Reference to the mesh database
66  const fvMesh& mesh_;
67 
68  //- Name of the MRF region
69  const word name_;
70 
71  //- Coefficients dictionary
72  dictionary coeffs_;
73 
74  //- MRF region active flag
75  bool active_;
76 
77  //- Name of cell zone
78  word cellZoneName_;
79 
80  //- Cell zone ID
81  label cellZoneID_;
82 
83  const wordRes excludedPatchNames_;
84 
85  labelList excludedPatchLabels_;
86 
87  //- Internal faces that are part of MRF
88  labelList internalFaces_;
89 
90  //- Outside faces (per patch) that move with the MRF
91  labelListList includedFaces_;
92 
93  //- Excluded faces (per patch) that do not move with the MRF
94  labelListList excludedFaces_;
95 
96  //- Origin of the axis
97  const vector origin_;
98 
99  //- Axis vector
100  vector axis_;
101 
102  //- Angular velocity (rad/sec)
103  scalar omega_;
104 
105  // Private Member Functions
106 
107  //- Divide faces in frame according to patch
108  void setMRFFaces();
109 
110  //- Make the given absolute mass/vol flux relative within the MRF region
111  template<class RhoFieldType>
112  void makeRelativeRhoFlux(
113  const RhoFieldType& rho,
114  surfaceScalarField& phi) const;
115 
116  //- Make the given absolute mass/vol flux relative within the MRF region
117  template<class RhoFieldType>
118  void makeRelativeRhoFlux(
119  const RhoFieldType& rho,
120  FieldField<fvsPatchField, scalar>& phi) const;
121 
122  //- Make the given absolute mass/vol flux relative within the MRF region
123  template<class RhoFieldType>
124  void makeRelativeRhoFlux(
125  const RhoFieldType& rho,
126  Field<scalar>& phi,
127  const label patchi) const;
128 
129  //- No copy construct
130  MRFZoneDF(const MRFZoneDF&) = delete;
131 
132  //- No copy assignment
133  void operator=(const MRFZoneDF&) = delete;
134 
135 public:
136  // Declare name of the class and its debug switch
137  ClassName("MRFZoneDF");
138 
139  // Constructors
140 
141  //- Construct from fvMesh
142  MRFZoneDF(
143  const word& name,
144  const fvMesh& mesh,
145  const dictionary& dict,
146  const word& cellZoneName = word::null);
147 
148  //- Return clone
149  autoPtr<MRFZoneDF> clone() const
150  {
151  NotImplemented;
152  return nullptr;
153  }
154 
155  // Member Functions
156 
157  //- Return const access to the MRF region name
158  inline const word& name() const;
159 
160  //- Return const access to the MRF active flag
161  inline bool active() const;
162 
163  //- Return the current Omega vector
164  vector Omega() const;
165 
166  //- Return the reference of rotation speed omega scalar
167  const scalar& getOmegaRef() const
168  {
169  return omega_;
170  }
171 
172  //- Update the mesh corresponding to given map
173  void updateMesh(const mapPolyMesh& mpm)
174  {
175  // Only updates face addressing
176  setMRFFaces();
177  }
178 
179  //- Add the Coriolis force contribution to the acceleration field
180  void addCoriolis(
181  const volVectorField& U,
182  volVectorField& ddtU) const;
183 
184  //- Make the given absolute velocity relative within the MRF region
185  void makeRelative(volVectorField& U) const;
186 
187  //- Make the given absolute flux relative within the MRF region
188  void makeRelative(surfaceScalarField& phi) const;
189 
190  //- Make the given absolute boundary flux relative
191  // within the MRF region
192  void makeRelative(FieldField<fvsPatchField, scalar>& phi) const;
193 
194  //- Make the given absolute patch flux relative
195  // within the MRF region
196  void makeRelative(Field<scalar>& phi, const label patchi) const;
197 
198  //- Make the given absolute mass-flux relative within the MRF region
199  void makeRelative(
200  const surfaceScalarField& rho,
201  surfaceScalarField& phi) const;
202 
203  //- Correct the boundary velocity for the rotation of the MRF region
204  void correctBoundaryVelocity(volVectorField& U) const;
205 
206  //- Update MRFZoneDF faces if the mesh topology changes
207  //void update();
208 
209  // I-O
210 
211  //- Write
212  //void writeData(Ostream& os) const;
213 
214  //- Read MRF dictionary
215  //bool read(const dictionary& dict);
216 };
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 } // End namespace Foam
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 #ifdef NoRepository
225 #include "MRFZoneTemplatesDF.C"
226 #endif
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 #include "MRFZoneIDF.H"
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 #endif
235 
236 // ************************************************************************* //
Foam::MRFZoneDF::name
const word & name() const
Definition: MRFZoneIDF.H:34
U
U
Definition: pEqnPimpleDyM.H:60
Foam::MRFZoneDF::active
bool active() const
Definition: MRFZoneIDF.H:39
Foam::MRFZoneDF::Omega
vector Omega() const
Definition: MRFZoneDF.C:262
Foam::MRFZoneDF
Definition: MRFZoneDF.H:61
Foam::MRFZoneDF::getOmegaRef
const scalar & getOmegaRef() const
Definition: MRFZoneDF.H:167
MRFZoneTemplatesDF.C
Foam::MRFZoneDF::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Definition: MRFZoneDF.C:356
Foam::MRFZoneDF::makeRelative
void makeRelative(volVectorField &U) const
Definition: MRFZoneDF.C:289
MRFZoneIDF.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::MRFZoneDF::addCoriolis
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Definition: MRFZoneDF.C:267
Foam::MRFZoneDF::updateMesh
void updateMesh(const mapPolyMesh &mpm)
Definition: MRFZoneDF.H:173
Foam::MRFZoneDF::clone
autoPtr< MRFZoneDF > clone() const
Definition: MRFZoneDF.H:149
Foam::MRFZoneDF::ClassName
ClassName("MRFZoneDF")
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8