BC

namespace Cubism::BC

Namespace for block field boundary conditions.

template<typename Lab>
class Absorbing : public Cubism::BC::Base<Lab>
#include <Absorbing.h>

Zeroth-Order absorbing BC.

Zeroth-Order absorbing boundary condition

Template Parameters
  • Lab: Type of FieldLab

template<typename Lab>
class Base
#include <Base.h>

Boundary condition base class.

Each boundary condition is applied for a specific dir < CUBISM_DIMENSION and corresponding side.

Template Parameters
  • Lab: Type of FieldLab

Subclassed by Cubism::BC::Absorbing< Lab >, Cubism::BC::Dirichlet< Lab >, Cubism::BC::Symmetry< Lab >

struct BoundaryInfo
#include <Base.h>

Boundary information meta data.

template<typename Lab>
class Dirichlet : public Cubism::BC::Base<Lab>
#include <Dirichlet.h>

Dirichlet BC.

Constant value Dirichlet boundary condition

Template Parameters
  • Lab: Type of FieldLab

template<typename Lab>
class Symmetry : public Cubism::BC::Base<Lab>
#include <Symmetry.h>

Symmetry BC (reflecting)

Symmetry/reflecting boundary condition

Template Parameters
  • Lab: Type of FieldLab

Block

namespace Cubism::Block

Namespace for block data types.

Typedefs

using CellField = Field<T, EntityType::Cell, DIM, State, Alloc>

Basic cell-centered data field.

Template Parameters
  • T: Data type (must be POD)

  • DIM: Field dimension

  • State: State type

  • Alloc: Allocator type

using NodeField = Field<T, EntityType::Node, DIM, State, Alloc>

Basic node-centered data field.

Template Parameters
  • T: Data type (must be POD)

  • DIM: Field dimension

  • State: State type

  • Alloc: Allocator type

using FaceField = Field<T, EntityType::Face, DIM, State, Alloc>

Basic face-centered data field.

Faces are stored individually for the dimensionality specified by CUBISM_DIMENSION at compile time. See the FaceContainer type for a container of size CUBISM_DIMENSION.

Template Parameters
  • T: Data type (must be POD)

  • DIM: Field dimension

  • State: State type

  • Alloc: Allocator type

using VectorField = TensorField<T, 1, Entity, DIM, State, Alloc>

Convenience type for vector fields.

Template Parameters
  • T: Field data type

  • Entity: Entity type

  • DIM: Field dimension

  • State: Field state type

  • Alloc: Memory allocator

template<typename T, Cubism::EntityType Entity, size_t DIM, typename BlockAlloc = AlignedBlockAllocator<T>>
class Data : public Cubism::Block::DataBase
#include <Data.h>

Generic block data.

Data type that manages memory allocation and data access for the specified index range spanned by the data.

Note

The memory allocation for a block may be larger than the minimum required data specified by the index range.

Template Parameters
  • T: Data type

  • Entity: Entity type

  • DIM: Data dimensionality

  • BlockAlloc: Allocator type

class DataBase
#include <Data.h>

Type-less base class for a data block.

For byte-level operations on the data only.

Subclassed by Cubism::Block::Data< T, Entity, DIM, BlockAlloc >, Cubism::Block::Data< T, Entity, DIM, Alloc< T > >, Cubism::Block::Data< TField::DataType, TField::EntityType, TField::IndexRangeType::Dim, Cubism::AlignedBlockAllocator< TField::DataType > >

template<typename TField>
class FaceContainer : public Cubism::Block::FieldContainer<TField>
#include <Field.h>

Container class for all faces in a CUBISM_DIMENSION-ional problem.

For CUBISM_DIMENSION in {1,2,3}, the face field for faces with normal in the X direction can be obtained with ff[0] or ff[Cubism::Dir::X] for example, where ff is of type FaceContainer.

Template Parameters
  • TField: Face field type (scalar or tensor)

template<typename T, Cubism::EntityType Entity, size_t DIM = CUBISM_DIMENSION, typename State = FieldState, template<typename> class Alloc = AlignedBlockAllocator>
class Field : public Cubism::Block::Data<T, Entity, DIM, Alloc<T>>
#include <Field.h>

Block scalar field base class.

