147c6ae99SBarry Smith /*
247c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel.
347c6ae99SBarry Smith */
447c6ae99SBarry Smith
5af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/
647c6ae99SBarry Smith
747c6ae99SBarry Smith /*
847c6ae99SBarry Smith Gets the natural number for each global number on the process.
947c6ae99SBarry Smith
10aa219208SBarry Smith Used by DMDAGetAO() and DMDAGlobalToNatural_Create()
1147c6ae99SBarry Smith */
DMDAGetNatural_Private(DM da,PetscInt * outNlocal,IS * isnatural)12d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNatural_Private(DM da, PetscInt *outNlocal, IS *isnatural)
13d71ae5a4SJacob Faibussowitsch {
14c73cfb54SMatthew G. Knepley PetscInt Nlocal, i, j, k, *lidx, lict = 0, dim = da->dim;
1547c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data;
1647c6ae99SBarry Smith
1747c6ae99SBarry Smith PetscFunctionBegin;
1847c6ae99SBarry Smith Nlocal = (dd->xe - dd->xs);
19c73cfb54SMatthew G. Knepley if (dim > 1) Nlocal *= (dd->ye - dd->ys);
20c73cfb54SMatthew G. Knepley if (dim > 2) Nlocal *= (dd->ze - dd->zs);
2147c6ae99SBarry Smith
229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nlocal, &lidx));
2347c6ae99SBarry Smith
24c73cfb54SMatthew G. Knepley if (dim == 1) {
2547c6ae99SBarry Smith for (i = dd->xs; i < dd->xe; i++) {
2647c6ae99SBarry Smith /* global number in natural ordering */
2747c6ae99SBarry Smith lidx[lict++] = i;
2847c6ae99SBarry Smith }
29c73cfb54SMatthew G. Knepley } else if (dim == 2) {
3047c6ae99SBarry Smith for (j = dd->ys; j < dd->ye; j++) {
3147c6ae99SBarry Smith for (i = dd->xs; i < dd->xe; i++) {
3247c6ae99SBarry Smith /* global number in natural ordering */
3347c6ae99SBarry Smith lidx[lict++] = i + j * dd->M * dd->w;
3447c6ae99SBarry Smith }
3547c6ae99SBarry Smith }
36c73cfb54SMatthew G. Knepley } else if (dim == 3) {
3747c6ae99SBarry Smith for (k = dd->zs; k < dd->ze; k++) {
3847c6ae99SBarry Smith for (j = dd->ys; j < dd->ye; j++) {
39ad540459SPierre Jolivet for (i = dd->xs; i < dd->xe; i++) lidx[lict++] = i + j * dd->M * dd->w + k * dd->M * dd->N * dd->w;
4047c6ae99SBarry Smith }
4147c6ae99SBarry Smith }
4247c6ae99SBarry Smith }
4347c6ae99SBarry Smith *outNlocal = Nlocal;
449566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)da), Nlocal, lidx, PETSC_OWN_POINTER, isnatural));
453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4647c6ae99SBarry Smith }
4747c6ae99SBarry Smith
48*5d83a8b1SBarry Smith /*@
4912b4a537SBarry Smith DMDASetAOType - Sets the type of application ordering to create with `DMDAGetAO()`, for a distributed array.
509db3d8bcSStefano Zampini
5120f4b53cSBarry Smith Collective
529db3d8bcSStefano Zampini
53d8d19677SJose E. Roman Input Parameters:
5472fd0fbdSBarry Smith + da - the `DMDA`
5512b4a537SBarry Smith - aotype - type of `AO`. `AOType` which can be `AOBASIC`, `AOADVANCED`, `AOMAPPING`, or `AOMEMORYSCALABLE`
569db3d8bcSStefano Zampini
579db3d8bcSStefano Zampini Level: intermediate
589db3d8bcSStefano Zampini
59dce8aebaSBarry Smith Note:
6012b4a537SBarry Smith It will generate an error if an `AO` has already been obtained with a call to `DMDAGetAO()` and the user sets a different `AOType`
619db3d8bcSStefano Zampini
6212b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate2d()`, `DMDAGetAO()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMLocalToGlobal()`
63db781477SPatrick Sanan `DMGlobalToLocalBegin()`, `DMGlobalToLocalEnd()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDAGetGlobalIndices()`, `DMDAGetOwnershipRanges()`,
6412b4a537SBarry Smith `AO`, `AOPetscToApplication()`, `AOApplicationToPetsc()`, `AOType`, `AOBASIC`, `AOADVANCED`, `AOMAPPING`, `AOMEMORYSCALABLE`
659db3d8bcSStefano Zampini @*/
DMDASetAOType(DM da,AOType aotype)66d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetAOType(DM da, AOType aotype)
67d71ae5a4SJacob Faibussowitsch {
689db3d8bcSStefano Zampini DM_DA *dd;
699db3d8bcSStefano Zampini PetscBool isdmda;
709db3d8bcSStefano Zampini
719db3d8bcSStefano Zampini PetscFunctionBegin;
72a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isdmda));
747a8be351SBarry Smith PetscCheck(isdmda, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Requires a DMDA as input");
759db3d8bcSStefano Zampini /* now we can safely dereference */
769db3d8bcSStefano Zampini dd = (DM_DA *)da->data;
779db3d8bcSStefano Zampini if (dd->ao) { /* check if the already computed AO has the same type as requested */
789db3d8bcSStefano Zampini PetscBool match;
799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dd->ao, aotype, &match));
807a8be351SBarry Smith PetscCheck(match, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot change AO type");
813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
829db3d8bcSStefano Zampini }
839566063dSJacob Faibussowitsch PetscCall(PetscFree(dd->aotype));
849566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(aotype, (char **)&dd->aotype));
853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
869db3d8bcSStefano Zampini }
879db3d8bcSStefano Zampini
8847c6ae99SBarry Smith /*@
89aa219208SBarry Smith DMDAGetAO - Gets the application ordering context for a distributed array.
9047c6ae99SBarry Smith
9120f4b53cSBarry Smith Collective
9247c6ae99SBarry Smith
9347c6ae99SBarry Smith Input Parameter:
9472fd0fbdSBarry Smith . da - the `DMDA`
9547c6ae99SBarry Smith
962fe279fdSBarry Smith Output Parameter:
97dce8aebaSBarry Smith . ao - the application ordering context for `DMDA`
9847c6ae99SBarry Smith
9947c6ae99SBarry Smith Level: intermediate
10047c6ae99SBarry Smith
10147c6ae99SBarry Smith Notes:
102dce8aebaSBarry Smith In this case, the `AO` maps to the natural grid ordering that would be used
103dce8aebaSBarry Smith for the `DMDA` if only 1 processor were employed (ordering most rapidly in the
10447c6ae99SBarry Smith x-direction, then y, then z). Multiple degrees of freedom are numbered
10547c6ae99SBarry Smith for each node (rather than 1 component for the whole grid, then the next
10647c6ae99SBarry Smith component, etc.)
10747c6ae99SBarry Smith
10812b4a537SBarry Smith Do NOT call `AODestroy()` on the `ao` returned by this function.
1097a612ce8SBarry Smith
11012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate2d()`, `DMDASetAOType()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMLocalToGlobal()`
111db781477SPatrick Sanan `DMGlobalToLocalBegin()`, `DMGlobalToLocalEnd()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDAGetOwnershipRanges()`,
112db781477SPatrick Sanan `AO`, `AOPetscToApplication()`, `AOApplicationToPetsc()`
11347c6ae99SBarry Smith @*/
DMDAGetAO(DM da,AO * ao)114d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetAO(DM da, AO *ao)
115d71ae5a4SJacob Faibussowitsch {
11663868d81SStefano Zampini DM_DA *dd;
11763868d81SStefano Zampini PetscBool isdmda;
11847c6ae99SBarry Smith
11947c6ae99SBarry Smith PetscFunctionBegin;
120a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1214f572ea9SToby Isaac PetscAssertPointer(ao, 2);
1229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isdmda));
1237a8be351SBarry Smith PetscCheck(isdmda, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Requires a DMDA as input");
12463868d81SStefano Zampini /* now we can safely dereference */
12563868d81SStefano Zampini dd = (DM_DA *)da->data;
12647c6ae99SBarry Smith
12747c6ae99SBarry Smith /*
12847c6ae99SBarry Smith Build the natural ordering to PETSc ordering mappings.
12947c6ae99SBarry Smith */
13047c6ae99SBarry Smith if (!dd->ao) {
13147c6ae99SBarry Smith IS ispetsc, isnatural;
13247c6ae99SBarry Smith PetscInt Nlocal;
13347c6ae99SBarry Smith
1349566063dSJacob Faibussowitsch PetscCall(DMDAGetNatural_Private(da, &Nlocal, &isnatural));
1359566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)da), Nlocal, dd->base, 1, &ispetsc));
1369566063dSJacob Faibussowitsch PetscCall(AOCreate(PetscObjectComm((PetscObject)da), &dd->ao));
1379566063dSJacob Faibussowitsch PetscCall(AOSetIS(dd->ao, isnatural, ispetsc));
1389566063dSJacob Faibussowitsch PetscCall(AOSetType(dd->ao, dd->aotype));
1399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispetsc));
1409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnatural));
14147c6ae99SBarry Smith }
14247c6ae99SBarry Smith *ao = dd->ao;
1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14447c6ae99SBarry Smith }
145