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