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