The Multiblock framework is accessed from a program via a set of routines. These routines are available in both C and Fortran90 versions. The C versions are all functions, and always return an error code of MBLK_SUCCESS or MBLK_FAILURE. The Fortran90 versions are all subroutines, and take an extra integer parameter ``err'' which will be set to MBLK_SUCCESS or MBLK_FAILURE.
int MBLK_Set_prefix(const char *prefix);
subroutine MBLK_Set_prefix(prefix,err)
character*, intent(in)::prefix
integer, intent(out)::err
This function is called to set the block filename prefix.
For example, if the input block files are named ``gridX00001.mblk''
and ``gridX00002.mblk'', the prefix is the string ``gridX''.
int MBLK_Set_nblocks(const int n);
subroutine MBLK_Set_nblocks(n,err)
integer, intent(in)::n
integer, intent(out)::err
This call is made to set the number of partitioned blocks to be used.
Each block is read from an input file and a separate driver
is spawned for each. The number of blocks determines the available
parallelism; so be sure to have at least as many blocks as processors.
We recommend using several times more blocks than processors, to ease
load balancing and allow adaptive overlap of computation and communication.
Be sure to set the number of blocks equal to the number of virtual processors (+vp command-line option).
int MBLK_Set_dim(const int n);
subroutine MBLK_Set_dim(n, err)
integer, intent(in)::n
integer, intent(out)::err
This call is made to set the number of spatial dimensions.
Only three dimensional computations are currently supported.
int MBLK_Get_nblocks(int* n);
subroutine MBLK_Get_nblocks(n,err)
integer,intent(out)::n
integer,intent(out)::err
Get the total number of blocks in the current computation.
Can only be called from the driver routine.
int MBLK_Get_myblock(int* m);
subroutine MBLK_Get_myblock(m,err)
integer,intent(out)::m
integer,intent(out)::err
Get the id of the current block, an integer from 0 to the number of blocks minus one.
Can only be called from the driver routine.
int MBLK_Get_blocksize(int* dims);
subroutine MBLK_Get_blocksize(dimsm,err)
integer,intent(out)::dims(3)
integer,intent(out)::err
Get the interior dimensions of the current block, in voxels.
The size of the array dims should be 3, and will be filled with
the , , and dimensions of the block.
Can only be called from the driver routine.
int MBLK_Get_nodelocs(const int* nodedim,double *nodelocs);
subroutine MBLK_Get_blocksize(nodedim,nodelocs,err)
integer,intent(in)::nodedims(3)
double precision,intent(out)::nodedims(3,nodedims(0),nodedims(1),nodedims(2))
integer,intent(out)::err
Get the locations of the nodes of the current block.
The 3-array nodedim should be the number of nodes you expect,
which must be exactly one more than the number of interior voxels.
|
|
You cannot obtain the locations of ghost nodes via this routine. To get the locations of ghost nodes, create a node-centered field containing the node locations and do an update field. Can only be called from the driver routine.
double MBLK_Timer(void);
function double precision :: MBLK_Timer()
Return the current wall clock time, in seconds. Resolution is machine-dependent, but is at worst 10ms.
void MBLK_Print_block(void);
subroutine MBLK_Print_block()
Print a debugging representation of the framework's information
about the current block.
void MBLK_Print(const char *str);
subroutine MBLK_Print(str)
character*, intent(in) :: str
Print the given string, prepended by the block id if called from the
driver. Works on all machines; unlike printf or
print *, which may not work on all parallel machines.
The Multiblock framework handles the exchange of boundary values between neighboring blocks. The basic mechanism to do this exchange is the field- numeric data items associated with each cell of a block. These items must be arranged in a regular 3D grid; but otherwise we make no assumptions about the meaning of a field.
You create a field once, with MBLK_Create_Field, then pass the resulting field ID to MBLK_Update_Field (which does the overlapping block communication) and/or MBLK_Reduce_Field (which applies a reduction over block values).
int MBLK_Create_Field(int *dimensions,int isVoxel,const int base_type,const int vec_len,const int offset,const int dist, int *fid);
subroutine MBLK_Create_Field(dimensions, isVoxel,base_type, vec_len, offset, dist, err)
integer, intent(in) :: dimensions, isVoxel, base_type, vec_len, offset, dist
integer, intent(out) :: fid, err
Creates and returns a Multiblock field ID, which can be passed to
MBLK_Update_Field and MBLK_Reduce_Field. Can only be called from
driver().
dimensions describes the size of the array the field is in. Dimensions is itself an array of size 3, giving the , , and sizes. The size should include the ghost regions- i.e., pass the actual allocated size of the array. isVoxel describes whether the data item is to be associated with a voxel (1, a volume-centered value) or the nodes (0, a node-centered value). base_type describes the type of each data item, one of:
vec_len describes the number of data items associated with each cell, an integer at least 1.
offset is the byte offset from the start of the array to the first interior cell's data items, a non-negative integer. This can be calculated using the offsetof() function; normally with offsetof(array(1,1,1),array(interiorX,interiorY,interiorZ)). Be sure to skip over any ghost regions.
dist is the byte offset from the first cell's data items to the second, a positive integer (normally the size of the data items). This can also be calculated using offsetof(); normally with offsetof(array(1,1,1),array(2,1,1)).
fid is the identifier for the field that is created by the function.
In the example below, we register a single double-precision value with
each voxel. The ghost region is 2 cells deep along all sides.
void MBLK_Update_field(const int fid,int ghostwidth, void *grid);
subroutine MBLK_Update_field(fid,ghostwidth, grid,err)
integer, intent(in) :: fid, ghostwidth
integer,intent(out) :: err
varies, intent(inout) :: grid
Update the values in the ghost regions specified when the field was created. This call sends this block's interior region out, and receives this block's boundary region from adjoining blocks.
Ghostwidth controls the thickness of the ghost region. To only exchange one cell on the boundary, pass 1. To exchange two cells, pass 2. To include diagonal regions, make the ghost width negative. A ghost width of zero would communicate no data.
|
MBLK_Update_field can only be called from driver, and to be useful, must be called from every block's driver routine.
MBLK_Update_field blocks till the field has been updated. After this routine returns, the given field will updated. If the update was successful MBLK_SUCCESS is returned and MBLK_FAILURE is returned in case of error.
void MBLK_Iupdate_field(const int fid,int ghostwidth, void *ingrid, void* outgrid);
subroutine MBLK_Iupdate_field(fid,ghostwidth, ingrid, outgrid,err)
integer, intent(in) :: fid, ghostwidth
integer,intent(out) :: err
varies,intent(in) :: ingrid
varies,intent(out) :: outgrid
Update the values in the ghost regions which were specified when the
field was created. For the example above the ghost regions will be
updated once for each step in the time loop.
MBLK_Iupdate_field can only be called from driver, and to be useful, must be called from every block's driver routine.
MBLK_Iupdate_field is a non blocking call similar to MPI_IRecv. After the routine returns the update may not yet be complete; and the outgrid may be in an inconsistent state. Before using the values the status of the update must be checked using MBLK_Test_update or MBLK_Wait_update.
There can be only one outstanding iupdate call in progress at any time.
int MBLK_Test_update(int *status);
subroutine MBLK_Test_update(status,err)
integer, intent(out) :: status,err
MBLK_Test_update is a call that is used in association with
MBLK_Iupdate_field from the driver sub routine. It tests whether
the preceeding iupdate has completed or not.
status is returned as MBLK_DONE if the update was completed or
MBLK_NOTDONE if the update is still pending.
Rather than looping if the update is still pending, call MBLK_Wait_update
to relinquish the CPU.
void MBLK_Wait_update(void);
subroutine MBLK_Wait_update()
MBLK_Wait_update call is a blocking call and is used in assoication
with MBLK_Iupdate_field call. It blocks until the update is completed.
void MBLK_Reduce_field(int fid,void *grid, void *out,int op);
subroutine MBLK_Reduce_field(fid,grid,outVal,op)
integer, intent(in) :: fid,op
varies, intent(in) :: grid
varies, intent(out) :: outVal
Combine a field from each block, according to op, across all blocks.
Only the interior values of the field will be combined; not the ghost cells.
After Reduce_Field returns, all blocks will have identical
values in outVal, which must be vec_len copies of base_type.
May only be called from driver, and to complete, must be called from every chunk's driver routine.
void MBLK_Reduce(int fid,void *inVal,void *outVal,int op);
subroutine MBLK_Reduce(fid,inVal,outVal,op)
integer, intent(in) :: fid,op
varies, intent(in) :: inVal
varies, intent(out) :: outVal
Combine a field from each block, acoording to op, across all blocks.
Fid is only used for the base_type and vec_len- offset and
dist are not used. After this call returns, all blocks will have
identical values in outVal. Op has the same values and meaning as
MBLK_Reduce_Field.
May only be called from driver, and to complete, must be called
from every blocks driver routine.
Most problems include some sort of boundary conditions. These conditions are normally applied in the ghost cells surrounding the actual computational domain. Examples of boundary conditions are imposed values, reflection walls, symmetry planes, inlets, and exits.
The Multiblock framework keeps track of where boundary conditions are to be applied. You register a subroutine that the framework will call to apply each type of external boundary condition.
int MBLK_Register_bc(const int bcnum, int ghostWidth, const MBLK_BcFn bcfn);
subroutine MBLK_Register_bc(bcnum, ghostwidth, bcfn, err)
integer,intent(in) :: bcnum, ghostWidth
integer,intent(out) :: err
subroutine :: bcfn
This is call is used to bind an external boundary condition subroutine,
written by you, to a boundary condition number.
MBLK_Register_bc should only be called from the driver.
When you ask the framework to apply boundary conditions, it will call this routine. The routine should be declared like:
param1 and param2 are not used by the framework- they are passed in unmodified from MBLK_Apply_bc and MBLK_Apply_bc_all. param1 and param2 typically contain the block data and dimensions.
start and end are 3-element arrays that give the ,, block locations where the boundary condition is to be applied. They are both inclusive and both relative to the block interior- you must shift them over your ghost cells. The C versions are 0-based (the first index is zero); the Fortran versions are 1-based (the first index is one).
For example, a Fortran subroutine to apply the constant value 1.0 across the boundary, with a 2-deep ghost region, would be:
int MBLK_Apply_bc(const int bcnum, void *param1,void *param2);
subroutine MBLK_Apply_bc(bcnum, param1,param2,err)
integer,intent(in)::bcnum
varies,intent(inout)::param1
varies,intent(inout)::param2
integer,intent(out)::err
MBLK_Apply_bc call is made to apply all boundry condition functions
of type bcnum to the block. param1 and param2 are passed unmodified
to the boundary condition function.
int MBLK_Apply_bc_all(void* param1, void* param2);
subroutine MBLK_Apply_bc_all(param1,param2, err)
integer,intent(out)::err
varies,intent(inout)::param1
varies,intent(inout)::param2
This call is same as MBLK_Apply_bc except it applies all
external boundary conditions to the block.
The CHARM++ runtime framework includes an automated, run-time load balancer, which will automatically monitor the performance of your parallel program. If needed, the load balancer can ``migrate'' mesh chunks from heavily-loaded processors to more lightly-loaded processors, improving the load balance and speeding up the program. For this to be useful, pass the +vpN argument with a larger number of blocks N than processors Because this is somewhat involved, you may refrain from calling MBLK_Migrate and migration will never take place.
The runtime system can automatically move your thread stack to the new processor, but you must write a PUP function to move any global or heap-allocated data to the new processor (global data is declared at file scope or static in C and COMMON in Fortran77; heap allocated data comes from C malloc, C++ new, or Fortran90 ALLOCATE). A PUP (Pack/UnPack) function performs both packing (converting heap data into a message) and unpacking (converting a message back into heap data). All your global and heap data must be collected into a single block (struct in C; user-defined TYPE in Fortran) so the PUP function can access it all.
Your PUP function will be passed a pointer to your heap data block and a special handle called a ``pupper'', which contains the network message to be sent. Your PUP function returns a pointer to your heap data block. In a PUP function, you pass all your heap data to routines named pup_type, where type is either a basic type (such as int, char, float, or double) or an array type (as before, but with a ``s'' suffix). Depending on the direction of packing, the pupper will either read from or write to the values you pass- normally, you shouldn't even know which. The only time you need to know the direction is when you are leaving a processor or just arriving. Correspondingly, the pupper passed to you may be deleting (indicating that you are leaving the processor, and should delete your heap storage after packing), unpacking (indicating you've just arrived on a processor, and should allocate your heap storage before unpacking), or neither (indicating the system is merely sizing a buffer, or checkpointing your values).
PUP functions are much easier to write than explain- a simple C heap block and the corresponding PUP function is:
This single PUP function can be used to copy the my_block data into a message buffer and free the old heap storage (deleting pupper); allocate storage on the new processor and copy the message data back (unpacking pupper); or save the heap data for debugging or checkpointing.
A Fortran block TYPE and corresponding PUP routine is as follows:
int MBLK_Register(void *block, MBLK_PupFn pup_ud, int* rid)
subroutine MBLK_Register(block,pup_ud, rid)
integer, intent(out)::rid
TYPE(varies), POINTER :: block
SUBROUTINE :: pup_ud
Associates the given data block and PUP function. Returns a block
ID, which can be passed to MBLK_Get_registered later. Can only be
called from driver. It returns MBLK_SUCESS if the call was successful
and MBLK_FAILURE in case of error. For the declarations above, you call
MBLK_Register as:
Note that Fortran blocks must be allocated on the stack in driver;
while C/C++ blocks may be allocated on the heap.
void MBLK_Migrate()
subroutine MBLK_Migrate()
Informs the load balancing system that you are ready to be
migrated, if needed. If the system decides to migrate you, the
PUP function passed to MBLK_Register will be called with a sizing
pupper, then a packing, deleting pupper. Your stack (and pupped
data) will then be sent to the destination machine, where your PUP
function will be called with an unpacking pupper. MBLK_Migrate
will then return, whereupon you should call MBLK_Get_registered to
get your unpacked data block. Can only be called from driver.
int MBLK_Get_Userdata(int n, void** block)
Return your unpacked userdata after migration- that is, the
return value of the unpacking call to your PUP function. Takes
the userdata ID returned by MBLK_Register. Can be called from
driver at any time.
Since Fortran blocks are always allocated on the stack, the system migrates them to the same location on the new processor, so no Get_registered call is needed from Fortran.
November 23, 2009
MBlock Homepage
Charm Homepage