MRFZoneTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 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>
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.internalField();
47  const vectorField& Sfi = Sf.internalField();
48  scalarField& phii = phi.internalField();
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.boundaryField());
58 }
59 
60 
61 template<class RhoFieldType>
63 (
64  const RhoFieldType& rho,
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>
102 (
103  const RhoFieldType& rho,
105 ) const
106 {
107  const surfaceVectorField& Cf = mesh_.Cf();
108  const surfaceVectorField& Sf = mesh_.Sf();
109 
110  const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
111 
112  const vectorField& Cfi = Cf.internalField();
113  const vectorField& Sfi = Sf.internalField();
114  scalarField& phii = phi.internalField();
115 
116  // Internal faces
117  forAll(internalFaces_, i)
118  {
119  label facei = internalFaces_[i];
120  phii[facei] += rho[facei]*(Omega ^ (Cfi[facei] - origin_)) & Sfi[facei];
121  }
122 
123  // Included patches
124  forAll(includedFaces_, patchi)
125  {
126  forAll(includedFaces_[patchi], i)
127  {
128  label patchFacei = includedFaces_[patchi][i];
129 
130  phi.boundaryField()[patchi][patchFacei] +=
131  rho.boundaryField()[patchi][patchFacei]
132  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
133  & Sf.boundaryField()[patchi][patchFacei];
134  }
135  }
136 
137  // Excluded patches
138  forAll(excludedFaces_, patchi)
139  {
140  forAll(excludedFaces_[patchi], i)
141  {
142  label patchFacei = excludedFaces_[patchi][i];
143 
144  phi.boundaryField()[patchi][patchFacei] +=
145  rho.boundaryField()[patchi][patchFacei]
146  * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
147  & Sf.boundaryField()[patchi][patchFacei];
148  }
149  }
150 }
151 
152 
153 // ************************************************************************* //
volFields.H
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
MRFZone.H
Foam::MRFZone::makeAbsoluteRhoFlux
void makeAbsoluteRhoFlux(const RhoFieldType &rho, surfaceScalarField &phi) const
Make the given relative mass/vol flux absolute within the MRF region.
Definition: MRFZoneTemplates.C:102
Foam::MRFZone::makeRelativeRhoFlux
void makeRelativeRhoFlux(const RhoFieldType &rho, surfaceScalarField &phi) const
Make the given absolute mass/vol flux relative within the MRF region.
Definition: MRFZoneTemplates.C:36
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
surfaceFields.H
Foam::surfaceFields.
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
fvMesh.H
rho
rho
Definition: pEqn.H:3
Foam::Vector< scalar >
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52