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