Actual source code: pmap.c
1: /*
2: This file contains routines for basic map object implementation.
3: */
5: #include <petsc/private/isimpl.h>
7: /*@
8: PetscLayoutCreate - Allocates `PetscLayout` object
10: Collective
12: Input Parameter:
13: . comm - the MPI communicator
15: Output Parameter:
16: . map - the new `PetscLayout`
18: Level: advanced
20: Notes:
21: Typical calling sequence
22: .vb
23: PetscLayoutCreate(MPI_Comm,PetscLayout *);
24: PetscLayoutSetBlockSize(PetscLayout,bs);
25: PetscLayoutSetSize(PetscLayout,N); // or PetscLayoutSetLocalSize(PetscLayout,n);
26: PetscLayoutSetUp(PetscLayout);
27: .ve
28: Alternatively,
29: .vb
30: PetscLayoutCreateFromSizes(comm,n,N,bs,&layout);
31: .ve
33: Optionally use any of the following
34: .vb
35: PetscLayoutGetSize(PetscLayout,PetscInt *);
36: PetscLayoutGetLocalSize(PetscLayout,PetscInt *);
37: PetscLayoutGetRange(PetscLayout,PetscInt *rstart,PetscInt *rend);
38: PetscLayoutGetRanges(PetscLayout,const PetscInt *range[]);
39: PetscLayoutDestroy(PetscLayout*);
40: .ve
42: The `PetscLayout` object and methods are intended to be used in the PETSc `Vec` and `Mat` implementations; it is often not needed in
43: user codes unless you really gain something in their use.
45: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutSetLocalSize()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutGetLocalSize()`,
46: `PetscLayout`, `PetscLayoutDestroy()`,
47: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`, `PetscLayoutSetUp()`,
48: `PetscLayoutCreateFromSizes()`
49: @*/
50: PetscErrorCode PetscLayoutCreate(MPI_Comm comm, PetscLayout *map)
51: {
52: PetscFunctionBegin;
53: PetscCall(PetscNew(map));
54: PetscCallMPI(MPI_Comm_size(comm, &(*map)->size));
55: (*map)->comm = comm;
56: (*map)->bs = -1;
57: (*map)->n = -1;
58: (*map)->N = -1;
59: (*map)->range = NULL;
60: (*map)->range_alloc = PETSC_TRUE;
61: (*map)->rstart = 0;
62: (*map)->rend = 0;
63: (*map)->setupcalled = PETSC_FALSE;
64: (*map)->oldn = -1;
65: (*map)->oldN = -1;
66: (*map)->oldbs = -1;
67: PetscFunctionReturn(PETSC_SUCCESS);
68: }
70: /*@
71: PetscLayoutCreateFromSizes - Allocates `PetscLayout` object and sets the layout sizes, and sets the layout up.
73: Collective
75: Input Parameters:
76: + comm - the MPI communicator
77: . n - the local size (or `PETSC_DECIDE`)
78: . N - the global size (or `PETSC_DECIDE`)
79: - bs - the block size (or `PETSC_DECIDE`)
81: Output Parameter:
82: . map - the new `PetscLayout`
84: Level: advanced
86: Note:
87: $ PetscLayoutCreateFromSizes(comm, n, N, bs, &layout);
88: is a shorthand for
89: .vb
90: PetscLayoutCreate(comm, &layout);
91: PetscLayoutSetLocalSize(layout, n);
92: PetscLayoutSetSize(layout, N);
93: PetscLayoutSetBlockSize(layout, bs);
94: PetscLayoutSetUp(layout);
95: .ve
97: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutGetLocalSize()`, `PetscLayout`, `PetscLayoutDestroy()`,
98: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`, `PetscLayoutSetUp()`, `PetscLayoutCreateFromRanges()`
99: @*/
100: PetscErrorCode PetscLayoutCreateFromSizes(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt bs, PetscLayout *map)
101: {
102: PetscFunctionBegin;
103: PetscCall(PetscLayoutCreate(comm, map));
104: PetscCall(PetscLayoutSetLocalSize(*map, n));
105: PetscCall(PetscLayoutSetSize(*map, N));
106: PetscCall(PetscLayoutSetBlockSize(*map, bs));
107: PetscCall(PetscLayoutSetUp(*map));
108: PetscFunctionReturn(PETSC_SUCCESS);
109: }
111: /*@
112: PetscLayoutDestroy - Frees a `PetscLayout` object and frees its range if that exists.
114: Collective
116: Input Parameter:
117: . map - the `PetscLayout`
119: Level: developer
121: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutSetLocalSize()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutGetLocalSize()`,
122: `PetscLayout`, `PetscLayoutCreate()`,
123: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`, `PetscLayoutSetUp()`
124: @*/
125: PetscErrorCode PetscLayoutDestroy(PetscLayout *map)
126: {
127: PetscFunctionBegin;
128: if (!*map) PetscFunctionReturn(PETSC_SUCCESS);
129: if (!(*map)->refcnt--) {
130: if ((*map)->range_alloc) PetscCall(PetscFree((*map)->range));
131: PetscCall(ISLocalToGlobalMappingDestroy(&(*map)->mapping));
132: PetscCall(PetscFree(*map));
133: }
134: *map = NULL;
135: PetscFunctionReturn(PETSC_SUCCESS);
136: }
138: /*@
139: PetscLayoutCreateFromRanges - Creates a new `PetscLayout` with the given ownership ranges and sets it up.
141: Collective
143: Input Parameters:
144: + comm - the MPI communicator
145: . range - the array of ownership ranges for each rank with length commsize+1
146: . mode - the copy mode for range
147: - bs - the block size (or `PETSC_DECIDE`)
149: Output Parameter:
150: . newmap - the new `PetscLayout`
152: Level: developer
154: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`,
155: `PetscLayoutGetLocalSize()`, `PetscLayout`, `PetscLayoutDestroy()`,
156: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`, `PetscLayoutSetUp()`, `PetscLayoutCreateFromSizes()`
157: @*/
158: PetscErrorCode PetscLayoutCreateFromRanges(MPI_Comm comm, const PetscInt range[], PetscCopyMode mode, PetscInt bs, PetscLayout *newmap)
159: {
160: PetscLayout map;
161: PetscMPIInt rank;
163: PetscFunctionBegin;
164: PetscCallMPI(MPI_Comm_rank(comm, &rank));
165: PetscCall(PetscLayoutCreate(comm, &map));
166: PetscCall(PetscLayoutSetBlockSize(map, bs));
167: switch (mode) {
168: case PETSC_COPY_VALUES:
169: PetscCall(PetscMalloc1(map->size + 1, &map->range));
170: PetscCall(PetscArraycpy(map->range, range, map->size + 1));
171: break;
172: case PETSC_USE_POINTER:
173: map->range_alloc = PETSC_FALSE;
174: break;
175: default:
176: map->range = (PetscInt *)range;
177: break;
178: }
179: map->rstart = map->range[rank];
180: map->rend = map->range[rank + 1];
181: map->n = map->rend - map->rstart;
182: map->N = map->range[map->size];
183: if (PetscDefined(USE_DEBUG)) { /* just check that n, N and bs are consistent */
184: PetscInt tmp;
185: PetscCall(MPIU_Allreduce(&map->n, &tmp, 1, MPIU_INT, MPI_SUM, map->comm));
186: PetscCheck(tmp == map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Sum of local lengths %" PetscInt_FMT " does not equal global length %" PetscInt_FMT ", my local length %" PetscInt_FMT ". The provided PetscLayout is wrong.", tmp, map->N, map->n);
187: if (map->bs > 1) PetscCheck(map->n % map->bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local size %" PetscInt_FMT " must be divisible by blocksize %" PetscInt_FMT, map->n, map->bs);
188: if (map->bs > 1) PetscCheck(map->N % map->bs == 0, map->comm, PETSC_ERR_PLIB, "Global size %" PetscInt_FMT " must be divisible by blocksize %" PetscInt_FMT, map->N, map->bs);
189: }
190: /* lock the layout */
191: map->setupcalled = PETSC_TRUE;
192: map->oldn = map->n;
193: map->oldN = map->N;
194: map->oldbs = map->bs;
195: *newmap = map;
196: PetscFunctionReturn(PETSC_SUCCESS);
197: }
199: /*@
200: PetscLayoutSetUp - given a map where you have set either the global or local
201: size sets up the map so that it may be used.
203: Collective
205: Input Parameter:
206: . map - pointer to the map
208: Level: developer
210: Notes:
211: Typical calling sequence
212: .vb
213: PetscLayoutCreate(MPI_Comm,PetscLayout *);
214: PetscLayoutSetBlockSize(PetscLayout,1);
215: PetscLayoutSetSize(PetscLayout,n) or PetscLayoutSetLocalSize(PetscLayout,N); or both
216: PetscLayoutSetUp(PetscLayout);
217: PetscLayoutGetSize(PetscLayout,PetscInt *);
218: .ve
220: If range exists, and local size is not set, everything gets computed from the range.
222: If the local size, global size are already set and range exists then this does nothing.
224: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutSetLocalSize()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutGetLocalSize()`,
225: `PetscLayout`, `PetscLayoutDestroy()`,
226: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`, `PetscLayoutCreate()`, `PetscSplitOwnership()`
227: @*/
228: PetscErrorCode PetscLayoutSetUp(PetscLayout map)
229: {
230: PetscMPIInt rank;
231: PetscInt p;
233: PetscFunctionBegin;
234: PetscCheck(!map->setupcalled || !(map->n != map->oldn || map->N != map->oldN), map->comm, PETSC_ERR_ARG_WRONGSTATE, "Layout is already setup with (local=%" PetscInt_FMT ",global=%" PetscInt_FMT "), cannot call setup again with (local=%" PetscInt_FMT ",global=%" PetscInt_FMT ")",
235: map->oldn, map->oldN, map->n, map->N);
236: if (map->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
238: if (map->n > 0 && map->bs > 1) PetscCheck(map->n % map->bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local size %" PetscInt_FMT " must be divisible by blocksize %" PetscInt_FMT, map->n, map->bs);
239: if (map->N > 0 && map->bs > 1) PetscCheck(map->N % map->bs == 0, map->comm, PETSC_ERR_PLIB, "Global size %" PetscInt_FMT " must be divisible by blocksize %" PetscInt_FMT, map->N, map->bs);
241: PetscCallMPI(MPI_Comm_rank(map->comm, &rank));
242: if (map->n > 0) map->n = map->n / PetscAbs(map->bs);
243: if (map->N > 0) map->N = map->N / PetscAbs(map->bs);
244: PetscCall(PetscSplitOwnership(map->comm, &map->n, &map->N));
245: map->n = map->n * PetscAbs(map->bs);
246: map->N = map->N * PetscAbs(map->bs);
247: if (!map->range) PetscCall(PetscMalloc1(map->size + 1, &map->range));
248: PetscCallMPI(MPI_Allgather(&map->n, 1, MPIU_INT, map->range + 1, 1, MPIU_INT, map->comm));
250: map->range[0] = 0;
251: for (p = 2; p <= map->size; p++) map->range[p] += map->range[p - 1];
253: map->rstart = map->range[rank];
254: map->rend = map->range[rank + 1];
256: /* lock the layout */
257: map->setupcalled = PETSC_TRUE;
258: map->oldn = map->n;
259: map->oldN = map->N;
260: map->oldbs = map->bs;
261: PetscFunctionReturn(PETSC_SUCCESS);
262: }
264: /*@
265: PetscLayoutDuplicate - creates a new `PetscLayout` with the same information as a given one. If the `PetscLayout` already exists it is destroyed first.
267: Collective
269: Input Parameter:
270: . in - input `PetscLayout` to be duplicated
272: Output Parameter:
273: . out - the copy
275: Level: developer
277: Note:
278: `PetscLayoutSetUp()` does not need to be called on the resulting `PetscLayout`
280: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutDestroy()`, `PetscLayoutSetUp()`, `PetscLayoutReference()`
281: @*/
282: PetscErrorCode PetscLayoutDuplicate(PetscLayout in, PetscLayout *out)
283: {
284: MPI_Comm comm = in->comm;
286: PetscFunctionBegin;
287: PetscCall(PetscLayoutDestroy(out));
288: PetscCall(PetscLayoutCreate(comm, out));
289: PetscCall(PetscMemcpy(*out, in, sizeof(struct _n_PetscLayout)));
290: if (in->range) {
291: PetscCall(PetscMalloc1((*out)->size + 1, &(*out)->range));
292: PetscCall(PetscArraycpy((*out)->range, in->range, (*out)->size + 1));
293: }
294: (*out)->refcnt = 0;
295: PetscFunctionReturn(PETSC_SUCCESS);
296: }
298: /*@
299: PetscLayoutReference - Causes a PETSc `Vec` or `Mat` to share a `PetscLayout` with one that already exists.
301: Collective
303: Input Parameter:
304: . in - input `PetscLayout` to be copied
306: Output Parameter:
307: . out - the reference location
309: Level: developer
311: Notes:
312: `PetscLayoutSetUp()` does not need to be called on the resulting `PetscLayout`
314: If the out location already contains a `PetscLayout` it is destroyed
316: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutDestroy()`, `PetscLayoutSetUp()`, `PetscLayoutDuplicate()`
317: @*/
318: PetscErrorCode PetscLayoutReference(PetscLayout in, PetscLayout *out)
319: {
320: PetscFunctionBegin;
321: in->refcnt++;
322: PetscCall(PetscLayoutDestroy(out));
323: *out = in;
324: PetscFunctionReturn(PETSC_SUCCESS);
325: }
327: /*@
328: PetscLayoutSetISLocalToGlobalMapping - sets a `ISLocalGlobalMapping` into a `PetscLayout`
330: Collective
332: Input Parameters:
333: + in - input `PetscLayout`
334: - ltog - the local to global mapping
336: Level: developer
338: Notes:
339: `PetscLayoutSetUp()` does not need to be called on the resulting `PetscLayout`
341: If the `PetscLayout` already contains a `ISLocalGlobalMapping` it is destroyed
343: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutDestroy()`, `PetscLayoutSetUp()`, `PetscLayoutDuplicate()`
344: @*/
345: PetscErrorCode PetscLayoutSetISLocalToGlobalMapping(PetscLayout in, ISLocalToGlobalMapping ltog)
346: {
347: PetscFunctionBegin;
348: if (ltog) {
349: PetscInt bs;
351: PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog, &bs));
352: PetscCheck(in->bs <= 0 || bs == 1 || in->bs == bs, in->comm, PETSC_ERR_PLIB, "Blocksize of layout %" PetscInt_FMT " must match that of mapping %" PetscInt_FMT " (or the latter must be 1)", in->bs, bs);
353: }
354: PetscCall(PetscObjectReference((PetscObject)ltog));
355: PetscCall(ISLocalToGlobalMappingDestroy(&in->mapping));
356: in->mapping = ltog;
357: PetscFunctionReturn(PETSC_SUCCESS);
358: }
360: /*@
361: PetscLayoutSetLocalSize - Sets the local size for a `PetscLayout` object.
363: Collective
365: Input Parameters:
366: + map - pointer to the map
367: - n - the local size
369: Level: developer
371: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutSetUp()`
372: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`
373: @*/
374: PetscErrorCode PetscLayoutSetLocalSize(PetscLayout map, PetscInt n)
375: {
376: PetscFunctionBegin;
377: PetscCheck(map->bs <= 1 || (n % map->bs) == 0, map->comm, PETSC_ERR_ARG_INCOMP, "Local size %" PetscInt_FMT " not compatible with block size %" PetscInt_FMT, n, map->bs);
378: map->n = n;
379: PetscFunctionReturn(PETSC_SUCCESS);
380: }
382: /*@
383: PetscLayoutGetLocalSize - Gets the local size for a `PetscLayout` object.
385: Not Collective
387: Input Parameter:
388: . map - pointer to the map
390: Output Parameter:
391: . n - the local size
393: Level: developer
395: Note:
396: Call this after the call to `PetscLayoutSetUp()`
398: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutSetUp()`
399: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`
400: @*/
401: PetscErrorCode PetscLayoutGetLocalSize(PetscLayout map, PetscInt *n)
402: {
403: PetscFunctionBegin;
404: *n = map->n;
405: PetscFunctionReturn(PETSC_SUCCESS);
406: }
408: /*@
409: PetscLayoutSetSize - Sets the global size for a `PetscLayout` object.
411: Logically Collective
413: Input Parameters:
414: + map - pointer to the map
415: - n - the global size
417: Level: developer
419: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutGetSize()`, `PetscLayoutSetUp()`
420: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`
421: @*/
422: PetscErrorCode PetscLayoutSetSize(PetscLayout map, PetscInt n)
423: {
424: PetscFunctionBegin;
425: map->N = n;
426: PetscFunctionReturn(PETSC_SUCCESS);
427: }
429: /*@
430: PetscLayoutGetSize - Gets the global size for a `PetscLayout` object.
432: Not Collective
434: Input Parameter:
435: . map - pointer to the map
437: Output Parameter:
438: . n - the global size
440: Level: developer
442: Note:
443: Call this after the call to `PetscLayoutSetUp()`
445: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutSetSize()`, `PetscLayoutSetUp()`
446: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetBlockSize()`
447: @*/
448: PetscErrorCode PetscLayoutGetSize(PetscLayout map, PetscInt *n)
449: {
450: PetscFunctionBegin;
451: *n = map->N;
452: PetscFunctionReturn(PETSC_SUCCESS);
453: }
455: /*@
456: PetscLayoutSetBlockSize - Sets the block size for a `PetscLayout` object.
458: Logically Collective
460: Input Parameters:
461: + map - pointer to the map
462: - bs - the size
464: Level: developer
466: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutGetBlockSize()`,
467: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutSetUp()`
468: @*/
469: PetscErrorCode PetscLayoutSetBlockSize(PetscLayout map, PetscInt bs)
470: {
471: PetscFunctionBegin;
472: if (bs < 0) PetscFunctionReturn(PETSC_SUCCESS);
473: PetscCheck(map->n <= 0 || (map->n % bs) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local size %" PetscInt_FMT " not compatible with block size %" PetscInt_FMT, map->n, bs);
474: if (map->mapping) {
475: PetscInt obs;
477: PetscCall(ISLocalToGlobalMappingGetBlockSize(map->mapping, &obs));
478: if (obs > 1) PetscCall(ISLocalToGlobalMappingSetBlockSize(map->mapping, bs));
479: }
480: map->bs = bs;
481: PetscFunctionReturn(PETSC_SUCCESS);
482: }
484: /*@
485: PetscLayoutGetBlockSize - Gets the block size for a `PetscLayout` object.
487: Not Collective
489: Input Parameter:
490: . map - pointer to the map
492: Output Parameter:
493: . bs - the size
495: Level: developer
497: Notes:
498: Call this after the call to `PetscLayoutSetUp()`
500: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutSetSize()`, `PetscLayoutSetUp()`
501: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutGetSize()`
502: @*/
503: PetscErrorCode PetscLayoutGetBlockSize(PetscLayout map, PetscInt *bs)
504: {
505: PetscFunctionBegin;
506: *bs = PetscAbs(map->bs);
507: PetscFunctionReturn(PETSC_SUCCESS);
508: }
510: /*@
511: PetscLayoutGetRange - gets the range of values owned by this process
513: Not Collective
515: Input Parameter:
516: . map - pointer to the map
518: Output Parameters:
519: + rstart - first index owned by this process
520: - rend - one more than the last index owned by this process
522: Level: developer
524: Note:
525: Call this after the call to `PetscLayoutSetUp()`
527: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutSetSize()`,
528: `PetscLayoutGetSize()`, `PetscLayoutGetRanges()`, `PetscLayoutSetBlockSize()`, `PetscLayoutSetUp()`
529: @*/
530: PetscErrorCode PetscLayoutGetRange(PetscLayout map, PetscInt *rstart, PetscInt *rend)
531: {
532: PetscFunctionBegin;
533: if (rstart) *rstart = map->rstart;
534: if (rend) *rend = map->rend;
535: PetscFunctionReturn(PETSC_SUCCESS);
536: }
538: /*@C
539: PetscLayoutGetRanges - gets the ranges of values owned by all processes
541: Not Collective
543: Input Parameter:
544: . map - pointer to the map
546: Output Parameter:
547: . range - start of each processors range of indices (the final entry is one more than the
548: last index on the last process). The length of the array is one more than the number of processes in the MPI
549: communicator owned by `map`
551: Level: developer
553: Note:
554: Call this after the call to `PetscLayoutSetUp()`
556: Fortran Notes:
557: In Fortran, use PetscLayoutGetRangesF90()
559: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutSetSize()`,
560: `PetscLayoutGetSize()`, `PetscLayoutGetRange()`, `PetscLayoutSetBlockSize()`, `PetscLayoutSetUp()`
561: @*/
562: PetscErrorCode PetscLayoutGetRanges(PetscLayout map, const PetscInt *range[])
563: {
564: PetscFunctionBegin;
565: *range = map->range;
566: PetscFunctionReturn(PETSC_SUCCESS);
567: }
569: /*@
570: PetscLayoutCompare - Compares two layouts
572: Not Collective
574: Input Parameters:
575: + mapa - pointer to the first map
576: - mapb - pointer to the second map
578: Output Parameter:
579: . congruent - `PETSC_TRUE` if the two layouts are congruent, `PETSC_FALSE` otherwise
581: Level: beginner
583: .seealso: [PetscLayout](sec_matlayout), `PetscLayoutCreate()`, `PetscLayoutSetLocalSize()`, `PetscLayoutGetLocalSize()`, `PetscLayoutGetBlockSize()`,
584: `PetscLayoutGetRange()`, `PetscLayoutGetRanges()`, `PetscLayoutSetSize()`, `PetscLayoutGetSize()`, `PetscLayoutSetUp()`
585: @*/
586: PetscErrorCode PetscLayoutCompare(PetscLayout mapa, PetscLayout mapb, PetscBool *congruent)
587: {
588: PetscFunctionBegin;
589: *congruent = PETSC_FALSE;
590: if (mapa->N == mapb->N && mapa->range && mapb->range && mapa->size == mapb->size) PetscCall(PetscArraycmp(mapa->range, mapb->range, mapa->size + 1, congruent));
591: PetscFunctionReturn(PETSC_SUCCESS);
592: }