1c4762a1bSJed Brown static char help[] = "Demonstrates constructing an application ordering.\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscsys.h>
4c4762a1bSJed Brown #include <petscao.h>
5c4762a1bSJed Brown #include <petscviewer.h>
6c4762a1bSJed Brown
main(int argc,char ** argv)7d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown PetscInt i, n = 5;
10c4762a1bSJed Brown PetscInt getpetsc[] = {0, 3, 4}, getapp[] = {2, 1, 9, 7};
11c4762a1bSJed Brown PetscInt getpetsc1[] = {0, 3, 4}, getapp1[] = {2, 1, 9, 7};
12c4762a1bSJed Brown PetscInt getpetsc2[] = {0, 3, 4}, getapp2[] = {2, 1, 9, 7};
13c4762a1bSJed Brown PetscInt getpetsc3[] = {0, 3, 4}, getapp3[] = {2, 1, 9, 7};
14c4762a1bSJed Brown PetscInt getpetsc4[] = {0, 3, 4}, getapp4[] = {2, 1, 9, 7};
15c4762a1bSJed Brown PetscMPIInt rank, size;
16c4762a1bSJed Brown IS ispetsc, isapp;
17c4762a1bSJed Brown AO ao;
18c4762a1bSJed Brown const PetscInt *app;
19c4762a1bSJed Brown
20327415f7SBarry Smith PetscFunctionBeginUser;
21c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
25c4762a1bSJed Brown
26c4762a1bSJed Brown /* create the index sets */
279566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, rank, size, &isapp));
289566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, n * rank, 1, &ispetsc)); /* natural numbering */
29c4762a1bSJed Brown
30c4762a1bSJed Brown /* create the application ordering */
319566063dSJacob Faibussowitsch PetscCall(AOCreateBasicIS(isapp, ispetsc, &ao));
329566063dSJacob Faibussowitsch PetscCall(AOView(ao, PETSC_VIEWER_STDOUT_WORLD));
33c4762a1bSJed Brown
349566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp));
359566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc));
36c4762a1bSJed Brown
379566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] 2,1,9,7 PetscToApplication %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, getapp[0], getapp[1], getapp[2], getapp[3]));
389566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] 0,3,4 ApplicationToPetsc %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, getpetsc[0], getpetsc[1], getpetsc[2]));
399566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
409566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao));
41c4762a1bSJed Brown
42c4762a1bSJed Brown /* test MemoryScalable ao */
43c4762a1bSJed Brown /*-------------------------*/
449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nTest AOCreateMemoryScalable: \n"));
459566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalableIS(isapp, ispetsc, &ao));
469566063dSJacob Faibussowitsch PetscCall(AOView(ao, PETSC_VIEWER_STDOUT_WORLD));
47c4762a1bSJed Brown
489566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp1));
499566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc1));
50c4762a1bSJed Brown
51*f5c5fea7SStefano Zampini /* Check accuracy */
52ad540459SPierre Jolivet for (i = 0; i < 4; i++) PetscCheck(getapp1[i] == getapp[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getapp1 %" PetscInt_FMT " != getapp %" PetscInt_FMT, getapp1[i], getapp[i]);
53ad540459SPierre Jolivet for (i = 0; i < 3; i++) PetscCheck(getpetsc1[i] == getpetsc[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getpetsc1 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT, getpetsc1[i], getpetsc[i]);
54c4762a1bSJed Brown
559566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao));
56c4762a1bSJed Brown
57c4762a1bSJed Brown /* test MemoryScalable ao: ispetsc = NULL */
58c4762a1bSJed Brown /*-----------------------------------------------*/
599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nTest AOCreateMemoryScalable with ispetsc=NULL:\n"));
609566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalableIS(isapp, NULL, &ao));
61c4762a1bSJed Brown
629566063dSJacob Faibussowitsch PetscCall(AOView(ao, PETSC_VIEWER_STDOUT_WORLD));
63c4762a1bSJed Brown
649566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp2));
659566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc2));
66c4762a1bSJed Brown
67*f5c5fea7SStefano Zampini /* Check accuracy */
68ad540459SPierre Jolivet for (i = 0; i < 4; i++) PetscCheck(getapp2[i] == getapp[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getapp2 %" PetscInt_FMT " != getapp %" PetscInt_FMT, getapp2[i], getapp[i]);
69ad540459SPierre Jolivet for (i = 0; i < 3; i++) PetscCheck(getpetsc2[i] == getpetsc[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getpetsc2 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT, getpetsc2[i], getpetsc[i]);
709566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao));
71c4762a1bSJed Brown
72c4762a1bSJed Brown /* test AOCreateMemoryScalable() ao: */
739566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isapp, &app));
749566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalable(PETSC_COMM_WORLD, n, app, NULL, &ao));
759566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isapp, &app));
76c4762a1bSJed Brown
779566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp4));
789566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc4));
79c4762a1bSJed Brown
80*f5c5fea7SStefano Zampini /* Check accuracy */
81ad540459SPierre Jolivet for (i = 0; i < 4; i++) PetscCheck(getapp4[i] == getapp[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getapp4 %" PetscInt_FMT " != getapp %" PetscInt_FMT, getapp4[i], getapp[i]);
82ad540459SPierre Jolivet for (i = 0; i < 3; i++) PetscCheck(getpetsc4[i] == getpetsc[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getpetsc4 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT, getpetsc4[i], getpetsc[i]);
839566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao));
84c4762a1bSJed Brown
85c4762a1bSJed Brown /* test general API */
86c4762a1bSJed Brown /*------------------*/
879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nTest general API: \n"));
889566063dSJacob Faibussowitsch PetscCall(AOCreate(PETSC_COMM_WORLD, &ao));
899566063dSJacob Faibussowitsch PetscCall(AOSetIS(ao, isapp, ispetsc));
909566063dSJacob Faibussowitsch PetscCall(AOSetType(ao, AOMEMORYSCALABLE));
919566063dSJacob Faibussowitsch PetscCall(AOSetFromOptions(ao));
92c4762a1bSJed Brown
93c4762a1bSJed Brown /* ispetsc and isapp are nolonger used. */
949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispetsc));
959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isapp));
96c4762a1bSJed Brown
979566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp3));
989566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc3));
99c4762a1bSJed Brown
1009566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] 2,1,9,7 PetscToApplication %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, getapp3[0], getapp3[1], getapp3[2], getapp3[3]));
1019566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] 0,3,4 ApplicationToPetsc %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, getpetsc3[0], getpetsc3[1], getpetsc3[2]));
1029566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
103c4762a1bSJed Brown
104*f5c5fea7SStefano Zampini /* Check accuracy */
105ad540459SPierre Jolivet for (i = 0; i < 4; i++) PetscCheck(getapp3[i] == getapp[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getapp3 %" PetscInt_FMT " != getapp %" PetscInt_FMT, getapp3[i], getapp[i]);
106ad540459SPierre Jolivet for (i = 0; i < 3; i++) PetscCheck(getpetsc3[i] == getpetsc[i], PETSC_COMM_SELF, PETSC_ERR_USER, "getpetsc3 %" PetscInt_FMT " != getpetsc %" PetscInt_FMT, getpetsc3[i], getpetsc[i]);
107c4762a1bSJed Brown
1089566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao));
1099566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
110b122ec5aSJacob Faibussowitsch return 0;
111c4762a1bSJed Brown }
112c4762a1bSJed Brown
113c4762a1bSJed Brown /*TEST
114c4762a1bSJed Brown
115c4762a1bSJed Brown test:
116c4762a1bSJed Brown
117c4762a1bSJed Brown test:
118c4762a1bSJed Brown suffix: 2
119c4762a1bSJed Brown nsize: 2
120c4762a1bSJed Brown
121c4762a1bSJed Brown test:
122c4762a1bSJed Brown suffix: 3
123c4762a1bSJed Brown nsize: 3
124c4762a1bSJed Brown
125c4762a1bSJed Brown test:
126c4762a1bSJed Brown suffix: 4
127c4762a1bSJed Brown nsize: 3
128c4762a1bSJed Brown args: -ao_type basic
129c4762a1bSJed Brown output_file: output/ex1_3.out
130c4762a1bSJed Brown
131c4762a1bSJed Brown TEST*/
132