Foam::ensightOutput::Detail Namespace Reference

Implementation details and output backends that would not normally be called directly by a user. More...

Functions

labelList getFaceSizes (const UList< face > &faces)
 Return sizes of faces in the list. More...
 
labelList getFaceSizes (const UIndirectList< face > &faces)
 Return sizes of faces in the list. More...
 
labelList getPolysNFaces (const polyMesh &mesh, const labelUList &addr)
 The number of faces per poly element. More...
 
labelList getPolysNPointsPerFace (const polyMesh &mesh, const labelUList &addr)
 The number of points for each face of the poly elements. More...
 
CompactListList< label > getPolysFacePoints (const polyMesh &mesh, const labelUList &addr, const labelList &pointMap)
 Generate 0-based point ids for each poly element face. More...
 
void writeLabelListList (ensightGeoFile &os, const labelUList &offsets, const labelUList &values, const label pointOffset)
 Write CompactListList<label> by components. More...
 
template<class LabelListListType >
void writeLabelListList (ensightGeoFile &os, const LabelListListType &listOfLists, const label pointOffset)
 Write a list of faces or cell shapes with one-entity per line. More...
 
template<template< typename > class FieldContainer, class Type >
void copyComponent (List< scalar > &cmptBuffer, const FieldContainer< Type > &input, const direction cmpt)
 
template<template< typename > class FieldContainer, class Type >
void writeFieldContent (ensightFile &os, const FieldContainer< Type > &fld, bool parallel)
 Write field content (component-wise) More...
 
template<template< typename > class FieldContainer>
bool writeCoordinates (ensightGeoFile &os, const label partId, const word &partName, const label nPoints, const FieldContainer< Foam::point > &fld, bool parallel)
 Write coordinates (component-wise) for the given part. More...
 
template<template< typename > class FieldContainer, class Type >
bool writeFieldComponents (ensightFile &os, const char *key, const FieldContainer< Type > &fld, bool parallel)
 Write field content (component-wise) for the given ensight element type. More...
 
template<class Type >
bool writeFaceSubField (ensightFile &os, const Field< Type > &fld, const ensightFaces &part, bool parallel)
 
template<class Type >
bool writeFaceLocalField (ensightFile &os, const Field< Type > &fld, const ensightFaces &part, bool parallel)
 
template<class Type >
label writeCloudFieldContent (ensightFile &os, const UList< Type > &fld, label count=0)
 Write cloud field data (serial) with rounding and newlines. More...
 

Detailed Description

Implementation details and output backends that would not normally be called directly by a user.

Function Documentation

◆ getFaceSizes() [1/2]

Foam::labelList getFaceSizes ( const UList< face > &  faces)

Return sizes of faces in the list.

Definition at line 40 of file ensightOutput.C.

References UList< T >::begin(), f(), and UList< T >::size().

Here is the call graph for this function:

◆ getFaceSizes() [2/2]

Foam::labelList getFaceSizes ( const UIndirectList< face > &  faces)

Return sizes of faces in the list.

Definition at line 59 of file ensightOutput.C.

References UList< T >::begin(), f(), and IndirectListBase< T, Addr >::size().

Here is the call graph for this function:

◆ getPolysNFaces()

Foam::labelList getPolysNFaces ( const polyMesh mesh,
const labelUList addr 
)

The number of faces per poly element.

const cellList& meshCells = mesh.cells();

Definition at line 78 of file ensightOutput.C.

References UList< T >::begin(), cellId, mesh, and UList< T >::size().

Here is the call graph for this function:

◆ getPolysNPointsPerFace()

Foam::labelList getPolysNPointsPerFace ( const polyMesh mesh,
const labelUList addr 
)

The number of points for each face of the poly elements.

const cellList& meshCells = mesh.cells();

Definition at line 102 of file ensightOutput.C.

References UList< T >::begin(), cellId, mesh, and UList< T >::size().

Here is the call graph for this function:

◆ getPolysFacePoints()

Foam::CompactListList< Foam::label > getPolysFacePoints ( const polyMesh mesh,
const labelUList addr,
const labelList pointMap 
)

Generate 0-based point ids for each poly element face.

The returned CompactListList is divided per output face

const cellList& meshCells = mesh.cells();

