xref: /petsc/src/mat/tests/ex89.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] ="Tests MatPtAP() for MPIMAIJ and MPIAIJ \n ";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscdmda.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown int main(int argc,char **argv)
6*c4762a1bSJed Brown {
7*c4762a1bSJed Brown   PetscErrorCode ierr;
8*c4762a1bSJed Brown   DM             coarsedm,finedm;
9*c4762a1bSJed Brown   PetscMPIInt    size,rank;
10*c4762a1bSJed Brown   PetscInt       M,N,Z,i,nrows;
11*c4762a1bSJed Brown   PetscScalar    one = 1.0;
12*c4762a1bSJed Brown   PetscReal      fill=2.0;
13*c4762a1bSJed Brown   Mat            A,P,C;
14*c4762a1bSJed Brown   PetscScalar    *array,none = -1.0,alpha;
15*c4762a1bSJed Brown   Vec            x,v1,v2,v3,v4;
16*c4762a1bSJed Brown   PetscReal      norm,norm_tmp,norm_tmp1,tol=100.*PETSC_MACHINE_EPSILON;
17*c4762a1bSJed Brown   PetscRandom    rdm;
18*c4762a1bSJed Brown   PetscBool      Test_3D=PETSC_FALSE,flg;
19*c4762a1bSJed Brown   const PetscInt *ia,*ja;
20*c4762a1bSJed Brown   PetscInt       dof;
21*c4762a1bSJed Brown   MPI_Comm       comm;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
24*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
25*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
27*c4762a1bSJed Brown   M = 10; N = 10; Z = 10;
28*c4762a1bSJed Brown   dof  = 10;
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown   /* Set up distributed array for fine grid */
35*c4762a1bSJed Brown   if (!Test_3D) {
36*c4762a1bSJed Brown     ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm);CHKERRQ(ierr);
37*c4762a1bSJed Brown   } else {
38*c4762a1bSJed Brown     ierr = DMDACreate3d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,Z,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,NULL,&coarsedm);CHKERRQ(ierr);
39*c4762a1bSJed Brown   }
40*c4762a1bSJed Brown   ierr = DMSetFromOptions(coarsedm);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = DMSetUp(coarsedm);CHKERRQ(ierr);
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   /* This makes sure the coarse DMDA has the same partition as the fine DMDA */
44*c4762a1bSJed Brown   ierr = DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm);CHKERRQ(ierr);
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown   /*------------------------------------------------------------*/
47*c4762a1bSJed Brown   ierr = DMSetMatType(finedm,MATAIJ);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = DMCreateMatrix(finedm,&A);CHKERRQ(ierr);
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   /* set val=one to A */
51*c4762a1bSJed Brown   if (size == 1) {
52*c4762a1bSJed Brown     ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
53*c4762a1bSJed Brown     if (flg) {
54*c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr);
55*c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
56*c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(A,&array);CHKERRQ(ierr);
57*c4762a1bSJed Brown     }
58*c4762a1bSJed Brown     ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
59*c4762a1bSJed Brown   } else {
60*c4762a1bSJed Brown     Mat AA,AB;
61*c4762a1bSJed Brown     ierr = MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL);CHKERRQ(ierr);
62*c4762a1bSJed Brown     ierr = MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
63*c4762a1bSJed Brown     if (flg) {
64*c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(AA,&array);CHKERRQ(ierr);
65*c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
66*c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(AA,&array);CHKERRQ(ierr);
67*c4762a1bSJed Brown     }
68*c4762a1bSJed Brown     ierr = MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
69*c4762a1bSJed Brown     ierr = MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
70*c4762a1bSJed Brown     if (flg) {
71*c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(AB,&array);CHKERRQ(ierr);
72*c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
73*c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(AB,&array);CHKERRQ(ierr);
74*c4762a1bSJed Brown     }
75*c4762a1bSJed Brown     ierr = MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
76*c4762a1bSJed Brown   }
77*c4762a1bSJed Brown   /* Create interpolation between the fine and coarse grids */
78*c4762a1bSJed Brown   ierr = DMCreateInterpolation(coarsedm,finedm,&P,NULL);CHKERRQ(ierr);
79*c4762a1bSJed Brown   /* Create vectors v1 and v2 that are compatible with A */
80*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&v1,NULL);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = VecDuplicate(v1,&v2);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = PetscRandomCreate(comm,&rdm);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown   /* Test P^T * A * P - MatPtAP() */
86*c4762a1bSJed Brown   /*------------------------------*/
87*c4762a1bSJed Brown   ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr);
88*c4762a1bSJed Brown   /* Test MAT_REUSE_MATRIX - reuse symbolic C */
89*c4762a1bSJed Brown   alpha=1.0;
90*c4762a1bSJed Brown   for (i=0; i<1; i++) {
91*c4762a1bSJed Brown     alpha -= 0.1;
92*c4762a1bSJed Brown     ierr   = MatScale(A,alpha);CHKERRQ(ierr);
93*c4762a1bSJed Brown     ierr   = MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr);
94*c4762a1bSJed Brown   }
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   /* Free intermediate data structures created for reuse of C=Pt*A*P */
97*c4762a1bSJed Brown   ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown   /* Create vector x that is compatible with P */
100*c4762a1bSJed Brown   ierr = MatCreateVecs(P,&x,NULL);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = VecDuplicate(x,&v3);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = VecDuplicate(x,&v4);CHKERRQ(ierr);
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown   norm = 0.0;
105*c4762a1bSJed Brown   for (i=0; i<10; i++) {
106*c4762a1bSJed Brown     ierr = VecSetRandom(x,rdm);CHKERRQ(ierr);
107*c4762a1bSJed Brown     ierr = MatMult(P,x,v1);CHKERRQ(ierr);
108*c4762a1bSJed Brown     ierr = MatMult(A,v1,v2);CHKERRQ(ierr);  /* v2 = A*P*x */
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown     ierr = MatMultTranspose(P,v2,v3);CHKERRQ(ierr); /* v3 = Pt*A*P*x */
111*c4762a1bSJed Brown     ierr = MatMult(C,x,v4);CHKERRQ(ierr);           /* v3 = C*x   */
112*c4762a1bSJed Brown     ierr = VecAXPY(v4,none,v3);CHKERRQ(ierr);
113*c4762a1bSJed Brown     ierr = VecNorm(v4,NORM_1,&norm_tmp);CHKERRQ(ierr);
114*c4762a1bSJed Brown     ierr = VecNorm(v3,NORM_1,&norm_tmp1);CHKERRQ(ierr);
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown     norm_tmp /= norm_tmp1;
117*c4762a1bSJed Brown     if (norm_tmp > norm) norm = norm_tmp;
118*c4762a1bSJed Brown   }
119*c4762a1bSJed Brown   if (norm >= tol && !rank) {
120*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatPtAP(), |v3 - v4|/|v3|: %g\n",(double)norm);CHKERRQ(ierr);
121*c4762a1bSJed Brown   }
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = MatDestroy(&P);CHKERRQ(ierr);
126*c4762a1bSJed Brown   ierr = VecDestroy(&v3);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = VecDestroy(&v4);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = VecDestroy(&v1);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = VecDestroy(&v2);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = DMDestroy(&finedm);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = DMDestroy(&coarsedm);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = PetscFinalize();
135*c4762a1bSJed Brown   return ierr;
136*c4762a1bSJed Brown }
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown /*TEST
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown    test:
141*c4762a1bSJed Brown       args: -M 10 -N 10 -Z 10
142*c4762a1bSJed Brown       output_file: output/ex89_1.out
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown    test:
145*c4762a1bSJed Brown       suffix: allatonce
146*c4762a1bSJed Brown       nsize: 4
147*c4762a1bSJed Brown       args: -M 10 -N 10 -Z 10 -matmaijptap_via allatonce
148*c4762a1bSJed Brown       output_file: output/ex89_1.out
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown    test:
151*c4762a1bSJed Brown       suffix: allatonce_merged
152*c4762a1bSJed Brown       nsize: 4
153*c4762a1bSJed Brown       args: -M 10 -M 5 -M 10 -matmaijptap_via allatonce_merged
154*c4762a1bSJed Brown       output_file: output/ex96_1.out
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown    test:
157*c4762a1bSJed Brown       suffix: allatonce_3D
158*c4762a1bSJed Brown       nsize: 4
159*c4762a1bSJed Brown       args: -M 10 -M 5 -M 10 -test_3D 1 -matmaijptap_via allatonce
160*c4762a1bSJed Brown       output_file: output/ex96_1.out
161*c4762a1bSJed Brown 
162*c4762a1bSJed Brown    test:
163*c4762a1bSJed Brown       suffix: allatonce_merged_3D
164*c4762a1bSJed Brown       nsize: 4
165*c4762a1bSJed Brown       args: -M 10 -M 5 -M 10 -test_3D 1 -matmaijptap_via allatonce_merged
166*c4762a1bSJed Brown       output_file: output/ex96_1.out
167*c4762a1bSJed Brown 
168*c4762a1bSJed Brown TEST*/
169