int HYPRE_StructVectorSetBoxValues( HYPRE_StructVector vector, int *ilower, int *iupper, double *values ) { hypre_Index new_ilower; hypre_Index new_iupper; hypre_Box *new_value_box; int d; int ierr = 0; hypre_ClearIndex(new_ilower); hypre_ClearIndex(new_iupper); for (d = 0; d < hypre_StructGridDim(hypre_StructVectorGrid(vector)); d++) { hypre_IndexD(new_ilower, d) = ilower[d]; hypre_IndexD(new_iupper, d) = iupper[d]; } new_value_box = hypre_BoxCreate(); hypre_BoxSetExtents(new_value_box, new_ilower, new_iupper); ierr = hypre_StructVectorSetBoxValues(vector, new_value_box, values, 0 ); hypre_BoxDestroy(new_value_box); return ierr; }
int HYPRE_StructMatrixSetBoxValues( HYPRE_StructMatrix matrix, int *ilower, int *iupper, int num_stencil_indices, int *stencil_indices, double *values ) { hypre_Index new_ilower; hypre_Index new_iupper; hypre_Box *new_value_box; int d; int ierr = 0; hypre_ClearIndex(new_ilower); hypre_ClearIndex(new_iupper); for (d = 0; d < hypre_StructGridDim(hypre_StructMatrixGrid(matrix)); d++) { hypre_IndexD(new_ilower, d) = ilower[d]; hypre_IndexD(new_iupper, d) = iupper[d]; } new_value_box = hypre_BoxCreate(); hypre_BoxSetExtents(new_value_box, new_ilower, new_iupper); ierr = hypre_StructMatrixSetBoxValues(matrix, new_value_box, num_stencil_indices, stencil_indices, values, 0); hypre_BoxDestroy(new_value_box); return (ierr); }
hypre_SStructSendInfoData * hypre_SStructSendInfo( hypre_StructGrid *fgrid, hypre_BoxManager *cboxman, hypre_Index rfactor ) { hypre_SStructSendInfoData *sendinfo_data; MPI_Comm comm= hypre_SStructVectorComm(fgrid); hypre_BoxArray *grid_boxes; hypre_Box *grid_box, cbox; hypre_Box *intersect_box, boxman_entry_box; hypre_BoxManEntry **boxman_entries; HYPRE_Int nboxman_entries; hypre_BoxArrayArray *send_boxes; HYPRE_Int **send_processes; HYPRE_Int **send_remote_boxnums; hypre_Index ilower, iupper, index; HYPRE_Int myproc, proc; HYPRE_Int cnt; HYPRE_Int i, j; hypre_ClearIndex(index); hypre_MPI_Comm_rank(comm, &myproc); sendinfo_data= hypre_CTAlloc(hypre_SStructSendInfoData, 1); /*------------------------------------------------------------------------ * Create the structured sendbox patterns. * * send_boxes are obtained by intersecting this proc's fgrid boxes * with cgrid's box_man. Intersecting BoxManEntries not on this proc * will give boxes that we will need to send data to- i.e., we scan * through the boxes of grid and find the processors that own a chunk * of it. *------------------------------------------------------------------------*/ intersect_box = hypre_CTAlloc(hypre_Box, 1); grid_boxes = hypre_StructGridBoxes(fgrid); send_boxes= hypre_BoxArrayArrayCreate(hypre_BoxArraySize(grid_boxes)); send_processes= hypre_CTAlloc(HYPRE_Int *, hypre_BoxArraySize(grid_boxes)); send_remote_boxnums= hypre_CTAlloc(HYPRE_Int *, hypre_BoxArraySize(grid_boxes)); hypre_ForBoxI(i, grid_boxes) { grid_box= hypre_BoxArrayBox(grid_boxes, i); /*--------------------------------------------------------------------- * Find the boxarray that must be sent. BoxManIntersect returns * the full extents of the boxes that intersect with the given box. * We further need to intersect each box in the list with the given * box to determine the actual box that needs to be sent. *---------------------------------------------------------------------*/ hypre_SStructIndexScaleF_C(hypre_BoxIMin(grid_box), index, rfactor, hypre_BoxIMin(&cbox)); hypre_SStructIndexScaleF_C(hypre_BoxIMax(grid_box), index, rfactor, hypre_BoxIMax(&cbox)); hypre_BoxManIntersect(cboxman, hypre_BoxIMin(&cbox), hypre_BoxIMax(&cbox), &boxman_entries, &nboxman_entries); cnt= 0; for (j= 0; j< nboxman_entries; j++) { hypre_SStructBoxManEntryGetProcess(boxman_entries[j], &proc); if (proc != myproc) { cnt++; } } send_processes[i] = hypre_CTAlloc(HYPRE_Int, cnt); send_remote_boxnums[i]= hypre_CTAlloc(HYPRE_Int, cnt); cnt= 0; for (j= 0; j< nboxman_entries; j++) { hypre_SStructBoxManEntryGetProcess(boxman_entries[j], &proc); /* determine the chunk of the boxman_entries[j] box that is needed */ hypre_BoxManEntryGetExtents(boxman_entries[j], ilower, iupper); hypre_BoxSetExtents(&boxman_entry_box, ilower, iupper); hypre_IntersectBoxes(&boxman_entry_box, &cbox, &boxman_entry_box); if (proc != myproc) { send_processes[i][cnt] = proc; hypre_SStructBoxManEntryGetBoxnum(boxman_entries[j], &send_remote_boxnums[i][cnt]); hypre_AppendBox(&boxman_entry_box, hypre_BoxArrayArrayBoxArray(send_boxes, i)); cnt++; } } hypre_TFree(boxman_entries); } /* hypre_ForBoxI(i, grid_boxes) */
/*-------------------------------------------------------------------------- * hypre_FacZeroCFSten: Zeroes the coarse stencil coefficients that reach * into an underlying coarsened refinement box. * Algo: For each cbox * { * 1) refine cbox and expand by one in each direction * 2) boxman_intersect with the fboxman * 3) loop over intersection boxes to see if stencil * reaches over. * } *--------------------------------------------------------------------------*/ HYPRE_Int hypre_FacZeroCFSten( hypre_SStructPMatrix *Af, hypre_SStructPMatrix *Ac, hypre_SStructGrid *grid, HYPRE_Int fine_part, hypre_Index rfactors ) { hypre_BoxManager *fboxman; hypre_BoxManEntry **boxman_entries; HYPRE_Int nboxman_entries; hypre_SStructPGrid *p_cgrid; hypre_Box fgrid_box; hypre_StructGrid *cgrid; hypre_BoxArray *cgrid_boxes; hypre_Box *cgrid_box; hypre_Box scaled_box; hypre_Box *shift_ibox; hypre_StructMatrix *smatrix; hypre_StructStencil *stencils; HYPRE_Int stencil_size; hypre_Index refine_factors, upper_shift; hypre_Index stride; hypre_Index stencil_shape; hypre_Index zero_index, ilower, iupper; HYPRE_Int nvars, var1, var2; HYPRE_Int ndim; hypre_Box *ac_dbox; HYPRE_Real *ac_ptr; hypre_Index loop_size; HYPRE_Int iac; HYPRE_Int ci, i, j; HYPRE_Int abs_shape; HYPRE_Int ierr = 0; p_cgrid = hypre_SStructPMatrixPGrid(Ac); nvars = hypre_SStructPMatrixNVars(Ac); ndim = hypre_SStructPGridNDim(p_cgrid); hypre_BoxInit(&fgrid_box, ndim); hypre_BoxInit(&scaled_box, ndim); hypre_ClearIndex(zero_index); hypre_ClearIndex(stride); hypre_ClearIndex(upper_shift); for (i= 0; i< ndim; i++) { stride[i]= 1; upper_shift[i]= rfactors[i]-1; } hypre_CopyIndex(rfactors, refine_factors); if (ndim < 3) { for (i= ndim; i< 3; i++) { refine_factors[i]= 1; } } for (var1= 0; var1< nvars; var1++) { cgrid= hypre_SStructPGridSGrid(hypre_SStructPMatrixPGrid(Ac), var1); cgrid_boxes= hypre_StructGridBoxes(cgrid); fboxman= hypre_SStructGridBoxManager(grid, fine_part, var1); /*------------------------------------------------------------------ * For each parent coarse box find all fboxes that may be connected * through a stencil entry- refine this box, expand it by one * in each direction, and boxman_intersect with fboxman *------------------------------------------------------------------*/ hypre_ForBoxI(ci, cgrid_boxes) { cgrid_box= hypre_BoxArrayBox(cgrid_boxes, ci); hypre_StructMapCoarseToFine(hypre_BoxIMin(cgrid_box), zero_index, refine_factors, hypre_BoxIMin(&scaled_box)); hypre_StructMapCoarseToFine(hypre_BoxIMax(cgrid_box), upper_shift, refine_factors, hypre_BoxIMax(&scaled_box)); hypre_SubtractIndexes(hypre_BoxIMin(&scaled_box), stride, 3, hypre_BoxIMin(&scaled_box)); hypre_AddIndexes(hypre_BoxIMax(&scaled_box), stride, 3, hypre_BoxIMax(&scaled_box)); hypre_BoxManIntersect(fboxman, hypre_BoxIMin(&scaled_box), hypre_BoxIMax(&scaled_box), &boxman_entries, &nboxman_entries); for (var2= 0; var2< nvars; var2++) { stencils= hypre_SStructPMatrixSStencil(Ac, var1, var2); if (stencils != NULL) { stencil_size= hypre_StructStencilSize(stencils); smatrix = hypre_SStructPMatrixSMatrix(Ac, var1, var2); ac_dbox = hypre_BoxArrayBox(hypre_StructMatrixDataSpace(smatrix), ci); /*--------------------------------------------------------- * Find the stencil coefficients that must be zeroed off. * Loop over all possible boxes. *---------------------------------------------------------*/ for (i= 0; i< stencil_size; i++) { hypre_CopyIndex(hypre_StructStencilElement(stencils, i), stencil_shape); AbsStencilShape(stencil_shape, abs_shape); if (abs_shape) /* non-centre stencils are zeroed */ { /* look for connecting fboxes that must be zeroed. */ for (j= 0; j< nboxman_entries; j++) { hypre_BoxManEntryGetExtents(boxman_entries[j], ilower, iupper); hypre_BoxSetExtents(&fgrid_box, ilower, iupper); shift_ibox= hypre_CF_StenBox(&fgrid_box, cgrid_box, stencil_shape, refine_factors, ndim); if ( hypre_BoxVolume(shift_ibox) ) { ac_ptr= hypre_StructMatrixExtractPointerByIndex(smatrix, ci, stencil_shape); hypre_BoxGetSize(shift_ibox, loop_size); hypre_BoxLoop1Begin(ndim, loop_size, ac_dbox, hypre_BoxIMin(shift_ibox), stride, iac); #ifdef HYPRE_USING_OPENMP #pragma omp parallel for private(HYPRE_BOX_PRIVATE,iac) HYPRE_SMP_SCHEDULE #endif hypre_BoxLoop1For(iac) { ac_ptr[iac] = 0.0; } hypre_BoxLoop1End(iac); } /* if ( hypre_BoxVolume(shift_ibox) ) */ hypre_BoxDestroy(shift_ibox); } /* for (j= 0; j< nboxman_entries; j++) */ } /* if (abs_shape) */ } /* for (i= 0; i< stencil_size; i++) */ } /* if (stencils != NULL) */ } /* for (var2= 0; var2< nvars; var2++) */ hypre_TFree(boxman_entries); } /* hypre_ForBoxI ci */
HYPRE_Int hypre_CreateCommInfoFromStencil( hypre_StructGrid *grid, hypre_StructStencil *stencil, hypre_CommInfo **comm_info_ptr ) { HYPRE_Int i,j,k, d, m, s; hypre_BoxArrayArray *send_boxes; hypre_BoxArrayArray *recv_boxes; HYPRE_Int **send_procs; HYPRE_Int **recv_procs; HYPRE_Int **send_rboxnums; HYPRE_Int **recv_rboxnums; hypre_BoxArrayArray *send_rboxes; hypre_BoxArrayArray *recv_rboxes; hypre_BoxArray *local_boxes; HYPRE_Int num_boxes; HYPRE_Int *local_ids; hypre_BoxManager *boxman; hypre_Index *stencil_shape; hypre_IndexRef stencil_offset; hypre_IndexRef pshift; hypre_Box *box; hypre_Box *hood_box; hypre_Box *grow_box; hypre_Box *extend_box; hypre_Box *int_box; hypre_Box *periodic_box; HYPRE_Int stencil_grid[3][3][3]; HYPRE_Int grow[3][2]; hypre_BoxManEntry **entries; hypre_BoxManEntry *entry; HYPRE_Int num_entries; hypre_BoxArray *neighbor_boxes = NULL; HYPRE_Int *neighbor_procs = NULL; HYPRE_Int *neighbor_ids = NULL; HYPRE_Int *neighbor_shifts = NULL; HYPRE_Int neighbor_count; HYPRE_Int neighbor_alloc; hypre_Index ilower, iupper; hypre_BoxArray *send_box_array; hypre_BoxArray *recv_box_array; hypre_BoxArray *send_rbox_array; hypre_BoxArray *recv_rbox_array; hypre_Box **cboxes; hypre_Box *cboxes_mem; HYPRE_Int *cboxes_neighbor_location; HYPRE_Int num_cboxes, cbox_alloc; HYPRE_Int istart[3], istop[3]; HYPRE_Int sgindex[3]; HYPRE_Int num_periods, loc, box_id, id, proc_id; HYPRE_Int myid; MPI_Comm comm; /*------------------------------------------------------ * Initializations *------------------------------------------------------*/ local_boxes = hypre_StructGridBoxes(grid); local_ids = hypre_StructGridIDs(grid); num_boxes = hypre_BoxArraySize(local_boxes); num_periods = hypre_StructGridNumPeriods(grid); boxman = hypre_StructGridBoxMan(grid); comm = hypre_StructGridComm(grid); hypre_MPI_Comm_rank(comm, &myid); for (k = 0; k < 3; k++) { for (j = 0; j < 3; j++) { for (i = 0; i < 3; i++) { stencil_grid[i][j][k] = 0; } } } /*------------------------------------------------------ * Compute the "grow" information from the stencil *------------------------------------------------------*/ stencil_shape = hypre_StructStencilShape(stencil); for (d = 0; d < 3; d++) { grow[d][0] = 0; grow[d][1] = 0; } for (s = 0; s < hypre_StructStencilSize(stencil); s++) { stencil_offset = stencil_shape[s]; for (d = 0; d < 3; d++) { m = stencil_offset[d]; istart[d] = 1; istop[d] = 1; if (m < 0) { istart[d] = 0; grow[d][0] = hypre_max(grow[d][0], -m); } else if (m > 0) { istop[d] = 2; grow[d][1] = hypre_max(grow[d][1], m); } } /* update stencil grid from the grow_stencil */ for (k = istart[2]; k <= istop[2]; k++) { for (j = istart[1]; j <= istop[1]; j++) { for (i = istart[0]; i <= istop[0]; i++) { stencil_grid[i][j][k] = 1; } } } } /*------------------------------------------------------ * Compute send/recv boxes and procs for each local box *------------------------------------------------------*/ /* initialize: for each local box, we create an array of send/recv info */ send_boxes = hypre_BoxArrayArrayCreate(num_boxes); recv_boxes = hypre_BoxArrayArrayCreate(num_boxes); send_procs = hypre_CTAlloc(HYPRE_Int *, num_boxes); recv_procs = hypre_CTAlloc(HYPRE_Int *, num_boxes); /* Remote boxnums and boxes describe data on the opposing processor, so some shifting of boxes is needed below for periodic neighbor boxes. Remote box info is also needed for receives to allow for reverse communication. */ send_rboxnums = hypre_CTAlloc(HYPRE_Int *, num_boxes); send_rboxes = hypre_BoxArrayArrayCreate(num_boxes); recv_rboxnums = hypre_CTAlloc(HYPRE_Int *, num_boxes); recv_rboxes = hypre_BoxArrayArrayCreate(num_boxes); grow_box = hypre_BoxCreate(); extend_box = hypre_BoxCreate(); int_box = hypre_BoxCreate(); periodic_box = hypre_BoxCreate(); /* storage we will use and keep track of the neighbors */ neighbor_alloc = 30; /* initial guess at max size */ neighbor_boxes = hypre_BoxArrayCreate(neighbor_alloc); neighbor_procs = hypre_CTAlloc(HYPRE_Int, neighbor_alloc); neighbor_ids = hypre_CTAlloc(HYPRE_Int, neighbor_alloc); neighbor_shifts = hypre_CTAlloc(HYPRE_Int, neighbor_alloc); /* storage we will use to collect all of the intersected boxes (the send and recv regions for box i (this may not be enough in the case of periodic boxes, so we will have to check) */ cbox_alloc = hypre_BoxManNEntries(boxman); cboxes_neighbor_location = hypre_CTAlloc(HYPRE_Int, cbox_alloc); cboxes = hypre_CTAlloc(hypre_Box *, cbox_alloc); cboxes_mem = hypre_CTAlloc(hypre_Box, cbox_alloc); /******* loop through each local box **************/ for (i = 0; i < num_boxes; i++) { /* get the box */ box = hypre_BoxArrayBox(local_boxes, i); /* box_id = local_ids[i]; the box id in the Box Manager is the box number, * and we use this to find out if a box has intersected with itself */ box_id = i; /* grow box local i according to the stencil*/ hypre_CopyBox(box, grow_box); for (d = 0; d < 3; d++) { hypre_BoxIMinD(grow_box, d) -= grow[d][0]; hypre_BoxIMaxD(grow_box, d) += grow[d][1]; } /* extend_box - to find the list of potential neighbors, we need to grow the local box a bit differently in case, for example, the stencil grows in one dimension [0] and not the other [1] */ hypre_CopyBox(box, extend_box); for (d = 0; d < 3; d++) { hypre_BoxIMinD(extend_box, d) -= hypre_max(grow[d][0],grow[d][1]); hypre_BoxIMaxD(extend_box, d) += hypre_max(grow[d][0],grow[d][1]); } /*------------------------------------------------ * Determine the neighbors of box i *------------------------------------------------*/ /* Do this by intersecting the extend box with the BoxManager. We must also check for periodic neighbors. */ neighbor_count = 0; hypre_BoxArraySetSize(neighbor_boxes, 0); /* shift the box by each period (k=0 is original box) */ for (k = 0; k < num_periods; k++) { hypre_CopyBox(extend_box, periodic_box); pshift = hypre_StructGridPShift(grid, k); hypre_BoxShiftPos(periodic_box, pshift); /* get the intersections */ hypre_BoxManIntersect(boxman, hypre_BoxIMin(periodic_box) , hypre_BoxIMax(periodic_box) , &entries , &num_entries); /* note: do we need to remove the intersection with our original box? no if periodic, yes if non-periodic (k=0) */ /* unpack entries (first check storage) */ if (neighbor_count + num_entries > neighbor_alloc) { neighbor_alloc = neighbor_count + num_entries + 5; neighbor_procs = hypre_TReAlloc(neighbor_procs, HYPRE_Int, neighbor_alloc); neighbor_ids = hypre_TReAlloc(neighbor_ids, HYPRE_Int, neighbor_alloc); neighbor_shifts = hypre_TReAlloc(neighbor_shifts, HYPRE_Int, neighbor_alloc); } /* check storage for the array */ hypre_BoxArraySetSize(neighbor_boxes, neighbor_count + num_entries); /* now unpack */ for (j = 0; j < num_entries; j++) { entry = entries[j]; proc_id = hypre_BoxManEntryProc(entry); id = hypre_BoxManEntryId(entry); /* don't keep box i in the non-periodic case*/ if (!k) { if((myid == proc_id) && (box_id == id)) { continue; } } hypre_BoxManEntryGetExtents(entry, ilower, iupper); hypre_BoxSetExtents(hypre_BoxArrayBox(neighbor_boxes, neighbor_count), ilower, iupper); /* shift the periodic boxes (needs to be the opposite of above) */ if (k) { hypre_BoxShiftNeg( hypre_BoxArrayBox(neighbor_boxes, neighbor_count), pshift); } neighbor_procs[neighbor_count] = proc_id; neighbor_ids[neighbor_count] = id; neighbor_shifts[neighbor_count] = k; neighbor_count++; } hypre_BoxArraySetSize(neighbor_boxes, neighbor_count); hypre_TFree(entries); } /* end of loop through periods k */ /* Now we have a list of all of the neighbors for box i! */ /* note: we don't want/need to remove duplicates - they should have different intersections (TO DO: put more thought into if there are ever any exceptions to this? - the intersection routine already eliminates duplicates - so what i mean is eliminating duplicates from multiple intersection calls in periodic case) */ /*------------------------------------------------ * Compute recv_box_array for box i *------------------------------------------------*/ /* check size of storage for cboxes */ /* let's make sure that we have enough storage in case each neighbor produces a send/recv region */ if (neighbor_count > cbox_alloc) { cbox_alloc = neighbor_count; cboxes_neighbor_location = hypre_TReAlloc(cboxes_neighbor_location, HYPRE_Int, cbox_alloc); cboxes = hypre_TReAlloc(cboxes, hypre_Box *, cbox_alloc); cboxes_mem = hypre_TReAlloc(cboxes_mem, hypre_Box, cbox_alloc); } /* Loop through each neighbor box. If the neighbor box intersects the grown box i (grown according to our stencil), then the intersection is a recv region. If the neighbor box was shifted to handle periodicity, we need to (positive) shift it back. */ num_cboxes = 0; for (k = 0; k < neighbor_count; k++) { hood_box = hypre_BoxArrayBox(neighbor_boxes, k); /* check the stencil grid to see if it makes sense to intersect */ for (d = 0; d < 3; d++) { sgindex[d] = 1; s = hypre_BoxIMinD(hood_box, d) - hypre_BoxIMaxD(box, d); if (s > 0) { sgindex[d] = 2; } s = hypre_BoxIMinD(box, d) - hypre_BoxIMaxD(hood_box, d); if (s > 0) { sgindex[d] = 0; } } /* it makes sense only if we have at least one non-zero entry */ if (stencil_grid[sgindex[0]][sgindex[1]][sgindex[2]]) { /* intersect - result is int_box */ hypre_IntersectBoxes(grow_box, hood_box, int_box); /* if we have a positive volume box, this is a recv region */ if (hypre_BoxVolume(int_box)) { /* keep track of which neighbor: k... */ cboxes_neighbor_location[num_cboxes] = k; cboxes[num_cboxes] = &cboxes_mem[num_cboxes]; /* keep the intersected box */ hypre_CopyBox(int_box, cboxes[num_cboxes]); num_cboxes++; } } } /* end of loop through each neighbor */ /* create recv_box_array and recv_procs for box i */ recv_box_array = hypre_BoxArrayArrayBoxArray(recv_boxes, i); hypre_BoxArraySetSize(recv_box_array, num_cboxes); recv_procs[i] = hypre_CTAlloc(HYPRE_Int, num_cboxes); recv_rboxnums[i] = hypre_CTAlloc(HYPRE_Int, num_cboxes); recv_rbox_array = hypre_BoxArrayArrayBoxArray(recv_rboxes, i); hypre_BoxArraySetSize(recv_rbox_array, num_cboxes); for (m = 0; m < num_cboxes; m++) { loc = cboxes_neighbor_location[m]; recv_procs[i][m] = neighbor_procs[loc]; recv_rboxnums[i][m] = neighbor_ids[loc]; hypre_CopyBox(cboxes[m], hypre_BoxArrayBox(recv_box_array, m)); /* if periodic, positive shift before copying to the rbox_array */ if (neighbor_shifts[loc]) /* periodic if shift != 0 */ { pshift = hypre_StructGridPShift(grid, neighbor_shifts[loc]); hypre_BoxShiftPos(cboxes[m], pshift); } hypre_CopyBox(cboxes[m], hypre_BoxArrayBox(recv_rbox_array, m)); cboxes[m] = NULL; } /*------------------------------------------------ * Compute send_box_array for box i *------------------------------------------------*/ /* Loop through each neighbor box. If the grown neighbor box intersects box i, then the intersection is a send region. If the neighbor box was shifted to handle periodicity, we need to (positive) shift it back. */ num_cboxes = 0; for (k = 0; k < neighbor_count; k++) { hood_box = hypre_BoxArrayBox(neighbor_boxes, k); /* check the stencil grid to see if it makes sense to intersect */ for (d = 0; d < 3; d++) { sgindex[d] = 1; s = hypre_BoxIMinD(box, d) - hypre_BoxIMaxD(hood_box, d); if (s > 0) { sgindex[d] = 2; } s = hypre_BoxIMinD(hood_box, d) - hypre_BoxIMaxD(box, d); if (s > 0) { sgindex[d] = 0; } } /* it makes sense only if we have at least one non-zero entry */ if (stencil_grid[sgindex[0]][sgindex[1]][sgindex[2]]) { /* grow the neighbor box and intersect */ hypre_CopyBox(hood_box, grow_box); for (d = 0; d < 3; d++) { hypre_BoxIMinD(grow_box, d) -= grow[d][0]; hypre_BoxIMaxD(grow_box, d) += grow[d][1]; } hypre_IntersectBoxes(box, grow_box, int_box); /* if we have a positive volume box, this is a send region */ if (hypre_BoxVolume(int_box)) { /* keep track of which neighbor: k... */ cboxes_neighbor_location[num_cboxes] = k; cboxes[num_cboxes] = &cboxes_mem[num_cboxes]; /* keep the intersected box */ hypre_CopyBox(int_box, cboxes[num_cboxes]); num_cboxes++; } } }/* end of loop through neighbors */ /* create send_box_array and send_procs for box i */ send_box_array = hypre_BoxArrayArrayBoxArray(send_boxes, i); hypre_BoxArraySetSize(send_box_array, num_cboxes); send_procs[i] = hypre_CTAlloc(HYPRE_Int, num_cboxes); send_rboxnums[i] = hypre_CTAlloc(HYPRE_Int, num_cboxes); send_rbox_array = hypre_BoxArrayArrayBoxArray(send_rboxes, i); hypre_BoxArraySetSize(send_rbox_array, num_cboxes); for (m = 0; m < num_cboxes; m++) { loc = cboxes_neighbor_location[m]; send_procs[i][m] = neighbor_procs[loc]; send_rboxnums[i][m] = neighbor_ids[loc]; hypre_CopyBox(cboxes[m], hypre_BoxArrayBox(send_box_array, m)); /* if periodic, positive shift before copying to the rbox_array */ if (neighbor_shifts[loc]) /* periodic if shift != 0 */ { pshift = hypre_StructGridPShift(grid, neighbor_shifts[loc]); hypre_BoxShiftPos(cboxes[m], pshift); } hypre_CopyBox(cboxes[m], hypre_BoxArrayBox(send_rbox_array, m)); cboxes[m] = NULL; } } /* end of loop through each local box */
HYPRE_Int hypre_FacZeroCData( void *fac_vdata, hypre_SStructMatrix *A ) { hypre_FACData *fac_data = fac_vdata; hypre_SStructGrid *grid; hypre_SStructPGrid *p_cgrid; hypre_StructGrid *cgrid; hypre_BoxArray *cgrid_boxes; hypre_Box *cgrid_box; hypre_BoxManager *fboxman; hypre_BoxManEntry **boxman_entries; HYPRE_Int nboxman_entries; hypre_Box scaled_box; hypre_Box intersect_box; hypre_SStructPMatrix *level_pmatrix; hypre_StructStencil *stencils; HYPRE_Int stencil_size; hypre_Index *refine_factors; hypre_Index temp_index; hypre_Index ilower, iupper; HYPRE_Int max_level = fac_data -> max_levels; HYPRE_Int *level_to_part = fac_data -> level_to_part; HYPRE_Int ndim = hypre_SStructMatrixNDim(A); HYPRE_Int part_crse = 0; HYPRE_Int part_fine = 1; HYPRE_Int level; HYPRE_Int nvars, var; HYPRE_Int ci, i, j, rem, intersect_size; double *values; HYPRE_Int ierr = 0; for (level= max_level; level> 0; level--) { level_pmatrix = hypre_SStructMatrixPMatrix(fac_data -> A_level[level], part_crse); grid = (fac_data -> grid_level[level]); refine_factors= &(fac_data -> refine_factors[level]); p_cgrid= hypre_SStructGridPGrid(grid, part_crse); nvars = hypre_SStructPGridNVars(p_cgrid); for (var= 0; var< nvars; var++) { stencils = hypre_SStructPMatrixSStencil(level_pmatrix, var, var); stencil_size= hypre_StructStencilSize(stencils); /*--------------------------------------------------------------------- * For each variable, find the underlying boxes for each coarse box. *---------------------------------------------------------------------*/ cgrid = hypre_SStructPGridSGrid(p_cgrid, var); cgrid_boxes = hypre_StructGridBoxes(cgrid); fboxman = hypre_SStructGridBoxManager(grid, part_fine, var); hypre_ForBoxI(ci, cgrid_boxes) { cgrid_box= hypre_BoxArrayBox(cgrid_boxes, ci); hypre_ClearIndex(temp_index); hypre_StructMapCoarseToFine(hypre_BoxIMin(cgrid_box), temp_index, *refine_factors, hypre_BoxIMin(&scaled_box)); for (i= 0; i< ndim; i++) { temp_index[i]= (*refine_factors)[i]-1; } hypre_StructMapCoarseToFine(hypre_BoxIMax(cgrid_box), temp_index, *refine_factors, hypre_BoxIMax(&scaled_box)); hypre_BoxManIntersect(fboxman, hypre_BoxIMin(&scaled_box), hypre_BoxIMax(&scaled_box), &boxman_entries, &nboxman_entries); for (i= 0; i< nboxman_entries; i++) { hypre_BoxManEntryGetExtents(boxman_entries[i], ilower, iupper); hypre_BoxSetExtents(&intersect_box, ilower, iupper); hypre_IntersectBoxes(&intersect_box, &scaled_box, &intersect_box); /* adjust the box so that it is divisible by refine_factors */ for (j= 0; j< ndim; j++) { rem= hypre_BoxIMin(&intersect_box)[j]%(*refine_factors)[j]; if (rem) { hypre_BoxIMin(&intersect_box)[j]+=(*refine_factors)[j] - rem; } } hypre_ClearIndex(temp_index); hypre_StructMapFineToCoarse(hypre_BoxIMin(&intersect_box), temp_index, *refine_factors, hypre_BoxIMin(&intersect_box)); hypre_StructMapFineToCoarse(hypre_BoxIMax(&intersect_box), temp_index, *refine_factors, hypre_BoxIMax(&intersect_box)); intersect_size= hypre_BoxVolume(&intersect_box); if (intersect_size > 0) { /*------------------------------------------------------------ * Coarse underlying box found. Now zero off. *------------------------------------------------------------*/ values= hypre_CTAlloc(double, intersect_size); for (j= 0; j< stencil_size; j++) { HYPRE_SStructMatrixSetBoxValues(fac_data -> A_level[level], part_crse, hypre_BoxIMin(&intersect_box), hypre_BoxIMax(&intersect_box), var, 1, &j, values); HYPRE_SStructMatrixSetBoxValues(A, level_to_part[level-1], hypre_BoxIMin(&intersect_box), hypre_BoxIMax(&intersect_box), var, 1, &j, values); } hypre_TFree(values); } /* if (intersect_size > 0) */ } /* for (i= 0; i< nboxman_entries; i++) */ hypre_TFree(boxman_entries); } /* hypre_ForBoxI(ci, cgrid_boxes) */ } /* for (var= 0; var< nvars; var++) */
int hypre_StructGridAssembleWithAP( hypre_StructGrid *grid ) { int ierr = 0; int tmp_i; int size, global_num_boxes, num_local_boxes; int i, j, d, k, index; int num_procs, myid; int *sendbuf8, *recvbuf8, *sendbuf2, *recvbuf2; int min_box_size, max_box_size; int global_min_box_size, global_max_box_size; int *ids; int max_regions, max_refinements, ologp; double gamma; hypre_Index min_index, max_index; int prune; hypre_Box *box; MPI_Comm comm = hypre_StructGridComm(grid); hypre_Box *bounding_box = hypre_StructGridBoundingBox(grid); hypre_BoxArray *local_boxes = hypre_StructGridBoxes(grid); int dim = hypre_StructGridDim(grid); hypre_BoxNeighbors *neighbors = hypre_StructGridNeighbors(grid); int max_distance = hypre_StructGridMaxDistance(grid); hypre_IndexRef periodic = hypre_StructGridPeriodic(grid); int *local_boxnums; double dbl_global_size, tmp_dbl; hypre_BoxArray *my_partition; int *part_ids, *part_boxnums; int *proc_array, proc_count, proc_alloc, count; int *tmp_proc_ids = NULL; int max_response_size; int *ap_proc_ids, *send_buf, *send_buf_starts; int *response_buf, *response_buf_starts; hypre_BoxArray *neighbor_boxes, *n_boxes_copy; int *neighbor_proc_ids, *neighbor_boxnums; int *order_index, *delete_array; int tmp_id, start, first_local; int grow, grow_array[6]; hypre_Box *grow_box; int *numghost; int ghostsize; hypre_Box *ghostbox; hypre_StructAssumedPart *assumed_part; hypre_DataExchangeResponse response_obj; int px = hypre_IndexX(periodic); int py = hypre_IndexY(periodic); int pz = hypre_IndexZ(periodic); int i_periodic = px ? 1 : 0; int j_periodic = py ? 1 : 0; int k_periodic = pz ? 1 : 0; int num_periods, multiple_ap, p; hypre_Box *result_box, *period_box; hypre_Index *pshifts; hypre_IndexRef pshift; #if NEIGH_PRINT double start_time, end_time; #endif /*--------------------------------------------- Step 1: Initializations -----------------------------------------------*/ prune = 1; /* default is to prune */ num_local_boxes = hypre_BoxArraySize(local_boxes); num_periods = (1+2*i_periodic) * (1+2*j_periodic) * (1+2*k_periodic); MPI_Comm_size(comm, &num_procs); MPI_Comm_rank(comm, &myid); /*--------------------------------------------- Step 2: Determine the global size, total number of boxes, and global bounding box. Also get the min and max box sizes since it is convenient to do so. -----------------------------------------------*/ if (neighbors == NULL) { /*these may not be needed - check later */ ids = hypre_TAlloc(int, num_local_boxes); /* for the vol and number of boxes */ sendbuf2 = hypre_CTAlloc(int, 2); recvbuf2 = hypre_CTAlloc(int, 2); size = 0; bounding_box = hypre_BoxCreate(); grow_box = hypre_BoxCreate(); if (num_local_boxes) { min_box_size = hypre_BoxVolume( hypre_BoxArrayBox(local_boxes, 0)); max_box_size = hypre_BoxVolume( hypre_BoxArrayBox(local_boxes, 0)); /* initialize min and max */ for (d=0; d<3; d++) { hypre_IndexD(min_index, d) = pow(2,30); hypre_IndexD(max_index, d) = -pow(2,30); } hypre_ForBoxI(i, local_boxes) { box = hypre_BoxArrayBox(local_boxes, i); /* get global size and number of boxes */ tmp_i = hypre_BoxVolume(box); size += tmp_i; min_box_size = hypre_min(min_box_size, tmp_i); max_box_size = hypre_max(max_box_size, tmp_i); /* set id */ ids[i] = i; /* 1/3/05 we need this for the case of holes in the domain. (I had commented it out on 12/04 - as I thought this was not necessary. */ /* zero volume boxes - still look at for getting the bounding box */ if (hypre_BoxVolume(box) == 0) /* zero volume boxes - still count */ { hypre_CopyBox(box, grow_box); for (d = 0; d < 3; d++) { if(!hypre_BoxSizeD(box, d)) { grow = (hypre_BoxIMinD(box, d) - hypre_BoxIMaxD(box, d) + 1)/2; grow_array[2*d] = grow; grow_array[2*d+1] = grow; } else { grow_array[2*d] = 0; grow_array[2*d+1] = 0; } } /* expand the box */ hypre_BoxExpand(grow_box, grow_array); box = grow_box; /*pointer copy*/ } /*now we have a vol > 0 box */ for (d = 0; d < dim; d++) /* for each dimension */ { hypre_IndexD(min_index, d) = hypre_min( hypre_IndexD(min_index, d), hypre_BoxIMinD(box, d)); hypre_IndexD(max_index, d) = hypre_max( hypre_IndexD(max_index, d), hypre_BoxIMaxD(box, d)); } }/*end for each box loop */ /* bounding box extents */ hypre_BoxSetExtents(bounding_box, min_index, max_index); }
/*-------------------------------------------------------------------------- * hypre_CFInterfaceExtents: Given a cgrid_box, a fgrid_box, and stencils, * find the extents of the C/F interface (interface nodes in the C box). * Boxes corresponding to stencil shifts are stored in the first stencil_size * boxes, and the union of these are appended to the end of the returned * box_array. *--------------------------------------------------------------------------*/ hypre_BoxArray * hypre_CFInterfaceExtents( hypre_Box *fgrid_box, hypre_Box *cgrid_box, hypre_StructStencil *stencils, hypre_Index rfactors ) { hypre_BoxArray *stencil_box_extents; hypre_BoxArray *union_boxes; hypre_Box *cfine_box; hypre_Box *box; hypre_Index stencil_shape, cstart, zero_index, neg_index; HYPRE_Int stencil_size; HYPRE_Int abs_stencil; HYPRE_Int ndim= hypre_StructStencilDim(stencils); HYPRE_Int i, j; hypre_ClearIndex(zero_index); hypre_ClearIndex(neg_index); for (i= 0; i< ndim; i++) { neg_index[i]= -1; } hypre_CopyIndex(hypre_BoxIMin(cgrid_box), cstart); stencil_size = hypre_StructStencilSize(stencils); stencil_box_extents= hypre_BoxArrayCreate(stencil_size); union_boxes = hypre_BoxArrayCreate(0); for (i= 0; i< stencil_size; i++) { hypre_CopyIndex(hypre_StructStencilElement(stencils, i), stencil_shape); AbsStencilShape(stencil_shape, abs_stencil); if (abs_stencil) /* only do if not the centre stencil */ { cfine_box= hypre_CF_StenBox(fgrid_box, cgrid_box, stencil_shape, rfactors, ndim); if ( hypre_BoxVolume(cfine_box) ) { hypre_AppendBox(cfine_box, union_boxes); hypre_CopyBox(cfine_box, hypre_BoxArrayBox(stencil_box_extents, i)); for (j= 0; j< ndim; j++) { hypre_BoxIMin(cfine_box)[j]-= cstart[j]; hypre_BoxIMax(cfine_box)[j]-= cstart[j]; } hypre_CopyBox(cfine_box, hypre_BoxArrayBox(stencil_box_extents, i)); } else { hypre_BoxSetExtents(hypre_BoxArrayBox(stencil_box_extents, i), zero_index, neg_index); } hypre_BoxDestroy(cfine_box); } else /* centre */ { hypre_BoxSetExtents(hypre_BoxArrayBox(stencil_box_extents, i), zero_index, neg_index); } } /*-------------------------------------------------------------------------- * Union the stencil_box_extents to get the full CF extents and append to * the end of the stencil_box_extents BoxArray. Then shift the unioned boxes * by cstart. *--------------------------------------------------------------------------*/ if (hypre_BoxArraySize(union_boxes) > 1) { hypre_UnionBoxes(union_boxes); } hypre_ForBoxI(i, union_boxes) { hypre_AppendBox(hypre_BoxArrayBox(union_boxes, i), stencil_box_extents); }