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-2016 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 template<class Type, class CombineOp>
37 (
38  const GeometricField<Type, fvsPatchField, surfaceMesh>& ssf,
39  const CombineOp& cop
40 )
41 {
42  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
43 
44  const fvMesh& mesh = ssf.mesh();
45 
46  tmp<volFieldType> tresult
47  (
48  new volFieldType
49  (
50  IOobject
51  (
52  "cellReduce(" + ssf.name() + ')',
53  ssf.instance(),
54  mesh,
55  IOobject::NO_READ,
56  IOobject::NO_WRITE
57  ),
58  mesh,
59  dimensioned<Type>("0", ssf.dimensions(), Zero),
60  extrapolatedCalculatedFvPatchField<Type>::typeName
61  )
62  );
63 
64  volFieldType& result = tresult.ref();
65 
66  const labelUList& own = mesh.owner();
67  const labelUList& nbr = mesh.neighbour();
68 
69  forAll(own, i)
70  {
71  label celli = own[i];
72  cop(result[celli], ssf[i]);
73  }
74  forAll(nbr, i)
75  {
76  label celli = nbr[i];
77  cop(result[celli], ssf[i]);
78  }
79 
80  result.correctBoundaryConditions();
81 
82  return tresult;
83 }
84 
85 
86 template<class Type, class CombineOp>
89 (
90  const tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>& tssf,
91  const CombineOp& cop
92 )
93 {
94  tmp<GeometricField<Type, fvPatchField, volMesh>> tvf(cellReduce(tssf, cop));
95 
96  tssf.clear();
97 
98  return tvf;
99 }
100 
101 
102 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
Construct a volume field from a surface field using a combine operator.
UList< label > labelUList
Definition: UList.H:64
dynamicFvMesh & mesh
static const zero Zero
Definition: zero.H:91
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop)
A class for managing temporary objects.
Definition: PtrList.H:54