Parameters
addrCell ids to write
pointMapPoint map to use

Definition at line 210 of file ensightOutput.C.

References cellId, f(), faceId(), mesh, nFaces(), nPoints, CompactListList< T >::offsets(), UList< T >::size(), and CompactListList< T >::values().

Here is the call graph for this function:

◆ writeLabelListList() [1/2]

void writeLabelListList ( ensightGeoFile os,
const labelUList offsets,
const labelUList values,
const label  pointOffset 
)

Write CompactListList<label> by components.

Definition at line 138 of file ensightOutput.C.

References os(), and UList< T >::size().

Here is the call graph for this function:

◆ writeLabelListList() [2/2]

void writeLabelListList ( ensightGeoFile os,
const LabelListListType &  listOfLists,
const label  pointOffset 
)

Write a list of faces or cell shapes with one-entity per line.

Definition at line 35 of file ensightOutputTemplates.C.

References forAll, and os().

Here is the call graph for this function:

◆ copyComponent()

void copyComponent ( List< scalar > &  cmptBuffer,
const FieldContainer< Type > &  input,
const direction  cmpt 
)

Copy specified field component into a scalar buffer works for various lists types. Must be adequately sized before calling

Definition at line 56 of file ensightOutputTemplates.C.

References UList< T >::begin(), Foam::component(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::nl, and UList< T >::size().

Referenced by writeFieldContent().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeFieldContent()

void writeFieldContent ( ensightFile os,
const FieldContainer< Type > &  fld,
bool  parallel 
)

Write field content (component-wise)

Parameters
parallelCollective write?

Definition at line 82 of file ensightOutputTemplates.C.

References UList< T >::cdata_bytes(), copyComponent(), UList< T >::data_bytes(), fld(), os(), procAddr(), DynamicList< T, SizeMin >::resize_nocopy(), and UList< T >::size_bytes().

Here is the call graph for this function:

◆ writeCoordinates()

bool writeCoordinates ( ensightGeoFile os,
const label  partId,
const word partName,
const label  nPoints,
const FieldContainer< Foam::point > &  fld,
bool  parallel 
)

Write coordinates (component-wise) for the given part.

Parameters
parallelCollective write?

Definition at line 152 of file ensightOutputTemplates.C.

References fld(), nPoints, and os().

Referenced by ensightOutputSurface::write(), ensightCells::write(), and ensightFaces::write().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeFieldComponents()

bool writeFieldComponents ( ensightFile os,
const char *  key,
const FieldContainer< Type > &  fld,
bool  parallel 
)

Write field content (component-wise) for the given ensight element type.

Parameters
parallelCollective write?

Definition at line 177 of file ensightOutputTemplates.C.

References fld(), os(), and reduce().

Referenced by ensightOutputSurface::writePointData(), and Foam::writeTrackField().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeFaceSubField()

bool writeFaceSubField ( ensightFile os,
const Field< Type > &  fld,
const ensightFaces part,
bool  parallel 
)

Write a sub-field of faces values as an indirect list, using the sub-list sizing information from ensightFaces

Parameters
parallelCollective write?

Definition at line 213 of file ensightOutputTemplates.C.

References fld(), ensightPart::index(), os(), ensightFaces::range(), reduce(), ensightFaces::size(), and ensightFaces::total().

Here is the call graph for this function:

◆ writeFaceLocalField()

bool writeFaceLocalField ( ensightFile os,
const Field< Type > &  fld,
const ensightFaces part,
bool  parallel 
)

Write a field of faces values as an indirect list, using the face order from ensightFaces

Parameters
parallelCollective write?

Definition at line 263 of file ensightOutputTemplates.C.

References Foam::exit(), ensightFaces::faceOrder(), Foam::FatalError, FatalErrorInFunction, fld(), ensightPart::index(), Foam::nl, os(), reduce(), ensightFaces::size(), UList< T >::size(), and ensightFaces::total().

Here is the call graph for this function:

◆ writeCloudFieldContent()

label writeCloudFieldContent ( ensightFile os,
const UList< Type > &  fld,
label  count = 0 
)

Write cloud field data (serial) with rounding and newlines.

Returns
the current output count
Parameters
countThe current output count