MRFZoneListDF.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 MRFZoneListDF_H
35 #define MRFZoneListDF_H
36 
37 #include "fvMesh.H"
38 #include "dictionary.H"
39 #include "fvMatricesFwd.H"
40 #include "MRFZoneDF.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // Forward declaration of friend functions and operators
48 class MRFZoneListDF;
49 Ostream& operator<<(Ostream& os, const MRFZoneListDF& models);
50 
51 /*---------------------------------------------------------------------------*\
52  Class MRFZoneListDF Declaration
53 \*---------------------------------------------------------------------------*/
54 
56  : public PtrList<MRFZoneDF>
57 {
58  // Private Member Functions
59 
60  //- No copy construct
61  MRFZoneListDF(const MRFZoneListDF&) = delete;
62 
63  //- No copy assignment
64  void operator=(const MRFZoneListDF&) = delete;
65 
66 protected:
67  // Protected data
68 
69  //- Reference to the mesh database
70  const fvMesh& mesh_;
71 
72 public:
73  //- Constructor
74  MRFZoneListDF(const fvMesh& mesh, const dictionary& dict);
75 
76  //- Destructor
78 
79  // Member Functions
80 
81  //- Return active status
82  bool active(const bool warn = false) const;
83 
84  //- Reset the source list
85  void reset(const dictionary& dict);
86 
87  //- Return the frame acceleration
88  tmp<volVectorField> DDt(
89  const volVectorField& U) const;
90 
91  //- Return the frame acceleration
92  tmp<volVectorField> DDt(
93  const volScalarField& rho,
94  const volVectorField& U) const;
95 
96  //- Make the given absolute velocity relative within the MRF region
97  void makeRelative(volVectorField& U) const;
98 
99  //- Make the given absolute flux relative within the MRF region
100  void makeRelative(surfaceScalarField& phi) const;
101 
102  //- Return the given absolute flux relative within the MRF region
103  tmp<surfaceScalarField> relative(
104  const tmp<surfaceScalarField>& phi) const;
105 
106  //- Return the given absolute boundary flux relative within
107  // the MRF region
108  tmp<FieldField<fvsPatchField, scalar>> relative(
109  const tmp<FieldField<fvsPatchField, scalar>>& tphi) const;
110 
111  //- Return the given absolute patch flux relative within
112  // the MRF region
113  tmp<Field<scalar>> relative(
114  const tmp<Field<scalar>>& tphi,
115  const label patchi) const;
116 
117  //- Make the given absolute mass-flux relative within the MRF region
118  void makeRelative(
119  const surfaceScalarField& rho,
120  surfaceScalarField& phi) const;
121 
122  //- Correct the boundary velocity for the rotation of the MRF region
123  void correctBoundaryVelocity(volVectorField& U) const;
124 
125  //- Update MRFZone faces if the mesh topology changes
126  //void update();
127 
128  //- Get the reference of the rotation speed omega from MRFZone
129  const scalar& getOmegaRef() const;
130 
131  // I-O
132 
133  //- Read dictionary
134  //bool read(const dictionary& dict);
135 
136  //- Write data to Ostream
137  //bool writeData(Ostream& os) const;
138 
139  //- Ostream operator
140  friend Ostream& operator<<(
141  Ostream& os,
142  const MRFZoneListDF& models);
143 };
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 } // End namespace Foam
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 #endif
152 
153 // ************************************************************************* //
Foam::MRFZoneListDF::DDt
tmp< volVectorField > DDt(const volVectorField &U) const
Definition: MRFZoneListDF.C:120
Foam::MRFZoneListDF::~MRFZoneListDF
~MRFZoneListDF()
Definition: MRFZoneListDF.C:53
U
U
Definition: pEqnPimpleDyM.H:60
MRFZoneDF.H
Foam::MRFZoneListDF::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Definition: MRFZoneListDF.C:164
Foam::MRFZoneListDF
Definition: MRFZoneListDF.H:55
Foam::operator<<
Ostream & operator<<(Ostream &os, const MRFZoneListDF &models)
Definition: MRFZoneListDF.C:248
Foam::MRFZoneListDF::operator<<
friend Ostream & operator<<(Ostream &os, const MRFZoneListDF &models)
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
Foam::MRFZoneListDF::reset
void reset(const dictionary &dict)
Definition: MRFZoneListDF.C:75
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::MRFZoneListDF::makeRelative
void makeRelative(volVectorField &U) const
Definition: MRFZoneListDF.C:148
Foam::MRFZoneListDF::getOmegaRef
const scalar & getOmegaRef() const
Definition: MRFZoneListDF.C:215
Foam::MRFZoneListDF::mesh_
const fvMesh & mesh_
Definition: MRFZoneListDF.H:70
Foam::MRFZoneListDF::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Definition: MRFZoneListDF.C:207
Foam::MRFZoneListDF::active
bool active(const bool warn=false) const
Definition: MRFZoneListDF.C:59
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8