Actual source code: ex104.c
petsc-3.8.4 2018-03-24
1: static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n";
2: /*
3: Example:
4: mpiexec -n <np> ./ex104 -mat_type elemental
5: */
7: #include <petscmat.h>
9: int main(int argc,char **argv)
10: {
11: Mat A,B,C,D;
12: PetscInt i,M=10,N=5,j,nrows,ncols,am,an,rstart,rend;
14: PetscRandom r;
15: PetscBool equal,iselemental;
16: PetscReal fill = 1.0;
17: IS isrows,iscols;
18: const PetscInt *rows,*cols;
19: PetscScalar *v,rval;
20: #if defined(PETSC_HAVE_ELEMENTAL)
21: PetscBool Test_MatMatMult=PETSC_TRUE;
22: #else
23: PetscBool Test_MatMatMult=PETSC_FALSE;
24: #endif
25: PetscMPIInt size;
27: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
28: MPI_Comm_size(PETSC_COMM_WORLD,&size);
30: PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);
31: PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);
32: MatCreate(PETSC_COMM_WORLD,&A);
33: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);
34: MatSetType(A,MATDENSE);
35: MatSetFromOptions(A);
36: MatSetUp(A);
37: PetscRandomCreate(PETSC_COMM_WORLD,&r);
38: PetscRandomSetFromOptions(r);
40: /* Set local matrix entries */
41: MatGetOwnershipIS(A,&isrows,&iscols);
42: ISGetLocalSize(isrows,&nrows);
43: ISGetIndices(isrows,&rows);
44: ISGetLocalSize(iscols,&ncols);
45: ISGetIndices(iscols,&cols);
46: PetscMalloc1(nrows*ncols,&v);
47: for (i=0; i<nrows; i++) {
48: for (j=0; j<ncols; j++) {
49: PetscRandomGetValue(r,&rval);
50: v[i*ncols+j] = rval;
51: }
52: }
53: MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
56: ISRestoreIndices(isrows,&rows);
57: ISRestoreIndices(iscols,&cols);
58: ISDestroy(&isrows);
59: ISDestroy(&iscols);
60: PetscRandomDestroy(&r);
62: /* Test MatMatMult() */
63: if (Test_MatMatMult) {
64: #if !defined(PETSC_HAVE_ELEMENTAL)
65: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL");
66: #endif
67: MatTranspose(A,MAT_INITIAL_MATRIX,&B); /* B = A^T */
68: MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C); /* C = B*A = A^T*A */
69: MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C);
71: /* Test B*A*x = C*x for n random vector x */
72: MatMatMultEqual(B,A,C,10,&equal);
73: if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x");
74: MatDestroy(&C);
76: MatMatMultSymbolic(B,A,fill,&C);
77: for (i=0; i<2; i++) {
78: /* Repeat the numeric product to test reuse of the previous symbolic product */
79: MatMatMultNumeric(B,A,C);
80:
81: MatMatMultEqual(B,A,C,10,&equal);
82: if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x");
83: }
84: MatDestroy(&C);
85: MatDestroy(&B);
86: }
88: /* Test MatTransposeMatMult() */
89: PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&iselemental);
90: if (!iselemental) {
91: MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D); /* D = A^T*A */
92: MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D);
93: /* MatView(D,PETSC_VIEWER_STDOUT_WORLD); */
94: MatTransposeMatMultEqual(A,A,D,10,&equal);
95: if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
96: MatDestroy(&D);
98: /* Test D*x = A^T*C*A*x, where C is in AIJ format */
99: MatGetLocalSize(A,&am,&an);
100: MatCreate(PETSC_COMM_WORLD,&C);
101: if (size == 1) {
102: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am);
103: } else {
104: MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE);
105: }
106: MatSetFromOptions(C);
107: MatSetUp(C);
108: MatGetOwnershipRange(C,&rstart,&rend);
109: v[0] = 1.0;
110: for (i=rstart; i<rend; i++) {
111: MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES);
112: }
113: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
114: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
116: /* B = C*A, D = A^T*B */
117: MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B);
118: MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D);
119: MatTransposeMatMultEqual(A,B,D,10,&equal);
120: if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x");
122: MatDestroy(&D);
123: MatDestroy(&C);
124: MatDestroy(&B);
125: }
127: MatDestroy(&A);
128: PetscFree(v);
129: PetscFinalize();
130: return ierr;
131: }