xref: /petsc/src/mat/tests/ex242.c (revision d24d420451aab1d024a728b511d7b61b40627436)
1*d24d4204SJose E. Roman 
2*d24d4204SJose E. Roman static char help[] = "Tests ScaLAPACK interface.\n\n";
3*d24d4204SJose E. Roman 
4*d24d4204SJose E. Roman #include <petscmat.h>
5*d24d4204SJose E. Roman 
6*d24d4204SJose E. Roman int main(int argc,char **args)
7*d24d4204SJose E. Roman {
8*d24d4204SJose E. Roman   Mat            Cdense,Caij,B,C,Ct;
9*d24d4204SJose E. Roman   Vec            d;
10*d24d4204SJose E. Roman   PetscInt       i,j,M = 5,N,mb = 2,nb,nrows,ncols,mloc,nloc;
11*d24d4204SJose E. Roman   const PetscInt *rows,*cols;
12*d24d4204SJose E. Roman   IS             isrows,iscols;
13*d24d4204SJose E. Roman   PetscErrorCode ierr;
14*d24d4204SJose E. Roman   PetscScalar    *v;
15*d24d4204SJose E. Roman   PetscMPIInt    rank;
16*d24d4204SJose E. Roman   PetscReal      Cnorm;
17*d24d4204SJose E. Roman   PetscBool      flg,mats_view=PETSC_FALSE;
18*d24d4204SJose E. Roman 
19*d24d4204SJose E. Roman   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20*d24d4204SJose E. Roman   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
21*d24d4204SJose E. Roman   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
22*d24d4204SJose E. Roman   N    = M;
23*d24d4204SJose E. Roman   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr);
24*d24d4204SJose E. Roman   ierr = PetscOptionsGetInt(NULL,NULL,"-mb",&mb,NULL);CHKERRQ(ierr);
25*d24d4204SJose E. Roman   nb   = mb;
26*d24d4204SJose E. Roman   ierr = PetscOptionsGetInt(NULL,NULL,"-nb",&nb,NULL);CHKERRQ(ierr);
27*d24d4204SJose E. Roman   ierr = PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);CHKERRQ(ierr);
28*d24d4204SJose E. Roman 
29*d24d4204SJose E. Roman   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
30*d24d4204SJose E. Roman   ierr = MatSetType(C,MATSCALAPACK);CHKERRQ(ierr);
31*d24d4204SJose E. Roman   mloc = PETSC_DECIDE;
32*d24d4204SJose E. Roman   ierr = PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M);CHKERRQ(ierr);
33*d24d4204SJose E. Roman   nloc = PETSC_DECIDE;
34*d24d4204SJose E. Roman   ierr = PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N);CHKERRQ(ierr);
35*d24d4204SJose E. Roman   ierr = MatSetSizes(C,mloc,nloc,M,N);CHKERRQ(ierr);
36*d24d4204SJose E. Roman   ierr = MatScaLAPACKSetBlockSizes(C,mb,nb);CHKERRQ(ierr);
37*d24d4204SJose E. Roman   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
38*d24d4204SJose E. Roman   ierr = MatSetUp(C);CHKERRQ(ierr);
39*d24d4204SJose E. Roman   /*ierr = MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&C);CHKERRQ(ierr); */
40*d24d4204SJose E. Roman 
41*d24d4204SJose E. Roman   ierr = MatGetOwnershipIS(C,&isrows,&iscols);CHKERRQ(ierr);
42*d24d4204SJose E. Roman   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
43*d24d4204SJose E. Roman   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
44*d24d4204SJose E. Roman   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
45*d24d4204SJose E. Roman   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
46*d24d4204SJose E. Roman   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
47*d24d4204SJose E. Roman   for (i=0;i<nrows;i++) {
48*d24d4204SJose E. Roman     for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(rows[i]+1+(cols[j]+1)*0.01);
49*d24d4204SJose E. Roman   }
50*d24d4204SJose E. Roman   ierr = MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
51*d24d4204SJose E. Roman   ierr = PetscFree(v);CHKERRQ(ierr);
52*d24d4204SJose E. Roman   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
53*d24d4204SJose E. Roman   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
54*d24d4204SJose E. Roman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55*d24d4204SJose E. Roman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56*d24d4204SJose E. Roman   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
57*d24d4204SJose E. Roman   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
58*d24d4204SJose E. Roman 
59*d24d4204SJose E. Roman   /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */
60*d24d4204SJose E. Roman   ierr = MatDuplicate(C,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
61*d24d4204SJose E. Roman   if (mats_view) {
62*d24d4204SJose E. Roman     ierr = PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n");CHKERRQ(ierr);
63*d24d4204SJose E. Roman     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
64*d24d4204SJose E. Roman   }
65*d24d4204SJose E. Roman   ierr = MatDestroy(&B);CHKERRQ(ierr);
66*d24d4204SJose E. Roman   ierr = MatConvert(C,MATDENSE,MAT_INITIAL_MATRIX,&Cdense);CHKERRQ(ierr);
67*d24d4204SJose E. Roman   ierr = MatMultEqual(C,Cdense,5,&flg);CHKERRQ(ierr);
68*d24d4204SJose E. Roman   if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Check fails: Cdense != C");
69*d24d4204SJose E. Roman 
70*d24d4204SJose E. Roman   /* Test MatNorm() */
71*d24d4204SJose E. Roman   ierr = MatNorm(C,NORM_1,&Cnorm);CHKERRQ(ierr);
72*d24d4204SJose E. Roman 
73*d24d4204SJose E. Roman   /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */
74*d24d4204SJose E. Roman   ierr = MatTranspose(C,MAT_INITIAL_MATRIX,&Ct);CHKERRQ(ierr);
75*d24d4204SJose E. Roman   ierr = MatConjugate(Ct);CHKERRQ(ierr);
76*d24d4204SJose E. Roman   if (mats_view) {
77*d24d4204SJose E. Roman     ierr = PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n");CHKERRQ(ierr);
78*d24d4204SJose E. Roman     ierr = MatView(Ct,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
79*d24d4204SJose E. Roman   }
80*d24d4204SJose E. Roman   ierr = MatZeroEntries(Ct);CHKERRQ(ierr);
81*d24d4204SJose E. Roman   if (M>N) { ierr = MatCreateVecs(C,&d,NULL);CHKERRQ(ierr); }
82*d24d4204SJose E. Roman   else { ierr = MatCreateVecs(C,NULL,&d);CHKERRQ(ierr); }
83*d24d4204SJose E. Roman   ierr = MatGetDiagonal(C,d);CHKERRQ(ierr);
84*d24d4204SJose E. Roman   if (mats_view) {
85*d24d4204SJose E. Roman     ierr = PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n");CHKERRQ(ierr);
86*d24d4204SJose E. Roman     ierr = VecView(d,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
87*d24d4204SJose E. Roman   }
88*d24d4204SJose E. Roman   if (M>N) {
89*d24d4204SJose E. Roman     ierr = MatDiagonalScale(C,NULL,d);CHKERRQ(ierr);
90*d24d4204SJose E. Roman   } else {
91*d24d4204SJose E. Roman     ierr = MatDiagonalScale(C,d,NULL);CHKERRQ(ierr);
92*d24d4204SJose E. Roman   }
93*d24d4204SJose E. Roman   if (mats_view) {
94*d24d4204SJose E. Roman     ierr = PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n");CHKERRQ(ierr);
95*d24d4204SJose E. Roman     ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
96*d24d4204SJose E. Roman   }
97*d24d4204SJose E. Roman 
98*d24d4204SJose E. Roman   /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */
99*d24d4204SJose E. Roman   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
100*d24d4204SJose E. Roman   ierr = MatSetType(B,MATSCALAPACK);CHKERRQ(ierr);
101*d24d4204SJose E. Roman   ierr = MatSetSizes(B,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
102*d24d4204SJose E. Roman   ierr = MatScaLAPACKSetBlockSizes(B,mb,nb);CHKERRQ(ierr);
103*d24d4204SJose E. Roman   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
104*d24d4204SJose E. Roman   ierr = MatSetUp(B);CHKERRQ(ierr);
105*d24d4204SJose E. Roman   /* ierr = MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&B);CHKERRQ(ierr); */
106*d24d4204SJose E. Roman   ierr = MatGetOwnershipIS(B,&isrows,&iscols);CHKERRQ(ierr);
107*d24d4204SJose E. Roman   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
108*d24d4204SJose E. Roman   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
109*d24d4204SJose E. Roman   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
110*d24d4204SJose E. Roman   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
111*d24d4204SJose E. Roman   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
112*d24d4204SJose E. Roman   for (i=0;i<nrows;i++) {
113*d24d4204SJose E. Roman     for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(1000*rows[i]+cols[j]);
114*d24d4204SJose E. Roman   }
115*d24d4204SJose E. Roman   ierr = MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
116*d24d4204SJose E. Roman   ierr = PetscFree(v);CHKERRQ(ierr);
117*d24d4204SJose E. Roman   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
118*d24d4204SJose E. Roman   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
119*d24d4204SJose E. Roman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
120*d24d4204SJose E. Roman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
121*d24d4204SJose E. Roman   if (mats_view) {
122*d24d4204SJose E. Roman     ierr = PetscPrintf(PETSC_COMM_WORLD,"B:\n");CHKERRQ(ierr);
123*d24d4204SJose E. Roman     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
124*d24d4204SJose E. Roman   }
125*d24d4204SJose E. Roman   ierr = MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
126*d24d4204SJose E. Roman   ierr = MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
127*d24d4204SJose E. Roman   ierr = MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr);
128*d24d4204SJose E. Roman   if (mats_view) {
129*d24d4204SJose E. Roman     ierr = PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n");CHKERRQ(ierr);
130*d24d4204SJose E. Roman     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
131*d24d4204SJose E. Roman   }
132*d24d4204SJose E. Roman   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
133*d24d4204SJose E. Roman   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
134*d24d4204SJose E. Roman   ierr = MatDestroy(&B);CHKERRQ(ierr);
135*d24d4204SJose E. Roman 
136*d24d4204SJose E. Roman   /* Test MatMatTransposeMult(): B = C*C^T */
137*d24d4204SJose E. Roman   ierr = MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
138*d24d4204SJose E. Roman   ierr = MatScale(C,2.0);CHKERRQ(ierr);
139*d24d4204SJose E. Roman   ierr = MatMatTransposeMult(C,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
140*d24d4204SJose E. Roman   ierr = MatMatTransposeMultEqual(C,C,B,10,&flg);CHKERRQ(ierr);
141*d24d4204SJose E. Roman   if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Chech fails: B != C*C^T");
142*d24d4204SJose E. Roman 
143*d24d4204SJose E. Roman   if (mats_view) {
144*d24d4204SJose E. Roman     ierr = PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n");CHKERRQ(ierr);
145*d24d4204SJose E. Roman     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
146*d24d4204SJose E. Roman   }
147*d24d4204SJose E. Roman 
148*d24d4204SJose E. Roman   /* Test MatMult() */
149*d24d4204SJose E. Roman   ierr = MatComputeOperator(C,MATAIJ,&Caij);CHKERRQ(ierr);
150*d24d4204SJose E. Roman   ierr = MatMultEqual(C,Caij,5,&flg);CHKERRQ(ierr);
151*d24d4204SJose E. Roman   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultEqual() fails");
152*d24d4204SJose E. Roman   ierr = MatMultTransposeEqual(C,Caij,5,&flg);CHKERRQ(ierr);
153*d24d4204SJose E. Roman   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeEqual() fails");
154*d24d4204SJose E. Roman 
155*d24d4204SJose E. Roman   /* Test MatMultAdd() and MatMultTransposeAddEqual() */
156*d24d4204SJose E. Roman   ierr = MatMultAddEqual(C,Caij,5,&flg);CHKERRQ(ierr);
157*d24d4204SJose E. Roman   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultAddEqual() fails");
158*d24d4204SJose E. Roman   ierr = MatMultTransposeAddEqual(C,Caij,5,&flg);CHKERRQ(ierr);
159*d24d4204SJose E. Roman   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeAddEqual() fails");
160*d24d4204SJose E. Roman 
161*d24d4204SJose E. Roman   /* Test MatMatMult() */
162*d24d4204SJose E. Roman   ierr = PetscOptionsHasName(NULL,NULL,"-test_matmatmult",&flg);CHKERRQ(ierr);
163*d24d4204SJose E. Roman   if (flg) {
164*d24d4204SJose E. Roman     Mat CC,CCaij;
165*d24d4204SJose E. Roman     ierr = MatMatMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CC);CHKERRQ(ierr);
166*d24d4204SJose E. Roman     ierr = MatMatMult(Caij,Caij,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CCaij);CHKERRQ(ierr);
167*d24d4204SJose E. Roman     ierr = MatMultEqual(CC,CCaij,5,&flg);CHKERRQ(ierr);
168*d24d4204SJose E. Roman     if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"CC != CCaij. MatMatMult() fails");
169*d24d4204SJose E. Roman     ierr = MatDestroy(&CCaij);CHKERRQ(ierr);
170*d24d4204SJose E. Roman     ierr = MatDestroy(&CC);CHKERRQ(ierr);
171*d24d4204SJose E. Roman   }
172*d24d4204SJose E. Roman 
173*d24d4204SJose E. Roman   ierr = MatDestroy(&Cdense);CHKERRQ(ierr);
174*d24d4204SJose E. Roman   ierr = MatDestroy(&Caij);CHKERRQ(ierr);
175*d24d4204SJose E. Roman   ierr = MatDestroy(&B);CHKERRQ(ierr);
176*d24d4204SJose E. Roman   ierr = MatDestroy(&C);CHKERRQ(ierr);
177*d24d4204SJose E. Roman   ierr = MatDestroy(&Ct);CHKERRQ(ierr);
178*d24d4204SJose E. Roman   ierr = VecDestroy(&d);CHKERRQ(ierr);
179*d24d4204SJose E. Roman   ierr = PetscFinalize();
180*d24d4204SJose E. Roman   return ierr;
181*d24d4204SJose E. Roman }
182*d24d4204SJose E. Roman 
183*d24d4204SJose E. Roman /*TEST
184*d24d4204SJose E. Roman 
185*d24d4204SJose E. Roman    build:
186*d24d4204SJose E. Roman       requires: scalapack
187*d24d4204SJose E. Roman 
188*d24d4204SJose E. Roman    test:
189*d24d4204SJose E. Roman       nsize: 2
190*d24d4204SJose E. Roman       args: -mb 5 -nb 5 -M 12 -N 10
191*d24d4204SJose E. Roman       requires: scalapack
192*d24d4204SJose E. Roman 
193*d24d4204SJose E. Roman    test:
194*d24d4204SJose E. Roman       suffix: 2
195*d24d4204SJose E. Roman       nsize: 6
196*d24d4204SJose E. Roman       args: -mb 8 -nb 6 -M 20 -N 50
197*d24d4204SJose E. Roman       requires: scalapack
198*d24d4204SJose E. Roman       output_file: output/ex242_1.out
199*d24d4204SJose E. Roman 
200*d24d4204SJose E. Roman    test:
201*d24d4204SJose E. Roman       suffix: 3
202*d24d4204SJose E. Roman       nsize: 3
203*d24d4204SJose E. Roman       args: -mb 2 -nb 2 -M 20 -N 20 -test_matmatmult
204*d24d4204SJose E. Roman       requires: scalapack
205*d24d4204SJose E. Roman       output_file: output/ex242_1.out
206*d24d4204SJose E. Roman 
207*d24d4204SJose E. Roman TEST*/
208