37template<
class Type,
class CombineOp>
41 const GeometricField<Type, fvsPatchField, surfaceMesh>& ssf,
46 typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
48 const fvMesh&
mesh = ssf.mesh();
50 tmp<volFieldType> tresult
56 "cellReduce(" + ssf.name() +
')',
63 dimensioned<Type>(
"initialValue", ssf.dimensions(), nullValue),
64 fvPatchFieldBase::extrapolatedCalculatedType()
68 volFieldType& result = tresult.ref();
75 const label celli = own[i];
76 cop(result[celli], ssf[i]);
80 const label celli = nbr[i];
81 cop(result[celli], ssf[i]);
86 const auto& pFaceCells =
mesh.boundary()[patchi].faceCells();
87 const auto& pssf = ssf.boundaryField()[patchi];
91 const label celli = pFaceCells[i];
92 cop(result[celli], pssf[i]);
96 result.correctBoundaryConditions();
102template<
class Type,
class CombineOp>
106 const tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>& tssf,
107 const CombineOp& cop,
108 const Type& nullValue
111 tmp<GeometricField<Type, fvPatchField, volMesh>> tvf
A class for managing temporary objects.
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.
#define forAll(list, i)
Loop across all elements in list.