xref: /petsc/src/mat/utils/multequal.c (revision 86a22c911aeec26a939927b0dea988da890dcc5b)
1*86a22c91SHong Zhang 
2*86a22c91SHong Zhang #include "src/mat/matimpl.h"  /*I   "petscmat.h"  I*/
3*86a22c91SHong Zhang 
4*86a22c91SHong Zhang #undef __FUNCT__
5*86a22c91SHong Zhang #define __FUNCT__ "MatMultEqual"
6*86a22c91SHong Zhang /*@
7*86a22c91SHong Zhang    MatMultEqual - Compares matrix-vector products of two matrices.
8*86a22c91SHong Zhang 
9*86a22c91SHong Zhang    Collective on Mat
10*86a22c91SHong Zhang 
11*86a22c91SHong Zhang    Input Parameters:
12*86a22c91SHong Zhang +  A - the first matrix
13*86a22c91SHong Zhang -  B - the second matrix
14*86a22c91SHong Zhang -  n - number of random vectors to be tested
15*86a22c91SHong Zhang 
16*86a22c91SHong Zhang    Output Parameter:
17*86a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
18*86a22c91SHong Zhang 
19*86a22c91SHong Zhang    Level: intermediate
20*86a22c91SHong Zhang 
21*86a22c91SHong Zhang    Concepts: matrices^equality between
22*86a22c91SHong Zhang @*/
23*86a22c91SHong Zhang PetscErrorCode MatMultEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
24*86a22c91SHong Zhang {
25*86a22c91SHong Zhang   PetscErrorCode ierr;
26*86a22c91SHong Zhang   Vec            x,s1,s2;
27*86a22c91SHong Zhang   PetscRandom    rctx;
28*86a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
29*86a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
30*86a22c91SHong Zhang 
31*86a22c91SHong Zhang   PetscFunctionBegin;
32*86a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
33*86a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
34*86a22c91SHong Zhang   if (am != bm || an != bn) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
35*86a22c91SHong Zhang   ierr = PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
36*86a22c91SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
37*86a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
38*86a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
39*86a22c91SHong Zhang   ierr = VecDuplicate(x,&s1);CHKERRQ(ierr);
40*86a22c91SHong Zhang   ierr = VecDuplicate(x,&s2);CHKERRQ(ierr);
41*86a22c91SHong Zhang 
42*86a22c91SHong Zhang   for (k=0; k<n; k++) {
43*86a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
44*86a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
45*86a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
46*86a22c91SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
47*86a22c91SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
48*86a22c91SHong Zhang     r1 -= r2;
49*86a22c91SHong Zhang     if (r1<-tol || r1>tol) {
50*86a22c91SHong Zhang       *flg = PETSC_FALSE;
51*86a22c91SHong Zhang     } else {
52*86a22c91SHong Zhang       *flg = PETSC_TRUE;
53*86a22c91SHong Zhang     }
54*86a22c91SHong Zhang   }
55*86a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
56*86a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
57*86a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
58*86a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
59*86a22c91SHong Zhang   PetscFunctionReturn(0);
60*86a22c91SHong Zhang }
61*86a22c91SHong Zhang 
62*86a22c91SHong Zhang #undef __FUNCT__
63*86a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual"
64*86a22c91SHong Zhang /*@
65*86a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
66*86a22c91SHong Zhang 
67*86a22c91SHong Zhang    Collective on Mat
68*86a22c91SHong Zhang 
69*86a22c91SHong Zhang    Input Parameters:
70*86a22c91SHong Zhang +  A - the first matrix
71*86a22c91SHong Zhang -  B - the second matrix
72*86a22c91SHong Zhang -  n - number of random vectors to be tested
73*86a22c91SHong Zhang 
74*86a22c91SHong Zhang    Output Parameter:
75*86a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
76*86a22c91SHong Zhang 
77*86a22c91SHong Zhang    Level: intermediate
78*86a22c91SHong Zhang 
79*86a22c91SHong Zhang    Concepts: matrices^equality between
80*86a22c91SHong Zhang @*/
81*86a22c91SHong Zhang PetscErrorCode MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
82*86a22c91SHong Zhang {
83*86a22c91SHong Zhang   PetscErrorCode ierr;
84*86a22c91SHong Zhang   Vec            x,y,s1,s2;
85*86a22c91SHong Zhang   PetscRandom    rctx;
86*86a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
87*86a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
88*86a22c91SHong Zhang 
89*86a22c91SHong Zhang   PetscFunctionBegin;
90*86a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
91*86a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
92*86a22c91SHong Zhang   if (am != bm || an != bn) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
93*86a22c91SHong Zhang   ierr = PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
94*86a22c91SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
95*86a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
96*86a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
97*86a22c91SHong Zhang   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
98*86a22c91SHong Zhang   ierr = VecDuplicate(x,&s1);CHKERRQ(ierr);
99*86a22c91SHong Zhang   ierr = VecDuplicate(x,&s2);CHKERRQ(ierr);
100*86a22c91SHong Zhang 
101*86a22c91SHong Zhang   for (k=0; k<n; k++) {
102*86a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
103*86a22c91SHong Zhang     ierr = VecSetRandom(rctx,y);CHKERRQ(ierr);
104*86a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
105*86a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
106*86a22c91SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
107*86a22c91SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
108*86a22c91SHong Zhang     r1 -= r2;
109*86a22c91SHong Zhang     if (r1<-tol || r1>tol) {
110*86a22c91SHong Zhang       *flg = PETSC_FALSE;
111*86a22c91SHong Zhang     } else {
112*86a22c91SHong Zhang       *flg = PETSC_TRUE;
113*86a22c91SHong Zhang     }
114*86a22c91SHong Zhang   }
115*86a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
116*86a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
117*86a22c91SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
118*86a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
119*86a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
120*86a22c91SHong Zhang   PetscFunctionReturn(0);
121*86a22c91SHong Zhang }
122