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