Actual source code: ex104.c
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;
13: PetscRandom r;
14: PetscBool equal, Aiselemental;
15: PetscReal fill = 1.0;
16: IS isrows, iscols;
17: const PetscInt *rows, *cols;
18: PetscScalar *v, rval;
19: #if defined(PETSC_HAVE_ELEMENTAL)
20: PetscBool Test_MatMatMult = PETSC_TRUE;
21: #else
22: PetscBool Test_MatMatMult = PETSC_FALSE;
23: #endif
24: PetscMPIInt size;
27: PetscInitialize(&argc, &argv, (char *)0, help);
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: PetscObjectTypeCompare((PetscObject)A, MATELEMENTAL, &Aiselemental);
64: /* Test MatCreateTranspose() and MatTranspose() */
65: MatCreateTranspose(A, &C);
66: MatTranspose(A, MAT_INITIAL_MATRIX, &B); /* B = A^T */
67: MatMultEqual(C, B, 10, &equal);
69: MatDestroy(&B);
71: MatDuplicate(A, MAT_COPY_VALUES, &B);
72: if (!Aiselemental) {
73: MatTranspose(B, MAT_INPLACE_MATRIX, &B);
74: MatMultEqual(C, B, 10, &equal);
76: }
77: MatDestroy(&B);
79: /* Test B = C*A for matrix type transpose and seqdense */
80: if (size == 1 && !Aiselemental) {
81: MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &B);
82: MatMatMultEqual(C, A, B, 10, &equal);
84: MatDestroy(&B);
85: }
86: MatDestroy(&C);
88: /* Test MatMatMult() */
89: if (Test_MatMatMult) {
90: #if !defined(PETSC_HAVE_ELEMENTAL)
92: #endif
93: MatTranspose(A, MAT_INITIAL_MATRIX, &B); /* B = A^T */
94: MatMatMult(B, A, MAT_INITIAL_MATRIX, fill, &C); /* C = B*A = A^T*A */
95: MatMatMult(B, A, MAT_REUSE_MATRIX, fill, &C);
96: MatMatMultEqual(B, A, C, 10, &equal);
99: /* Test MatDuplicate for matrix product */
100: MatDuplicate(C, MAT_COPY_VALUES, &D);
102: MatDestroy(&D);
103: MatDestroy(&C);
104: MatDestroy(&B);
105: }
107: /* Test MatTransposeMatMult() */
108: if (!Aiselemental) {
109: MatTransposeMatMult(A, A, MAT_INITIAL_MATRIX, fill, &D); /* D = A^T*A */
110: MatTransposeMatMult(A, A, MAT_REUSE_MATRIX, fill, &D);
111: MatTransposeMatMultEqual(A, A, D, 10, &equal);
114: /* Test MatDuplicate for matrix product */
115: MatDuplicate(D, MAT_COPY_VALUES, &C);
116: MatDestroy(&C);
117: MatDestroy(&D);
119: /* Test D*x = A^T*C*A*x, where C is in AIJ format */
120: MatGetLocalSize(A, &am, &an);
121: MatCreate(PETSC_COMM_WORLD, &C);
122: if (size == 1) {
123: MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, am, am);
124: } else {
125: MatSetSizes(C, am, am, PETSC_DECIDE, PETSC_DECIDE);
126: }
127: MatSetFromOptions(C);
128: MatSetUp(C);
129: MatGetOwnershipRange(C, &rstart, &rend);
130: v[0] = 1.0;
131: for (i = rstart; i < rend; i++) MatSetValues(C, 1, &i, 1, &i, v, INSERT_VALUES);
132: MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY);
133: MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY);
135: /* B = C*A, D = A^T*B */
136: MatMatMult(C, A, MAT_INITIAL_MATRIX, 1.0, &B);
137: MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, fill, &D);
138: MatTransposeMatMultEqual(A, B, D, 10, &equal);
141: MatDestroy(&D);
142: MatDestroy(&C);
143: MatDestroy(&B);
144: }
146: /* Test MatMatTransposeMult() */
147: if (!Aiselemental) {
148: PetscReal diff, scale;
149: PetscInt am, an, aM, aN;
151: MatGetLocalSize(A, &am, &an);
152: MatGetSize(A, &aM, &aN);
153: MatCreateDense(PetscObjectComm((PetscObject)A), PETSC_DECIDE, an, aM + 10, aN, NULL, &B);
154: MatSetRandom(B, NULL);
155: MatMatTransposeMult(A, B, MAT_INITIAL_MATRIX, fill, &D); /* D = A*A^T */
157: /* Test MatDuplicate for matrix product */
158: MatDuplicate(D, MAT_COPY_VALUES, &C);
160: MatMatTransposeMult(A, B, MAT_REUSE_MATRIX, fill, &D);
161: MatAXPY(C, -1., D, SAME_NONZERO_PATTERN);
163: MatNorm(C, NORM_FROBENIUS, &diff);
164: MatNorm(D, NORM_FROBENIUS, &scale);
166: MatDestroy(&C);
168: MatMatTransposeMultEqual(A, B, D, 10, &equal);
170: MatDestroy(&D);
171: MatDestroy(&B);
172: }
174: MatDestroy(&A);
175: PetscFree(v);
176: PetscFinalize();
177: return 0;
178: }
180: /*TEST
182: test:
183: output_file: output/ex104.out
185: test:
186: suffix: 2
187: nsize: 2
188: output_file: output/ex104.out
190: test:
191: suffix: 3
192: nsize: 4
193: output_file: output/ex104.out
194: args: -M 23 -N 31
196: test:
197: suffix: 4
198: nsize: 4
199: output_file: output/ex104.out
200: args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic
202: test:
203: suffix: 5
204: nsize: 4
205: output_file: output/ex104.out
206: args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv
208: test:
209: suffix: 6
210: args: -mat_type elemental
211: requires: elemental
212: output_file: output/ex104.out
214: test:
215: suffix: 7
216: nsize: 2
217: args: -mat_type elemental
218: requires: elemental
219: output_file: output/ex104.out
221: TEST*/