xref: /petsc/src/vec/is/is/impls/block/block.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
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     PetscCheckFalse(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) {
175     PetscCall(ISView_Binary(is,viewer));
176   }
177   PetscFunctionReturn(0);
178 }
179 
180 static PetscErrorCode ISSort_Block(IS is)
181 {
182   IS_Block       *sub = (IS_Block*)is->data;
183   PetscInt       bs, n;
184 
185   PetscFunctionBegin;
186   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
187   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
188   PetscCall(PetscIntSortSemiOrdered(n/bs,sub->idx));
189   PetscFunctionReturn(0);
190 }
191 
192 static PetscErrorCode ISSortRemoveDups_Block(IS is)
193 {
194   IS_Block       *sub = (IS_Block*)is->data;
195   PetscInt       bs, n, nb;
196   PetscBool      sorted;
197 
198   PetscFunctionBegin;
199   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
200   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
201   nb   = n/bs;
202   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sorted));
203   if (sorted) {
204     PetscCall(PetscSortedRemoveDupsInt(&nb,sub->idx));
205   } else {
206     PetscCall(PetscSortRemoveDupsInt(&nb,sub->idx));
207   }
208   PetscCall(PetscLayoutDestroy(&is->map));
209   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), nb*bs, PETSC_DECIDE, bs, &is->map));
210   PetscFunctionReturn(0);
211 }
212 
213 static PetscErrorCode ISSorted_Block(IS is,PetscBool  *flg)
214 {
215   PetscFunctionBegin;
216   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,flg));
217   PetscFunctionReturn(0);
218 }
219 
220 static PetscErrorCode ISSortedLocal_Block(IS is,PetscBool *flg)
221 {
222   IS_Block       *sub = (IS_Block*)is->data;
223   PetscInt       n, bs, i, *idx;
224 
225   PetscFunctionBegin;
226   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
227   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
228   n   /= bs;
229   idx  = sub->idx;
230   for (i = 1; i < n; i++) if (idx[i] < idx[i - 1]) break;
231   if (i < n) *flg = PETSC_FALSE;
232   else       *flg = PETSC_TRUE;
233   PetscFunctionReturn(0);
234 }
235 
236 static PetscErrorCode ISUniqueLocal_Block(IS is,PetscBool *flg)
237 {
238   IS_Block       *sub = (IS_Block*)is->data;
239   PetscInt       n, bs, i, *idx, *idxcopy = NULL;
240   PetscBool      sortedLocal;
241 
242   PetscFunctionBegin;
243   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
244   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
245   n   /= bs;
246   idx  = sub->idx;
247   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal));
248   if (!sortedLocal) {
249     PetscCall(PetscMalloc1(n, &idxcopy));
250     PetscCall(PetscArraycpy(idxcopy, idx, n));
251     PetscCall(PetscIntSortSemiOrdered(n, idxcopy));
252     idx = idxcopy;
253   }
254   for (i = 1; i < n; i++) if (idx[i] == idx[i - 1]) break;
255   if (i < n) *flg = PETSC_FALSE;
256   else       *flg = PETSC_TRUE;
257   PetscCall(PetscFree(idxcopy));
258   PetscFunctionReturn(0);
259 }
260 
261 static PetscErrorCode ISPermutationLocal_Block(IS is,PetscBool *flg)
262 {
263   IS_Block       *sub = (IS_Block*)is->data;
264   PetscInt       n, bs, i, *idx, *idxcopy = NULL;
265   PetscBool      sortedLocal;
266 
267   PetscFunctionBegin;
268   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
269   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
270   n   /= bs;
271   idx  = sub->idx;
272   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal));
273   if (!sortedLocal) {
274     PetscCall(PetscMalloc1(n, &idxcopy));
275     PetscCall(PetscArraycpy(idxcopy, idx, n));
276     PetscCall(PetscIntSortSemiOrdered(n, idxcopy));
277     idx = idxcopy;
278   }
279   for (i = 0; i < n; i++) if (idx[i] != i) break;
280   if (i < n) *flg = PETSC_FALSE;
281   else       *flg = PETSC_TRUE;
282   PetscCall(PetscFree(idxcopy));
283   PetscFunctionReturn(0);
284 }
285 
286 static PetscErrorCode ISIntervalLocal_Block(IS is,PetscBool *flg)
287 {
288   IS_Block       *sub = (IS_Block*)is->data;
289   PetscInt       n, bs, i, *idx;
290 
291   PetscFunctionBegin;
292   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
293   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
294   n   /= bs;
295   idx  = sub->idx;
296   for (i = 1; i < n; i++) if (idx[i] != idx[i - 1] + 1) break;
297   if (i < n) *flg = PETSC_FALSE;
298   else       *flg = PETSC_TRUE;
299   PetscFunctionReturn(0);
300 }
301 
302 static PetscErrorCode ISDuplicate_Block(IS is,IS *newIS)
303 {
304   IS_Block       *sub = (IS_Block*)is->data;
305   PetscInt        bs, n;
306 
307   PetscFunctionBegin;
308   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
309   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
310   n   /= bs;
311   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)is),bs,n,sub->idx,PETSC_COPY_VALUES,newIS));
312   PetscFunctionReturn(0);
313 }
314 
315 static PetscErrorCode ISCopy_Block(IS is,IS isy)
316 {
317   IS_Block       *is_block = (IS_Block*)is->data,*isy_block = (IS_Block*)isy->data;
318   PetscInt       bs, n, N, bsy, ny, Ny;
319 
320   PetscFunctionBegin;
321   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
322   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
323   PetscCall(PetscLayoutGetSize(is->map, &N));
324   PetscCall(PetscLayoutGetBlockSize(isy->map, &bsy));
325   PetscCall(PetscLayoutGetLocalSize(isy->map, &ny));
326   PetscCall(PetscLayoutGetSize(isy->map, &Ny));
327   PetscCheckFalse(n != ny || N != Ny || bs != bsy,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Index sets incompatible");
328   PetscCall(PetscArraycpy(isy_block->idx,is_block->idx,(n/bs)));
329   PetscFunctionReturn(0);
330 }
331 
332 static PetscErrorCode ISOnComm_Block(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
333 {
334   IS_Block       *sub = (IS_Block*)is->data;
335   PetscInt       bs, n;
336 
337   PetscFunctionBegin;
338   PetscCheckFalse(mode == PETSC_OWN_POINTER,comm,PETSC_ERR_ARG_WRONG,"Cannot use PETSC_OWN_POINTER");
339   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
340   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
341   PetscCall(ISCreateBlock(comm,bs,n/bs,sub->idx,mode,newis));
342   PetscFunctionReturn(0);
343 }
344 
345 static PetscErrorCode ISSetBlockSize_Block(IS is,PetscInt bs)
346 {
347   PetscFunctionBegin;
348   PetscCheckFalse(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);
349   PetscCall(PetscLayoutSetBlockSize(is->map, bs));
350   PetscFunctionReturn(0);
351 }
352 
353 static PetscErrorCode ISToGeneral_Block(IS inis)
354 {
355   IS_Block       *sub   = (IS_Block*)inis->data;
356   PetscInt       bs,n;
357   const PetscInt *idx;
358 
359   PetscFunctionBegin;
360   PetscCall(ISGetBlockSize(inis,&bs));
361   PetscCall(ISGetLocalSize(inis,&n));
362   PetscCall(ISGetIndices(inis,&idx));
363   if (bs == 1) {
364     PetscCopyMode mode = sub->allocated ? PETSC_OWN_POINTER : PETSC_USE_POINTER;
365     sub->allocated = PETSC_FALSE; /* prevent deallocation when changing the subtype*/
366     PetscCall(ISSetType(inis,ISGENERAL));
367     PetscCall(ISGeneralSetIndices(inis,n,idx,mode));
368   } else {
369     PetscCall(ISSetType(inis,ISGENERAL));
370     PetscCall(ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER));
371   }
372   PetscFunctionReturn(0);
373 }
374 
375 static struct _ISOps myops = { ISGetIndices_Block,
376                                ISRestoreIndices_Block,
377                                ISInvertPermutation_Block,
378                                ISSort_Block,
379                                ISSortRemoveDups_Block,
380                                ISSorted_Block,
381                                ISDuplicate_Block,
382                                ISDestroy_Block,
383                                ISView_Block,
384                                ISLoad_Default,
385                                ISCopy_Block,
386                                ISToGeneral_Block,
387                                ISOnComm_Block,
388                                ISSetBlockSize_Block,
389                                NULL,
390                                ISLocate_Block,
391                                /* we can have specialized local routines for determining properties,
392                                 * but unless the block size is the same on each process (which is not guaranteed at
393                                 * the moment), then trying to do something specialized for global properties is too
394                                 * complicated */
395                                ISSortedLocal_Block,
396                                NULL,
397                                ISUniqueLocal_Block,
398                                NULL,
399                                ISPermutationLocal_Block,
400                                NULL,
401                                ISIntervalLocal_Block,
402                                NULL};
403 
404 /*@
405    ISBlockSetIndices - Set integers representing blocks of indices in an index set.
406 
407    Collective on IS
408 
409    Input Parameters:
410 +  is - the index set
411 .  bs - number of elements in each block
412 .   n - the length of the index set (the number of blocks)
413 .  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
414 -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported
415 
416    Notes:
417    When the communicator is not MPI_COMM_SELF, the operations on the
418    index sets, IS, are NOT conceptually the same as MPI_Group operations.
419    The index sets are then distributed sets of indices and thus certain operations
420    on them are collective.
421 
422    Example:
423    If you wish to index the values {0,1,4,5}, then use
424    a block size of 2 and idx of {0,2}.
425 
426    Level: beginner
427 
428 .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather(), ISCreateBlock(), ISBLOCK, ISGeneralSetIndices()
429 @*/
430 PetscErrorCode  ISBlockSetIndices(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
431 {
432   PetscFunctionBegin;
433   PetscCall(ISClearInfoCache(is,PETSC_FALSE));
434   PetscCall(PetscUseMethod(is,"ISBlockSetIndices_C",(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode),(is,bs,n,idx,mode)));
435   PetscFunctionReturn(0);
436 }
437 
438 static PetscErrorCode  ISBlockSetIndices_Block(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
439 {
440   PetscInt       i,min,max;
441   IS_Block       *sub = (IS_Block*)is->data;
442   PetscLayout    map;
443 
444   PetscFunctionBegin;
445   PetscCheckFalse(bs < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1");
446   PetscCheckFalse(n < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
447   if (n) PetscValidIntPointer(idx,4);
448 
449   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is),n*bs,is->map->N,bs,&map));
450   PetscCall(PetscLayoutDestroy(&is->map));
451   is->map = map;
452 
453   if (sub->allocated) PetscCall(PetscFree(sub->idx));
454   if (mode == PETSC_COPY_VALUES) {
455     PetscCall(PetscMalloc1(n,&sub->idx));
456     PetscCall(PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt)));
457     PetscCall(PetscArraycpy(sub->idx,idx,n));
458     sub->allocated = PETSC_TRUE;
459   } else if (mode == PETSC_OWN_POINTER) {
460     sub->idx = (PetscInt*) idx;
461     PetscCall(PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt)));
462     sub->allocated = PETSC_TRUE;
463   } else if (mode == PETSC_USE_POINTER) {
464     sub->idx = (PetscInt*) idx;
465     sub->allocated = PETSC_FALSE;
466   }
467 
468   if (n) {
469     min = max = idx[0];
470     for (i=1; i<n; i++) {
471       if (idx[i] < min) min = idx[i];
472       if (idx[i] > max) max = idx[i];
473     }
474     is->min = bs*min;
475     is->max = bs*max+bs-1;
476   } else {
477     is->min = PETSC_MAX_INT;
478     is->max = PETSC_MIN_INT;
479   }
480   PetscFunctionReturn(0);
481 }
482 
483 /*@
484    ISCreateBlock - Creates a data structure for an index set containing
485    a list of integers. Each integer represents a fixed block size set of indices.
486 
487    Collective
488 
489    Input Parameters:
490 +  comm - the MPI communicator
491 .  bs - number of elements in each block
492 .  n - the length of the index set (the number of blocks)
493 .  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
494 -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported in this routine
495 
496    Output Parameter:
497 .  is - the new index set
498 
499    Notes:
500    When the communicator is not MPI_COMM_SELF, the operations on the
501    index sets, IS, are NOT conceptually the same as MPI_Group operations.
502    The index sets are then distributed sets of indices and thus certain operations
503    on them are collective.
504 
505    Example:
506    If you wish to index the values {0,1,6,7}, then use
507    a block size of 2 and idx of {0,3}.
508 
509    Level: beginner
510 
511 .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather(), ISBlockSetIndices(), ISBLOCK, ISGENERAL
512 @*/
513 PetscErrorCode  ISCreateBlock(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode,IS *is)
514 {
515   PetscFunctionBegin;
516   PetscValidPointer(is,6);
517   PetscCheckFalse(bs < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1");
518   PetscCheckFalse(n < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
519   if (n) PetscValidIntPointer(idx,4);
520 
521   PetscCall(ISCreate(comm,is));
522   PetscCall(ISSetType(*is,ISBLOCK));
523   PetscCall(ISBlockSetIndices(*is,bs,n,idx,mode));
524   PetscFunctionReturn(0);
525 }
526 
527 static PetscErrorCode  ISBlockGetIndices_Block(IS is,const PetscInt *idx[])
528 {
529   IS_Block *sub = (IS_Block*)is->data;
530 
531   PetscFunctionBegin;
532   *idx = sub->idx;
533   PetscFunctionReturn(0);
534 }
535 
536 static PetscErrorCode  ISBlockRestoreIndices_Block(IS is,const PetscInt *idx[])
537 {
538   PetscFunctionBegin;
539   PetscFunctionReturn(0);
540 }
541 
542 /*@C
543    ISBlockGetIndices - Gets the indices associated with each block.
544 
545    Not Collective
546 
547    Input Parameter:
548 .  is - the index set
549 
550    Output Parameter:
551 .  idx - the integer indices, one for each block and count of block not indices
552 
553    Level: intermediate
554 
555 .seealso: ISGetIndices(), ISBlockRestoreIndices(), ISBLOCK, ISBlockSetIndices(), ISCreateBlock()
556 @*/
557 PetscErrorCode  ISBlockGetIndices(IS is,const PetscInt *idx[])
558 {
559   PetscFunctionBegin;
560   PetscCall(PetscUseMethod(is,"ISBlockGetIndices_C",(IS,const PetscInt*[]),(is,idx)));
561   PetscFunctionReturn(0);
562 }
563 
564 /*@C
565    ISBlockRestoreIndices - Restores the indices associated with each block.
566 
567    Not Collective
568 
569    Input Parameter:
570 .  is - the index set
571 
572    Output Parameter:
573 .  idx - the integer indices
574 
575    Level: intermediate
576 
577 .seealso: ISRestoreIndices(), ISBlockGetIndices()
578 @*/
579 PetscErrorCode  ISBlockRestoreIndices(IS is,const PetscInt *idx[])
580 {
581   PetscFunctionBegin;
582   PetscCall(PetscUseMethod(is,"ISBlockRestoreIndices_C",(IS,const PetscInt*[]),(is,idx)));
583   PetscFunctionReturn(0);
584 }
585 
586 /*@
587    ISBlockGetLocalSize - Returns the local number of blocks in the index set.
588 
589    Not Collective
590 
591    Input Parameter:
592 .  is - the index set
593 
594    Output Parameter:
595 .  size - the local number of blocks
596 
597    Level: intermediate
598 
599 .seealso: ISGetBlockSize(), ISBlockGetSize(), ISGetSize(), ISCreateBlock(), ISBLOCK
600 @*/
601 PetscErrorCode  ISBlockGetLocalSize(IS is,PetscInt *size)
602 {
603   PetscFunctionBegin;
604   PetscCall(PetscUseMethod(is,"ISBlockGetLocalSize_C",(IS,PetscInt*),(is,size)));
605   PetscFunctionReturn(0);
606 }
607 
608 static PetscErrorCode  ISBlockGetLocalSize_Block(IS is,PetscInt *size)
609 {
610   PetscInt       bs, n;
611 
612   PetscFunctionBegin;
613   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
614   PetscCall(PetscLayoutGetLocalSize(is->map, &n));
615   *size = n/bs;
616   PetscFunctionReturn(0);
617 }
618 
619 /*@
620    ISBlockGetSize - Returns the global number of blocks in the index set.
621 
622    Not Collective
623 
624    Input Parameter:
625 .  is - the index set
626 
627    Output Parameter:
628 .  size - the global number of blocks
629 
630    Level: intermediate
631 
632 .seealso: ISGetBlockSize(), ISBlockGetLocalSize(), ISGetSize(), ISCreateBlock(), ISBLOCK
633 @*/
634 PetscErrorCode  ISBlockGetSize(IS is,PetscInt *size)
635 {
636   PetscFunctionBegin;
637   PetscCall(PetscUseMethod(is,"ISBlockGetSize_C",(IS,PetscInt*),(is,size)));
638   PetscFunctionReturn(0);
639 }
640 
641 static PetscErrorCode  ISBlockGetSize_Block(IS is,PetscInt *size)
642 {
643   PetscInt       bs, N;
644 
645   PetscFunctionBegin;
646   PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
647   PetscCall(PetscLayoutGetSize(is->map, &N));
648   *size = N/bs;
649   PetscFunctionReturn(0);
650 }
651 
652 PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is)
653 {
654   IS_Block       *sub;
655 
656   PetscFunctionBegin;
657   PetscCall(PetscNewLog(is,&sub));
658   is->data = (void *) sub;
659   PetscCall(PetscMemcpy(is->ops,&myops,sizeof(myops)));
660   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",ISBlockSetIndices_Block));
661   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",ISBlockGetIndices_Block));
662   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",ISBlockRestoreIndices_Block));
663   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",ISBlockGetSize_Block));
664   PetscCall(PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",ISBlockGetLocalSize_Block));
665   PetscFunctionReturn(0);
666 }
667