1c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIAIJWithSplitArrays().\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **argv) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown Mat A,B; 8c4762a1bSJed Brown PetscInt i,j,column,*ooj; 9c4762a1bSJed Brown PetscInt *di,*dj,*oi,*oj,nd; 10c4762a1bSJed Brown const PetscInt *garray; 11c4762a1bSJed Brown PetscScalar *oa,*da; 12c4762a1bSJed Brown PetscScalar value; 13c4762a1bSJed Brown PetscRandom rctx; 14c4762a1bSJed Brown PetscBool equal,done; 15c4762a1bSJed Brown Mat AA,AB; 16c4762a1bSJed Brown PetscMPIInt size,rank; 17c4762a1bSJed Brown 18*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 202c71b3e2SJacob Faibussowitsch PetscCheckFalse(size == 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Must run with 2 or more processes"); 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Create a mpiaij matrix for checking */ 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,5,5,PETSC_DECIDE,PETSC_DECIDE,0,NULL,0,NULL,&A)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rctx)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown for (i=5*rank; i<5*rank+5; i++) { 31c4762a1bSJed Brown for (j=0; j<5*size; j++) { 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rctx,&value)); 33c4762a1bSJed Brown column = (PetscInt) (5*size*PetscRealPart(value)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rctx,&value)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&column,&value,INSERT_VALUES)); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown } 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 40c4762a1bSJed Brown 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJGetSeqAIJ(A,&AA,&AB,&garray)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&di,(const PetscInt**)&dj,&done)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AA,&da)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&oi,(const PetscInt**)&oj,&done)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AB,&oa)); 46c4762a1bSJed Brown 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(oi[5],&ooj)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(ooj,oj,oi[5])); 49c4762a1bSJed Brown /* modify the column entries in the non-diagonal portion back to global numbering */ 50c4762a1bSJed Brown for (i=0; i<oi[5]; i++) { 51c4762a1bSJed Brown ooj[i] = garray[ooj[i]]; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIAIJWithSplitArrays(PETSC_COMM_WORLD,5,5,PETSC_DETERMINE,PETSC_DETERMINE,di,dj,da,oi,ooj,oa,&B)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(A,B,&equal)); 57c4762a1bSJed Brown 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&di,(const PetscInt**)&dj,&done)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AA,&da)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&oi,(const PetscInt**)&oj,&done)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AB,&oa)); 62c4762a1bSJed Brown 6328b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Likely a bug in MatCreateMPIAIJWithSplitArrays()"); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Free spaces */ 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ooj)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(oj)); 71*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 72*b122ec5aSJacob Faibussowitsch return 0; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /*TEST 76c4762a1bSJed Brown 77c4762a1bSJed Brown test: 78c4762a1bSJed Brown nsize: 3 79c4762a1bSJed Brown 80c4762a1bSJed Brown TEST*/ 81