예제 #1
0
 hypre_ForBoxI(i, boxes)
    {
       cbox_array = hypre_BoxArrayArrayBoxArray(dept_boxes, i);
       hypre_BoxArraySetSize(cbox_array, 1);
       cbox = hypre_BoxArrayBox(cbox_array, 0);
       hypre_CopyBox(hypre_BoxArrayBox(boxes, i), cbox);
    }
예제 #2
0
int
hypre_SubtractBoxArrays( hypre_BoxArray *box_array1,
                         hypre_BoxArray *box_array2,
                         hypre_BoxArray *tmp_box_array )
{
   int ierr = 0;
              
   hypre_BoxArray *diff_boxes     = box_array1;
   hypre_BoxArray *new_diff_boxes = tmp_box_array;
   hypre_BoxArray  box_array;
   hypre_Box      *box1;
   hypre_Box      *box2;
   int             i, k;

   hypre_ForBoxI(i, box_array2)
      {
         box2 = hypre_BoxArrayBox(box_array2, i);

         /* compute new_diff_boxes = (diff_boxes - box2) */
         hypre_BoxArraySetSize(new_diff_boxes, 0);
         hypre_ForBoxI(k, diff_boxes)
            {
               box1 = hypre_BoxArrayBox(diff_boxes, k);
               hypre_SubtractBoxes(box1, box2, new_diff_boxes);
            }
예제 #3
0
HYPRE_Int
hypre_BoxBoundaryIntersect( hypre_Box *box,
                            hypre_StructGrid *grid,
                            HYPRE_Int d,
                            HYPRE_Int dir,
                            hypre_BoxArray *boundary )
{
   HYPRE_Int           ndim = hypre_BoxNDim(box);
   hypre_BoxManager   *boxman;
   hypre_BoxManEntry **entries;
   hypre_BoxArray     *int_boxes, *tmp_boxes;
   hypre_Box          *bbox, *ibox;
   HYPRE_Int           nentries, i;

   /* set bbox to the box surface of interest */
   hypre_BoxArraySetSize(boundary, 1);
   bbox = hypre_BoxArrayBox(boundary, 0);
   hypre_CopyBox(box, bbox);
   if (dir > 0)
   {
      hypre_BoxIMinD(bbox, d) = hypre_BoxIMaxD(bbox, d);
   }
   else if (dir < 0)
   {
      hypre_BoxIMaxD(bbox, d) = hypre_BoxIMinD(bbox, d);
   }

   /* temporarily shift bbox in direction dir and intersect with the grid */
   hypre_BoxIMinD(bbox, d) += dir;
   hypre_BoxIMaxD(bbox, d) += dir;
   boxman = hypre_StructGridBoxMan(grid);
   hypre_BoxManIntersect(boxman, hypre_BoxIMin(bbox), hypre_BoxIMax(bbox),
                         &entries, &nentries);
   hypre_BoxIMinD(bbox, d) -= dir;
   hypre_BoxIMaxD(bbox, d) -= dir;

   /* shift intersected boxes in direction -dir and subtract from bbox */
   int_boxes  = hypre_BoxArrayCreate(nentries, ndim);
   tmp_boxes  = hypre_BoxArrayCreate(0, ndim);
   for (i = 0; i < nentries; i++)
   {
      ibox = hypre_BoxArrayBox(int_boxes, i);
      hypre_BoxManEntryGetExtents(
         entries[i], hypre_BoxIMin(ibox), hypre_BoxIMax(ibox));
      hypre_BoxIMinD(ibox, d) -= dir;
      hypre_BoxIMaxD(ibox, d) -= dir;
   }
   hypre_SubtractBoxArrays(boundary, int_boxes, tmp_boxes);

   hypre_BoxArrayDestroy(int_boxes);
   hypre_BoxArrayDestroy(tmp_boxes);
   hypre_TFree(entries);

   return hypre_error_flag;
}
예제 #4
0
   hypre_ForBoxI(i, boxes)
      {
         cbox_array = hypre_BoxArrayArrayBoxArray(indt_boxes, i);
         hypre_BoxArraySetSize(cbox_array, 1);
         cbox = hypre_BoxArrayBox(cbox_array, 0);
         hypre_CopyBox(hypre_BoxArrayBox(boxes, i), cbox);

         for (d = 0; d < 3; d++)
         {
            if ( (border[d][0]) )
            {
               hypre_BoxIMinD(cbox, d) += border[d][0];
            }
            if ( (border[d][1]) )
            {
               hypre_BoxIMaxD(cbox, d) -= border[d][1];
            }
         }
      }
예제 #5
0
int
hypre_CreateComputeInfo( hypre_StructGrid      *grid,
                         hypre_StructStencil   *stencil,
                         hypre_ComputeInfo    **compute_info_ptr )
{
   int                      ierr = 0;

   hypre_CommInfo          *comm_info;
   hypre_BoxArrayArray     *indt_boxes;
   hypre_BoxArrayArray     *dept_boxes;

   hypre_BoxArray          *boxes;

   hypre_BoxArray          *cbox_array;
   hypre_Box               *cbox;

   int                      i;

#ifdef HYPRE_OVERLAP_COMM_COMP
   hypre_Box               *rembox;
   hypre_Index             *stencil_shape;
   int                      border[3][2] = {{0, 0}, {0, 0}, {0, 0}};
   int                      cbox_array_size;
   int                      s, d;
#endif

   /*------------------------------------------------------
    * Extract needed grid info
    *------------------------------------------------------*/

   boxes = hypre_StructGridBoxes(grid);

   /*------------------------------------------------------
    * Get communication info
    *------------------------------------------------------*/

   hypre_CreateCommInfoFromStencil(grid, stencil, &comm_info);

#ifdef HYPRE_OVERLAP_COMM_COMP

   /*------------------------------------------------------
    * Compute border info
    *------------------------------------------------------*/

   stencil_shape = hypre_StructStencilShape(stencil);
   for (s = 0; s < hypre_StructStencilSize(stencil); s++)
   {
      for (d = 0; d < 3; d++)
      {
         i = hypre_IndexD(stencil_shape[s], d);
         if (i < 0)
         {
            border[d][0] = hypre_max(border[d][0], -i);
         }
         else if (i > 0)
         {
            border[d][1] = hypre_max(border[d][1], i);
         }
      }
   }

   /*------------------------------------------------------
    * Set up the dependent boxes
    *------------------------------------------------------*/

   dept_boxes = hypre_BoxArrayArrayCreate(hypre_BoxArraySize(boxes));

   rembox = hypre_BoxCreate();
   hypre_ForBoxI(i, boxes)
      {
         cbox_array = hypre_BoxArrayArrayBoxArray(dept_boxes, i);
         hypre_BoxArraySetSize(cbox_array, 6);

         hypre_CopyBox(hypre_BoxArrayBox(boxes, i), rembox);
         cbox_array_size = 0;
         for (d = 0; d < 3; d++)
         {
            if ( (hypre_BoxVolume(rembox)) && (border[d][0]) )
            {
               cbox = hypre_BoxArrayBox(cbox_array, cbox_array_size);
               hypre_CopyBox(rembox, cbox);
               hypre_BoxIMaxD(cbox, d) =
                  hypre_BoxIMinD(cbox, d) + border[d][0] - 1;
               hypre_BoxIMinD(rembox, d) =
                  hypre_BoxIMinD(cbox, d) + border[d][0];
               cbox_array_size++;
            }
            if ( (hypre_BoxVolume(rembox)) && (border[d][1]) )
            {
               cbox = hypre_BoxArrayBox(cbox_array, cbox_array_size);
               hypre_CopyBox(rembox, cbox);
               hypre_BoxIMinD(cbox, d) =
                  hypre_BoxIMaxD(cbox, d) - border[d][1] + 1;
               hypre_BoxIMaxD(rembox, d) =
                  hypre_BoxIMaxD(cbox, d) - border[d][1];
               cbox_array_size++;
            }
         }
         hypre_BoxArraySetSize(cbox_array, cbox_array_size);
      }
예제 #6
0
int
hypre_SubtractBoxes( hypre_Box      *box1,
                     hypre_Box      *box2,
                     hypre_BoxArray *box_array )
{
   int ierr = 0;
              
   hypre_Box  *box;
   hypre_Box  *rembox;
   int         d, size;

   /*------------------------------------------------------
    * Set the box array size to the maximum possible,
    * plus one, to have space for the remainder box.
    *------------------------------------------------------*/

   size = hypre_BoxArraySize(box_array);
   hypre_BoxArraySetSize(box_array, (size + 7));

   /*------------------------------------------------------
    * Subtract the boxes by cutting box1 in x, y, then z
    *------------------------------------------------------*/

   rembox = hypre_BoxArrayBox(box_array, (size + 6));
   hypre_CopyBox(box1, rembox);

   for (d = 0; d < 3; d++)
   {
      /* if the boxes do not intersect, the subtraction is trivial */
      if ( (hypre_BoxIMinD(box2, d) > hypre_BoxIMaxD(rembox, d)) ||
           (hypre_BoxIMaxD(box2, d) < hypre_BoxIMinD(rembox, d)) )
      {
         size = hypre_BoxArraySize(box_array) - 7;
         hypre_CopyBox(box1, hypre_BoxArrayBox(box_array, size));
         size++;
         break;
      }

      /* update the box array */
      else
      {
         if ( hypre_BoxIMinD(box2, d) > hypre_BoxIMinD(rembox, d) )
         {
            box = hypre_BoxArrayBox(box_array, size);
            hypre_CopyBox(rembox, box);
            hypre_BoxIMaxD(box, d) = hypre_BoxIMinD(box2, d) - 1;
            hypre_BoxIMinD(rembox, d) = hypre_BoxIMinD(box2, d);
            if ( hypre_BoxVolume(box)>0 ) size++;
         }
         if ( hypre_BoxIMaxD(box2, d) < hypre_BoxIMaxD(rembox, d) )
         {
            box = hypre_BoxArrayBox(box_array, size);
            hypre_CopyBox(rembox, box);
            hypre_BoxIMinD(box, d) = hypre_BoxIMaxD(box2, d) + 1;
            hypre_BoxIMaxD(rembox, d) = hypre_BoxIMaxD(box2, d);
            if ( hypre_BoxVolume(box)>0 ) size++;
         }
      }
   }
   hypre_BoxArraySetSize(box_array, size);

   return ierr;
}
예제 #7
0
HYPRE_Int
hypre_StructCoarsen( hypre_StructGrid  *fgrid,
                     hypre_Index        index,
                     hypre_Index        stride,
                     HYPRE_Int          prune,
                     hypre_StructGrid **cgrid_ptr )
{
   hypre_StructGrid *cgrid;

   MPI_Comm          comm;
   HYPRE_Int         ndim;

   hypre_BoxArray   *my_boxes;

   hypre_Index       periodic;
   hypre_Index       ilower, iupper;

   hypre_Box        *box;
   hypre_Box        *new_box;
   hypre_Box        *bounding_box;

   HYPRE_Int         i, j, myid, count;
   HYPRE_Int         info_size, max_nentries;
   HYPRE_Int         num_entries;
   HYPRE_Int        *fids, *cids;
   hypre_Index       new_dist;
   hypre_IndexRef    max_distance;
   HYPRE_Int         proc, id;
   HYPRE_Int         coarsen_factor, known;
   HYPRE_Int         num, last_proc;
#if 0
   hypre_StructAssumedPart *fap = NULL, *cap = NULL;
#endif
   hypre_BoxManager   *fboxman, *cboxman;

   hypre_BoxManEntry *entries;
   hypre_BoxManEntry  *entry;
     
   void               *entry_info = NULL;
 
#if TIME_DEBUG  
   HYPRE_Int tindex;
   char new_title[80];
   hypre_sprintf(new_title,"Coarsen.%d",s_coarsen_num);
   tindex = hypre_InitializeTiming(new_title);
   s_coarsen_num++;

   hypre_BeginTiming(tindex);
#endif

   hypre_SetIndex(ilower, 0);
   hypre_SetIndex(iupper, 0);

   /* get relevant information from the fine grid */
   fids = hypre_StructGridIDs(fgrid);
   fboxman = hypre_StructGridBoxMan(fgrid);
   comm  = hypre_StructGridComm(fgrid);
   ndim  = hypre_StructGridNDim(fgrid);
   max_distance = hypre_StructGridMaxDistance(fgrid);
   
   /* initial */
   hypre_MPI_Comm_rank(comm, &myid );

   /* create new coarse grid */
   hypre_StructGridCreate(comm, ndim, &cgrid);

   /* coarsen my boxes and create the coarse grid ids (same as fgrid) */
   my_boxes = hypre_BoxArrayDuplicate(hypre_StructGridBoxes(fgrid));
   cids = hypre_TAlloc(HYPRE_Int,  hypre_BoxArraySize(my_boxes));
   for (i = 0; i < hypre_BoxArraySize(my_boxes); i++)
   {
      box = hypre_BoxArrayBox(my_boxes, i);
      hypre_StructCoarsenBox(box, index, stride);
      cids[i] = fids[i];
   }
   
   /* prune? */
   /* zero volume boxes are needed when forming P and P^T */ 
   if (prune)
   {
      count = 0;    
      hypre_ForBoxI(i, my_boxes)
      {
         box = hypre_BoxArrayBox(my_boxes, i);
         if (hypre_BoxVolume(box))
         {
            hypre_CopyBox(box, hypre_BoxArrayBox(my_boxes, count));
            cids[count] = cids[i];
            count++;
         }
      }
      hypre_BoxArraySetSize(my_boxes, count);
   }
예제 #8
0
int
hypre_PointRelaxSetup( void               *relax_vdata,
                       hypre_StructMatrix *A,
                       hypre_StructVector *b,
                       hypre_StructVector *x           )
{
   hypre_PointRelaxData *relax_data = (hypre_PointRelaxData *)relax_vdata;

   int                    num_pointsets    = (relax_data -> num_pointsets);
   int                   *pointset_sizes   = (relax_data -> pointset_sizes);
   hypre_Index           *pointset_strides = (relax_data -> pointset_strides);
   hypre_Index          **pointset_indices = (relax_data -> pointset_indices);
   hypre_StructVector    *t;
   int                    diag_rank;
   hypre_ComputePkg     **compute_pkgs;

   hypre_Index            unit_stride;
   hypre_Index            diag_index;
   hypre_IndexRef         stride;
   hypre_IndexRef         index;
                       
   hypre_StructGrid      *grid;
   hypre_StructStencil   *stencil;
                       
   hypre_BoxArrayArray   *send_boxes;
   hypre_BoxArrayArray   *recv_boxes;
   int                  **send_processes;
   int                  **recv_processes;
   hypre_BoxArrayArray   *indt_boxes;
   hypre_BoxArrayArray   *dept_boxes;

   hypre_BoxArrayArray   *orig_indt_boxes;
   hypre_BoxArrayArray   *orig_dept_boxes;
   hypre_BoxArrayArray   *box_aa;
   hypre_BoxArray        *box_a;
   hypre_Box             *box;
   int                    box_aa_size;
   int                    box_a_size;
   hypre_BoxArrayArray   *new_box_aa;
   hypre_BoxArray        *new_box_a;
   hypre_Box             *new_box;

   double                 scale;
   int                    frac;

   int                    i, j, k, p, m, compute_i;
   int                    ierr = 0;
                       
   /*----------------------------------------------------------
    * Set up the temp vector
    *----------------------------------------------------------*/

   if ((relax_data -> t) == NULL)
   {
      t = hypre_StructVectorCreate(hypre_StructVectorComm(b),
                                   hypre_StructVectorGrid(b));
      hypre_StructVectorSetNumGhost(t, hypre_StructVectorNumGhost(b));
      hypre_StructVectorInitialize(t);
      hypre_StructVectorAssemble(t);
      (relax_data -> t) = t;
   }

   /*----------------------------------------------------------
    * Find the matrix diagonal
    *----------------------------------------------------------*/

   grid    = hypre_StructMatrixGrid(A);
   stencil = hypre_StructMatrixStencil(A);

   hypre_SetIndex(diag_index, 0, 0, 0);
   diag_rank = hypre_StructStencilElementRank(stencil, diag_index);

   /*----------------------------------------------------------
    * Set up the compute packages
    *----------------------------------------------------------*/

   hypre_SetIndex(unit_stride, 1, 1, 1);

   compute_pkgs = hypre_CTAlloc(hypre_ComputePkg *, num_pointsets);

   for (p = 0; p < num_pointsets; p++)
   {
      hypre_CreateComputeInfo(grid, stencil,
                           &send_boxes, &recv_boxes,
                           &send_processes, &recv_processes,
                           &orig_indt_boxes, &orig_dept_boxes);

      stride = pointset_strides[p];

      for (compute_i = 0; compute_i < 2; compute_i++)
      {
         switch(compute_i)
         {
            case 0:
            box_aa = orig_indt_boxes;
            break;

            case 1:
            box_aa = orig_dept_boxes;
            break;
         }
         box_aa_size = hypre_BoxArrayArraySize(box_aa);
         new_box_aa = hypre_BoxArrayArrayCreate(box_aa_size);

         for (i = 0; i < box_aa_size; i++)
         {
            box_a = hypre_BoxArrayArrayBoxArray(box_aa, i);
            box_a_size = hypre_BoxArraySize(box_a);
            new_box_a = hypre_BoxArrayArrayBoxArray(new_box_aa, i);
            hypre_BoxArraySetSize(new_box_a, box_a_size * pointset_sizes[p]);

            k = 0;
            for (m = 0; m < pointset_sizes[p]; m++)
            {
               index  = pointset_indices[p][m];

               for (j = 0; j < box_a_size; j++)
               {
                  box = hypre_BoxArrayBox(box_a, j);
                  new_box = hypre_BoxArrayBox(new_box_a, k);
                  
                  hypre_CopyBox(box, new_box);
                  hypre_ProjectBox(new_box, index, stride);
                  
                  k++;
               }
            }
         }

         switch(compute_i)
         {
            case 0:
            indt_boxes = new_box_aa;
            break;

            case 1:
            dept_boxes = new_box_aa;
            break;
         }
      }

      hypre_ComputePkgCreate(send_boxes, recv_boxes,
                             unit_stride, unit_stride,
                             send_processes, recv_processes,
                             indt_boxes, dept_boxes,
                             stride, grid,
                             hypre_StructVectorDataSpace(x), 1,
                             &compute_pkgs[p]);

      hypre_BoxArrayArrayDestroy(orig_indt_boxes);
      hypre_BoxArrayArrayDestroy(orig_dept_boxes);
   }

   /*----------------------------------------------------------
    * Set up the relax data structure
    *----------------------------------------------------------*/

   (relax_data -> A) = hypre_StructMatrixRef(A);
   (relax_data -> x) = hypre_StructVectorRef(x);
   (relax_data -> b) = hypre_StructVectorRef(b);
   (relax_data -> diag_rank)    = diag_rank;
   (relax_data -> compute_pkgs) = compute_pkgs;

   /*-----------------------------------------------------
    * Compute flops
    *-----------------------------------------------------*/

   scale = 0.0;
   for (p = 0; p < num_pointsets; p++)
   {
      stride = pointset_strides[p];
      frac   = hypre_IndexX(stride);
      frac  *= hypre_IndexY(stride);
      frac  *= hypre_IndexZ(stride);
      scale += (pointset_sizes[p] / frac);
   }
   (relax_data -> flops) = scale * (hypre_StructMatrixGlobalSize(A) +
                                    hypre_StructVectorGlobalSize(x));

   return ierr;
}
예제 #9
0
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 */
예제 #10
0
HYPRE_Int
hypre_GeneralBoxBoundaryIntersect( hypre_Box *box,
                            hypre_StructGrid *grid,
                            hypre_Index stencil_element,
                            hypre_BoxArray *boundary )
{
   hypre_BoxManager   *boxman;
   hypre_BoxManEntry **entries;
   hypre_BoxArray     *int_boxes, *tmp_boxes;
   hypre_Box          *bbox, *ibox;
   HYPRE_Int           nentries, i, j;
   HYPRE_Int          *dd;
   HYPRE_Int           ndim;

   ndim = hypre_StructGridNDim(grid);
   dd = hypre_CTAlloc(HYPRE_Int, ndim);

   for (i=0; i < ndim; i++)
     dd[i] = hypre_IndexD(stencil_element, i);

   /* set bbox to the box surface of interest */
   hypre_BoxArraySetSize(boundary, 1);
   bbox = hypre_BoxArrayBox(boundary, 0);
   hypre_CopyBox(box, bbox);

   /* temporarily shift bbox in direction dir and intersect with the grid */
   for (i=0; i < ndim; i++)
   {
      hypre_BoxIMinD(bbox, i) += dd[i];
      hypre_BoxIMaxD(bbox, i) += dd[i];
   }

   boxman = hypre_StructGridBoxMan(grid);
   hypre_BoxManIntersect(boxman, hypre_BoxIMin(bbox), hypre_BoxIMax(bbox),
                         &entries, &nentries);
   for (i=0; i < ndim; i++)
   {
      hypre_BoxIMinD(bbox, i) -= dd[i];
      hypre_BoxIMaxD(bbox, i) -= dd[i];
   }

   /* shift intersected boxes in direction -dir and subtract from bbox */
   int_boxes  = hypre_BoxArrayCreate(nentries, ndim);
   tmp_boxes  = hypre_BoxArrayCreate(0, ndim);
   for (i = 0; i < nentries; i++)
   {
      ibox = hypre_BoxArrayBox(int_boxes, i);
      hypre_BoxManEntryGetExtents(
         entries[i], hypre_BoxIMin(ibox), hypre_BoxIMax(ibox));
      for (j=0; j < ndim; j++)
      {
         hypre_BoxIMinD(ibox, j) -= dd[j];
         hypre_BoxIMaxD(ibox, j) -= dd[j];
      }
   }
   hypre_SubtractBoxArrays(boundary, int_boxes, tmp_boxes);

   hypre_BoxArrayDestroy(int_boxes);
   hypre_BoxArrayDestroy(tmp_boxes);
   hypre_TFree(entries);
   hypre_TFree(dd);

   return hypre_error_flag;
}