xref: /petsc/src/mat/tests/ex111.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Tests sequential and parallel MatMatMatMult() and MatPtAP(). Modified from ex96.c \n\
3c4762a1bSJed Brown   -Mx <xg>, where <xg> = number of coarse grid points in the x-direction\n\
4c4762a1bSJed Brown   -My <yg>, where <yg> = number of coarse grid points in the y-direction\n\
5c4762a1bSJed Brown   -Mz <zg>, where <zg> = number of coarse grid points in the z-direction\n\
6c4762a1bSJed Brown   -Npx <npx>, where <npx> = number of processors in the x-direction\n\
7c4762a1bSJed Brown   -Npy <npy>, where <npy> = number of processors in the y-direction\n\
8c4762a1bSJed Brown   -Npz <npz>, where <npz> = number of processors in the z-direction\n\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown     Example of usage: mpiexec -n 3 ./ex41 -Mx 10 -My 10 -Mz 10
12c4762a1bSJed Brown */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscdm.h>
15c4762a1bSJed Brown #include <petscdmda.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /* User-defined application contexts */
18c4762a1bSJed Brown typedef struct {
19c4762a1bSJed Brown   PetscInt mx,my,mz;            /* number grid points in x, y and z direction */
20c4762a1bSJed Brown   Vec      localX,localF;       /* local vectors with ghost region */
21c4762a1bSJed Brown   DM       da;
22c4762a1bSJed Brown   Vec      x,b,r;               /* global vectors */
23c4762a1bSJed Brown   Mat      J;                   /* Jacobian on grid */
24c4762a1bSJed Brown } GridCtx;
25c4762a1bSJed Brown typedef struct {
26c4762a1bSJed Brown   GridCtx  fine;
27c4762a1bSJed Brown   GridCtx  coarse;
28c4762a1bSJed Brown   PetscInt ratio;
29c4762a1bSJed Brown   Mat      Ii;                  /* interpolation from coarse to fine */
30c4762a1bSJed Brown } AppCtx;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown #define COARSE_LEVEL 0
33c4762a1bSJed Brown #define FINE_LEVEL   1
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown       Mm_ratio - ration of grid lines between fine and coarse grids.
37c4762a1bSJed Brown */
38c4762a1bSJed Brown int main(int argc,char **argv)
39c4762a1bSJed Brown {
40c4762a1bSJed Brown   AppCtx         user;
41c4762a1bSJed Brown   PetscMPIInt    size,rank;
42c4762a1bSJed Brown   PetscInt       m,n,M,N,i,nrows;
43c4762a1bSJed Brown   PetscScalar    one = 1.0;
44c4762a1bSJed Brown   PetscReal      fill=2.0;
45c20d7725SJed Brown   Mat            A,P,R,C,PtAP,D;
46c4762a1bSJed Brown   PetscScalar    *array;
47c4762a1bSJed Brown   PetscRandom    rdm;
48c4762a1bSJed Brown   PetscBool      Test_3D=PETSC_FALSE,flg;
49c4762a1bSJed Brown   const PetscInt *ia,*ja;
50c4762a1bSJed Brown 
51*327415f7SBarry Smith   PetscFunctionBeginUser;
529566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Get size of fine grids and coarse grids */
57c4762a1bSJed Brown   user.ratio     = 2;
58c4762a1bSJed Brown   user.coarse.mx = 4; user.coarse.my = 4; user.coarse.mz = 4;
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL));
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL));
64c4762a1bSJed Brown   if (user.coarse.mz) Test_3D = PETSC_TRUE;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   user.fine.mx = user.ratio*(user.coarse.mx-1)+1;
67c4762a1bSJed Brown   user.fine.my = user.ratio*(user.coarse.my-1)+1;
68c4762a1bSJed Brown   user.fine.mz = user.ratio*(user.coarse.mz-1)+1;
69c4762a1bSJed Brown 
70dd400576SPatrick Sanan   if (rank == 0) {
71c4762a1bSJed Brown     if (!Test_3D) {
729566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"coarse grids: %" PetscInt_FMT " %" PetscInt_FMT "; fine grids: %" PetscInt_FMT " %" PetscInt_FMT "\n",user.coarse.mx,user.coarse.my,user.fine.mx,user.fine.my));
73c4762a1bSJed Brown     } else {
749566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"coarse grids: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "; fine grids: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",user.coarse.mx,user.coarse.my,user.coarse.mz,user.fine.mx,user.fine.my,user.fine.mz));
75c4762a1bSJed Brown     }
76c4762a1bSJed Brown   }
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* Set up distributed array for fine grid */
79c4762a1bSJed Brown   if (!Test_3D) {
809566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.fine.mx,user.fine.my,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.fine.da));
81c4762a1bSJed Brown   } else {
82d0609cedSBarry Smith     PetscCall(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.fine.mx,user.fine.my,user.fine.mz,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,
83d0609cedSBarry Smith                            1,1,NULL,NULL,NULL,&user.fine.da));
84c4762a1bSJed Brown   }
859566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.fine.da));
869566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.fine.da));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Create and set A at fine grids */
899566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(user.fine.da,MATAIJ));
909566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.fine.da,&A));
919566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
929566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&M,&N));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* set val=one to A (replace with random values!) */
959566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
969566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
97c4762a1bSJed Brown   if (size == 1) {
989566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
99c4762a1bSJed Brown     if (flg) {
1009566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(A,&array));
101c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
1029566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(A,&array));
103c4762a1bSJed Brown     }
1049566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
105c4762a1bSJed Brown   } else {
106c4762a1bSJed Brown     Mat AA,AB;
1079566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL));
1089566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
109c4762a1bSJed Brown     if (flg) {
1109566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AA,&array));
111c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
1129566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AA,&array));
113c4762a1bSJed Brown     }
1149566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
1159566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
116c4762a1bSJed Brown     if (flg) {
1179566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AB,&array));
118c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
1199566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AB,&array));
120c4762a1bSJed Brown     }
1219566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown   /* Set up distributed array for coarse grid */
124c4762a1bSJed Brown   if (!Test_3D) {
1259566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.coarse.mx,user.coarse.my,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.coarse.da));
126c4762a1bSJed Brown   } else {
1279566063dSJacob Faibussowitsch     PetscCall(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.coarse.mx,user.coarse.my,user.coarse.mz,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&user.coarse.da));
128c4762a1bSJed Brown   }
1299566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.coarse.da));
1309566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.coarse.da));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* Create interpolation between the fine and coarse grids */
1339566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* Get R = P^T */
1369566063dSJacob Faibussowitsch   PetscCall(MatTranspose(P,MAT_INITIAL_MATRIX,&R));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* C = R*A*P */
139c20d7725SJed Brown   /* Developer's API */
1409566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(R,A,P,&D));
1419566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(D,MATPRODUCT_ABC));
1429566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(D));
1439566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(D));
1449566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(D));
1459566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(D)); /* Test reuse symbolic D */
146c20d7725SJed Brown 
147c20d7725SJed Brown   /* User's API */
148c20d7725SJed Brown   { /* Test MatMatMatMult_Basic() */
149c20d7725SJed Brown     Mat Adense,Cdense;
1509566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense));
1519566063dSJacob Faibussowitsch     PetscCall(MatMatMatMult(R,Adense,P,MAT_INITIAL_MATRIX,fill,&Cdense));
1529566063dSJacob Faibussowitsch     PetscCall(MatMatMatMult(R,Adense,P,MAT_REUSE_MATRIX,fill,&Cdense));
153c20d7725SJed Brown 
1549566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(D,Cdense,10,&flg));
15528b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"D*v != Cdense*v");
1569566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Adense));
1579566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cdense));
158c20d7725SJed Brown   }
159c20d7725SJed Brown 
1609566063dSJacob Faibussowitsch   PetscCall(MatMatMatMult(R,A,P,MAT_INITIAL_MATRIX,fill,&C));
1619566063dSJacob Faibussowitsch   PetscCall(MatMatMatMult(R,A,P,MAT_REUSE_MATRIX,fill,&C));
1629566063dSJacob Faibussowitsch   PetscCall(MatProductClear(C));
163c4762a1bSJed Brown 
164c20d7725SJed Brown   /* Test D == C */
1659566063dSJacob Faibussowitsch   PetscCall(MatEqual(D,C,&flg));
16628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"D != C");
167c20d7725SJed Brown 
168c4762a1bSJed Brown   /* Test C == PtAP */
1699566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&PtAP));
1709566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&PtAP));
1719566063dSJacob Faibussowitsch   PetscCall(MatEqual(C,PtAP,&flg));
17228b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"C != PtAP");
1739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&PtAP));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* Clean up */
1769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1779566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1789566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.fine.da));
1799566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.coarse.da));
1809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
1819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&R));
1829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
1849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
185b122ec5aSJacob Faibussowitsch   return 0;
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /*TEST
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    test:
191c4762a1bSJed Brown 
192c4762a1bSJed Brown    test:
193c4762a1bSJed Brown       suffix: 2
194c4762a1bSJed Brown       nsize: 2
195c4762a1bSJed Brown       args: -matmatmatmult_via scalable
196c4762a1bSJed Brown 
197c4762a1bSJed Brown    test:
198c4762a1bSJed Brown       suffix: 3
199c4762a1bSJed Brown       nsize: 2
200c4762a1bSJed Brown       args: -matmatmatmult_via nonscalable
201c4762a1bSJed Brown       output_file: output/ex111_1.out
202c4762a1bSJed Brown 
203c4762a1bSJed Brown TEST*/
204