Actual source code: multequal.c
1: #define PETSCMAT_DLL
3: #include private/matimpl.h
7: /*@
8: MatMultEqual - Compares matrix-vector products of two matrices.
10: Collective on Mat
12: Input Parameters:
13: + A - the first matrix
14: - B - the second matrix
15: - n - number of random vectors to be tested
17: Output Parameter:
18: . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
20: Level: intermediate
22: Concepts: matrices^equality between
23: @*/
24: PetscErrorCode MatMultEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
25: {
27: Vec x,s1,s2;
28: PetscRandom rctx;
29: PetscReal r1,r2,tol=1.e-10;
30: PetscInt am,an,bm,bn,k;
31: PetscScalar none = -1.0;
36: MatGetLocalSize(A,&am,&an);
37: MatGetLocalSize(B,&bm,&bn);
38: if (am != bm || an != bn) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
41: #if defined(PETSC_USE_SCALAR_SINGLE)
42: tol = 1.e-5;
43: #endif
44: PetscRandomCreate(((PetscObject)A)->comm,&rctx);
45: PetscRandomSetFromOptions(rctx);
46: VecCreate(((PetscObject)A)->comm,&x);
47: VecSetSizes(x,an,PETSC_DECIDE);
48: VecSetFromOptions(x);
49:
50: VecCreate(((PetscObject)A)->comm,&s1);
51: VecSetSizes(s1,am,PETSC_DECIDE);
52: VecSetFromOptions(s1);
53: VecDuplicate(s1,&s2);
54:
55: *flg = PETSC_TRUE;
56: for (k=0; k<n; k++) {
57: VecSetRandom(x,rctx);
58: MatMult(A,x,s1);
59: MatMult(B,x,s2);
60: VecNorm(s2,NORM_INFINITY,&r2);
61: if (r2 < tol){
62: VecNorm(s1,NORM_INFINITY,&r1);
63: } else {
64: VecAXPY(s2,none,s1);
65: VecNorm(s2,NORM_INFINITY,&r1);
66: r1 /= r2;
67: }
68: if (r1 > tol) {
69: *flg = PETSC_FALSE;
70: PetscInfo2(A,"Error: %D-th MatMult() %G\n",k,r1);
71: break;
72: }
73: }
74: PetscRandomDestroy(rctx);
75: VecDestroy(x);
76: VecDestroy(s1);
77: VecDestroy(s2);
78: return(0);
79: }
83: /*@
84: MatMultAddEqual - Compares matrix-vector products of two matrices.
86: Collective on Mat
88: Input Parameters:
89: + A - the first matrix
90: - B - the second matrix
91: - n - number of random vectors to be tested
93: Output Parameter:
94: . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
96: Level: intermediate
98: Concepts: matrices^equality between
99: @*/
100: PetscErrorCode MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
101: {
103: Vec x,y,s1,s2;
104: PetscRandom rctx;
105: PetscReal r1,r2,tol=1.e-10;
106: PetscInt am,an,bm,bn,k;
107: PetscScalar none = -1.0;
110: MatGetLocalSize(A,&am,&an);
111: MatGetLocalSize(B,&bm,&bn);
112: if (am != bm || an != bn) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
114: PetscRandomCreate(((PetscObject)A)->comm,&rctx);
115: PetscRandomSetFromOptions(rctx);
116: VecCreate(((PetscObject)A)->comm,&x);
117: VecSetSizes(x,an,PETSC_DECIDE);
118: VecSetFromOptions(x);
120: VecCreate(((PetscObject)A)->comm,&s1);
121: VecSetSizes(s1,am,PETSC_DECIDE);
122: VecSetFromOptions(s1);
123: VecDuplicate(s1,&s2);
124: VecDuplicate(s1,&y);
125:
126: *flg = PETSC_TRUE;
127: for (k=0; k<n; k++) {
128: VecSetRandom(x,rctx);
129: VecSetRandom(y,rctx);
130: MatMultAdd(A,x,y,s1);
131: MatMultAdd(B,x,y,s2);
132: VecNorm(s2,NORM_INFINITY,&r2);
133: if (r2 < tol){
134: VecNorm(s1,NORM_INFINITY,&r1);
135: } else {
136: VecAXPY(s2,none,s1);
137: VecNorm(s2,NORM_INFINITY,&r1);
138: r1 /= r2;
139: }
140: if (r1 > tol) {
141: *flg = PETSC_FALSE;
142: PetscInfo2(A,"Error: %d-th MatMultAdd() %G\n",k,r1);
143: break;
144: }
145: }
146: PetscRandomDestroy(rctx);
147: VecDestroy(x);
148: VecDestroy(y);
149: VecDestroy(s1);
150: VecDestroy(s2);
151: return(0);
152: }
156: /*@
157: MatMultTransposeEqual - Compares matrix-vector products of two matrices.
159: Collective on Mat
161: Input Parameters:
162: + A - the first matrix
163: - B - the second matrix
164: - n - number of random vectors to be tested
166: Output Parameter:
167: . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
169: Level: intermediate
171: Concepts: matrices^equality between
172: @*/
173: PetscErrorCode MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
174: {
176: Vec x,s1,s2;
177: PetscRandom rctx;
178: PetscReal r1,r2,tol=1.e-10;
179: PetscInt am,an,bm,bn,k;
180: PetscScalar none = -1.0;
183: MatGetLocalSize(A,&am,&an);
184: MatGetLocalSize(B,&bm,&bn);
185: if (am != bm || an != bn) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
187: PetscRandomCreate(((PetscObject)A)->comm,&rctx);
188: PetscRandomSetFromOptions(rctx);
189: VecCreate(((PetscObject)A)->comm,&x);
190: VecSetSizes(x,am,PETSC_DECIDE);
191: VecSetFromOptions(x);
192:
193: VecCreate(((PetscObject)A)->comm,&s1);
194: VecSetSizes(s1,an,PETSC_DECIDE);
195: VecSetFromOptions(s1);
196: VecDuplicate(s1,&s2);
197:
198: *flg = PETSC_TRUE;
199: for (k=0; k<n; k++) {
200: VecSetRandom(x,rctx);
201: MatMultTranspose(A,x,s1);
202: MatMultTranspose(B,x,s2);
203: VecNorm(s2,NORM_INFINITY,&r2);
204: if (r2 < tol){
205: VecNorm(s1,NORM_INFINITY,&r1);
206: } else {
207: VecAXPY(s2,none,s1);
208: VecNorm(s2,NORM_INFINITY,&r1);
209: r1 /= r2;
210: }
211: if (r1 > tol) {
212: *flg = PETSC_FALSE;
213: PetscInfo2(A,"Error: %d-th MatMultTranspose() %G\n",k,r1);
214: break;
215: }
216: }
217: PetscRandomDestroy(rctx);
218: VecDestroy(x);
219: VecDestroy(s1);
220: VecDestroy(s2);
221: return(0);
222: }
226: /*@
227: MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
229: Collective on Mat
231: Input Parameters:
232: + A - the first matrix
233: - B - the second matrix
234: - n - number of random vectors to be tested
236: Output Parameter:
237: . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
239: Level: intermediate
241: Concepts: matrices^equality between
242: @*/
243: PetscErrorCode MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
244: {
246: Vec x,y,s1,s2;
247: PetscRandom rctx;
248: PetscReal r1,r2,tol=1.e-10;
249: PetscInt am,an,bm,bn,k;
250: PetscScalar none = -1.0;
253: MatGetLocalSize(A,&am,&an);
254: MatGetLocalSize(B,&bm,&bn);
255: if (am != bm || an != bn) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
257: PetscRandomCreate(((PetscObject)A)->comm,&rctx);
258: PetscRandomSetFromOptions(rctx);
259: VecCreate(((PetscObject)A)->comm,&x);
260: VecSetSizes(x,am,PETSC_DECIDE);
261: VecSetFromOptions(x);
263: VecCreate(((PetscObject)A)->comm,&s1);
264: VecSetSizes(s1,an,PETSC_DECIDE);
265: VecSetFromOptions(s1);
266: VecDuplicate(s1,&s2);
267: VecDuplicate(s1,&y);
268:
269: *flg = PETSC_TRUE;
270: for (k=0; k<n; k++) {
271: VecSetRandom(x,rctx);
272: VecSetRandom(y,rctx);
273: MatMultTransposeAdd(A,x,y,s1);
274: MatMultTransposeAdd(B,x,y,s2);
275: VecNorm(s2,NORM_INFINITY,&r2);
276: if (r2 < tol){
277: VecNorm(s1,NORM_INFINITY,&r1);
278: } else {
279: VecAXPY(s2,none,s1);
280: VecNorm(s2,NORM_INFINITY,&r1);
281: r1 /= r2;
282: }
283: if (r1 > tol) {
284: *flg = PETSC_FALSE;
285: PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %G\n",k,r1);
286: break;
287: }
288: }
289: PetscRandomDestroy(rctx);
290: VecDestroy(x);
291: VecDestroy(y);
292: VecDestroy(s1);
293: VecDestroy(s2);
294: return(0);
295: }