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