MRFZoneTemplatesDF.C
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 #include "MRFZoneDF.H"
35 #include "fvMesh.H"
36 #include "volFields.H"
37 #include "surfaceFields.H"
38 #include "fvMatrices.H"
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 template<class RhoFieldType>
43 void Foam::MRFZoneDF::makeRelativeRhoFlux(
44  const RhoFieldType& rho,
45  surfaceScalarField& phi) const
46 {
47  if (!active_)
48  {
49  return;
50  }
51 
52  const surfaceVectorField& Cf = mesh_.Cf();
53  const surfaceVectorField& Sf = mesh_.Sf();
54 
55  const vector Omega = omega_ * axis_;
56 
57  const vectorField& Cfi = Cf;
58  const vectorField& Sfi = Sf;
59  scalarField& phii = phi.primitiveFieldRef();
60 
61  // Internal faces
62  forAll(internalFaces_, i)
63  {
64  label facei = internalFaces_[i];
65  phii[facei] -= rho[facei] * (Omega ^ (Cfi[facei] - origin_)) & Sfi[facei];
66  }
67 
68  makeRelativeRhoFlux(rho.boundaryField(), phi.boundaryFieldRef());
69 }
70 
71 template<class RhoFieldType>
72 void Foam::MRFZoneDF::makeRelativeRhoFlux(
73  const RhoFieldType& rho,
74  FieldField<fvsPatchField, scalar>& phi) const
75 {
76  if (!active_)
77  {
78  return;
79  }
80 
81  const surfaceVectorField& Cf = mesh_.Cf();
82  const surfaceVectorField& Sf = mesh_.Sf();
83 
84  const vector Omega = omega_ * axis_;
85 
86  // Included patches
87  forAll(includedFaces_, patchi)
88  {
89  forAll(includedFaces_[patchi], i)
90  {
91  label patchFacei = includedFaces_[patchi][i];
92 
93  phi[patchi][patchFacei] = 0.0;
94  }
95  }
96 
97  // Excluded patches
98  forAll(excludedFaces_, patchi)
99  {
100  forAll(excludedFaces_[patchi], i)
101  {
102  label patchFacei = excludedFaces_[patchi][i];
103 
104  phi[patchi][patchFacei] -=
105  rho[patchi][patchFacei]
106  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
107  & Sf.boundaryField()[patchi][patchFacei];
108  }
109  }
110 }
111 
112 template<class RhoFieldType>
113 void Foam::MRFZoneDF::makeRelativeRhoFlux(
114  const RhoFieldType& rho,
115  Field<scalar>& phi,
116  const label patchi) const
117 {
118  /*
119  if (!active_)
120  {
121  return;
122  }
123 
124  const surfaceVectorField& Cf = mesh_.Cf();
125  const surfaceVectorField& Sf = mesh_.Sf();
126 
127  const vector Omega = omega_ * axis_;
128 
129  // Included patches
130  forAll(includedFaces_[patchi], i)
131  {
132  label patchFacei = includedFaces_[patchi][i];
133 
134  phi[patchFacei] = 0.0;
135  }
136 
137  // Excluded patches
138  forAll(excludedFaces_[patchi], i)
139  {
140  label patchFacei = excludedFaces_[patchi][i];
141 
142  phi[patchFacei] -=
143  rho[patchFacei]
144  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
145  & Sf.boundaryField()[patchi][patchFacei];
146  }
147  */
148 }
149 
150 // ************************************************************************* //
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
MRFZoneDF.H
Foam::MRFZoneDF::Omega
vector Omega() const
Definition: MRFZoneDF.C:262
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8