xref: /petsc/src/vec/is/is/impls/block/block.c (revision 95b2e421c83884142534a24efe8ff5b5a32bde75)
1 
2 /*
3      Provides the functions for index sets (IS) defined by a list of integers.
4    These are for blocks of data, each block is indicated with a single integer.
5 */
6 #include <petsc/private/isimpl.h>               /*I  "petscis.h"     I*/
7 #include <petscviewer.h>
8 
9 typedef struct {
10   PetscBool sorted;      /* are the blocks sorted? */
11   PetscBool allocated;   /* did we allocate the index array ourselves? */
12   PetscInt  *idx;
13 } IS_Block;
14 
15 static PetscErrorCode ISDestroy_Block(IS is)
16 {
17   IS_Block       *sub = (IS_Block*)is->data;
18 
19   PetscFunctionBegin;
20   if (sub->allocated) PetscCall(PetscFree(sub->idx));
21   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",NULL));
22   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",NULL));
23   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",NULL));
24   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",NULL));
25   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",NULL));
26   PetscCall(PetscFree(is->data));
27   PetscFunctionReturn(0);
28 }
29 
30 static PetscErrorCode ISLocate_Block(IS is,PetscInt key,PetscInt *location)
31 {
32   IS_Block       *sub = (IS_Block*)is->data;
33   PetscInt       numIdx, i, bs, bkey, mkey;
34   PetscBool      sorted;
35 
36   PetscFunctionBegin;
37   PetscCall(PetscLayoutGetBlockSize(is->map,&bs));
38   PetscCall(PetscLayoutGetSize(is->map,&numIdx));
39   numIdx /= bs;
40   bkey    = key / bs;
41   mkey    = key % bs;
42   if (mkey < 0) {
43     bkey--;
44     mkey += bs;
45   }
46   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sorted));
47   if (sorted) {
48     PetscCall(PetscFindInt(bkey,numIdx,sub->idx,location));
49   } else {
50     const PetscInt *idx = sub->idx;
51 
52     *location = -1;
53     for (i = 0; i < numIdx; i++) {
54       if (idx[i] == bkey) {
55         *location = i;
56         break;
57       }
58     }
59   }
60   if (*location >= 0) {
61     *location = *location * bs + mkey;
62   }
63   PetscFunctionReturn(0);
64 }
65 
66 static PetscErrorCode ISGetIndices_Block(IS in,const PetscInt *idx[])
67 {
68   IS_Block       *sub = (IS_Block*)in->data;
69   PetscInt       i,j,k,bs,n,*ii,*jj;
70 
71   PetscFunctionBegin;
72   PetscCall(PetscLayoutGetBlockSize(in->map, &bs));
73   PetscCall(PetscLayoutGetLocalSize(in->map, &n));
74   n   /= bs;
75   if (bs == 1) *idx = sub->idx;
76   else {
77     if (n) {
78       PetscCall(PetscMalloc1(bs*n,&jj));
79       *idx = jj;
80       k    = 0;
81       ii   = sub->idx;
82       for (i=0; i<n; i++)
83         for (j=0; j<bs; j++)
84           jj[k++] = bs*ii[i] + j;
85     } else {
86       /* do not malloc for zero size because F90Array1dCreate() inside ISRestoreArrayF90() does not keep array when zero length array */
87       *idx = NULL;
88     }
89   }
90   PetscFunctionReturn(0);
91 }
92 
93 static PetscErrorCode ISRestoreIndices_Block(IS is,const PetscInt *idx[])
94 {
95   IS_Block       *sub = (IS_Block*)is->data;
96   PetscInt       bs;
97 
98   PetscFunctionBegin;
99   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
100   if (bs != 1) {
101     PetscCall(PetscFree(*(void**)idx));
102   } else {
103     /* F90Array1dCreate() inside ISRestoreArrayF90() does not keep array when zero length array */
104     PetscCheck(is->map->n <= 0 || *idx == sub->idx,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must restore with value from ISGetIndices()");
105   }
106   PetscFunctionReturn(0);
107 }
108 
109 static PetscErrorCode ISInvertPermutation_Block(IS is,PetscInt nlocal,IS *isout)
110 {
111   IS_Block       *sub = (IS_Block*)is->data;
112   PetscInt       i,*ii,bs,n,*idx = sub->idx;
113   PetscMPIInt    size;
114 
115   PetscFunctionBegin;
116   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is),&size));
117   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
118   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
119   n   /= bs;
120   if (size == 1) {
121     PetscCall(PetscMalloc1(n,&ii));
122     for (i=0; i<n; i++) ii[idx[i]] = i;
123     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,n,ii,PETSC_OWN_POINTER,isout));
124     PetscCall(ISSetPermutation(*isout));
125   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No inversion written yet for block IS");
126   PetscFunctionReturn(0);
127 }
128 
129 static PetscErrorCode ISView_Block(IS is, PetscViewer viewer)
130 {
131   IS_Block       *sub = (IS_Block*)is->data;
132   PetscInt       i,bs,n,*idx = sub->idx;
133   PetscBool      iascii,ibinary;
134 
135   PetscFunctionBegin;
136   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
137   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
138   n   /= bs;
139   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
140   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&ibinary));
141   if (iascii) {
142     PetscViewerFormat fmt;
143 
144     PetscCall(PetscViewerGetFormat(viewer,&fmt));
145     if (fmt == PETSC_VIEWER_ASCII_MATLAB) {
146       IS             ist;
147       const char     *name;
148       const PetscInt *idx;
149       PetscInt       n;
150 
151       PetscCall(PetscObjectGetName((PetscObject)is,&name));
152       PetscCall(ISGetLocalSize(is,&n));
153       PetscCall(ISGetIndices(is,&idx));
154       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idx,PETSC_USE_POINTER,&ist));
155       PetscCall(PetscObjectSetName((PetscObject)ist,name));
156       PetscCall(ISView(ist,viewer));
157       PetscCall(ISDestroy(&ist));
158       PetscCall(ISRestoreIndices(is,&idx));
159     } else {
160       PetscBool isperm;
161 
162       PetscCall(ISGetInfo(is,IS_PERMUTATION,IS_GLOBAL,PETSC_FALSE,&isperm));
163       if (isperm) PetscCall(PetscViewerASCIIPrintf(viewer,"Block Index set is permutation\n"));
164       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
165       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Block size %" PetscInt_FMT "\n",bs));
166       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Number of block indices in set %" PetscInt_FMT "\n",n));
167       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"The first indices of each block are\n"));
168       for (i=0; i<n; i++) {
169         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Block %" PetscInt_FMT " Index %" PetscInt_FMT "\n",i,idx[i]));
170       }
171       PetscCall(PetscViewerFlush(viewer));
172       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
173     }
174   } else if (ibinary) PetscCall(ISView_Binary(is,viewer));
175   PetscFunctionReturn(0);
176 }
177 
178 static PetscErrorCode ISSort_Block(IS is)
179 {
180   IS_Block       *sub = (IS_Block*)is->data;
181   PetscInt       bs, n;
182 
183   PetscFunctionBegin;
184   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
185   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
186   PetscCall(PetscIntSortSemiOrdered(n/bs,sub->idx));
187   PetscFunctionReturn(0);
188 }
189 
190 static PetscErrorCode ISSortRemoveDups_Block(IS is)
191 {
192   IS_Block       *sub = (IS_Block*)is->data;
193   PetscInt       bs, n, nb;
194   PetscBool      sorted;
195 
196   PetscFunctionBegin;
197   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
198   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
199   nb   = n/bs;
200   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sorted));
201   if (sorted) {
202     PetscCall(PetscSortedRemoveDupsInt(&nb,sub->idx));
203   } else {
204     PetscCall(PetscSortRemoveDupsInt(&nb,sub->idx));
205   }
206   PetscCall(PetscLayoutDestroy(&is->map));
207   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), nb*bs, PETSC_DECIDE, bs, &is->map));
208   PetscFunctionReturn(0);
209 }
210 
211 static PetscErrorCode ISSorted_Block(IS is,PetscBool  *flg)
212 {
213   PetscFunctionBegin;
214   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,flg));
215   PetscFunctionReturn(0);
216 }
217 
218 static PetscErrorCode ISSortedLocal_Block(IS is,PetscBool *flg)
219 {
220   IS_Block       *sub = (IS_Block*)is->data;
221   PetscInt       n, bs, i, *idx;
222 
223   PetscFunctionBegin;
224   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
225   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
226   n   /= bs;
227   idx  = sub->idx;
228   for (i = 1; i < n; i++) if (idx[i] < idx[i - 1]) break;
229   if (i < n) *flg = PETSC_FALSE;
230   else       *flg = PETSC_TRUE;
231   PetscFunctionReturn(0);
232 }
233 
234 static PetscErrorCode ISUniqueLocal_Block(IS is,PetscBool *flg)
235 {
236   IS_Block       *sub = (IS_Block*)is->data;
237   PetscInt       n, bs, i, *idx, *idxcopy = NULL;
238   PetscBool      sortedLocal;
239 
240   PetscFunctionBegin;
241   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
242   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
243   n   /= bs;
244   idx  = sub->idx;
245   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal));
246   if (!sortedLocal) {
247     PetscCall(PetscMalloc1(n, &idxcopy));
248     PetscCall(PetscArraycpy(idxcopy, idx, n));
249     PetscCall(PetscIntSortSemiOrdered(n, idxcopy));
250     idx = idxcopy;
251   }
252   for (i = 1; i < n; i++) if (idx[i] == idx[i - 1]) break;
253   if (i < n) *flg = PETSC_FALSE;
254   else       *flg = PETSC_TRUE;
255   PetscCall(PetscFree(idxcopy));
256   PetscFunctionReturn(0);
257 }
258 
259 static PetscErrorCode ISPermutationLocal_Block(IS is,PetscBool *flg)
260 {
261   IS_Block       *sub = (IS_Block*)is->data;
262   PetscInt       n, bs, i, *idx, *idxcopy = NULL;
263   PetscBool      sortedLocal;
264 
265   PetscFunctionBegin;
266   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
267   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
268   n   /= bs;
269   idx  = sub->idx;
270   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal));
271   if (!sortedLocal) {
272     PetscCall(PetscMalloc1(n, &idxcopy));
273     PetscCall(PetscArraycpy(idxcopy, idx, n));
274     PetscCall(PetscIntSortSemiOrdered(n, idxcopy));
275     idx = idxcopy;
276   }
277   for (i = 0; i < n; i++) if (idx[i] != i) break;
278   if (i < n) *flg = PETSC_FALSE;
279   else       *flg = PETSC_TRUE;
280   PetscCall(PetscFree(idxcopy));
281   PetscFunctionReturn(0);
282 }
283 
284 static PetscErrorCode ISIntervalLocal_Block(IS is,PetscBool *flg)
285 {
286   IS_Block       *sub = (IS_Block*)is->data;
287   PetscInt       n, bs, i, *idx;
288 
289   PetscFunctionBegin;
290   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
291   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
292   n   /= bs;
293   idx  = sub->idx;
294   for (i = 1; i < n; i++) if (idx[i] != idx[i - 1] + 1) break;
295   if (i < n) *flg = PETSC_FALSE;
296   else       *flg = PETSC_TRUE;
297   PetscFunctionReturn(0);
298 }
299 
300 static PetscErrorCode ISDuplicate_Block(IS is,IS *newIS)
301 {
302   IS_Block       *sub = (IS_Block*)is->data;
303   PetscInt        bs, n;
304 
305   PetscFunctionBegin;
306   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
307   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
308   n   /= bs;
309   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)is),bs,n,sub->idx,PETSC_COPY_VALUES,newIS));
310   PetscFunctionReturn(0);
311 }
312 
313 static PetscErrorCode ISCopy_Block(IS is,IS isy)
314 {
315   IS_Block       *is_block = (IS_Block*)is->data,*isy_block = (IS_Block*)isy->data;
316   PetscInt       bs, n;
317 
318   PetscFunctionBegin;
319   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
320   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
321   PetscCall(PetscArraycpy(isy_block->idx, is_block->idx, n/bs));
322   PetscFunctionReturn(0);
323 }
324 
325 static PetscErrorCode ISOnComm_Block(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
326 {
327   IS_Block       *sub = (IS_Block*)is->data;
328   PetscInt       bs, n;
329 
330   PetscFunctionBegin;
331   PetscCheck(mode != PETSC_OWN_POINTER,comm,PETSC_ERR_ARG_WRONG,"Cannot use PETSC_OWN_POINTER");
332   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
333   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
334   PetscCall(ISCreateBlock(comm,bs,n/bs,sub->idx,mode,newis));
335   PetscFunctionReturn(0);
336 }
337 
338 static PetscErrorCode ISSetBlockSize_Block(IS is,PetscInt bs)
339 {
340   PetscFunctionBegin;
341   PetscCheck(is->map->bs <= 0 || bs == is->map->bs,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot change blocksize %" PetscInt_FMT " (to %" PetscInt_FMT ") if ISType is ISBLOCK",is->map->bs,bs);
342   PetscCall(PetscLayoutSetBlockSize(is->map, bs));
343   PetscFunctionReturn(0);
344 }
345 
346 static PetscErrorCode ISToGeneral_Block(IS inis)
347 {
348   IS_Block       *sub   = (IS_Block*)inis->data;
349   PetscInt       bs,n;
350   const PetscInt *idx;
351 
352   PetscFunctionBegin;
353   PetscCall(ISGetBlockSize(inis,&bs));
354   PetscCall(ISGetLocalSize(inis,&n));
355   PetscCall(ISGetIndices(inis,&idx));
356   if (bs == 1) {
357     PetscCopyMode mode = sub->allocated ? PETSC_OWN_POINTER : PETSC_USE_POINTER;
358     sub->allocated = PETSC_FALSE; /* prevent deallocation when changing the subtype*/
359     PetscCall(ISSetType(inis,ISGENERAL));
360     PetscCall(ISGeneralSetIndices(inis,n,idx,mode));
361   } else {
362     PetscCall(ISSetType(inis,ISGENERAL));
363     PetscCall(ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER));
364   }
365   PetscFunctionReturn(0);
366 }
367 
368 static struct _ISOps myops = { ISGetIndices_Block,
369                                ISRestoreIndices_Block,
370                                ISInvertPermutation_Block,
371                                ISSort_Block,
372                                ISSortRemoveDups_Block,
373                                ISSorted_Block,
374                                ISDuplicate_Block,
375                                ISDestroy_Block,
376                                ISView_Block,
377                                ISLoad_Default,
378                                ISCopy_Block,
379                                ISToGeneral_Block,
380                                ISOnComm_Block,
381                                ISSetBlockSize_Block,
382                                NULL,
383                                ISLocate_Block,
384                                /* we can have specialized local routines for determining properties,
385                                 * but unless the block size is the same on each process (which is not guaranteed at
386                                 * the moment), then trying to do something specialized for global properties is too
387                                 * complicated */
388                                ISSortedLocal_Block,
389                                NULL,
390                                ISUniqueLocal_Block,
391                                NULL,
392                                ISPermutationLocal_Block,
393                                NULL,
394                                ISIntervalLocal_Block,
395                                NULL};
396 
397 /*@
398    ISBlockSetIndices - Set integers representing blocks of indices in an index set.
399 
400    Collective on IS
401 
402    Input Parameters:
403 +  is - the index set
404 .  bs - number of elements in each block
405 .   n - the length of the index set (the number of blocks)
406 .  idx - the list of integers, one for each block, the integers contain the index of the first index of each block divided by the block size
407 -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported
408 
409    Notes:
410    When the communicator is not MPI_COMM_SELF, the operations on the
411    index sets, IS, are NOT conceptually the same as MPI_Group operations.
412    The index sets are then distributed sets of indices and thus certain operations
413    on them are collective.
414 
415    Example:
416    If you wish to index the values {0,1,4,5}, then use
417    a block size of 2 and idx of {0,2}.
418 
419    Level: beginner
420 
421 .seealso: `ISCreateStride()`, `ISCreateGeneral()`, `ISAllGather()`, `ISCreateBlock()`, `ISBLOCK`, `ISGeneralSetIndices()`
422 @*/
423 PetscErrorCode  ISBlockSetIndices(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
424 {
425   PetscFunctionBegin;
426   PetscCall(ISClearInfoCache(is,PETSC_FALSE));
427   PetscUseMethod(is,"ISBlockSetIndices_C",(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode),(is,bs,n,idx,mode));
428   PetscFunctionReturn(0);
429 }
430 
431 static PetscErrorCode  ISBlockSetIndices_Block(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
432 {
433   PetscInt       i,min,max;
434   IS_Block       *sub = (IS_Block*)is->data;
435   PetscLayout    map;
436 
437   PetscFunctionBegin;
438   PetscCheck(bs >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1");
439   PetscCheck(n >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
440   if (n) PetscValidIntPointer(idx,4);
441 
442   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is),n*bs,is->map->N,bs,&map));
443   PetscCall(PetscLayoutDestroy(&is->map));
444   is->map = map;
445 
446   if (sub->allocated) PetscCall(PetscFree(sub->idx));
447   if (mode == PETSC_COPY_VALUES) {
448     PetscCall(PetscMalloc1(n,&sub->idx));
449     PetscCall(PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt)));
450     PetscCall(PetscArraycpy(sub->idx,idx,n));
451     sub->allocated = PETSC_TRUE;
452   } else if (mode == PETSC_OWN_POINTER) {
453     sub->idx = (PetscInt*) idx;
454     PetscCall(PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt)));
455     sub->allocated = PETSC_TRUE;
456   } else if (mode == PETSC_USE_POINTER) {
457     sub->idx = (PetscInt*) idx;
458     sub->allocated = PETSC_FALSE;
459   }
460 
461   if (n) {
462     min = max = idx[0];
463     for (i=1; i<n; i++) {
464       if (idx[i] < min) min = idx[i];
465       if (idx[i] > max) max = idx[i];
466     }
467     is->min = bs*min;
468     is->max = bs*max+bs-1;
469   } else {
470     is->min = PETSC_MAX_INT;
471     is->max = PETSC_MIN_INT;
472   }
473   PetscFunctionReturn(0);
474 }
475 
476 /*@
477    ISCreateBlock - Creates a data structure for an index set containing
478    a list of integers. Each integer represents a fixed block size set of indices.
479 
480    Collective
481 
482    Input Parameters:
483 +  comm - the MPI communicator
484 .  bs - number of elements in each block
485 .  n - the length of the index set (the number of blocks)
486 .  idx - the list of integers, one for each block, the integers contain the index of the first entry of each block divided by the block size
487 -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported in this routine
488 
489    Output Parameter:
490 .  is - the new index set
491 
492    Notes:
493    When the communicator is not MPI_COMM_SELF, the operations on the
494    index sets, IS, are NOT conceptually the same as MPI_Group operations.
495    The index sets are then distributed sets of indices and thus certain operations
496    on them are collective.
497 
498    Example:
499    If you wish to index the values {0,1,6,7}, then use
500    a block size of 2 and idx of {0,3}.
501 
502    Level: beginner
503 
504 .seealso: `ISCreateStride()`, `ISCreateGeneral()`, `ISAllGather()`, `ISBlockSetIndices()`, `ISBLOCK`, `ISGENERAL`
505 @*/
506 PetscErrorCode  ISCreateBlock(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode,IS *is)
507 {
508   PetscFunctionBegin;
509   PetscValidPointer(is,6);
510   PetscCheck(bs >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1");
511   PetscCheck(n >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
512   if (n) PetscValidIntPointer(idx,4);
513 
514   PetscCall(ISCreate(comm,is));
515   PetscCall(ISSetType(*is,ISBLOCK));
516   PetscCall(ISBlockSetIndices(*is,bs,n,idx,mode));
517   PetscFunctionReturn(0);
518 }
519 
520 static PetscErrorCode  ISBlockGetIndices_Block(IS is,const PetscInt *idx[])
521 {
522   IS_Block *sub = (IS_Block*)is->data;
523 
524   PetscFunctionBegin;
525   *idx = sub->idx;
526   PetscFunctionReturn(0);
527 }
528 
529 static PetscErrorCode  ISBlockRestoreIndices_Block(IS is,const PetscInt *idx[])
530 {
531   PetscFunctionBegin;
532   PetscFunctionReturn(0);
533 }
534 
535 /*@C
536    ISBlockGetIndices - Gets the indices associated with each block.
537 
538    Not Collective
539 
540    Input Parameter:
541 .  is - the index set
542 
543    Output Parameter:
544 .  idx - the integer indices, one for each block and count of block not indices
545 
546    Level: intermediate
547 
548 .seealso: `ISGetIndices()`, `ISBlockRestoreIndices()`, `ISBLOCK`, `ISBlockSetIndices()`, `ISCreateBlock()`
549 @*/
550 PetscErrorCode  ISBlockGetIndices(IS is,const PetscInt *idx[])
551 {
552   PetscFunctionBegin;
553   PetscUseMethod(is,"ISBlockGetIndices_C",(IS,const PetscInt*[]),(is,idx));
554   PetscFunctionReturn(0);
555 }
556 
557 /*@C
558    ISBlockRestoreIndices - Restores the indices associated with each block.
559 
560    Not Collective
561 
562    Input Parameter:
563 .  is - the index set
564 
565    Output Parameter:
566 .  idx - the integer indices
567 
568    Level: intermediate
569 
570 .seealso: `ISRestoreIndices()`, `ISBlockGetIndices()`
571 @*/
572 PetscErrorCode  ISBlockRestoreIndices(IS is,const PetscInt *idx[])
573 {
574   PetscFunctionBegin;
575   PetscUseMethod(is,"ISBlockRestoreIndices_C",(IS,const PetscInt*[]),(is,idx));
576   PetscFunctionReturn(0);
577 }
578 
579 /*@
580    ISBlockGetLocalSize - Returns the local number of blocks in the index set.
581 
582    Not Collective
583 
584    Input Parameter:
585 .  is - the index set
586 
587    Output Parameter:
588 .  size - the local number of blocks
589 
590    Level: intermediate
591 
592 .seealso: `ISGetBlockSize()`, `ISBlockGetSize()`, `ISGetSize()`, `ISCreateBlock()`, `ISBLOCK`
593 @*/
594 PetscErrorCode  ISBlockGetLocalSize(IS is,PetscInt *size)
595 {
596   PetscFunctionBegin;
597   PetscUseMethod(is,"ISBlockGetLocalSize_C",(IS,PetscInt*),(is,size));
598   PetscFunctionReturn(0);
599 }
600 
601 static PetscErrorCode  ISBlockGetLocalSize_Block(IS is,PetscInt *size)
602 {
603   PetscInt       bs, n;
604 
605   PetscFunctionBegin;
606   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
607   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
608   *size = n/bs;
609   PetscFunctionReturn(0);
610 }
611 
612 /*@
613    ISBlockGetSize - Returns the global number of blocks in the index set.
614 
615    Not Collective
616 
617    Input Parameter:
618 .  is - the index set
619 
620    Output Parameter:
621 .  size - the global number of blocks
622 
623    Level: intermediate
624 
625 .seealso: `ISGetBlockSize()`, `ISBlockGetLocalSize()`, `ISGetSize()`, `ISCreateBlock()`, `ISBLOCK`
626 @*/
627 PetscErrorCode  ISBlockGetSize(IS is,PetscInt *size)
628 {
629   PetscFunctionBegin;
630   PetscUseMethod(is,"ISBlockGetSize_C",(IS,PetscInt*),(is,size));
631   PetscFunctionReturn(0);
632 }
633 
634 static PetscErrorCode  ISBlockGetSize_Block(IS is,PetscInt *size)
635 {
636   PetscInt       bs, N;
637 
638   PetscFunctionBegin;
639   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
640   PetscCall(PetscLayoutGetSize(is->map, &N));
641   *size = N/bs;
642   PetscFunctionReturn(0);
643 }
644 
645 PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is)
646 {
647   IS_Block       *sub;
648 
649   PetscFunctionBegin;
650   PetscCall(PetscNewLog(is,&sub));
651   is->data = (void *) sub;
652   PetscCall(PetscMemcpy(is->ops,&myops,sizeof(myops)));
653   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",ISBlockSetIndices_Block));
654   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",ISBlockGetIndices_Block));
655   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",ISBlockRestoreIndices_Block));
656   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",ISBlockGetSize_Block));
657   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",ISBlockGetLocalSize_Block));
658   PetscFunctionReturn(0);
659 }
660