31 template<
class Type,
class GeoMesh>
36 const Field<Type>& completeField,
37 const bool preserveCouples,
38 const bool preserveProcessorOnly
41 tmp<FieldField<GeoMesh::template PatchField, Type>> tbf
43 new FieldField<GeoMesh::template PatchField, Type>
48 FieldField<GeoMesh::template PatchField, Type>& bf = tbf.ref();
57 !preserveProcessorOnly
104 template<
class Type,
class GeoMesh>
109 const FieldField<GeoMesh::template PatchField, Type>& bField,
110 const bool preserveCouples
113 tmp<FieldField<GeoMesh::template PatchField, Type>> tbf
115 new FieldField<GeoMesh::template PatchField, Type>
120 FieldField<GeoMesh::template PatchField, Type>& bf = tbf.ref();
163 template<
class Type,
class GeoMesh>
170 const bool preserveCouples
179 slicedBoundaryField(
mesh, completeField, preserveCouples)
192 template<
class Type,
class GeoMesh>
200 const bool preserveCouples,
201 const bool preserveProcessorOnly
215 preserveProcessorOnly
229 template<
class Type,
class GeoMesh>
234 const bool preserveCouples
243 slicedBoundaryField(gf.
mesh(), gf.boundaryField(), preserveCouples)
253 template<
class Type,
class GeoMesh>
265 slicedBoundaryField(gf.
mesh(), gf.boundaryField(), true)
273 template<
class Type,
class GeoMesh>
286 template<
class Type,
class GeoMesh>
297 template<
class Type,
class GeoMesh>
303 label completeSize = GeoMesh::size(
mesh);
314 = this->primitiveField();
331 ) = this->boundaryField()[
patchi];
346 return tCompleteField;
350 template<
class Type,
class GeoMesh>
#define forAll(list, i)
Loop across all elements in list.
static const DimensionedField< Type, GeoMesh, PrimitiveField > & null()
Return a null DimensionedField.
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
SubField< Type > subField
Declare type of subField.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
void correctBoundaryConditions()
Correct boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Specialisation of GeometricField which holds slices of given complete fields in a form that they act ...
tmp< SlicedGeometricField< Type, GeoMesh > > clone() const
Clone.
SlicedGeometricField(const IOobject &, const Mesh &, const dimensionSet &, const Field< Type > &completeField, const bool preserveCouples=true)
Construct from components and field to slice.
tmp< Field< Type > > splice() const
Splice the sliced field and return the complete field.
~SlicedGeometricField()
Destructor.
void correctBoundaryConditions()
Correct boundary field.
Pre-declare related SubField type.
void shallowCopy(const UList< T > &)
Copy the pointer held by the given UList.
label size() const
Return the number of elements in the UPtrList.
Dimension set for the base types.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.