1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Demonstrates constructing an application ordering.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscsys.h> 5c4762a1bSJed Brown #include <petscao.h> 6c4762a1bSJed Brown #include <petscviewer.h> 7c4762a1bSJed Brown 8*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 9*d71ae5a4SJacob Faibussowitsch { 10c4762a1bSJed Brown PetscInt i, n = 5; 11c4762a1bSJed Brown PetscInt getpetsc[] = {0, 3, 4}, getapp[] = {2, 1, 9, 7}; 12c4762a1bSJed Brown PetscInt getpetsc1[] = {0, 3, 4}, getapp1[] = {2, 1, 9, 7}; 13c4762a1bSJed Brown PetscInt getpetsc2[] = {0, 3, 4}, getapp2[] = {2, 1, 9, 7}; 14c4762a1bSJed Brown PetscInt getpetsc3[] = {0, 3, 4}, getapp3[] = {2, 1, 9, 7}; 15c4762a1bSJed Brown PetscInt getpetsc4[] = {0, 3, 4}, getapp4[] = {2, 1, 9, 7}; 16c4762a1bSJed Brown PetscMPIInt rank, size; 17c4762a1bSJed Brown IS ispetsc, isapp; 18c4762a1bSJed Brown AO ao; 19c4762a1bSJed Brown const PetscInt *app; 20c4762a1bSJed Brown 21327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* create the index sets */ 289566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, rank, size, &isapp)); 299566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, n * rank, 1, &ispetsc)); /* natural numbering */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* create the application ordering */ 329566063dSJacob Faibussowitsch PetscCall(AOCreateBasicIS(isapp, ispetsc, &ao)); 339566063dSJacob Faibussowitsch PetscCall(AOView(ao, PETSC_VIEWER_STDOUT_WORLD)); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp)); 369566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc)); 37c4762a1bSJed Brown 389566063dSJacob 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])); 399566063dSJacob 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])); 409566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 419566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* test MemoryScalable ao */ 44c4762a1bSJed Brown /*-------------------------*/ 459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nTest AOCreateMemoryScalable: \n")); 469566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalableIS(isapp, ispetsc, &ao)); 479566063dSJacob Faibussowitsch PetscCall(AOView(ao, PETSC_VIEWER_STDOUT_WORLD)); 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp1)); 509566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc1)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Check accuracy */; 53ad540459SPierre 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]); 54ad540459SPierre 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]); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* test MemoryScalable ao: ispetsc = NULL */ 59c4762a1bSJed Brown /*-----------------------------------------------*/ 609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nTest AOCreateMemoryScalable with ispetsc=NULL:\n")); 619566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalableIS(isapp, NULL, &ao)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(AOView(ao, PETSC_VIEWER_STDOUT_WORLD)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp2)); 669566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc2)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Check accuracy */; 69ad540459SPierre 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]); 70ad540459SPierre 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]); 719566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* test AOCreateMemoryScalable() ao: */ 749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isapp, &app)); 759566063dSJacob Faibussowitsch PetscCall(AOCreateMemoryScalable(PETSC_COMM_WORLD, n, app, NULL, &ao)); 769566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isapp, &app)); 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp4)); 799566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc4)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Check accuracy */; 82ad540459SPierre 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]); 83ad540459SPierre 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]); 849566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* test general API */ 87c4762a1bSJed Brown /*------------------*/ 889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nTest general API: \n")); 899566063dSJacob Faibussowitsch PetscCall(AOCreate(PETSC_COMM_WORLD, &ao)); 909566063dSJacob Faibussowitsch PetscCall(AOSetIS(ao, isapp, ispetsc)); 919566063dSJacob Faibussowitsch PetscCall(AOSetType(ao, AOMEMORYSCALABLE)); 929566063dSJacob Faibussowitsch PetscCall(AOSetFromOptions(ao)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* ispetsc and isapp are nolonger used. */ 959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispetsc)); 969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isapp)); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(AOPetscToApplication(ao, 4, getapp3)); 999566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, 3, getpetsc3)); 100c4762a1bSJed Brown 1019566063dSJacob 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])); 1029566063dSJacob 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])); 1039566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Check accuracy */; 106ad540459SPierre 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]); 107ad540459SPierre 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]); 108c4762a1bSJed Brown 1099566063dSJacob Faibussowitsch PetscCall(AODestroy(&ao)); 1109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 111b122ec5aSJacob Faibussowitsch return 0; 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /*TEST 115c4762a1bSJed Brown 116c4762a1bSJed Brown test: 117c4762a1bSJed Brown 118c4762a1bSJed Brown test: 119c4762a1bSJed Brown suffix: 2 120c4762a1bSJed Brown nsize: 2 121c4762a1bSJed Brown 122c4762a1bSJed Brown test: 123c4762a1bSJed Brown suffix: 3 124c4762a1bSJed Brown nsize: 3 125c4762a1bSJed Brown 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown suffix: 4 128c4762a1bSJed Brown nsize: 3 129c4762a1bSJed Brown args: -ao_type basic 130c4762a1bSJed Brown output_file: output/ex1_3.out 131c4762a1bSJed Brown 132c4762a1bSJed Brown TEST*/ 133