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.
January 17, 2008
MBlock Homepage
Charm Homepage