Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fvcCellReduce.H"
30#include "fvMesh.H"
31#include "volFields.H"
32#include "surfaceFields.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37template<class Type, class CombineOp>
40(
41 const GeometricField<Type, fvsPatchField, surfaceMesh>& ssf,
42 const CombineOp& cop,
43 const Type& nullValue
44)
45{
46 typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
47
48 const fvMesh& mesh = ssf.mesh();
49
50 tmp<volFieldType> tresult
51 (
52 new volFieldType
53 (
54 IOobject
55 (
56 "cellReduce(" + ssf.name() + ')',
57 ssf.instance(),
58 mesh,
59 IOobject::NO_READ,
60 IOobject::NO_WRITE
61 ),
62 mesh,
63 dimensioned<Type>("initialValue", ssf.dimensions(), nullValue),
64 fvPatchFieldBase::extrapolatedCalculatedType()
65 )
66 );
67
68 volFieldType& result = tresult.ref();
69
70 const labelUList& own = mesh.owner();
71 const labelUList& nbr = mesh.neighbour();
72
73 forAll(own, i)
74 {
75 const label celli = own[i];
76 cop(result[celli], ssf[i]);
77 }
78 forAll(nbr, i)
79 {
80 const label celli = nbr[i];
81 cop(result[celli], ssf[i]);
82 }
83
84 forAll(mesh.boundary(), patchi)
85 {
86 const auto& pFaceCells = mesh.boundary()[patchi].faceCells();
87 const auto& pssf = ssf.boundaryField()[patchi];
88
89 forAll(pssf, i)
90 {
91 const label celli = pFaceCells[i];
92 cop(result[celli], pssf[i]);
93 }
94 }
95
96 result.correctBoundaryConditions();
97
98 return tresult;
99}
100
101
102template<class Type, class CombineOp>
105(
106 const tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>& tssf,
107 const CombineOp& cop,
108 const Type& nullValue
109)
110{
111 tmp<GeometricField<Type, fvPatchField, volMesh>> tvf
112 (
113 cellReduce
114 (
115 tssf,
116 cop,
117 nullValue
118 )
119 );
120
121 tssf.clear();
122
123 return tvf;
124}
125
126
127// ************************************************************************* //
A class for managing temporary objects.
Definition tmp.H:75
dynamicFvMesh & mesh
Construct a volume field from a surface field using a combine operator.
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
UList< label > labelUList
A UList of labels.
Definition UList.H:75
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.