/* -------------------------------------- *\ LAPI Active Message Completion Handler \* -------------------------------------- */ void DDI_AM_Compl_hndlr(lapi_handle_t *handl,void *param) { size_t nbytes; char *stack = (char *) param; DDI_Patch *patch = (DDI_Patch *) param; double *buffer = NULL; MAX_DEBUG((stdout,"%s: Entering LAPI Active Message Completion Handler.\n",DDI_Id())) nbytes = sizeof(DDI_Patch); nbytes += (nbytes % sizeof(double)); stack += nbytes; buffer = (double *) stack; nbytes += patch->ds_buffer_size; switch(patch->oper) { case DDI_GET: DDI_Get_lapi_server(patch,buffer); break; case DDI_PUT: DDI_Put_lapi_server(patch,buffer); break; case DDI_ACC: DDI_Acc_lapi_server(patch,buffer); break; default: fprintf(stdout,"%s: unknown operation in data server handler.\n",DDI_Id()); Fatal_error(911); } if(patch->oper == DDI_ACC) DDI_Fence_release(patch->handle); DDI_Memory_heap_free(param,nbytes); MAX_DEBUG((stdout,"%s: Leaving LAPI AM Completion Handler.\n",DDI_Id())) }
void SMP_destroy(int handle) { SMP_Data *node = (SMP_Data *) gv(smp_base_data); SMP_Data *prev = (SMP_Data *) gv(smp_base_data); DEBUG_OUT(LVL3,(stdout,"%s: smp_destroy - %i.\n",DDI_Id(),handle)) if(node == NULL) { fprintf(stdout,"%s: Warning no SMP arrays to destroy.\n",DDI_Id()); Fatal_error(911); } while(node) { if(node->handle == handle) { # if defined USE_SYSV MAX_DEBUG((stdout,"%s: detaching from shm addr %x, removing semid=%i.\n", DDI_Id(),node->addr,node->semid)) shmdt(node->addr); DDI_Sem_remove(node->semid); # else free(node->addr); # endif if(node == gv(smp_base_data)) gv(smp_base_data) = (SMP_Data *) node->next; else prev->next = node->next; free(node); return; } prev = (SMP_Data *) node; node = (SMP_Data *) node->next; } fprintf(stdout,"%s: Shared Memory handle %i was not found.\n",DDI_Id(),handle); Fatal_error(911); }
int SMP_create(size_t size) { int semflg = 0600; DDI_Comm *comm = (DDI_Comm *) Comm_find(DDI_COMM_WORLD); /* hardcoded atm */ SMP_Data *new_data = (SMP_Data *) Malloc(sizeof(SMP_Data)); SMP_Data *end_data = (SMP_Data *) SMP_find_end(); STD_DEBUG((stdout,"%s: Entered DDI_SMP_Create.\n",DDI_Id())) if(end_data) end_data->next = (void *) new_data; else gv(smp_base_data) = (SMP_Data *) new_data; # if defined USE_SYSV Comm_sync_smp(comm); if(comm->me_local == 0) { new_data->handle = gv(smp_data_id)++; new_data->shmid = gv(shmid) = Shmget(IPC_PRIVATE,size,SHM_R|SHM_W); new_data->semid = Semget(IPC_PRIVATE,1,semflg); new_data->size = size; new_data->next = NULL; } Comm_bcast_smp(new_data,sizeof(SMP_Data),0,comm); new_data->addr = Shmat(new_data->shmid,0,0); MAX_DEBUG((stdout,"%s: SMP memory [%i] shmid=%i, semid=%i, addr=%x.\n", DDI_Id(),new_data->handle,new_data->shmid,new_data->semid,new_data->addr)) Comm_sync_smp(comm); if(comm->me_local == 0) { Shmctl(new_data->shmid,IPC_RMID,NULL); gv(shmid) = 0; } Comm_sync_smp(comm); # else new_data->handle = gv(smp_data_id)++; new_data->size = size; new_data->next = NULL; new_data->addr = Malloc(size); /* MWS: May 2010 It appears above that systems without SysV memory are expected to allocate Process-replicated memory instead of Node-replicated, and get on with it. If these are duplicated, at full size, as it appears, that's likely devastating for the system total memory usage. The parallel CCSD(T) on IBM Blue Gene/P got into a deadlock, but other systems with sockets or MPI seem to work if allowed to proceed. At this time, we kill off just the BG here... */ # ifdef IBMBG fprintf(stdout,"DDI compiled w/o SysV operating system support.\n"); fprintf(stdout,"IBM/BG parallel CCSD(T) cannot run w/o SysV.\n"); Fatal_error(911); # endif # endif return new_data->handle; }
/* ---------------------------- *\ FORTRAN Wrapper for DDI_Init \* ---------------------------- */ void F77_Init() { int i,j,lenmax=256,argc=iargc_(); char **argv = NULL; char arg[256]; STD_DEBUG((stdout," DDI: Entering F77 DDI_Init.\n")) /* -------------------------- *\ Get command line arguments \* -------------------------- */ if(argc) { argc++; argv = Malloc(argc*sizeof(char*)); for(i=0; i<argc; i++) { for(j=0; j<256; j++) arg[j]=' '; # if defined CRAY getarg_(&i,arg,&lenmax); # else getarg_(&i,arg,lenmax); # endif for(j=0; j<256 && arg[j] != ' '; j++); arg[j] = 0; argv[i] = (char *) strdup(arg); } } MAX_DEBUG((stdout," DDI: Calling DDI_Init.\n")) /* -------------- *\ Initialize DDI \* -------------- */ DDI_Init(argc,argv); }
/* -------------------------------------------------------------- *\ DDI_GetAccP(handle,patch,buff) ============================ [IN] handle - Handle of the distributed array to be accessed. [IN] patch - structure containing ilo, ihi, jlo, jhi, etc. [IN] buff - Data segment to be operated on. \* -------------------------------------------------------------- */ void DDI_GetAccP(int handle,DDI_Patch *patch,void *buff) { /* --------------- *\ Local Variables \* --------------- */ char ack=57; int i,np,me,nn,my,remote_id,nsubp; int ranks[MAX_NODES]; DDI_Patch subp[MAX_NODES]; char *working_buffer = (char *) buff; # if defined DDI_LAPI DDI_Patch *local_patch = NULL; lapi_cntr_t cntr[MAX_NODES]; # endif STD_DEBUG((stdout,"%s: Entering DDI_GetAccP.\n",DDI_Id())) /* -------------------- *\ Process OR Node Rank \* -------------------- */ DDI_NProc(&np,&me); DDI_NNode(&nn,&my); /* ------------------------------------- *\ Ensure the patch has the correct info \* ------------------------------------- */ patch->oper = DDI_GETACC; patch->handle = handle; /* ---------------------------------- *\ Check calling arguments for errors \* ---------------------------------- */ # if defined DDI_CHECK_ARGS if(handle < 0 || handle >= gv(ndda)) { fprintf(stdout,"%s: Invalid handle [%i] in DDI_GetAcc.\n",DDI_Id(),handle); Fatal_error(911); } if(patch->ilo > patch->ihi || patch->ilo < 0 || patch->ihi >= gv(nrow)[handle]) { fprintf(stdout,"%s: Invalid row dimensions during DDI_GetAcc => ilo=%i ihi=%i.\n",DDI_Id(),patch->ilo,patch->ihi); Fatal_error(911); } if(patch->jlo > patch->jhi || patch->jlo < 0 || patch->jhi >= gv(ncol)[handle]) { fprintf(stdout,"%s: Invalid colum dimensions during DDI_GetAcc => jlo=%i jhi=%i.\n",DDI_Id(),patch->jlo,patch->jhi); Fatal_error(911); } # endif /* ------------------------------ *\ Log some simple profiling info \* ------------------------------ */ # if defined DDI_COUNTERS gv(acc_profile).ncalls++; gv(acc_profile).nbytes += DDI_Patch_sizeof(patch); # endif /* ------------------------------------------------------- *\ Determine where the pieces of the requested patch exist \* ------------------------------------------------------- */ DDI_Subpatch(handle,patch,&nsubp,ranks,subp); MAX_DEBUG((stdout,"%s: %i subpatches.\n",DDI_Id(),nsubp)) /* ------------------------------------------------------------------- *\ Send data requests for all non-local pieces of the requested patch. Operate immediately to GetAcc a local portion of the patch. \* ------------------------------------------------------------------- */ for(i=0; i<nsubp; i++) { ULTRA_DEBUG((stdout,"%s: GetAccumulating subpatch %i.\n",DDI_Id(),i)) /* ------------------------------------------------------------- *\ Using SysV, take advantage of shared-memory for a local patch \* ------------------------------------------------------------- */ # if defined USE_SYSV /* ------------------------------------------------ *\ Determine if the ith patch is local to 'my' node \* ------------------------------------------------ */ if(ranks[i] == my) { MAX_DEBUG((stdout,"%s: Subpatch %i is local.\n",DDI_Id(),i)) /* ---------------------------------------------------- *\ Using LAPI, perform the local Getacc after all the data requests have been sent ==> maximize concurrency. \* ---------------------------------------------------- */ # if defined DDI_LAPI local_patch = &subp[i]; local_patch->cp_buffer_addr = working_buffer; # else /* --------------------------------------------- *\ Otherwise, perform the local Getacc immediately. \* --------------------------------------------- */ DDI_GetAcc_local(&subp[i],working_buffer); # endif /* ------------------------------------------------------- *\ Move the working buffer to the next patch and continue. \* ------------------------------------------------------- */ working_buffer += subp[i].size; continue; } # endif /* --------------------------------- *\ If the current patch is NOT local \* --------------------------------- */ remote_id = ranks[i]; /* ----------------------------------------------- *\ Using LAPI, then include some extra information \* ----------------------------------------------- */ # if defined DDI_LAPI subp[i].cp_lapi_id = gv(lapi_map)[me]; subp[i].cp_lapi_cntr = (void *) &cntr[i]; subp[i].cp_buffer_addr = (void *) working_buffer; LAPI_Setcntr(gv(lapi_hnd),&cntr[i],0); ULTRA_DEBUG((stdout,"%s: cp_lapi_id=%i.\n",DDI_Id(),gv(lapi_map)[me])) ULTRA_DEBUG((stdout,"%s: cp_lapi_cntr=%x.\n",DDI_Id(),&cntr[i])) ULTRA_DEBUG((stdout,"%s: cp_buffer_addr=%x.\n",DDI_Id(),working_buffer)) # endif /* -------------------------------- *\ Send data request for subpatch i \* -------------------------------- */ MAX_DEBUG((stdout,"%s: Sending data request to node %i.\n",DDI_Id(),remote_id)) DDI_Send_request(&subp[i],&remote_id,NULL); MAX_DEBUG((stdout,"%s: data request sent to global process %i.\n",DDI_Id(),remote_id)) /* ------------------------------------------------------------ *\ Receive an acknowledgement that the data server has raised a fence that will protect the distributed array from get or put access until all accumulates have finished. This block- ing receive ensures that the current process executing this accumulate can *NOT* finish, until the fence has been raised \* ------------------------------------------------------------ */ # if !defined DDI_LAPI # if defined USE_SYSV MAX_DEBUG((stdout,"%s: Receiving remote fence ACK.\n",DDI_Id())) DDI_Recv(&ack,1,remote_id); # endif /* ---------------------------- *\ Recv subpatch from remote_id \* ---------------------------- */ MAX_DEBUG((stdout,"%s: Sending subpatch %i to %i.\n",DDI_Id(),i,remote_id)) DDI_Send(working_buffer,subp[i].size,remote_id); DDI_Recv(working_buffer,subp[i].size,remote_id); # endif /* ------------ *\ Shift buffer \* ------------ */ working_buffer += subp[i].size; } /* ----------------------------------------------------------- *\ Using LAPI, perform the local Getaccumulate (if needed) as the remote processes are getting the data to Getaccumulate on the target processes. Then wait for all the data to be copied out of the buffer before returning. \* ----------------------------------------------------------- */ # if defined DDI_LAPI /* ------------------------------------ *\ GetAccumulating local patch (if exists) \* ------------------------------------ */ if(local_patch) DDI_GetAcc_local(local_patch,local_patch->cp_buffer_addr); /* ---------------------------------------------------------- *\ Wait for all remote LAPI_Gets to finish copying local data \* ---------------------------------------------------------- */ for(i=0; i<nsubp; i++) { if(subp[i].cp_lapi_cntr) { ULTRA_DEBUG((stdout,"%s: Wait for subpatch %i to be copied.\n",DDI_Id(),i)) LAPI_Waitcntr(gv(lapi_hnd),&cntr[i],3,NULL); ULTRA_DEBUG((stdout,"%s: Subpatch %i copy completed.\n",DDI_Id(),i)) } }