36 template<
class Type,
class CombineOp>
40 const GeometricField<Type, fvsPatchField, surfaceMesh>& ssf,
45 typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
47 const fvMesh&
mesh = ssf.mesh();
49 tmp<volFieldType> tresult
55 "cellReduce(" + ssf.name() +
')',
62 dimensioned<Type>(
"initialValue", ssf.dimensions(), nullValue),
67 volFieldType& result = tresult.ref();
75 cop(result[celli], ssf[i]);
80 cop(result[celli], ssf[i]);
83 result.correctBoundaryConditions();
89 template<
class Type,
class CombineOp>
93 const tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>& tssf,
98 tmp<GeometricField<Type, fvPatchField, volMesh>> tvf
Construct a volume field from a surface field using a combine operator.
UList< label > labelUList
A UList of labels.
#define forAll(list, i)
Loop across all elements in list.
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.
const word extrapolatedCalculatedType
A combined zero-gradient and calculated patch field type.