xref: /petsc/src/mat/tests/ex104.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n";
2*c4762a1bSJed Brown /*
3*c4762a1bSJed Brown  Example:
4*c4762a1bSJed Brown    mpiexec -n <np> ./ex104 -mat_type elemental
5*c4762a1bSJed Brown */
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown #include <petscmat.h>
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown int main(int argc,char **argv)
10*c4762a1bSJed Brown {
11*c4762a1bSJed Brown   Mat            A,B,C,D;
12*c4762a1bSJed Brown   PetscInt       i,M=10,N=5,j,nrows,ncols,am,an,rstart,rend;
13*c4762a1bSJed Brown   PetscErrorCode ierr;
14*c4762a1bSJed Brown   PetscRandom    r;
15*c4762a1bSJed Brown   PetscBool      equal,iselemental;
16*c4762a1bSJed Brown   PetscReal      fill = 1.0;
17*c4762a1bSJed Brown   IS             isrows,iscols;
18*c4762a1bSJed Brown   const PetscInt *rows,*cols;
19*c4762a1bSJed Brown   PetscScalar    *v,rval;
20*c4762a1bSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
21*c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_TRUE;
22*c4762a1bSJed Brown #else
23*c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_FALSE;
24*c4762a1bSJed Brown #endif
25*c4762a1bSJed Brown   PetscMPIInt    size;
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
28*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   /* Set local matrix entries */
41*c4762a1bSJed Brown   ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
47*c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
48*c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
49*c4762a1bSJed Brown       ierr         = PetscRandomGetValue(r,&rval);CHKERRQ(ierr);
50*c4762a1bSJed Brown       v[i*ncols+j] = rval;
51*c4762a1bSJed Brown     }
52*c4762a1bSJed Brown   }
53*c4762a1bSJed Brown   ierr = MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   /* Test MatTranspose() */
63*c4762a1bSJed Brown   ierr = MatCreateTranspose(A,&C);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */
65*c4762a1bSJed Brown   ierr = MatMultEqual(C,B,10,&equal);CHKERRQ(ierr);
66*c4762a1bSJed Brown   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T");
67*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_REUSE_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */
68*c4762a1bSJed Brown   ierr = MatMultEqual(C,B,10,&equal);CHKERRQ(ierr);
69*c4762a1bSJed Brown   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T");
70*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = MatTranspose(B,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = MatMultEqual(C,B,10,&equal);CHKERRQ(ierr);
74*c4762a1bSJed Brown   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T");
75*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   /* Test MatMatMult() */
79*c4762a1bSJed Brown   if (Test_MatMatMult) {
80*c4762a1bSJed Brown #if !defined(PETSC_HAVE_ELEMENTAL)
81*c4762a1bSJed Brown     if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL");
82*c4762a1bSJed Brown #endif
83*c4762a1bSJed Brown     ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */
84*c4762a1bSJed Brown     ierr = MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); /* C = B*A = A^T*A */
85*c4762a1bSJed Brown     ierr = MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
88*c4762a1bSJed Brown     ierr = MatDuplicate(C,MAT_COPY_VALUES,&D);CHKERRQ(ierr);
89*c4762a1bSJed Brown     ierr = MatDestroy(&D);CHKERRQ(ierr);
90*c4762a1bSJed Brown 
91*c4762a1bSJed Brown     /* Test B*A*x = C*x for n random vector x */
92*c4762a1bSJed Brown     ierr = MatMatMultEqual(B,A,C,10,&equal);CHKERRQ(ierr);
93*c4762a1bSJed Brown     if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x");
94*c4762a1bSJed Brown     ierr = MatDestroy(&C);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown     ierr = MatMatMultSymbolic(B,A,fill,&C);CHKERRQ(ierr);
97*c4762a1bSJed Brown     for (i=0; i<2; i++) {
98*c4762a1bSJed Brown       /* Repeat the numeric product to test reuse of the previous symbolic product */
99*c4762a1bSJed Brown       ierr = MatMatMultNumeric(B,A,C);CHKERRQ(ierr);
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown       ierr = MatMatMultEqual(B,A,C,10,&equal);CHKERRQ(ierr);
102*c4762a1bSJed Brown       if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x");
103*c4762a1bSJed Brown     }
104*c4762a1bSJed Brown     ierr = MatDestroy(&C);CHKERRQ(ierr);
105*c4762a1bSJed Brown     ierr = MatDestroy(&B);CHKERRQ(ierr);
106*c4762a1bSJed Brown   }
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown   /* Test MatTransposeMatMult() */
109*c4762a1bSJed Brown   ierr = PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&iselemental);CHKERRQ(ierr);
110*c4762a1bSJed Brown   if (!iselemental) {
111*c4762a1bSJed Brown     ierr = MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); /* D = A^T*A */
112*c4762a1bSJed Brown     ierr = MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr);
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
115*c4762a1bSJed Brown     ierr = MatDuplicate(D,MAT_COPY_VALUES,&C);CHKERRQ(ierr);
116*c4762a1bSJed Brown     ierr = MatDestroy(&C);CHKERRQ(ierr);
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown     ierr = MatTransposeMatMultEqual(A,A,D,10,&equal);CHKERRQ(ierr);
119*c4762a1bSJed Brown     if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
120*c4762a1bSJed Brown     ierr = MatDestroy(&D);CHKERRQ(ierr);
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown     /* Test D*x = A^T*C*A*x, where C is in AIJ format */
123*c4762a1bSJed Brown     ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
124*c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
125*c4762a1bSJed Brown     if (size == 1) {
126*c4762a1bSJed Brown       ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am);CHKERRQ(ierr);
127*c4762a1bSJed Brown     } else {
128*c4762a1bSJed Brown       ierr = MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
129*c4762a1bSJed Brown     }
130*c4762a1bSJed Brown     ierr = MatSetFromOptions(C);CHKERRQ(ierr);
131*c4762a1bSJed Brown     ierr = MatSetUp(C);CHKERRQ(ierr);
132*c4762a1bSJed Brown     ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
133*c4762a1bSJed Brown     v[0] = 1.0;
134*c4762a1bSJed Brown     for (i=rstart; i<rend; i++) {
135*c4762a1bSJed Brown       ierr = MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES);CHKERRQ(ierr);
136*c4762a1bSJed Brown     }
137*c4762a1bSJed Brown     ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
138*c4762a1bSJed Brown     ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown     /* B = C*A, D = A^T*B */
141*c4762a1bSJed Brown     ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
142*c4762a1bSJed Brown     ierr = MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr);
143*c4762a1bSJed Brown     ierr = MatTransposeMatMultEqual(A,B,D,10,&equal);CHKERRQ(ierr);
144*c4762a1bSJed Brown     if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x");
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown     ierr = MatDestroy(&D);CHKERRQ(ierr);
147*c4762a1bSJed Brown     ierr = MatDestroy(&C);CHKERRQ(ierr);
148*c4762a1bSJed Brown     ierr = MatDestroy(&B);CHKERRQ(ierr);
149*c4762a1bSJed Brown   }
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown   /* Test MatMatTransposeMult() */
152*c4762a1bSJed Brown   if (!iselemental) {
153*c4762a1bSJed Brown     PetscReal diff, scale;
154*c4762a1bSJed Brown     PetscInt  am, an, aM, aN;
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown     ierr = MatGetLocalSize(A, &am, &an);CHKERRQ(ierr);
157*c4762a1bSJed Brown     ierr = MatGetSize(A, &aM, &aN);CHKERRQ(ierr);
158*c4762a1bSJed Brown     ierr = MatCreateDense(PetscObjectComm((PetscObject)A),PETSC_DECIDE, an, aM + 10, aN, NULL, &B);CHKERRQ(ierr);
159*c4762a1bSJed Brown     ierr = MatSetRandom(B, NULL);CHKERRQ(ierr);
160*c4762a1bSJed Brown     ierr = MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161*c4762a1bSJed Brown     ierr = MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162*c4762a1bSJed Brown     ierr = MatMatTransposeMult(A,B,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); /* D = A*A^T */
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
165*c4762a1bSJed Brown     ierr = MatDuplicate(D,MAT_COPY_VALUES,&C);CHKERRQ(ierr);
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown     ierr = MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr);
168*c4762a1bSJed Brown     ierr = MatAXPY(C, -1., D, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown     ierr = MatNorm(C, NORM_FROBENIUS, &diff);CHKERRQ(ierr);
171*c4762a1bSJed Brown     ierr = MatNorm(D, NORM_FROBENIUS, &scale);CHKERRQ(ierr);
172*c4762a1bSJed Brown     if (diff > PETSC_SMALL * scale) SETERRQ(PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX");
173*c4762a1bSJed Brown     ierr = MatDestroy(&C);CHKERRQ(ierr);
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown     ierr = MatMatTransposeMultEqual(A,B,D,10,&equal);CHKERRQ(ierr);
176*c4762a1bSJed Brown     if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
177*c4762a1bSJed Brown     ierr = MatDestroy(&D);CHKERRQ(ierr);
178*c4762a1bSJed Brown     ierr = MatDestroy(&B);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   }
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = PetscFree(v);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = PetscFinalize();
185*c4762a1bSJed Brown   return ierr;
186*c4762a1bSJed Brown }
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown /*TEST
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown     test:
191*c4762a1bSJed Brown       output_file: output/ex104.out
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown     test:
194*c4762a1bSJed Brown       suffix: 2
195*c4762a1bSJed Brown       nsize: 2
196*c4762a1bSJed Brown       output_file: output/ex104.out
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown     test:
199*c4762a1bSJed Brown       suffix: 3
200*c4762a1bSJed Brown       nsize: 4
201*c4762a1bSJed Brown       output_file: output/ex104.out
202*c4762a1bSJed Brown       args: -M 23 -N 31
203*c4762a1bSJed Brown 
204*c4762a1bSJed Brown     test:
205*c4762a1bSJed Brown       suffix: 4
206*c4762a1bSJed Brown       nsize: 4
207*c4762a1bSJed Brown       output_file: output/ex104.out
208*c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic
209*c4762a1bSJed Brown 
210*c4762a1bSJed Brown     test:
211*c4762a1bSJed Brown       suffix: 5
212*c4762a1bSJed Brown       nsize: 4
213*c4762a1bSJed Brown       output_file: output/ex104.out
214*c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv
215*c4762a1bSJed Brown 
216*c4762a1bSJed Brown TEST*/
217