MRFZoneTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "MRFZone.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "fvMatrices.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class RhoFieldType>
35 void Foam::MRFZone::makeRelativeRhoFlux
36 (
37  const RhoFieldType& rho,
39 ) const
40 {
41  if (!active_)
42  {
43  return;
44  }
45 
46  const surfaceVectorField& Cf = mesh_.Cf();
47  const surfaceVectorField& Sf = mesh_.Sf();
48 
49  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
50 
51  const vectorField& Cfi = Cf;
52  const vectorField& Sfi = Sf;
53  scalarField& phii = phi.primitiveFieldRef();
54 
55  // Internal faces
56  forAll(internalFaces_, i)
57  {
58  label facei = internalFaces_[i];
59  phii[facei] -= rho[facei]*(Omega ^ (Cfi[facei] - origin_)) & Sfi[facei];
60  }
61 
62  makeRelativeRhoFlux(rho.boundaryField(), phi.boundaryFieldRef());
63 }
64 
65 
66 template<class RhoFieldType>
67 void Foam::MRFZone::makeRelativeRhoFlux
68 (
69  const RhoFieldType& rho,
70  FieldField<fvsPatchField, scalar>& phi
71 ) const
72 {
73  if (!active_)
74  {
75  return;
76  }
77 
78  const surfaceVectorField& Cf = mesh_.Cf();
79  const surfaceVectorField& Sf = mesh_.Sf();
80 
81  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
82 
83  // Included patches
84  forAll(includedFaces_, patchi)
85  {
86  forAll(includedFaces_[patchi], i)
87  {
88  label patchFacei = includedFaces_[patchi][i];
89 
90  phi[patchi][patchFacei] = 0.0;
91  }
92  }
93 
94  // Excluded patches
95  forAll(excludedFaces_, patchi)
96  {
97  forAll(excludedFaces_[patchi], i)
98  {
99  label patchFacei = excludedFaces_[patchi][i];
100 
101  phi[patchi][patchFacei] -=
102  rho[patchi][patchFacei]
103  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
104  & Sf.boundaryField()[patchi][patchFacei];
105  }
106  }
107 }
108 
109 
110 template<class RhoFieldType>
111 void Foam::MRFZone::makeRelativeRhoFlux
112 (
113  const RhoFieldType& rho,
114  Field<scalar>& phi,
115  const label patchi
116 ) const
117 {
118  if (!active_)
119  {
120  return;
121  }
122 
123  const surfaceVectorField& Cf = mesh_.Cf();
124  const surfaceVectorField& Sf = mesh_.Sf();
125 
126  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
127 
128  // Included patches
129  forAll(includedFaces_[patchi], i)
130  {
131  label patchFacei = includedFaces_[patchi][i];
132 
133  phi[patchFacei] = 0.0;
134  }
135 
136  // Excluded patches
137  forAll(excludedFaces_[patchi], i)
138  {
139  label patchFacei = excludedFaces_[patchi][i];
140 
141  phi[patchFacei] -=
142  rho[patchFacei]
143  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
144  & Sf.boundaryField()[patchi][patchFacei];
145  }
146 }
147 
148 
149 template<class RhoFieldType>
150 void Foam::MRFZone::makeAbsoluteRhoFlux
151 (
152  const RhoFieldType& rho,
153  surfaceScalarField& phi
154 ) const
155 {
156  if (!active_)
157  {
158  return;
159  }
160 
161  const surfaceVectorField& Cf = mesh_.Cf();
162  const surfaceVectorField& Sf = mesh_.Sf();
163 
164  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
165 
166  const vectorField& Cfi = Cf;
167  const vectorField& Sfi = Sf;
168  scalarField& phii = phi.primitiveFieldRef();
169 
170  // Internal faces
171  forAll(internalFaces_, i)
172  {
173  label facei = internalFaces_[i];
174  phii[facei] += rho[facei]*(Omega ^ (Cfi[facei] - origin_)) & Sfi[facei];
175  }
176 
177  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
178 
179 
180  // Included patches
181  forAll(includedFaces_, patchi)
182  {
183  forAll(includedFaces_[patchi], i)
184  {
185  label patchFacei = includedFaces_[patchi][i];
186 
187  phibf[patchi][patchFacei] +=
188  rho.boundaryField()[patchi][patchFacei]
189  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
190  & Sf.boundaryField()[patchi][patchFacei];
191  }
192  }
193 
194  // Excluded patches
195  forAll(excludedFaces_, patchi)
196  {
197  forAll(excludedFaces_[patchi], i)
198  {
199  label patchFacei = excludedFaces_[patchi][i];
200 
201  phibf[patchi][patchFacei] +=
202  rho.boundaryField()[patchi][patchFacei]
203  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
204  & Sf.boundaryField()[patchi][patchFacei];
205  }
206  }
207 }
208 
209 
210 template<class Type>
212 (
214 ) const
215 {
216  if (!active_)
217  {
218  return;
219  }
220 
221  Field<Type>& phii = phi.primitiveFieldRef();
222 
223  forAll(internalFaces_, i)
224  {
225  phii[internalFaces_[i]] = Zero;
226  }
227 
229  phi.boundaryFieldRef();
230 
231  forAll(includedFaces_, patchi)
232  {
233  forAll(includedFaces_[patchi], i)
234  {
235  phibf[patchi][includedFaces_[patchi][i]] = Zero;
236  }
237  }
238 
239  forAll(excludedFaces_, patchi)
240  {
241  forAll(excludedFaces_[patchi], i)
242  {
243  phibf[patchi][excludedFaces_[patchi][i]] = Zero;
244  }
245  }
246 }
247 
248 
249 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const surfaceVectorField & Sf() const
Return cell face area vectors.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Generic GeometricField class.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
Pre-declare SubField and related Field type.
Definition: Field.H:56
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void zero(GeometricField< Type, fvsPatchField, surfaceMesh > &phi) const
Zero the MRF region of the given field.
static const zero Zero
Definition: zero.H:97
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
label patchi
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:29
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Field< vector > vectorField
Specialisation of Field<T> for vector.
mesh Sf()
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField