char * hypre_ReAlloc( char *ptr, int size ) { #ifdef HYPRE_USE_UMALLOC if (ptr == NULL) { ptr = hypre_MAlloc(size); } else if (size == 0) { hypre_Free(ptr); } else { int threadid = hypre_GetThreadID(); ptr = _urealloc_(ptr, size); } #else ptr = realloc(ptr, size); #endif #if 1 if ((ptr == NULL) && (size > 0)) { hypre_OutOfMemory(size); } #endif return ptr; }
char * hypre_SharedMAlloc( int size ) { char *ptr; int unthreaded = pthread_equal(initial_thread, pthread_self()); int I_call_malloc = unthreaded || pthread_equal(hypre_thread[0],pthread_self()); if (I_call_malloc) { global_alloc_ptr = hypre_MAlloc( size ); } hypre_barrier(&talloc_mtx, unthreaded); ptr = global_alloc_ptr; hypre_barrier(&talloc_mtx, unthreaded); return ptr; }
HYPRE_Int hypre_DataExchangeList(HYPRE_Int num_contacts, HYPRE_Int *contact_proc_list, void *contact_send_buf, HYPRE_Int *contact_send_buf_starts, HYPRE_Int contact_obj_size, HYPRE_Int response_obj_size, hypre_DataExchangeResponse *response_obj, HYPRE_Int max_response_size, HYPRE_Int rnum, MPI_Comm comm, void **p_response_recv_buf, HYPRE_Int **p_response_recv_buf_starts) { /*------------------------------------------- * parameters: * * num_contacts = how many procs to contact * contact_proc_list = list of processors to contact * contact_send_buf = array of data to send * contact_send_buf_starts = index for contact_send_buf corresponding to * contact_proc_list * contact_obj_size = sizeof() one obj in contact list * response_obj_size = sizeof() one obj in response_recv_buf * response_obj = this will give us the function we need to * fill the reponse as well as * any data we might need to accomplish that * max_response_size = max size of a single response expected (do NOT * need to be an absolute upper bound) * rnum = two consequentive exchanges should have different * rnums. Alternate rnum = 1 * and rnum=2 - these flags will be even (so odd * numbered tags could be used in calling code) * p_response_recv_buf = where to receive the reponses - will be allocated * in this function * p_response_recv_buf_starts = index of p_response_buf corresponding to * contact_buf_list - will be allocated here *-------------------------------------------*/ HYPRE_Int num_procs, myid; HYPRE_Int i; HYPRE_Int terminate, responses_complete; HYPRE_Int children_complete; HYPRE_Int contact_flag; HYPRE_Int proc; HYPRE_Int contact_size; HYPRE_Int size, post_size, copy_size; HYPRE_Int total_size, count; void *start_ptr = NULL, *index_ptr=NULL; HYPRE_Int *int_ptr=NULL; void *response_recv_buf = NULL; void *send_response_buf = NULL; HYPRE_Int *response_recv_buf_starts = NULL; void *initial_recv_buf = NULL; void *recv_contact_buf = NULL; HYPRE_Int recv_contact_buf_size = 0; HYPRE_Int response_message_size = 0; HYPRE_Int overhead; HYPRE_Int max_response_size_bytes; HYPRE_Int max_response_total_bytes; void **post_array = NULL; /*this must be set to null or realloc will crash */ HYPRE_Int post_array_storage = 0; HYPRE_Int post_array_size = 0; HYPRE_Int num_post_recvs =0; void **contact_ptrs = NULL, **response_ptrs=NULL, **post_ptrs=NULL; hypre_BinaryTree tree; hypre_MPI_Request *response_requests, *contact_requests; hypre_MPI_Status *response_statuses, *contact_statuses; hypre_MPI_Request *post_send_requests = NULL, *post_recv_requests = NULL; hypre_MPI_Status *post_send_statuses = NULL, *post_recv_statuses = NULL; hypre_MPI_Request *term_requests, term_request1, request_parent; hypre_MPI_Status *term_statuses, term_status1, status_parent; hypre_MPI_Status status, fill_status; const HYPRE_Int contact_tag = 1000*rnum; const HYPRE_Int response_tag = 1002*rnum; const HYPRE_Int term_tag = 1004*rnum; const HYPRE_Int post_tag = 1006*rnum; hypre_MPI_Comm_size(comm, &num_procs ); hypre_MPI_Comm_rank(comm, &myid ); /* ---------initializations ----------------*/ /* if the response_obj_size or contact_obj_size is 0, set to sizeof(HYPRE_Int) */ if (!response_obj_size) response_obj_size = sizeof(HYPRE_Int); if (!contact_obj_size) contact_obj_size = sizeof(HYPRE_Int); max_response_size_bytes = max_response_size*response_obj_size; /* pre-allocate the max space for responding to contacts */ overhead = ceil((HYPRE_Real) sizeof(HYPRE_Int)/response_obj_size); /*for appending an integer*/ max_response_total_bytes = (max_response_size+overhead)*response_obj_size; response_obj->send_response_overhead = overhead; response_obj->send_response_storage = max_response_size; /*send_response_buf = hypre_MAlloc(max_response_total_bytes);*/ send_response_buf = hypre_CAlloc(max_response_size+overhead, response_obj_size); /*allocate space for inital recv array for the responses - give each processor size max_response_size */ initial_recv_buf = hypre_MAlloc(max_response_total_bytes*num_contacts); response_recv_buf_starts = hypre_CTAlloc(HYPRE_Int, num_contacts+1); contact_ptrs = hypre_TAlloc( void *, num_contacts); response_ptrs = hypre_TAlloc(void *, num_contacts); /*-------------SEND CONTACTS AND POST RECVS FOR RESPONSES---*/ for (i=0; i<= num_contacts; i++) { response_recv_buf_starts[i] = i*(max_response_size+overhead); } /* Send "contact" messages to the list of processors and pre-post receives to wait for their response*/ responses_complete = 1; if (num_contacts > 0 ) { responses_complete = 0; response_requests = hypre_CTAlloc(hypre_MPI_Request, num_contacts); response_statuses = hypre_CTAlloc(hypre_MPI_Status, num_contacts); contact_requests = hypre_CTAlloc(hypre_MPI_Request, num_contacts); contact_statuses = hypre_CTAlloc(hypre_MPI_Status, num_contacts); /* post receives - could be confirmation or data*/ /* the size to post is max_response_total_bytes*/ for (i=0; i< num_contacts; i++) { /* response_ptrs[i] = initial_recv_buf + i*max_response_total_bytes ; */ response_ptrs[i] = (void *)((char *) initial_recv_buf + i*max_response_total_bytes) ; hypre_MPI_Irecv(response_ptrs[i], max_response_total_bytes, hypre_MPI_BYTE, contact_proc_list[i], response_tag, comm, &response_requests[i]); } /* send out contact messages */ start_ptr = contact_send_buf; for (i=0; i< num_contacts; i++) { contact_ptrs[i] = start_ptr; size = contact_send_buf_starts[i+1] - contact_send_buf_starts[i] ; hypre_MPI_Isend(contact_ptrs[i], size*contact_obj_size, hypre_MPI_BYTE, contact_proc_list[i], contact_tag, comm, &contact_requests[i]); /* start_ptr += (size*contact_obj_size); */ start_ptr = (void *) ((char *) start_ptr + (size*contact_obj_size)); } } /*------------BINARY TREE-----------------------*/ /*Now let's find out our binary tree information and initialize for the termination check sweep */ terminate = 1; /*indicates whether we can stop probing for contact */ children_complete = 1;/*indicates whether we have recv. term messages from our children*/ if (num_procs > 1) { hypre_CreateBinaryTree(myid, num_procs, &tree); /* we will get a message from all of our children when they have received responses for all of their contacts. So post receives now */ term_requests = hypre_CTAlloc(hypre_MPI_Request, tree.num_child); term_statuses = hypre_CTAlloc(hypre_MPI_Status, tree.num_child); for (i=0; i< tree.num_child; i++) { hypre_MPI_Irecv(NULL, 0, HYPRE_MPI_INT, tree.child_id[i], term_tag, comm, &term_requests[i]); } terminate = 0; children_complete = 0; } else if (num_procs ==1 && num_contacts > 0 ) /* added 11/08 */ { terminate = 0; } /*---------PROBE LOOP-----------------------------------------*/ /*Look for incoming contact messages - don't know how many I will get!*/ while (!terminate) { /* did I receive any contact messages? */ hypre_MPI_Iprobe(hypre_MPI_ANY_SOURCE, contact_tag, comm, &contact_flag, &status); while (contact_flag) { /* received contacts - from who and what do we do ?*/ proc = status.hypre_MPI_SOURCE; hypre_MPI_Get_count(&status, hypre_MPI_BYTE, &contact_size); contact_size = contact_size/contact_obj_size; /*---------------FILL RESPONSE ------------------------*/ /*first receive the contact buffer - then call a function to determine how to populate the send buffer for the reponse*/ /* do we have enough space to recv it? */ if(contact_size > recv_contact_buf_size) { recv_contact_buf = hypre_ReAlloc(recv_contact_buf, contact_obj_size*contact_size); recv_contact_buf_size = contact_size; } /* this must be blocking - can't fill recv without the buffer*/ hypre_MPI_Recv(recv_contact_buf, contact_size*contact_obj_size, hypre_MPI_BYTE, proc, contact_tag, comm, &fill_status); response_obj->fill_response(recv_contact_buf, contact_size, proc, response_obj, comm, &send_response_buf, &response_message_size ); /* we need to append the size of the send obj */ /* first we copy out any part that may be needed to send later so we don't overwrite */ post_size = response_message_size - max_response_size; if (post_size > 0) /*we will need to send the extra information later */ { /*hypre_printf("myid = %d, post_size = %d\n", myid, post_size);*/ if (post_array_size == post_array_storage) { /* allocate room for more posts - add 20*/ post_array_storage += 20; post_array = hypre_TReAlloc(post_array, void *, post_array_storage); post_send_requests = hypre_TReAlloc(post_send_requests, hypre_MPI_Request, post_array_storage); } /* allocate space for the data this post only*/ /* this should not happen often (unless a poor max_size has been chosen) - so we will allocate space for the data as needed */ size = post_size*response_obj_size; post_array[post_array_size] = hypre_MAlloc(size); /* index_ptr = send_response_buf + max_response_size_bytes */; index_ptr = (void *) ((char *) send_response_buf + max_response_size_bytes); memcpy(post_array[post_array_size], index_ptr, size); /*now post any part of the message that is too long with a non-blocking send and a different tag */ hypre_MPI_Isend(post_array[post_array_size], size, hypre_MPI_BYTE, proc, post_tag, /*hypre_MPI_COMM_WORLD, */ comm, &post_send_requests[post_array_size]); post_array_size++; } /*now append the size information into the overhead storage */ /* index_ptr = send_response_buf + max_response_size_bytes; */ index_ptr = (void *) ((char *) send_response_buf + max_response_size_bytes); memcpy(index_ptr, &response_message_size, sizeof(HYPRE_Int)); /*send the block of data that includes the overhead */ /* this is a blocking send - the recv has already been posted */ hypre_MPI_Send(send_response_buf, max_response_total_bytes, hypre_MPI_BYTE, proc, response_tag, comm); /*--------------------------------------------------------------*/ /* look for any more contact messages*/ hypre_MPI_Iprobe(hypre_MPI_ANY_SOURCE, contact_tag, comm, &contact_flag, &status); } /* no more contact messages waiting - either (1) check to see if we have received all of our response messages (2) participate in termination (check for messages from children) (3) participate in termination sweep (check for message from parent) */ if (!responses_complete) { hypre_MPI_Testall(num_contacts, response_requests, &responses_complete, response_statuses); if (responses_complete && num_procs == 1) terminate = 1; /*added 11/08 */ } else if(!children_complete) /* have all of our children received all of their response messages?*/ { hypre_MPI_Testall(tree.num_child, term_requests, &children_complete, term_statuses); /* if we have gotten term messages from all of our children, send a term message to our parent. Then post a receive to hear back from parent */ if (children_complete & (myid > 0)) /*root does not have a parent*/ { hypre_MPI_Isend(NULL, 0, HYPRE_MPI_INT, tree.parent_id, term_tag, comm, &request_parent); hypre_MPI_Irecv(NULL, 0, HYPRE_MPI_INT, tree.parent_id, term_tag, comm, &term_request1); } } else /*have we gotten a term message from our parent? */ { if (myid == 0) /* root doesn't have a parent */ { terminate = 1; } else { hypre_MPI_Test(&term_request1, &terminate, &term_status1); } if (terminate) /*tell children to terminate */ { if (myid > 0 ) hypre_MPI_Wait(&request_parent, &status_parent); for (i=0; i< tree.num_child; i++) { /*a blocking send - recv has been posted already*/ hypre_MPI_Send(NULL, 0, HYPRE_MPI_INT, tree.child_id[i], term_tag, comm); } } } }
HYPRE_Int HYPRE_SStructSplitSetup( HYPRE_SStructSolver solver, HYPRE_SStructMatrix A, HYPRE_SStructVector b, HYPRE_SStructVector x ) { hypre_SStructVector *y; HYPRE_Int nparts; HYPRE_Int *nvars; void ****smatvec_data; HYPRE_Int (***ssolver_solve)(); HYPRE_Int (***ssolver_destroy)(); void ***ssolver_data; HYPRE_Int ssolver = (solver -> ssolver); MPI_Comm comm; hypre_SStructGrid *grid; hypre_SStructPMatrix *pA; hypre_SStructPVector *px; hypre_SStructPVector *py; hypre_StructMatrix *sA; hypre_StructVector *sx; hypre_StructVector *sy; HYPRE_StructMatrix sAH; HYPRE_StructVector sxH; HYPRE_StructVector syH; HYPRE_Int (*ssolve)(); HYPRE_Int (*sdestroy)(); void *sdata; HYPRE_Int part, vi, vj; comm = hypre_SStructVectorComm(b); grid = hypre_SStructVectorGrid(b); HYPRE_SStructVectorCreate(comm, grid, &y); HYPRE_SStructVectorInitialize(y); HYPRE_SStructVectorAssemble(y); nparts = hypre_SStructMatrixNParts(A); nvars = hypre_TAlloc(HYPRE_Int, nparts); smatvec_data = hypre_TAlloc(void ***, nparts); ssolver_solve = (HYPRE_Int (***)()) hypre_MAlloc((sizeof(HYPRE_Int (**)()) * nparts)); ssolver_destroy = (HYPRE_Int (***)()) hypre_MAlloc((sizeof(HYPRE_Int (**)()) * nparts)); ssolver_data = hypre_TAlloc(void **, nparts); for (part = 0; part < nparts; part++) { pA = hypre_SStructMatrixPMatrix(A, part); px = hypre_SStructVectorPVector(x, part); py = hypre_SStructVectorPVector(y, part); nvars[part] = hypre_SStructPMatrixNVars(pA); smatvec_data[part] = hypre_TAlloc(void **, nvars[part]); ssolver_solve[part] = (HYPRE_Int (**)()) hypre_MAlloc((sizeof(HYPRE_Int (*)()) * nvars[part])); ssolver_destroy[part] = (HYPRE_Int (**)()) hypre_MAlloc((sizeof(HYPRE_Int (*)()) * nvars[part])); ssolver_data[part] = hypre_TAlloc(void *, nvars[part]); for (vi = 0; vi < nvars[part]; vi++) { smatvec_data[part][vi] = hypre_TAlloc(void *, nvars[part]); for (vj = 0; vj < nvars[part]; vj++) { sA = hypre_SStructPMatrixSMatrix(pA, vi, vj); sx = hypre_SStructPVectorSVector(px, vj); smatvec_data[part][vi][vj] = NULL; if (sA != NULL) { smatvec_data[part][vi][vj] = hypre_StructMatvecCreate(); hypre_StructMatvecSetup(smatvec_data[part][vi][vj], sA, sx); } } sA = hypre_SStructPMatrixSMatrix(pA, vi, vi); sx = hypre_SStructPVectorSVector(px, vi); sy = hypre_SStructPVectorSVector(py, vi); sAH = (HYPRE_StructMatrix) sA; sxH = (HYPRE_StructVector) sx; syH = (HYPRE_StructVector) sy; switch(ssolver) { default: /* If no solver is matched, use Jacobi, but throw and error */ if (ssolver != HYPRE_Jacobi) { hypre_error(HYPRE_ERROR_GENERIC); } /* don't break */ case HYPRE_Jacobi: HYPRE_StructJacobiCreate(comm, (HYPRE_StructSolver *)&sdata); HYPRE_StructJacobiSetMaxIter(sdata, 1); HYPRE_StructJacobiSetTol(sdata, 0.0); if (solver -> zero_guess) { HYPRE_StructJacobiSetZeroGuess(sdata); } HYPRE_StructJacobiSetup(sdata, sAH, syH, sxH); ssolve = HYPRE_StructJacobiSolve; sdestroy = HYPRE_StructJacobiDestroy; break; case HYPRE_SMG: HYPRE_StructSMGCreate(comm, (HYPRE_StructSolver *)&sdata); HYPRE_StructSMGSetMemoryUse(sdata, 0); HYPRE_StructSMGSetMaxIter(sdata, 1); HYPRE_StructSMGSetTol(sdata, 0.0); if (solver -> zero_guess) { HYPRE_StructSMGSetZeroGuess(sdata); } HYPRE_StructSMGSetNumPreRelax(sdata, 1); HYPRE_StructSMGSetNumPostRelax(sdata, 1); HYPRE_StructSMGSetLogging(sdata, 0); HYPRE_StructSMGSetPrintLevel(sdata, 0); HYPRE_StructSMGSetup(sdata, sAH, syH, sxH); ssolve = HYPRE_StructSMGSolve; sdestroy = HYPRE_StructSMGDestroy; break; case HYPRE_PFMG: HYPRE_StructPFMGCreate(comm, (HYPRE_StructSolver *)&sdata); HYPRE_StructPFMGSetMaxIter(sdata, 1); HYPRE_StructPFMGSetTol(sdata, 0.0); if (solver -> zero_guess) { HYPRE_StructPFMGSetZeroGuess(sdata); } HYPRE_StructPFMGSetRelaxType(sdata, 1); HYPRE_StructPFMGSetNumPreRelax(sdata, 1); HYPRE_StructPFMGSetNumPostRelax(sdata, 1); HYPRE_StructPFMGSetLogging(sdata, 0); HYPRE_StructPFMGSetPrintLevel(sdata, 0); HYPRE_StructPFMGSetup(sdata, sAH, syH, sxH); ssolve = HYPRE_StructPFMGSolve; sdestroy = HYPRE_StructPFMGDestroy; break; } ssolver_solve[part][vi] = ssolve; ssolver_destroy[part][vi] = sdestroy; ssolver_data[part][vi] = sdata; } } (solver -> y) = y; (solver -> nparts) = nparts; (solver -> nvars) = nvars; (solver -> smatvec_data) = smatvec_data; (solver -> ssolver_solve) = ssolver_solve; (solver -> ssolver_destroy) = ssolver_destroy; (solver -> ssolver_data) = ssolver_data; if ((solver -> tol) > 0.0) { hypre_SStructMatvecCreate(&(solver -> matvec_data)); hypre_SStructMatvecSetup((solver -> matvec_data), A, x); } return hypre_error_flag; }