1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] ="Tests MatCreateMPIAIJWithArrays() abd MatUpdateMPIAIJWithArrays()\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown /* 7c4762a1bSJed Brown * This is an extremely simple example to test MatUpdateMPIAIJWithArrays() 8c4762a1bSJed Brown * 9c4762a1bSJed Brown * A = 10c4762a1bSJed Brown 11c4762a1bSJed Brown 1 2 0 3 0 0 12c4762a1bSJed Brown 0 4 5 0 0 6 13c4762a1bSJed Brown 7 0 8 0 9 0 14c4762a1bSJed Brown 0 10 11 12 0 13 15c4762a1bSJed Brown 0 14 15 0 0 16 16c4762a1bSJed Brown 17 0 0 0 0 18 17c4762a1bSJed Brown * 18c4762a1bSJed Brown * */ 19c4762a1bSJed Brown 20c4762a1bSJed Brown int main(int argc,char **argv) 21c4762a1bSJed Brown { 22c4762a1bSJed Brown Mat A; 23c4762a1bSJed Brown PetscInt i[3][3] = {{0, 3, 6},{0, 3, 7},{0, 3, 5}}; 24c4762a1bSJed Brown PetscInt j[3][7] = {{0, 1, 3, 1, 2, 5, -1},{0, 2, 4, 1, 2, 3, 5},{1, 2, 5, 0, 5, -1, -1}}; 25c4762a1bSJed Brown PetscScalar a[3][7] = {{1, 2, 3, 4, 5, 6, -1}, {7, 8, 9, 10, 11, 12, 13},{14, 15, 16, 17, 18, -1, -1}}; 26c4762a1bSJed Brown PetscScalar anew[3][7] = {{10, 20, 30, 40, 50, 60, -1}, {70, 80, 90, 100, 110, 120, 130},{140, 150, 160, 170, 180, -1, -1}}; 27c4762a1bSJed Brown MPI_Comm comm; 28c4762a1bSJed Brown PetscMPIInt rank,size; 29c4762a1bSJed Brown 30*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 31c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 325f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm,&rank)); 335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 342c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 3,comm,PETSC_ERR_ARG_INCOMP,"You have to use three MPI processes to run this example "); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIAIJWithArrays(comm,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],a[rank],&A)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,NULL)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatUpdateMPIAIJWithArrays(A,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],anew[rank])); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,NULL)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatUpdateMPIAIJWithArrays(A,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],a[rank])); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,NULL)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 42*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 43*b122ec5aSJacob Faibussowitsch return 0; 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown /*TEST 47c4762a1bSJed Brown test: 48c4762a1bSJed Brown nsize: 3 49c4762a1bSJed Brown 50c4762a1bSJed Brown TEST*/ 51