/*@ BSdepend.h - This file defines all the message-passing macros for the BlockSolve95 package. System Description: All the macros can MPI routines. Note that there is blocking version and a nonblocking version based on whether NO_BLOCKING_SEND is defined. @*/ #ifndef __BSdependh #define __BSdependh /* include MPI */ #include "mpi.h" #define MPI_Aint int #include #if defined(PARCH_freebsd) #include #endif #if defined(PARCH_sun4) && !defined(__cplusplus) && defined(_Gnu_) extern int fprintf(FILE*,const char*,...); extern void sscanf(char *,...); extern int printf(const char *,...); extern void *malloc(long unsigned int ); extern void free(void *); extern int srand(unsigned int); extern int rand(); #endif #if defined(__cplusplus) extern "C" { extern void *malloc(long unsigned int ); extern void free(void *); extern int abs(int); extern double fabs(double); extern void srand(int); extern int rand(); extern double drand48(); extern void srand48(long); extern int exit(int); } #endif #define DEBUG_ALL 1 #ifdef DEBUG_ALL #define DEBUG_ACT printf("ERROR: Code %d occured at Line %d in File %s\n",__BSERROR_STATUS,__LINE__,__FILE__); #define DEBUG_ACTC(msg) printf("ERROR: Code %d occured at Line %d in File %s: %s\n",__BSERROR_STATUS,__LINE__,__FILE__,msg); #define SETERR(code_num) {__BSERROR_STATUS=code_num;DEBUG_ACT;} #define SETERRC(code_num,str) {__BSERROR_STATUS=code_num;DEBUG_ACTC(str);} #define GETERR (__BSERROR_STATUS) #define CHKERR(code_num) {if (GETERR) {DEBUG_ACT;return;}} #define CHKERRN(code_num) {if (GETERR) {DEBUG_ACT;return 0;}} #else #define SETERR(code_num) #define SETERRC(code_num,str) #define CHKERR(a) #define CHKERRN(a) #endif #ifdef TRMALLOC #define MALLOC(c) TrMalloc(c,__LINE__,__FILE__) #define FREE(c) TrFree(c,__LINE__,__FILE__) #define CHECK_MEMORY() \ { \ int ierr99; \ ierr99 = TrDump(stderr); \ } #else #define MALLOC(c) malloc(c) #define FREE(c) free(c) #define CHECK_MEMORY() #endif #define RECVSYNCNOMEM(type,buf,length,data_type,p_info,mpistat) \ MPI_Recv(buf,length,data_type,MPI_ANY_SOURCE,type,p_info->procset,&mpistat) #define RECVSYNCUNSZ(type,buf,length,data_type,p_info,mpistat) \ { \ MPI_Aint size99; \ MPI_Status mpistat99; \ MPI_Probe(MPI_ANY_SOURCE,type,p_info->procset,&mpistat99); \ MPI_Get_count(&mpistat99,data_type,&length); \ MPI_Type_size(data_type,&size99); \ MY_MALLOC(buf,(void *),size99*length,1); \ MPI_Recv(buf,length,data_type,mpistat99.MPI_SOURCE, \ mpistat99.MPI_TAG,p_info->procset,&mpistat); \ } #define RECVSYNCUNSZN(type,buf,length,data_type,p_info,mpistat) \ { \ MPI_Aint size99; \ MPI_Status mpistat99; \ MPI_Probe(MPI_ANY_SOURCE,type,p_info->procset,&mpistat99); \ MPI_Get_count(&mpistat99,data_type,&length); \ MPI_Type_size(data_type,&size99); \ MY_MALLOCN(buf,(void *),size99*length,1); \ MPI_Recv(buf,length,data_type,mpistat99.MPI_SOURCE, \ mpistat99.MPI_TAG,p_info->procset,&mpistat); \ } #define MSGFREERECV(msg) MY_FREE(msg) #define RECVASYNCNOMEMFORCE(type,buf,length,data_type,msg_id,p_info) \ MPI_Irecv(buf,length,data_type,MPI_ANY_SOURCE,type,p_info->procset,&(msg_id)) #define SENDASYNCNOMEMFORCE(type,buf,size,to_proc,data_type,msg_id,p_info) \ MPI_Irsend(buf,size,data_type,to_proc,type,p_info->procset,&(msg_id)) #define SENDASYNCNOMEM(type,buf,size,to_proc,data_type,msg_id,p_info) \ MPI_Isend((void *)buf,size,data_type,to_proc,type,p_info->procset,(MPI_Request *)(&(msg_id))) #define SENDSYNCNOMEM(type,buf,size,to_proc,data_type,p_info) \ MPI_Send((void *)buf,size,data_type,to_proc,type,(p_info->procset)) #define SENDWAITNOMEM(type,buf,size,to_proc,data_type,msg_id) \ { \ MPI_Status stat99; \ MPI_Wait(&(msg_id),&stat99); \ } #define GISUM(sum_vec,vec_len,work_vec,procset) \ { \ int i99, *iptr1 = (int *) sum_vec, *iptr2 = (int *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (iptr2)[i99] = (iptr1)[i99]; \ } \ MPI_Allreduce(iptr2,iptr1,vec_len,MPI_INT,MPI_SUM,procset); \ } #define GIOR(sum_vec,vec_len,work_vec,procset) \ { \ int i99, *iptr1 = (int *) sum_vec, *iptr2 = (int *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (iptr2)[i99] = (iptr1)[i99]; \ } \ MPI_Allreduce(iptr2,iptr1,vec_len,MPI_INT,MPI_LOR,procset); \ } #define GIMIN(sum_vec,vec_len,work_vec,procset) \ { \ int i99, *iptr1 = (int *) sum_vec, *iptr2 = (int *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (iptr2)[i99] = (iptr1)[i99]; \ } \ MPI_Allreduce(iptr2,iptr1,vec_len,MPI_INT,MPI_MIN,procset); \ } #define GIMAX(sum_vec,vec_len,work_vec,procset) \ { \ int i99, *iptr1 = (int *) sum_vec, *iptr2 = (int *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (iptr2)[i99] = (iptr1)[i99]; \ } \ MPI_Allreduce(iptr2,iptr1,vec_len,MPI_INT,MPI_MAX,procset); \ } #define GDMAX(sum_vec,vec_len,work_vec,procset) \ { \ int i99; \ double *dptr1 = (double *) sum_vec, *dptr2 = (double *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (dptr2)[i99] = (dptr1)[i99]; \ } \ MPI_Allreduce(dptr2,dptr1,vec_len,MPI_DOUBLE,MPI_MAX,procset); \ } #define GDMIN(sum_vec,vec_len,work_vec,procset) \ { \ int i99; \ double *dptr1 = (double *) sum_vec, *dptr2 = (double *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (dptr2)[i99] = (dptr1)[i99]; \ } \ MPI_Allreduce(dptr2,dptr1,vec_len,MPI_DOUBLE,MPI_MIN,procset); \ } #define GDSUM(sum_vec,vec_len,work_vec,procset) \ { \ int i99; \ double *dptr1 = (double *) sum_vec, *dptr2 = (double *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (dptr2)[i99] = (dptr1)[i99]; \ } \ MPI_Allreduce(dptr2,dptr1,vec_len,MPI_DOUBLE,MPI_SUM,procset); \ } #define GFSUM(sum_vec,vec_len,work_vec,procset) \ { \ int i99; \ float *fptr1 = (float *) sum_vec, *fptr2 = (float *) work_vec; \ for (i99=0;i99<(vec_len);i99++) { \ (fptr2)[i99] = (fptr1)[i99]; \ } \ MPI_Allreduce(fptr2,fptr1,vec_len,MPI_FLOAT,MPI_SUM,procset); \ } #define GSYNC(procset) MPI_Barrier(procset) #define PSISROOT(procinfo) ((procinfo->my_id == 0) ? 1 : 0) typedef MPI_Comm ProcSet; /* the function PSNbrTree returns the parent, left child, or right child */ /* processor number is the processors are organized conceptually as a tree */ /* op_code is one of the following */ #define PS_PARENT 0 #define PS_LCHILD 1 #define PS_RCHILD 2 #define PSNbrTree(op_code,nbr_id,procset) \ { \ int my_id99, np99; \ MPI_Comm_rank(procset,&my_id99); \ MPI_Comm_size(procset,&np99); \ switch(op_code) { \ case PS_PARENT: \ if (my_id99 == 0) { \ nbr_id = -1; \ } else { \ nbr_id = (my_id99-1) / 2; \ } \ break; \ case PS_LCHILD: \ nbr_id = (2*my_id99) + 1; \ break; \ case PS_RCHILD: \ nbr_id = (2*my_id99) + 2; \ break; \ } \ if (nbr_id >= np99) nbr_id = -1; \ } /* PICall(procedure,argc,argv) gets things started */ #define PICall(procedure,argc,argv) \ { \ MPI_Init(&argc,&argv); \ procedure(argc,argv); \ MPI_Finalize(); \ } /* set up the right "Fortran naming" definitions for use in BSsparse.h */ /* somehow the name of your architecture must be identified */ /* the following setup works for the *vast* majority of systems */ /* this affects only the linking of the blas and lapack routines called */ /* FORTRANCAPS: Names are uppercase, no trailing underscore FORTRANUNDERSCORE: Names are lowercase, trailing underscore */ #if defined(PARCH_cray) || defined(PARCH_NCUBE) || defined(PARCH_t3d) #define FORTRANCAPS #elif !defined(PARCH_rs6000) && !defined(PARCH_NeXT) && !defined(PARCH_hpux) #define FORTRANUNDERSCORE #endif /* ********************************************************************** */ /* this is a set of macros that retrofit the code for machines on which blocking sends are a bad idea. These macros make all sends asynchronous by allocating buffers for the synchronous sends and cleaning up those buffers as we go */ /* To turn this code "on" simply define NO_BLOCKING_SEND */ /* For now, the only architecture that we know needs this is the rs6000 */ /* and, perhaps, the HPs */ #if defined(PARCH_rs6000) || defined(PARCH_hpux) || defined(PARCH_t3d) || defined(PARCH_IRIX) #define NO_BLOCKING_SEND 1 #endif typedef struct __BSmsg_list { int msg_type; char *msg_buf; int msg_len; int msg_to; MPI_Datatype msg_data_type; MPI_Request msg_id; struct __BSmsg_list *next; } BSmsg_list; #ifdef NO_BLOCKING_SEND #define MY_SEND_SYNC(Mmsg_list,Mmsg_type,Mmsg,Mmsg_len,Mmsg_to,Mmsg_data_type,Mp_info) \ { \ BSmsg_list *node_99; \ int i99; \ MPI_Aint size99; \ char *tmsg_ptr99; \ MY_MALLOC(node_99,(BSmsg_list *),sizeof(BSmsg_list),1); \ node_99->next = Mmsg_list; \ Mmsg_list = node_99; \ MPI_Type_size(Mmsg_data_type,&size99); \ if (Mmsg_len == 0) { \ MY_MALLOC(node_99->msg_buf,(char *),size99,2); \ } else { \ MY_MALLOC(node_99->msg_buf,(char *),Mmsg_len*size99,3); \ } \ node_99->msg_type = Mmsg_type; \ node_99->msg_len = Mmsg_len; \ node_99->msg_to = Mmsg_to; \ node_99->msg_data_type = Mmsg_data_type; \ tmsg_ptr99 = (char *) Mmsg; \ for (i99=0;i99msg_buf[i99] = tmsg_ptr99[i99]; \ } \ SENDASYNCNOMEM(Mmsg_type,node_99->msg_buf,Mmsg_len,Mmsg_to,Mmsg_data_type,\ node_99->msg_id,Mp_info); \ CHECK_SEND_LIST(Mmsg_list); \ } #else #define MY_SEND_SYNC(Mmsg_list,Mmsg_type,Mmsg,Mmsg_len,Mmsg_to,Mmsg_data_type,Mp_info) \ { \ SENDSYNCNOMEM(Mmsg_type,Mmsg,Mmsg_len,Mmsg_to,Mmsg_data_type,Mp_info); \ } #endif #define MCHECK_SEND_LIST(Mmsg_list) \ { \ BSmsg_list *node_ptr_99, *prev_node_ptr_99, *tnode_ptr_99; \ int fin99; \ MPI_Status stat99; \ node_ptr_99 = Mmsg_list; \ prev_node_ptr_99 = Mmsg_list; \ while (node_ptr_99 != NULL) { \ MPI_Test(&(node_ptr_99->msg_id),&fin99,&stat99); \ if (fin99) { \ tnode_ptr_99 = node_ptr_99; \ if (node_ptr_99 == Mmsg_list) { \ Mmsg_list = Mmsg_list->next; \ prev_node_ptr_99 = Mmsg_list; \ node_ptr_99 = Mmsg_list; \ } else { \ prev_node_ptr_99->next = node_ptr_99->next; \ node_ptr_99 = node_ptr_99->next; \ } \ MY_FREE(tnode_ptr_99->msg_buf); \ MY_FREE(tnode_ptr_99); \ } else { \ prev_node_ptr_99 = node_ptr_99; \ node_ptr_99 = node_ptr_99->next; \ } \ } \ } #ifdef NO_BLOCKING_SEND #define CHECK_SEND_LIST(Mmsg_list) MCHECK_SEND_LIST(Mmsg_list) #else #define CHECK_SEND_LIST(Mmsg_list) #endif #define MFINISH_SEND_LIST(Mmsg_list) \ { \ BSmsg_list *tnode99; \ while(Mmsg_list != NULL) { \ tnode99 = Mmsg_list; \ Mmsg_list = Mmsg_list->next; \ SENDWAITNOMEM(tnode99->msg_type,tnode99->msg_buf,tnode99->msg_len, \ tnode99->msg_to,tnode99->msg_data_type,tnode99->msg_id); \ MY_FREE(tnode99->msg_buf); \ MY_FREE(tnode99); \ } \ } #ifdef NO_BLOCKING_SEND #define FINISH_SEND_LIST(Mmsg_list) MFINISH_SEND_LIST(Mmsg_list) #else #define FINISH_SEND_LIST(Mmsg_list) #endif #endif