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; 28*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 29c4762a1bSJed Brown 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutCreate(comm, &layout)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetSize(layout, n_global)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetLocalSize(layout, PETSC_DECIDE)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetUp(layout)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutGetLocalSize(layout, &local_size)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutGetRange(layout, &start, &end)); 36c4762a1bSJed Brown 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(local_size,&app_indices)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(local_size,&petsc_indices)); 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 */ 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(comm, local_size, &app_indices[0], PETSC_COPY_VALUES, &app_is)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(comm, local_size, &petsc_indices[0], PETSC_COPY_VALUES, &petsc_is)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOCreate(comm, &app2petsc)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOSetIS(app2petsc, app_is, petsc_is)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOSetType(app2petsc, AOMEMORYSCALABLE)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOSetFromOptions(app2petsc)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&app_is)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&petsc_is)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOView(app2petsc, PETSC_VIEWER_STDOUT_WORLD)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Test AOApplicationToPetsc */ 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n_loc,&ia)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n_loc,&ia0)); 59dd400576SPatrick 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 } 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(ia0,ia,n_loc)); 79c4762a1bSJed Brown 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOApplicationToPetsc(app2petsc, n_loc, ia)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown for (i=0; i<n_loc; ++i) { 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"proc = %d : %" PetscInt_FMT " -> %" PetscInt_FMT " \n", rank, ia0[i], ia[i])); 84c4762a1bSJed Brown } 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(AODestroy(&app2petsc)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutDestroy(&layout)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(app_indices)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(petsc_indices)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ia)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ia0)); 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