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