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  const Type& nullValue
41 )
42 {
43  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
44 
45  const fvMesh& mesh = ssf.mesh();
46 
47  tmp<volFieldType> tresult
48  (
49  new volFieldType
50  (
51  IOobject
52  (
53  "cellReduce(" + ssf.name() + ')',
54  ssf.instance(),
55  mesh,
56  IOobject::NO_READ,
57  IOobject::NO_WRITE
58  ),
59  mesh,
60  dimensioned<Type>("initialValue", ssf.dimensions(), nullValue),
61  extrapolatedCalculatedFvPatchField<Type>::typeName
62  )
63  );
64 
65  volFieldType& result = tresult.ref();
66 
67  const labelUList& own = mesh.owner();
68  const labelUList& nbr = mesh.neighbour();
69 
70  forAll(own, i)
71  {
72  label celli = own[i];
73  cop(result[celli], ssf[i]);
74  }
75  forAll(nbr, i)
76  {
77  label celli = nbr[i];
78  cop(result[celli], ssf[i]);
79  }
80 
81  result.correctBoundaryConditions();
82 
83  return tresult;
84 }
85 
86 
87 template<class Type, class CombineOp>
90 (
91  const tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>& tssf,
92  const CombineOp& cop,
93  const Type& nullValue
94 )
95 {
96  tmp<GeometricField<Type, fvPatchField, volMesh>> tvf
97  (
99  (
100  tssf,
101  cop,
102  nullValue
103  )
104  );
105 
106  tssf.clear();
107 
108  return tvf;
109 }
110 
111 
112 // ************************************************************************* //
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
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
A class for managing temporary objects.
Definition: PtrList.H:53