74template<
class Type>
class Field;
182 inline explicit
Field(const label len);
185 inline
Field(const label len, const Type& val);
203 inline explicit
Field(const
UList<Type>& list);
216 template<
int SizeMin>
236 const
UList<Type>& mapF,
252 const
UList<Type>& mapF,
254 const
bool applyFlip = true
260 const
UList<Type>& mapF,
262 const Type& defaultValue,
263 const
bool applyFlip = true
271 const
UList<Type>& defaultValues,
272 const
bool applyFlip = true
280 const
bool applyFlip = true
289 const Type& defaultValue,
290 const
bool applyFlip = true
299 const
UList<Type>& defaultValues,
300 const
bool applyFlip = true
341 template<
class Type2>
412 const bool applyFlip =
true
420 const bool applyFlip =
true
427 const bool applyFlip =
true
496 void clamp_range(
const Type& lower,
const Type& upper);
502 template<
class VSForm>
503 VSForm
block(
const label start)
const;
544 template<
int SizeMin>
553 template<
class Form,
class Cmpt, direction nCmpt>
580 friend Ostream& operator<< <Type>
583 friend Ostream& operator<< <Type>
593 expr.evaluate(*
this);
600 expr.evaluate(*
this);
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Expression templates for List.
static bool unifiedGeometricField
GeometricField with extra capacity for flattened boundary fields. Uses opt-switch "unifiedGeometricFi...
static bool allowConstructFromLargerSize
Permit read construct from a larger size.
static int localBoundaryConsistency() noexcept
Get flag for local boundary consistency checks.
static void warnLocalBoundaryConsistencyCompat(const dictionary &)
Warn about keyword changes for local boundary consistency checks.
constexpr FieldBase() noexcept
Default construct.
static scalar localBoundaryTolerance_
Tolerance for local boundary field consistency checks. Uses opt-switch "localBoundaryConsistency::tol...
static int localBoundaryConsistency(int val) noexcept
Set flag for local boundary consistency checks.
static int localBoundaryConsistency_
Local boundary field consistency checks. Uses opt-switch "localBoundaryConsistency".
static const char *const typeName
Typename for Field.
Abstract base class to hold the Field mapping addressing and weights.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
void clamp_range(const MinMax< Type > &range)
Clamp field values (in-place) to the specified range.
void operator=(Foam::zero)
Assign entries to zero.
static tmp< Field< Type > > NewCalculatedType(const Field< Type2 > &f)
Return a pointer to a new calculatedFvPatchFieldField created on freestore without setting patchField...
void map(const tmp< Field< Type > > &tmapF, const FieldMapper &map, const bool applyFlip=true)
Map from the given tmp field.
void operator+=(const tmp< Field< Type > > &)
void operator-=(const Type &)
static autoPtr< Field< T > > New(Istream &is)
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors).
void map(const UList< Type > &mapF, const FieldMapper &map, const bool applyFlip=true)
Map from the given field.
void operator/=(const scalar &)
void map(const tmp< Field< Type > > &tmapF, const labelListList &mapAddressing, const scalarListList &weights)
Interpolative map from the given tmp field.
void operator=(const Field< Type > &)
Copy assignment.
void operator/=(const tmp< Field< scalar > > &)
void operator=(Field< Type > &&rhs)
Move assignment.
void clamp_range(const Type &lower, const Type &upper)
Clamp field values (in-place) to the specified range.
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
void operator=(const SubField< Type > &rhs)
void operator*=(const scalar &)
void operator+=(const Type &)
void operator=(const UList< Type > &rhs)
void operator+=(const UList< Type > &)
void operator-=(const UList< Type > &)
void operator=(DynamicList< Type, SizeMin > &&rhs)
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place).
void operator=(const tmp< Field< Type > > &)
SubField< Type > slice(const labelRange &range)
Return SubField slice (non-const access) - with range checking.
void operator*=(const tmp< Field< scalar > > &)
pTraits< T >::cmptType cmptType
void clamp_max(const UList< Type > &upper)
Impose upper (ceiling) clamp on the field values (in-place).
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
static const Field< Type > & null() noexcept
Return a null Field (reference to a nullObject). Behaves like an empty Field.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
void clamp_min(const UList< Type > &lower)
Impose lower (floor) clamp on the field values (in-place).
constexpr Field() noexcept
Default construct.
void operator=(const Expression::ListExpression< E > &expr)
Assign values from expression.
void operator=(List< Type > &&rhs)
void map(const UList< Type > &mapF, const labelListList &mapAddressing, const scalarListList &weights)
Interpolative map from the given field.
void operator=(const VectorSpace< Form, Cmpt, nCmpt > &)
void rmap(const tmp< Field< Type > > &tmapF, const labelUList &mapAddressing, const UList< scalar > &weights)
Interpolative reverse map from the given tmp field.
void operator*=(const UList< scalar > &)
void replace(const direction, const cmptType &)
Replace a component field of the field.
void assign(const entry &e, const label len)
Assign from a primitive dictionary entry with the following behaviour:
void replace(const direction, const tmp< Field< cmptType > > &)
Replace a component field of the field.
tmp< Field< T > > clone() const
bool assign(const word &key, const dictionary &dict, const label len, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ)
Lookup a primitive dictionary entry by (literal) name and assign its contents to this (behaviour as d...
void operator/=(const UList< scalar > &)
void operator-=(const tmp< Field< Type > > &)
void negate()
Inplace negate this field (negative).
void map(const tmp< Field< Type > > &tmapF, const labelUList &mapAddressing)
1 to 1 map from the given tmp field
void clamp_max(const Type &upper)
Impose upper (ceiling) clamp on the field values (in-place).
void map(const UList< T > &mapF, const labelUList &mapAddressing)
void operator=(const Type &val)
Assign entries to the given value.
void operator=(const IndirectListBase< Type, Addr > &rhs)
Copy assign from IndirectList.
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing, const UList< scalar > &weights)
Interpolative reverse map from the given field.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
const SubField< Type > slice(const labelRange &range) const
Return SubField slice (const access) - with range checking.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
void rmap(const tmp< Field< Type > > &tmapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given tmp field
VSForm block(const label start) const
SubField< Type > slice(const label pos, label len=-1)
Return SubField slice (non-const access) - no range checking.
const SubField< Type > slice(const label pos, label len=-1) const
Return SubField slice (const access) - no range checking.
Field(const Expression::ListExpression< E > &expr)
Construct from value expression.
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
readOption
Enumeration defining read preferences.
@ MUST_READ
Reading required.
Base for lists with indirect addressing, templated on the list contents type and the addressing type....
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A min/max value pair with additional methods. In addition to conveniently storing values,...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
SubField is a Field obtained as a section of another Field, without its own allocation....
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
auto expr() const
Wrap value as expression.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A keyword and a list of tokens is an 'entry'.
A range or interval of labels defined by a start and a size.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
A traits class, which is primarily used for primitives and vector-space.
constexpr refCount() noexcept
Default construct, initializing count to 0.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Direction is an 8-bit unsigned integer type used to represent Cartesian directions,...
OBJstream os(runTime.globalPath()/outputName)
List< scalarList > scalarListList
List of scalarList.
dimensionedScalar pos(const dimensionedScalar &ds)
const T & NullObjectRef() noexcept
Const reference (of type T) to the nullObject.
List< labelList > labelListList
List of labelList.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
UList< label > labelUList
A UList of labels.
Various functors for unary and binary operations. Can be used for parallel combine-reduce operations ...
A non-counting (dummy) refCount.