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 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  (
50  (
51  "cellReduce(" + ssf.name() + ')',
52  mesh,
53  dimensioned<Type>
54  (
55  ssf.dimensions(),
56  nullValue
57  ),
58  extrapolatedCalculatedFvPatchField<Type>::typeName
59  )
60  );
61 
62  volFieldType& result = tresult.ref();
63 
64  const labelUList& own = mesh.owner();
65  const labelUList& nbr = mesh.neighbour();
66 
67  forAll(own, i)
68  {
69  label celli = own[i];
70  cop(result[celli], ssf[i]);
71  }
72  forAll(nbr, i)
73  {
74  label celli = nbr[i];
75  cop(result[celli], ssf[i]);
76  }
77 
78  result.correctBoundaryConditions();
79 
80  return tresult;
81 }
82 
83 
84 template<class Type, class CombineOp>
87 (
88  const tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>& tssf,
89  const CombineOp& cop,
90  const Type& nullValue
91 )
92 {
93  tmp<GeometricField<Type, fvPatchField, volMesh>> tvf
94  (
96  (
97  tssf(),
98  cop,
99  nullValue
100  )
101  );
102 
103  tssf.clear();
104 
105  return tvf;
106 }
107 
108 
109 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Construct a volume field from a surface field using a combine operator.
fvMesh & mesh
UList< label > labelUList
Definition: UList.H:65
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