Generic block scalar field type used by Grid classes to compose a certain topology of block fields.

Template Parameters
  • T: Field data type

  • Entity: Entity type

  • DIM: Field dimension

  • State: Field state type

  • Alloc: Memory allocator

template<typename TField>
class FieldContainer
#include <Field.h>

Field container type.

This is an actively managed field container. Unlike std::vector, the destructors of the container elements are called upon destruction of the container. An incomplete container contains nullptr for some of its components.

Template Parameters
  • TField: Type of field

Subclassed by Cubism::Block::FaceContainer< TField >

template<typename TField>
class FieldLab : public Cubism::Block::Data<TField::DataType, TField::EntityType, TField::IndexRangeType::Dim, Cubism::AlignedBlockAllocator<TField::DataType>>
#include <FieldLab.h>

Field laboratory.

A FieldLab is an extended data structure to include ghost cells for a given stencil. Loading a lab takes care of loading the ghost cells from neighboring block fields and applies boundary conditions if present. The default is periodic if no boundary conditions are specified otherwise.

Template Parameters
  • TField: Field type to map to the lab

template<typename FieldType, size_t DIM>
struct FieldLabLoader
template<typename FieldType>
struct FieldLabLoader<FieldType, 3>
struct FieldState
#include <Field.h>

Default meta data (state) of a block field.

Minimal (default) state of a field is empty. Custom field state types my add additional state to describe meta data of a field. A field state type must define copy semantics. Components in TensorField types and FaceContainer types share one instance of a state.

template<typename T, size_t RANK, Cubism::EntityType Entity, size_t DIM = CUBISM_DIMENSION, typename State = FieldState, template<typename> class Alloc = AlignedBlockAllocator>
struct FieldTypeFactory
#include <Field.h>

Field type factory for tensor fields.

Template Parameters
  • T: Field data type

  • RANK: Tensor rank

  • Entity: Entity type

  • DIM: Field dimension

  • State: Field state type

  • Alloc: Memory allocator

template<typename T, size_t DIM, typename State, template< typename > class Alloc> Face, DIM, State, Alloc >
#include <Field.h>

Field type factory for face scalar fields.

Template Parameters
  • T: Field data type

  • DIM: Field dimension

  • State: Field state type

  • Alloc: Memory allocator

template<typename T, Cubism::EntityType Entity, size_t DIM, typename State, template<typename> class Alloc>
struct FieldTypeFactory<T, 0, Entity, DIM, State, Alloc>
#include <Field.h>

Field type factory for scalar fields.

Template Parameters
  • T: Field data type

  • Entity: Entity type

  • DIM: Field dimension

  • State: Field state type

  • Alloc: Memory allocator

template<typename T, size_t RANK, size_t DIM, typename State, template< typename > class Alloc> Face, DIM, State, Alloc >
#include <Field.h>

Field type factory for face tensor fields.

Template Parameters
  • T: Field data type

  • RANK: Tensor rank

  • DIM: Field dimension

  • State: Field state type

  • Alloc: Memory allocator

template<typename TField>
class FieldView : public TField
#include <Field.h>

Field view type.

Provides a view (shallow copy) for scalar fields, tensor fields or face field containers. The corresponding field interface is inherited.

Note

A view type never owns memory.

Template Parameters
  • TField: Field type

template<typename FContainer, Cubism::FieldClass Class, size_t RANK>
class PeriodicIndexFunctor
template<typename FContainer, Cubism::FieldClass Class, size_t RANK>
struct ScalarFieldMap
template<typename FContainer, Cubism::FieldClass Class>
struct ScalarFieldMap<FContainer, Class, 0>
template<typename FContainer> FaceContainer, 0 >
template<typename FContainer, size_t RANK> FaceContainer, RANK >
template<typename T, size_t RANK, Cubism::EntityType Entity, size_t DIM = CUBISM_DIMENSION, typename State = FieldState, template<typename> class Alloc = AlignedBlockAllocator>
class TensorField : public Cubism::Block::FieldContainer<Field<T, Entity, DIM, State, Alloc>>
#include <Field.h>

Generic tensor field.

Template Parameters
  • T: Field data type

  • RANK: Tensor rank

  • Entity: Entity type

  • DIM: Field dimension

  • State: Field state type

  • Alloc: Memory allocator

Core

namespace Cubism::Core

Namespace for Cubism core components.

Typedefs

using Index = std::ptrdiff_t
using MultiIndex = typename IndexRange<DIM>::PointType

Alias for a multi-dimensional index.

Template Parameters
  • DIM: Dimension of associated index space

Functions

template<typename T, size_t DIM>
void swap(Vector<T, DIM> &va, Vector<T, DIM> &vb) noexcept

Non-STL swap function for vectors.

template<typename T, size_t DIM>
std::istream &operator>>(std::istream &stream, Vector<T, DIM> &vec)

Input stream for Vector<T, DIM> types.

template<typename T, size_t DIM>
std::ostream &operator<<(std::ostream &stream, const Vector<T, DIM> &vec)

Output stream for Vector<T, DIM> types.

template<size_t DIM>
class EntityIterator : public Cubism::Core::MultiIndexIterator<DIM>
template<size_t DIM>
struct IndexConverter
#include <Index.h>

Converts indices from and to one dimensional and multi-dimensional, respectively.

Template Parameters
  • DIM: Dimension of the index space

template<>
struct IndexConverter<1>
template<>
struct IndexConverter<2>
template<>
struct IndexConverter<3>
template<size_t DIM>
class IndexRange : public Cubism::Core::Range<Index, DIM>
#include <Index.h>

Rectangular index range.

Defines a simple consecutive index space.

Template Parameters
  • DIM: Dimension of the index space

template<size_t DIM>
struct MIIForward
template<>
struct MIIForward<1>
template<>
struct MIIForward<2>
template<>
struct MIIForward<3>
template<size_t DIM>
class MultiIndexIterator

Subclassed by Cubism::Core::EntityIterator< DIM >

template<typename T, size_t DIM>
class Range
#include <Range.h>

Rectangular range.

Template Parameters
  • T: Data type that describes coordinates in the range

  • DIM: Dimension

template<size_t DIM>
class Stencil
#include <Stencil.h>

Describes a stencil.

Template Parameters
  • DIM: Stencil dimensionality

template<typename T, size_t DIM>
class Vector
#include <Vector.h>

Vector class with support for common operations.

Wraps around std::array. The data type T must follow the rules of aggregate initialization. The vector dimension DIM should be low-dimensional when used for automatic variables on the stack.

Template Parameters
  • T: Data type

  • DIM: Dimension

Grid

namespace Cubism::Grid

Namespace for grid data types composed of block data types.

template<typename T, size_t RANK, Cubism::EntityType Entity, typename State, typename Mesh>
struct BlockFieldAssembler
#include <BlockFieldAssembler.h>

Block field assembler for an externally allocated region of memory.

Template Parameters
  • T: Field data type

  • RANK: Tensor rank

  • Entity: Entity type

  • State: Field state type

  • Mesh: Mesh type of global and block meshes

template<typename T, typename Mesh, Cubism::EntityType Entity = Cubism::EntityType::Cell, size_t RANK = 0, typename UserState = Block::FieldState, template<typename> class Alloc = AlignedBlockAllocator>
class Cartesian
#include <Cartesian.h>

Cartesian block (tensor) field.

Cartesian topology composed of block Field.h for the specified entity type. As opposed to an individual block Field.h, this class manages a structure of arrays (SoA) memory layout for all the blocks in the Cartesian topology instead of just individual blocks. See the CartesianMPI.h grid section for a distributed variant of this class. The field state can be extended with the UserState extension. The UserState type must be trivially copyable.

Template Parameters
  • T: Field data type

  • Mesh: Mesh type to be associated with fields

  • Entity: Entity type

  • RANK: Rank of (tensor) fields

  • UserState: Type for field state user extension

  • Alloc: Allocator for field data

Subclassed by Cubism::Grid::CartesianMPI< T, Mesh, Entity, RANK, UserState, Alloc >

template<typename T, typename Mesh, Cubism::EntityType Entity = Cubism::EntityType::Cell, size_t RANK = 0, typename UserState = Block::FieldState, template<typename> class Alloc = AlignedBlockAllocator>
class CartesianMPI : public Cubism::Grid::Cartesian<T, Mesh, Entity, RANK, UserState, Alloc>
#include <CartesianMPI.h>

Cartesian MPI block (tensor) field.

Cartesian topology composed of block Field.h for the specified entity type. As opposed to an individual block Field.h, this class manages a structure of arrays (SoA) memory layout for all the blocks in the rank local Cartesian topology instead of just individual blocks. See the Cartesian.h grid section for a non-distributed variant of this class as well as the UserState extension.

Template Parameters
  • T: Field data type

  • Mesh: Mesh type to be associated with fields

  • Entity: Entity type

  • RANK: Rank of (tensor) fields

  • Alloc: Allocator for field data

IO

namespace Cubism::IO

Namespace for input/output operations.

Input/Output routines are functions that can be used with a block Field.h object or block field compounds described in Grid. The components of this namespace form the contents of libCubismIO.a. Consequently an application must link to -lCubismIO to make use of I/O group.

Functions

template<typename FileDataType, typename Grid, typename Mesh, typename Dir = size_t>
void CartesianWriteHDF(const std::string &fname, const std::string &aname, const Grid &grid, const Mesh &mesh, const double time, const Dir face_dir = 0, const bool create_xdmf = true)

Write Cartesian grid data to HDF file.

Write the data carried by grid to an HDF5 container file. The data that is written to the file is specified by the index space described in mesh.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Mesh: Mesh type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Output full filename without file extension

  • aname: Name of quantity in grid

  • grid: Input grid

  • mesh: Input mesh corresponding to the extracted data

  • time: Current time

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

  • create_xdmf: Flag for XDMF wrapper

template<typename FileDataType, typename Grid, typename Dir = size_t>
void CartesianWriteHDF(const std::string &fname, const std::string &aname, const Grid &grid, const double time, const Dir face_dir = 0, const bool create_xdmf = true)

Write Cartesian grid data to HDF file.

Convenience wrapper to dump a full grid to an HDF container file.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Output full filename without file extension

  • aname: Name of quantity in grid

  • grid: Input grid

  • time: Current time

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

  • create_xdmf: Flag for XDMF wrapper

template<typename FileDataType, typename Grid, typename Mesh, typename Dir = size_t>
void CartesianReadHDF(const std::string &fname, Grid &grid, const Mesh &mesh, const Dir face_dir = 0)

Read Cartesian grid data from HDF file.

Read the data of an HDF5 container file into grid. The data that is read from the file is specified by the index space described in mesh.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Mesh: Mesh type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Input full filename without file extension

  • grid: Grid populated with file data

  • mesh: Grid (sub)mesh

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

template<typename FileDataType, typename Grid, typename Dir = size_t>
void CartesianReadHDF(const std::string &fname, Grid &grid, const Dir face_dir = 0)

Read Cartesian grid data from HDF file.

Convenience wrapper to read a full grid from an HDF container file.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Input full filename without file extension

  • grid: Grid populated with file data

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

template<typename FileDataType, typename Grid, typename Mesh, typename Dir = size_t>
void CartesianMPIWriteHDF(const std::string &fname, const std::string &aname, const Grid &grid, const Mesh &mesh, const double time, const Dir face_dir = 0, const bool create_xdmf = true)

Write Cartesian MPI grid data to HDF file.

Write the data carried by the MPI grid to an HDF5 container file. The data that is written to the file is specified by the index space described in mesh.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Mesh: Mesh type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Output full filename without file extension

  • aname: Name of quantity in grid

  • grid: Input grid

  • mesh: Input mesh corresponding to the extracted data

  • time: Current time

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

  • create_xdmf: Flag for XDMF wrapper

template<typename FileDataType, typename Grid, typename Dir = size_t>
void CartesianMPIWriteHDF(const std::string &fname, const std::string &aname, const Grid &grid, const double time, const Dir face_dir = 0, const bool create_xdmf = true)

Write Cartesian MPI grid data to HDF file.

Convenience wrapper to dump a full MPI grid to an HDF container file.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Output full filename without file extension

  • aname: Name of quantity in grid

  • grid: Input grid

  • time: Current time

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

  • create_xdmf: Flag for XDMF wrapper

template<typename FileDataType, typename Grid, typename Mesh, typename Dir = size_t>
void CartesianMPIReadHDF(const std::string &fname, Grid &grid, const Mesh &mesh, const Dir face_dir = 0)

Read Cartesian MPI grid data from HDF file.

Read the data of an HDF5 container file into the MPI grid. The data that is read from the file is specified by the index space described in mesh.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Mesh: Mesh type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Input full filename without file extension

  • grid: Grid populated with file data

  • mesh: Grid (sub)mesh

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

template<typename FileDataType, typename Grid, typename Dir = size_t>
void CartesianMPIReadHDF(const std::string &fname, Grid &grid, const Dir face_dir = 0)

Read Cartesian grid data from HDF file.

Convenience wrapper to read a full MPI grid from an HDF container file.

Template Parameters
  • FileDataType: HDF file data type

  • Grid: Grid type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Input full filename without file extension

  • grid: Grid populated with file data

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

template<typename Field, typename Buffer>
void Field2AOS(const Field &f, const typename Field::IndexRangeType &r, Buffer *buf, const size_t dface = 0)

Write field data into AoS buffer.

Copy the data from a structure of arrays (SoA) field into an array of structures (AoS) buffer for I/O operation. This is a low-level function which can be used in high-level I/O interfaces. The index range r may describe a sub-region of the index range spanned by the field f. The size of the output buffer buf is determined by the index range r and Field::NComponents.

Template Parameters
  • Field: Field type

  • Buffer: Data type of AoS buffer

Parameters
  • f: Input field

  • r: Index space for the copy (describes the memory region of buf)

  • buf: Output buffer

  • dface: Face direction (relevant for Cubism::EntityType::Face)

template<typename Field, typename Buffer>
void AOS2Field(const Buffer *buf, const typename Field::IndexRangeType &r, Field &f, const size_t dface = 0)

Read AoS buffer into field data.

Copy the data from an array of structures (AoS) buffer into a structure of arrays field. This is a low-level function which can be used in high-level I/O interfaces. The index range r may describe a sub-region of the index range spanned by the field f. The size of the input buffer buf is determined by the index range r and Field::NComponents.

Template Parameters
  • Field: Field type

  • Buffer: Data type of AoS buffer

Parameters
  • buf: Input buffer

  • r: Index space for the read (describes the memory region of buf)

  • f: Output field

  • dface: Face direction (relevant for Cubism::EntityType::Face)

template<typename FileDataType, typename Field, typename Mesh, typename Dir = size_t>
void FieldWriteHDF(const std::string &fname, const std::string &aname, const Field &field, const Mesh &mesh, const double time, const Dir face_dir = 0, const bool create_xdmf = true)

Write field data to HDF file.

Write the data carried by field to an HDF5 container file. The data that is written to the file is specified by the index space described in mesh.

Note

The variables field and mesh are only related by the index space they span. If there is no common intersection, no data will be written to the file. If the mesh spans a larger index space and full or partially contains the index space spanned by field then the index space spanned by mesh is clipped to that spanned by field.

Example:

#include "Cubism/Block/Field.h"
#include "Cubism/Core/Index.h"
#include "Cubism/Mesh/StructuredUniform.h"
#include "Cubism/IO/FieldHDF.h"

using Mesh = Cubism::Mesh::StructuredUniform<double, 3>;
using MeshIntegrity = typename Mesh::MeshIntegrity;
using PointType = typename Mesh::PointType;
using IRange = typename Mesh::IndexRangeType;
using MIndex = typename Mesh::MultiIndex;

using CellField = Cubism::Block::CellField<double, Mesh::Dim>;

int main(void)
{
    // mesh in [0,1]^3 with 16^3 cells
    const PointType end(1);
    const MIndex cells(16);
    Mesh m(end, cells, MeshIntegrity::FullMesh);

    // cell field spans an index space [6,10)^3
    const IRange cell_range(6, 10);
    CellField cf(cell_range);

    // write field to file (since the mesh fully contains cell_range, the
    // file contains 64 doubles at the coordinates defined in m using the
    // indices in cell_range
    Cubism::IO::FieldWriteHDF<double>(
        "filename", "attributename", cf, m, 0.0);
    return 0;
}

