fvcCellReduce.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) 2013-2022 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 SurfaceField<Type>& ssf,
39  const CombineOp& cop,
40  const Type& nullValue
41 )
42 {
43  const fvMesh& mesh = ssf.mesh();
44 
45  tmp<VolField<Type>> tresult
46  (
48  (
49  "cellReduce(" + ssf.name() + ')',
50  mesh,
51  dimensioned<Type>
52  (
53  ssf.dimensions(),
54  nullValue
55  ),
56  extrapolatedCalculatedFvPatchField<Type>::typeName
57  )
58  );
59 
60  VolField<Type>& result = tresult.ref();
61 
62  const labelUList& own = mesh.owner();
63  const labelUList& nbr = mesh.neighbour();
64 
65  forAll(own, i)
66  {
67  label celli = own[i];
68  cop(result[celli], ssf[i]);
69  }
70  forAll(nbr, i)
71  {
72  label celli = nbr[i];
73  cop(result[celli], ssf[i]);
74  }
75 
76  result.correctBoundaryConditions();
77 
78  return tresult;
79 }
80 
81 
82 template<class Type, class CombineOp>
85 (
86  const tmp<SurfaceField<Type>>& tssf,
87  const CombineOp& cop,
88  const Type& nullValue
89 )
90 {
91  tmp<VolField<Type>> tvf
92  (
94  (
95  tssf(),
96  cop,
97  nullValue
98  )
99  );
100 
101  tssf.clear();
102 
103  return tvf;
104 }
105 
106 
107 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A class for managing temporary objects.
Definition: tmp.H:55
Construct a volume field from a surface field using a combine operator.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > cellReduce(const SurfaceField< Type > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
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
UList< label > labelUList
Definition: UList.H:65
Foam::surfaceFields.