fml  0.1-0
Fused Matrix Library
fml::grid Class Reference

2-dimensional MPI process grid. More...

#include <grid.hh>

Public Member Functions

 grid ()
 Create a new grid object. Does not initialize any BLACS or MPI data.
 
 grid (const gridshape gridtype)
 Create a new grid object by initializing a new BLACS process grid. More...
 
void set (const int blacs_context)
 Create a grid object from an existing BLACS process grid. More...
 
void exit ()
 Exits the BLACS grid, but does not shutdown BLACS/MPI.
 
void finalize (const bool mpi_continue=false)
 Shuts down BLACS, and optionally MPI. More...
 
void printf (const int row, const int col, const char *fmt,...) const
 Helper wrapper around the C standard I/O 'printf()' function. Conceptually similar to guarding a normal 'printf()' function with a check for 'row==myrow() && col==mycol()'. More...
 
void info () const
 Print some brief information about the BLACS grid. The printing is done by row 0 and col 0.
 
bool rank0 () const
 Check if the executing process is rank 0, i.e., if the process row and column are 0.
 
bool ingrid () const
 Check if the executing process is in the grid, i.e., if neither the process row nor column are -1.
 
void barrier (const char scope='A') const
 Execute a barrier across the specified scope of the BLACS grid. More...
 
bool valid_grid () const
 Is the BLACS grid valid?
 
void send (const int m, const int n, const int *x, const int rdest=0, const int cdest=0) const
 Point-to-point send. Should be matched by a corresponding 'recv' call. More...
 
void send (const int m, const int n, const float *x, const int rdest=0, const int cdest=0) const
 
void send (const int m, const int n, const double *x, const int rdest=0, const int cdest=0) const
 
void send (const int m, const int n, const int ldx, const int *x, const int rdest=0, const int cdest=0) const
 
void send (const int m, const int n, const int ldx, const float *x, const int rdest=0, const int cdest=0) const
 
void send (const int m, const int n, const int ldx, const double *x, const int rdest=0, const int cdest=0) const
 
void recv (const int m, const int n, int *x, const int rsrc=0, const int csrc=0) const
 Point-to-point receive. Should be matched by a corresponding 'send' call. More...
 
void recv (const int m, const int n, float *x, const int rsrc=0, const int csrc=0) const
 
void recv (const int m, const int n, double *x, const int rsrc=0, const int csrc=0) const
 
void recv (const int m, const int n, const int ldx, int *x, const int rsrc=0, const int csrc=0) const
 
void recv (const int m, const int n, const int ldx, float *x, const int rsrc=0, const int csrc=0) const
 
void recv (const int m, const int n, const int ldx, double *x, const int rsrc=0, const int csrc=0) const
 
void allreduce (const int m, const int n, int *x, const char scope='A', const blacsops op=BLACS_SUM) const
 Sum reduce operation across all processes in the grid. More...
 
void allreduce (const int m, const int n, float *x, const char scope='A', const blacsops op=BLACS_SUM) const
 
void allreduce (const int m, const int n, double *x, const char scope='A', const blacsops op=BLACS_SUM) const
 
void reduce (const int m, const int n, int *x, const char scope='A', const blacsops op=BLACS_SUM, const int rdest=0, const int cdest=0) const
 Sum reduce operation. More...
 
void reduce (const int m, const int n, float *x, const char scope='A', const blacsops op=BLACS_SUM, const int rdest=0, const int cdest=0) const
 
void reduce (const int m, const int n, double *x, const char scope='A', const blacsops op=BLACS_SUM, const int rdest=0, const int cdest=0) const
 
void bcast (const int m, const int n, int *x, const char scope='A', const int rsrc=0, const int csrc=0) const
 Broadcast. More...
 
void bcast (const int m, const int n, float *x, const char scope='A', const int rsrc=0, const int csrc=0) const
 
void bcast (const int m, const int n, double *x, const char scope='A', const int rsrc=0, const int csrc=0) const
 
int ictxt () const
 
int nprocs () const
 The total number of processes bound to the BLACS context.
 
int nprow () const
 The number of processes rows in the BLACS context.
 
int npcol () const
 The number of processes columns in the BLACS context.
 
int myrow () const
 The process row (0-based index) of the calling process.
 
int mycol () const
 The process column (0-based index) of the calling process.
 

Protected Attributes

