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-2020 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  const surfaceVectorField& Cf = mesh_.Cf();
42  const surfaceVectorField& Sf = mesh_.Sf();
43 
44  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
45 
46  const vectorField& Cfi = Cf;
47  const vectorField& Sfi = Sf;
48  scalarField& phii = phi.primitiveFieldRef();
49 
50  // Internal faces
51  forAll(internalFaces_, i)
52  {
53  label facei = internalFaces_[i];
54  phii[facei] -= rho[facei]*(Omega ^ (Cfi[facei] - origin_)) & Sfi[facei];
55  }
56 
57  makeRelativeRhoFlux(rho.boundaryField(), phi.boundaryFieldRef());
58 }
59 
60 
61 template<class RhoFieldType>
62 void Foam::MRFZone::makeRelativeRhoFlux
63 (
64  const RhoFieldType& rho,
65  FieldField<fvsPatchField, scalar>& phi
66 ) const
67 {
68  const surfaceVectorField& Cf = mesh_.Cf();
69  const surfaceVectorField& Sf = mesh_.Sf();
70 
71  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
72 
73  // Included patches
74  forAll(includedFaces_, patchi)
75  {
76  forAll(includedFaces_[patchi], i)
77  {
78  label patchFacei = includedFaces_[patchi][i];
79 
80  phi[patchi][patchFacei] = 0.0;
81  }
82  }
83 
84  // Excluded patches
85  forAll(excludedFaces_, patchi)
86  {
87  forAll(excludedFaces_[patchi], i)
88  {
89  label patchFacei = excludedFaces_[patchi][i];
90 
91  phi[patchi][patchFacei] -=
92  rho[patchi][patchFacei]
93  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
94  & Sf.boundaryField()[patchi][patchFacei];
95  }
96  }
97 }
98 
99 
100 template<class RhoFieldType>
101 void Foam::MRFZone::makeRelativeRhoFlux
102 (
103  const RhoFieldType& rho,
104  Field<scalar>& phi,
105  const label patchi
106 ) const
107 {
108  const surfaceVectorField& Cf = mesh_.Cf();
109  const surfaceVectorField& Sf = mesh_.Sf();
110 
111  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
112 
113  // Included patches
114  forAll(includedFaces_[patchi], i)
115  {
116  label patchFacei = includedFaces_[patchi][i];
117 
118  phi[patchFacei] = 0.0;
119  }
120 
121  // Excluded patches
122  forAll(excludedFaces_[patchi], i)
123  {
124  label patchFacei = excludedFaces_[patchi][i];
125 
126  phi[patchFacei] -=
127  rho[patchFacei]
128  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
129  & Sf.boundaryField()[patchi][patchFacei];
130  }
131 }
132 
133 
134 template<class RhoFieldType>
135 void Foam::MRFZone::makeAbsoluteRhoFlux
136 (
137  const RhoFieldType& rho,
138  surfaceScalarField& phi
139 ) const
140 {
141  const surfaceVectorField& Cf = mesh_.Cf();
142  const surfaceVectorField& Sf = mesh_.Sf();
143 
144  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
145 
146  const vectorField& Cfi = Cf;
147  const vectorField& Sfi = Sf;
148  scalarField& phii = phi.primitiveFieldRef();
149 
150  // Internal faces
151  forAll(internalFaces_, i)
152  {
153  label facei = internalFaces_[i];
154  phii[facei] += rho[facei]*(Omega ^ (Cfi[facei] - origin_)) & Sfi[facei];
155  }
156 
157  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
158 
159 
160  // Included patches
161  forAll(includedFaces_, patchi)
162  {
163  forAll(includedFaces_[patchi], i)
164  {
165  label patchFacei = includedFaces_[patchi][i];
166 
167  phibf[patchi][patchFacei] +=
168  rho.boundaryField()[patchi][patchFacei]
169  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
170  & Sf.boundaryField()[patchi][patchFacei];
171  }
172  }
173 
174  // Excluded patches
175  forAll(excludedFaces_, patchi)
176  {
177  forAll(excludedFaces_[patchi], i)
178  {
179  label patchFacei = excludedFaces_[patchi][i];
180 
181  phibf[patchi][patchFacei] +=
182  rho.boundaryField()[patchi][patchFacei]
183  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
184  & Sf.boundaryField()[patchi][patchFacei];
185  }
186  }
187 }
188 
189 
190 template<class Type>
192 (
194 ) const
195 {
196  Field<Type>& phii = phi.primitiveFieldRef();
197 
198  forAll(internalFaces_, i)
199  {
200  phii[internalFaces_[i]] = Zero;
201  }
202 
204  phi.boundaryFieldRef();
205 
206  forAll(includedFaces_, patchi)
207  {
208  forAll(includedFaces_[patchi], i)
209  {
210  phibf[patchi][includedFaces_[patchi][i]] = Zero;
211  }
212  }
213 
214  forAll(excludedFaces_, patchi)
215  {
216  forAll(excludedFaces_[patchi], i)
217  {
218  phibf[patchi][excludedFaces_[patchi][i]] = Zero;
219  }
220  }
221 }
222 
223 
224 // ************************************************************************* //
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
virtual Type value(const scalar x) const =0
Return value as a function of scalar x.
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.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField