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