Actual source code: ex96.c
2: static char help[] ="Tests sequential and parallel MatMatMult() and MatPtAP()\n\
3: -Mx <xg>, where <xg> = number of coarse grid points in the x-direction\n\
4: -My <yg>, where <yg> = number of coarse grid points in the y-direction\n\
5: -Mz <zg>, where <zg> = number of coarse grid points in the z-direction\n\
6: -Npx <npx>, where <npx> = number of processors in the x-direction\n\
7: -Npy <npy>, where <npy> = number of processors in the y-direction\n\
8: -Npz <npz>, where <npz> = number of processors in the z-direction\n\n";
10: /* This test is modified from ~src/ksp/examples/tests/ex19.c. */
12: #include petscksp.h
13: #include petscda.h
14: #include petscmg.h
15: #include src/mat/impls/aij/seq/aij.h
16: #include src/mat/impls/aij/mpi/mpiaij.h
18: /* User-defined application contexts */
19: typedef struct {
20: PetscInt mx,my,mz; /* number grid points in x, y and z direction */
21: Vec localX,localF; /* local vectors with ghost region */
22: DA da;
23: Vec x,b,r; /* global vectors */
24: Mat J; /* Jacobian on grid */
25: } GridCtx;
26: typedef struct {
27: GridCtx fine;
28: GridCtx coarse;
29: KSP ksp_coarse;
30: PetscInt ratio;
31: Mat I; /* interpolation from coarse to fine */
32: } AppCtx;
34: #define COARSE_LEVEL 0
35: #define FINE_LEVEL 1
37: /*
38: Mm_ratio - ration of grid lines between fine and coarse grids.
39: */
42: int main(int argc,char **argv)
43: {
45: AppCtx user;
46: PetscInt Nx=PETSC_DECIDE,Ny=PETSC_DECIDE,Nz=PETSC_DECIDE;
47: PetscMPIInt size,rank;
48: PetscInt m,n,M,N,i,nrows,*ia,*ja;
49: PetscScalar one = 1.0;
50: PetscReal fill=2.0;
51: Mat A,A_tmp,P,C;
52: PetscScalar *array,none = -1.0,alpha;
53: PetscTruth flg;
54: Vec x,v1,v2,v3,v4;
55: PetscReal norm,norm_tmp,norm_tmp1,tol=1.e-12;
56: PetscRandom rdm;
57: PetscTruth Test_MatMatMult=PETSC_TRUE,Test_MatPtAP=PETSC_TRUE,Test_3D=PETSC_FALSE;
59: PetscInitialize(&argc,&argv,PETSC_NULL,help);
60: PetscOptionsGetReal(PETSC_NULL,"-tol",&tol,PETSC_NULL);
62: user.ratio = 2;
63: user.coarse.mx = 2; user.coarse.my = 2; user.coarse.mz = 0;
64: PetscOptionsGetInt(PETSC_NULL,"-Mx",&user.coarse.mx,PETSC_NULL);
65: PetscOptionsGetInt(PETSC_NULL,"-My",&user.coarse.my,PETSC_NULL);
66: PetscOptionsGetInt(PETSC_NULL,"-Mz",&user.coarse.mz,PETSC_NULL);
67: PetscOptionsGetInt(PETSC_NULL,"-ratio",&user.ratio,PETSC_NULL);
68: if (user.coarse.mz) Test_3D = PETSC_TRUE;
70: user.fine.mx = user.ratio*(user.coarse.mx-1)+1;
71: user.fine.my = user.ratio*(user.coarse.my-1)+1;
72: user.fine.mz = user.ratio*(user.coarse.mz-1)+1;
73: MPI_Comm_size(PETSC_COMM_WORLD,&size);
74: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
75: PetscOptionsGetInt(PETSC_NULL,"-Npx",&Nx,PETSC_NULL);
76: PetscOptionsGetInt(PETSC_NULL,"-Npy",&Ny,PETSC_NULL);
77: PetscOptionsGetInt(PETSC_NULL,"-Npz",&Nz,PETSC_NULL);
79: /* Set up distributed array for fine grid */
80: if (!Test_3D){
81: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.fine.mx,
82: user.fine.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.fine.da);
83: } else {
84: DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,
85: user.fine.mx,user.fine.my,user.fine.mz,Nx,Ny,Nz,
86: 1,1,PETSC_NULL,PETSC_NULL,PETSC_NULL,&user.fine.da);
87: }
89: DAGetMatrix(user.fine.da,MATAIJ,&A);
90: MatGetLocalSize(A,&m,&n);
91: MatGetSize(A,&M,&N);
92: /* set val=one to A */
93: if (size == 1){
94: MatGetRowIJ(A,0,PETSC_FALSE,&nrows,&ia,&ja,&flg);
95: if (flg){
96: MatGetArray(A,&array);
97: for (i=0; i<ia[nrows]; i++) array[i] = one;
98: MatRestoreArray(A,&array);
99: }
100: MatRestoreRowIJ(A,0,PETSC_FALSE,&nrows,&ia,&ja,&flg);
101: } else {
102: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data;
103: Mat_SeqAIJ *a=(Mat_SeqAIJ*)(aij->A)->data, *b=(Mat_SeqAIJ*)(aij->B)->data;
104: /* A_part */
105: for (i=0; i<a->i[m]; i++) a->a[i] = one;
106: /* B_part */
107: for (i=0; i<b->i[m]; i++) b->a[i] = one;
108:
109: }
110: /* MatView(A, PETSC_VIEWER_STDOUT_WORLD); */
112: /* Set up distributed array for coarse grid */
113: if (!Test_3D){
114: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.coarse.mx,
115: user.coarse.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.coarse.da);
116: } else {
117: DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,
118: user.coarse.mx,user.coarse.my,user.coarse.mz,Nx,Ny,Nz,
119: 1,1,PETSC_NULL,PETSC_NULL,PETSC_NULL,&user.coarse.da);
120: }
122: /* Create interpolation between the levels */
123: DAGetInterpolation(user.coarse.da,user.fine.da,&P,PETSC_NULL);
124:
125: MatGetLocalSize(P,&m,&n);
126: MatGetSize(P,&M,&N);
128: /* Create vectors v1 and v2 that are compatible with A */
129: VecCreate(PETSC_COMM_WORLD,&v1);
130: MatGetLocalSize(A,&m,PETSC_NULL);
131: VecSetSizes(v1,m,PETSC_DECIDE);
132: VecSetFromOptions(v1);
133: VecDuplicate(v1,&v2);
134: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rdm);
136: /* Test MatMatMult(): C = A*P */
137: /*----------------------------*/
138: if (Test_MatMatMult){
139: MatDuplicate(A,MAT_COPY_VALUES,&A_tmp);
140: MatMatMult(A_tmp,P,MAT_INITIAL_MATRIX,fill,&C);
141:
142: /* Test MAT_REUSE_MATRIX - reuse symbolic C */
143: alpha=1.0;
144: for (i=0; i<2; i++){
145: alpha -=0.1;
146: MatScale(&alpha,A_tmp);
147: MatMatMult(A_tmp,P,MAT_REUSE_MATRIX,fill,&C);
148: }
149: /*
150: if (rank == 0) PetscPrintf(PETSC_COMM_SELF, " \nA*P: \n");
151: MatView(C, PETSC_VIEWER_STDOUT_WORLD);
152: */
154: /* Create vector x that is compatible with P */
155: VecCreate(PETSC_COMM_WORLD,&x);
156: MatGetLocalSize(P,PETSC_NULL,&n);
157: VecSetSizes(x,n,PETSC_DECIDE);
158: VecSetFromOptions(x);
160: norm = 0.0;
161: for (i=0; i<10; i++) {
162: VecSetRandom(rdm,x);
163: MatMult(P,x,v1);
164: MatMult(A_tmp,v1,v2); /* v2 = A*P*x */
165: MatMult(C,x,v1); /* v1 = C*x */
166: VecAXPY(&none,v2,v1);
167: VecNorm(v1,NORM_1,&norm_tmp);
168: VecNorm(v2,NORM_1,&norm_tmp1);
169: norm_tmp /= norm_tmp1;
170: if (norm_tmp > norm) norm = norm_tmp;
171: }
172: if (norm >= tol && !rank) {
173: PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult(), |v1 - v2|/|v2|: %g\n",norm);
174: }
175:
176: VecDestroy(x);
177: MatDestroy(C);
178: MatDestroy(A_tmp);
179: }
181: /* Test P^T * A * P - MatPtAP() */
182: /*------------------------------*/
183: if (Test_MatPtAP){
184: MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C);
185: MatGetLocalSize(C,&m,&n);
186:
187: /* Test MAT_REUSE_MATRIX - reuse symbolic C */
188: alpha=1.0;
189: for (i=0; i<1; i++){
190: alpha -=0.1;
191: MatScale(&alpha,A);
192: MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C);
193: }
195: /* Create vector x that is compatible with P */
196: VecCreate(PETSC_COMM_WORLD,&x);
197: MatGetLocalSize(P,&m,&n);
198: VecSetSizes(x,n,PETSC_DECIDE);
199: VecSetFromOptions(x);
200:
201: VecCreate(PETSC_COMM_WORLD,&v3);
202: VecSetSizes(v3,n,PETSC_DECIDE);
203: VecSetFromOptions(v3);
204: VecDuplicate(v3,&v4);
206: norm = 0.0;
207: for (i=0; i<10; i++) {
208: VecSetRandom(rdm,x);
209: MatMult(P,x,v1);
210: MatMult(A,v1,v2); /* v2 = A*P*x */
212: MatMultTranspose(P,v2,v3); /* v3 = Pt*A*P*x */
213: MatMult(C,x,v4); /* v3 = C*x */
214: VecAXPY(&none,v3,v4);
215: VecNorm(v4,NORM_1,&norm_tmp);
216: VecNorm(v3,NORM_1,&norm_tmp1);
217: norm_tmp /= norm_tmp1;
218: if (norm_tmp > norm) norm = norm_tmp;
219: }
220: if (norm >= tol && !rank) {
221: PetscPrintf(PETSC_COMM_SELF,"Error: MatPtAP(), |v3 - v4|/|v3|: %g\n",norm);
222: }
223:
224: MatDestroy(C);
225: VecDestroy(v3);
226: VecDestroy(v4);
227: VecDestroy(x);
228:
229: }
231: /* Clean up */
232: MatDestroy(A);
233: PetscRandomDestroy(rdm);
234: VecDestroy(v1);
235: VecDestroy(v2);
236: DADestroy(user.fine.da);
237: DADestroy(user.coarse.da);
238: MatDestroy(P);
240: PetscFinalize();
242: return 0;
243: }