int _ictxt
 
int _nprocs
 
int _nprow
 
int _npcol
 
int _myrow
 
int _mycol
 

Detailed Description

2-dimensional MPI process grid.

Constructor & Destructor Documentation

◆ grid()

fml::grid::grid ( const gridshape  gridtype)
inline

Create a new grid object by initializing a new BLACS process grid.

Parameters
[in]gridtypeShould be one of PROC_GRID_SQUARE, PROC_GRID_WIDE, or PROC_GRID_TALL.

Exceptions
If 'gridtype' is not one of PROC_GRID_SQUARE, PROC_GRID_WIDE, or PROC_GRID_TALL, the method will throw a 'runtime_error' exception.

Member Function Documentation

◆ allreduce()

void fml::grid::allreduce ( const int  m,
const int  n,
int *  x,
const char  scope = 'A',
const blacsops  op = BLACS_SUM 
) const
inline

Sum reduce operation across all processes in the grid.

Parameters
[in]m,nDimensions (number of rows/cols) of the data 'x'.
[in,out]xThe data to reduce.
scopeThe scope of the operation. For just rows use 'R', just columns use 'C', and for all processes use 'A'.

◆ barrier()

void fml::grid::barrier ( const char  scope = 'A') const
inline

Execute a barrier across the specified scope of the BLACS grid.

Parameters
scopeThe scope of the operation. For just rows use 'R', just columns use 'C', and for all processes use 'A'.

◆ bcast()

void fml::grid::bcast ( const int  m,
const int  n,
int *  x,
const char  scope = 'A',
const int  rsrc = 0,
const int  csrc = 0 
) const
inline

Broadcast.

Parameters
[in]m,nDimensions (number of rows/cols) of the data 'x'.
[in,out]xThe data to broadcast.
scopeThe scope of the operation. For just rows use 'R', just columns use 'C', and for all processes use 'A'.
[in]rsrc,csrcThe process row/column of the BLACS grid that is broadcasting.

◆ finalize()

void fml::grid::finalize ( const bool  mpi_continue = false)
inline

Shuts down BLACS, and optionally MPI.

Parameters
mpi_continueShould MPI continue, i.e., not be shut down too?

◆ ictxt()

int fml::grid::ictxt ( ) const
inline

The BLACS integer context.

◆ printf()

void fml::grid::printf ( const int  row,
const int  col,
const char *  fmt,
  ... 
) const
inline

Helper wrapper around the C standard I/O 'printf()' function. Conceptually similar to guarding a normal 'printf()' function with a check for 'row==myrow() && col==mycol()'.

Parameters
[in]row,colThe process row/column that should do the printing.
[in]fmtThe printf format string.
[in]...additional arguments to printf.

◆ recv()

void fml::grid::recv ( const int  m,
const int  n,
int *  x,
const int  rsrc = 0,
const int  csrc = 0 
) const
inline

Point-to-point receive. Should be matched by a corresponding 'send' call.

Parameters
[in]m,nDimensions (number of rows/cols) of the data 'x'.
[in]xThe data to receive.
[in]rsrc,csrcThe row/col source in the BLACS grid.

◆ reduce()

void fml::grid::reduce ( const int  m,
const int  n,
int *  x,
const char  scope = 'A',
const blacsops  op = BLACS_SUM,
const int  rdest = 0,
const int  cdest = 0 
) const
inline

Sum reduce operation.

Parameters
[in]m,nDimensions (number of rows/cols) of the data 'x'.
[in,out]xThe data to reduce.
scopeThe scope of the operation. For just rows use 'R', just columns use 'C', and for all processes use 'A'.
[in]rdest,cdestThe row/column of the BLACS grid to receive the final answer.

◆ send()

void fml::grid::send ( const int  m,
const int  n,
const int *  x,
const int  rdest = 0,
const int  cdest = 0 
) const
inline

Point-to-point send. Should be matched by a corresponding 'recv' call.

Parameters
[in]m,nDimensions (number of rows/cols) of the data x.
[in]ldxLeading dimension of matrix x.
[in]xThe data to send.
[in]rdest,cdestThe row/col destination in the BLACS grid.

◆ set()

void fml::grid::set ( const int  blacs_context)
inline

Create a grid object from an existing BLACS process grid.

Parameters
blacs_contextThe BLACS integer context number.

The documentation for this class was generated from the following file: