1c4762a1bSJed Brown /* 2c4762a1bSJed Brown Created by Huy Vo on 12/3/18. 3c4762a1bSJed Brown */ 4c4762a1bSJed Brown static char help[] = "Test memory scalable AO.\n\n"; 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include<petsc.h> 7c4762a1bSJed Brown #include<petscvec.h> 8c4762a1bSJed Brown #include<petscmat.h> 9c4762a1bSJed Brown #include<petscao.h> 10c4762a1bSJed Brown 11c4762a1bSJed Brown int main(int argc, char *argv[]) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown PetscInt ierr; 14c4762a1bSJed Brown PetscInt n_global = 16; 15c4762a1bSJed Brown MPI_Comm comm; 16c4762a1bSJed Brown PetscLayout layout; 17c4762a1bSJed Brown PetscInt local_size; 18c4762a1bSJed Brown PetscInt start, end; 19c4762a1bSJed Brown PetscMPIInt rank; 20c4762a1bSJed Brown PetscInt *app_indices,*petsc_indices,*ia,*ia0; 21c4762a1bSJed Brown PetscInt i; 22c4762a1bSJed Brown AO app2petsc; 23c4762a1bSJed Brown IS app_is, petsc_is; 24c4762a1bSJed Brown const PetscInt n_loc = 8; 25c4762a1bSJed Brown 26c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, (char *) 0, help); if (ierr) return ierr; 27c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 28ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 29c4762a1bSJed Brown 30c4762a1bSJed Brown ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscLayoutSetSize(layout, n_global);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscLayoutSetLocalSize(layout, PETSC_DECIDE);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = PetscLayoutGetLocalSize(layout, &local_size);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = PetscLayoutGetRange(layout, &start, &end);CHKERRQ(ierr); 36c4762a1bSJed Brown 37c4762a1bSJed Brown ierr = PetscMalloc1(local_size,&app_indices);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscMalloc1(local_size,&petsc_indices);CHKERRQ(ierr); 39c4762a1bSJed Brown /* Add values for local indices for usual states */ 40c4762a1bSJed Brown for (i = 0; i < local_size; ++i) { 41c4762a1bSJed Brown app_indices[i] = start + i; 42c4762a1bSJed Brown petsc_indices[i] = end -1 - i; 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Create the AO object that maps from lexicographic ordering to Petsc Vec ordering */ 46c4762a1bSJed Brown ierr = ISCreateGeneral(comm, local_size, &app_indices[0], PETSC_COPY_VALUES, &app_is);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = ISCreateGeneral(comm, local_size, &petsc_indices[0], PETSC_COPY_VALUES, &petsc_is);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = AOCreate(comm, &app2petsc);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = AOSetIS(app2petsc, app_is, petsc_is);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = AOSetType(app2petsc, AOMEMORYSCALABLE);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = AOSetFromOptions(app2petsc);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = ISDestroy(&app_is);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = ISDestroy(&petsc_is);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = AOView(app2petsc, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Test AOApplicationToPetsc */ 57c4762a1bSJed Brown ierr = PetscMalloc1(n_loc,&ia);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscMalloc1(n_loc,&ia0);CHKERRQ(ierr); 59*dd400576SPatrick Sanan if (rank == 0) { 60c4762a1bSJed Brown ia[0] = 0; 61c4762a1bSJed Brown ia[1] = -1; 62c4762a1bSJed Brown ia[2] = 1; 63c4762a1bSJed Brown ia[3] = 2; 64c4762a1bSJed Brown ia[4] = -1; 65c4762a1bSJed Brown ia[5] = 4; 66c4762a1bSJed Brown ia[6] = 5; 67c4762a1bSJed Brown ia[7] = 6; 68c4762a1bSJed Brown } else { 69c4762a1bSJed Brown ia[0] = -1; 70c4762a1bSJed Brown ia[1] = 8; 71c4762a1bSJed Brown ia[2] = 9; 72c4762a1bSJed Brown ia[3] = 10; 73c4762a1bSJed Brown ia[4] = -1; 74c4762a1bSJed Brown ia[5] = 12; 75c4762a1bSJed Brown ia[6] = 13; 76c4762a1bSJed Brown ia[7] = 14; 77c4762a1bSJed Brown } 78c4762a1bSJed Brown ierr = PetscArraycpy(ia0,ia,n_loc);CHKERRQ(ierr); 79c4762a1bSJed Brown 80c4762a1bSJed Brown ierr = AOApplicationToPetsc(app2petsc, n_loc, ia);CHKERRQ(ierr); 81c4762a1bSJed Brown 82c4762a1bSJed Brown for (i=0; i<n_loc; ++i) { 83c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"proc = %d : %D -> %D \n", rank, ia0[i], ia[i]);CHKERRQ(ierr); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = AODestroy(&app2petsc);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = PetscFree(app_indices);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = PetscFree(petsc_indices);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = PetscFree(ia);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = PetscFree(ia0);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = PetscFinalize(); 93c4762a1bSJed Brown return ierr; 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /*TEST 97c4762a1bSJed Brown 98c4762a1bSJed Brown test: 99c4762a1bSJed Brown nsize: 2 100c4762a1bSJed Brown 101c4762a1bSJed Brown TEST*/ 102