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