Template Parameters
  • FileDataType: HDF file data type

  • Field: Field type

  • Mesh: Mesh type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Output full filename without file extension

  • aname: Name of quantity in field

  • field: Input field

  • mesh: Input mesh

  • time: Current time

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

  • create_xdmf: Flag for XDMF wrapper

template<typename FileDataType, typename Field, typename Mesh, typename Dir = size_t>
void FieldReadHDF(const std::string &fname, Field &field, const Mesh &mesh, const Dir face_dir = 0)

Read field data from HDF file.

Read the data of an HDF5 container file into field. The data that is read from the file is specified by the index space described in mesh.

Template Parameters
  • FileDataType: HDF file data type

  • Field: Field type

  • Mesh: Mesh type

  • Dir: Special type that defines a cast to size_t

Parameters
  • fname: Input full filename without file extension

  • field: Field populated with file data

  • mesh: Field (sub)mesh

  • face_dir: Face direction (relevant for Cubism::EntityType::Face)

template<Cubism::FieldClass Class>
struct AOSDriver
template<> FaceContainer >
template<> Scalar >
template<> Tensor >
template<typename FileDataType, typename Mesh, Cubism::MeshClass Class>
struct HDFDriver
#include <HDFDriver.h>

HDF read/write interface.

Template Parameters
  • FileDataType: File data taype

  • Mesh: Mesh type

  • Class: Mesh class

template<typename FileDataType, typename Mesh, Cubism::MeshClass Class>
struct HDFDriverMPI
#include <HDFDriver.h>

HDF MPI read/write interface.

Template Parameters
  • FileDataType: File data taype

  • Mesh: Mesh type

  • Class: Mesh class

template<typename DataType, Cubism::MeshClass Class>
struct XDMFDriver
template<typename DataType> Uniform >

Mesh

namespace Cubism::Mesh

Namespace for mesh data types.

template<typename TReal, size_t DIM>
class StructuredBase
#include <StructuredBase.h>

Structured mesh base class.

Defines the (pure virtual) interface for a structured mesh.

Template Parameters
  • TReal: Float type for mesh entities

  • DIM: Mesh dimension

Subclassed by Cubism::Mesh::StructuredUniform< TReal, DIM >

template<typename TReal, size_t DIM>
class StructuredUniform : public Cubism::Mesh::StructuredBase<TReal, DIM>
#include <StructuredUniform.h>

Structured uniform mesh.

Template Parameters
  • TReal: Float type for mesh entities

  • DIM: Mesh dimension

Util

namespace Cubism::Util

Namespace for input/output operations.

The members of this namespace are optional utilities that are not required to implement a working application. Its components form the content of libCubismUtil.a. An application using these utilities must link to -lCubismUtil.

Functions

std::ostream &operator<<(std::ostream &os, const INIParser &p)

Parameters
  • os: Output stream

  • p: INI parser instance

Variables

const char *CubismVersion

Cubism version at build time. Falls back to VERSION for a release.

Command:

git describe --long --dirty --broken

Example: v1.1.8-5-g824e676-dirty

The string starts with the version number, followed by the number of commits ahead of tagged commit (5 in this case), followed by g and the short hash (SHA-1) of HEAD, followed by the state of the working tree (if any).

class Histogram : public Cubism::Util::Sampler
#include <Histogram.h>

MPI profiling using histograms.

Collects samples for a profiled quantity of interest on individual ranks. Can be used to detect inhomogeneities among MPI ranks.

class INIParser
#include <INIParser.h>

INI config file parser.

Simple configuration file parser for use on the application level. CubismNova does not depend on this parser. The .ini file format is explained in more detail on Wikipedia for example. The python configparser module supports this format for convenient simulation case setup using python.

class Profiler : public Cubism::Util::Sampler
#include <Profiler.h>

Runtime profiler.

Used to collect runtime samples for a code section that is enclosed by the push() and pop() methods. Used for profiling.

class Sampler
#include <Sampler.h>

Sample collector.

Used to collect time samples (by default) for a code section that is enclosed by the seedSample() and collectSample() methods. Used for profiling.

Subclassed by Cubism::Util::Histogram, Cubism::Util::Profiler

class Timer
#include <Timer.h>

Simple timer class.