fvcCellReduce.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) 2013-2015 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 "fvcCellReduce.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fvc
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type, class CombineOp>
45 tmp<GeometricField<Type, fvPatchField, volMesh> > cellReduce
46 (
48  const CombineOp& cop,
49  const Type& nullValue
50 )
51 {
53 
54  const fvMesh& mesh = ssf.mesh();
55 
56  tmp<volFieldType> tresult
57  (
58  new volFieldType
59  (
60  IOobject
61  (
62  "cellReduce(" + ssf.name() + ')',
63  ssf.instance(),
64  mesh,
67  ),
68  mesh,
69  dimensioned<Type>("initialValue", ssf.dimensions(), nullValue),
71  )
72  );
73 
74  volFieldType& result = tresult();
75 
76  const labelUList& own = mesh.owner();
77  const labelUList& nbr = mesh.neighbour();
78 
79  // Internal field
80  const Field<Type>& iField = ssf.internalField();
81  forAll(iField, faceI)
82  {
83  label cellOwn = own[faceI];
84  cop(result[cellOwn], iField[faceI]);
85 
86  label cellNbr = nbr[faceI];
87  cop(result[cellNbr], iField[faceI]);
88  }
89 
90  // Boundary field
91  forAll(ssf.boundaryField(), patchI)
92  {
93  const fvsPatchField<Type>& pf = ssf.boundaryField()[patchI];
94  const label start = pf.patch().start();
95 
96  forAll(pf, i)
97  {
98  label faceI = start + i;
99  label cellI = own[faceI];
100  cop(result[cellI], pf[i]);
101  }
102  }
103 
104  result.correctBoundaryConditions();
105 
106  return tresult;
107 }
108 
109 
110 template<class Type, class CombineOp>
112 (
114  const CombineOp& cop,
115  const Type& nullValue
116 )
117 {
119  tvf(cellReduce(cop, tssf, nullValue));
120 
121  tssf.clear();
122  return tvf;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 } // End namespace fvc
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
132 } // End namespace Foam
133 
134 // ************************************************************************* //
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fvc::cellReduce
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf, const CombineOp &cop, const Type &nullValue)
Definition: fvcCellReduce.C:46
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvsPatchField< Type >
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
surfaceFields.H
Foam::surfaceFields.
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:63
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< Type >
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned< Type >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:278
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
fvcCellReduce.H
Construct a volume field from a surface field using a combine operator.
zeroGradientFvPatchFields.H