xref: /petsc/src/mat/tests/ex110.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1 static char help[] = "Testing MatCreateMPIAIJWithSplitArrays().\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc,char **argv)
6 {
7   Mat               A,B;
8   PetscInt          i,j,column,*ooj;
9   PetscInt          *di,*dj,*oi,*oj,nd;
10   const PetscInt    *garray;
11   PetscScalar       *oa,*da;
12   PetscScalar       value;
13   PetscRandom       rctx;
14   PetscErrorCode    ierr;
15   PetscBool         equal,done;
16   Mat               AA,AB;
17   PetscMPIInt       size,rank;
18 
19   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
20   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
21   PetscCheckFalse(size == 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Must run with 2 or more processes");
22   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23 
24   /* Create a mpiaij matrix for checking */
25   CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,5,5,PETSC_DECIDE,PETSC_DECIDE,0,NULL,0,NULL,&A));
26   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE));
27   CHKERRQ(MatSetUp(A));
28   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
29   CHKERRQ(PetscRandomSetFromOptions(rctx));
30 
31   for (i=5*rank; i<5*rank+5; i++) {
32     for (j=0; j<5*size; j++) {
33       CHKERRQ(PetscRandomGetValue(rctx,&value));
34       column = (PetscInt) (5*size*PetscRealPart(value));
35       CHKERRQ(PetscRandomGetValue(rctx,&value));
36       CHKERRQ(MatSetValues(A,1,&i,1,&column,&value,INSERT_VALUES));
37     }
38   }
39   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
40   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
41 
42   CHKERRQ(MatMPIAIJGetSeqAIJ(A,&AA,&AB,&garray));
43   CHKERRQ(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&di,(const PetscInt**)&dj,&done));
44   CHKERRQ(MatSeqAIJGetArray(AA,&da));
45   CHKERRQ(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&oi,(const PetscInt**)&oj,&done));
46   CHKERRQ(MatSeqAIJGetArray(AB,&oa));
47 
48   CHKERRQ(PetscMalloc1(oi[5],&ooj));
49   CHKERRQ(PetscArraycpy(ooj,oj,oi[5]));
50   /* modify the column entries in the non-diagonal portion back to global numbering */
51   for (i=0; i<oi[5]; i++) {
52     ooj[i] = garray[ooj[i]];
53   }
54 
55   CHKERRQ(MatCreateMPIAIJWithSplitArrays(PETSC_COMM_WORLD,5,5,PETSC_DETERMINE,PETSC_DETERMINE,di,dj,da,oi,ooj,oa,&B));
56   CHKERRQ(MatSetUp(B));
57   CHKERRQ(MatEqual(A,B,&equal));
58 
59   CHKERRQ(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&di,(const PetscInt**)&dj,&done));
60   CHKERRQ(MatSeqAIJRestoreArray(AA,&da));
61   CHKERRQ(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&oi,(const PetscInt**)&oj,&done));
62   CHKERRQ(MatSeqAIJRestoreArray(AB,&oa));
63 
64   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Likely a bug in MatCreateMPIAIJWithSplitArrays()");
65 
66   /* Free spaces */
67   CHKERRQ(PetscFree(ooj));
68   CHKERRQ(PetscRandomDestroy(&rctx));
69   CHKERRQ(MatDestroy(&A));
70   CHKERRQ(MatDestroy(&B));
71   CHKERRQ(PetscFree(oj));
72   ierr = PetscFinalize();
73   return ierr;
74 }
75 
76 /*TEST
77 
78    test:
79       nsize: 3
80 
81 TEST*/
82