xref: /petsc/src/mat/tests/ex110.c (revision be096a4675ce3b20dd7f9c195e6046e69d9955cf)
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 
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20*be096a46SBarry Smith   PetscCheck(size > 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 2 or more processes");
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* Create a mpiaij matrix for checking */
249566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,5,5,PETSC_DECIDE,PETSC_DECIDE,0,NULL,0,NULL,&A));
259566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE));
269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
279566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   for (i=5*rank; i<5*rank+5; i++) {
31c4762a1bSJed Brown     for (j=0; j<5*size; j++) {
329566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rctx,&value));
33c4762a1bSJed Brown       column = (PetscInt) (5*size*PetscRealPart(value));
349566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rctx,&value));
359566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,1,&column,&value,INSERT_VALUES));
36c4762a1bSJed Brown     }
37c4762a1bSJed Brown   }
389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A,&AA,&AB,&garray));
429566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&di,(const PetscInt**)&dj,&done));
439566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArray(AA,&da));
449566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&oi,(const PetscInt**)&oj,&done));
459566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArray(AB,&oa));
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(oi[5],&ooj));
489566063dSJacob Faibussowitsch   PetscCall(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 
549566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAIJWithSplitArrays(PETSC_COMM_WORLD,5,5,PETSC_DETERMINE,PETSC_DETERMINE,di,dj,da,oi,ooj,oa,&B));
559566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
569566063dSJacob Faibussowitsch   PetscCall(MatEqual(A,B,&equal));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&di,(const PetscInt**)&dj,&done));
599566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArray(AA,&da));
609566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nd,(const PetscInt**)&oi,(const PetscInt**)&oj,&done));
619566063dSJacob Faibussowitsch   PetscCall(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 */
669566063dSJacob Faibussowitsch   PetscCall(PetscFree(ooj));
679566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
709566063dSJacob Faibussowitsch   PetscCall(PetscFree(oj));
719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
72b122ec5aSJacob 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