xref: /petsc/src/mat/utils/multequal.c (revision 63879571fdc812aab9a52703ff2eccb2f9c1b75e)
186a22c91SHong Zhang 
286a22c91SHong Zhang #include "src/mat/matimpl.h"  /*I   "petscmat.h"  I*/
386a22c91SHong Zhang 
486a22c91SHong Zhang #undef __FUNCT__
586a22c91SHong Zhang #define __FUNCT__ "MatMultEqual"
686a22c91SHong Zhang /*@
786a22c91SHong Zhang    MatMultEqual - Compares matrix-vector products of two matrices.
886a22c91SHong Zhang 
986a22c91SHong Zhang    Collective on Mat
1086a22c91SHong Zhang 
1186a22c91SHong Zhang    Input Parameters:
1286a22c91SHong Zhang +  A - the first matrix
1386a22c91SHong Zhang -  B - the second matrix
1486a22c91SHong Zhang -  n - number of random vectors to be tested
1586a22c91SHong Zhang 
1686a22c91SHong Zhang    Output Parameter:
1786a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
1886a22c91SHong Zhang 
1986a22c91SHong Zhang    Level: intermediate
2086a22c91SHong Zhang 
2186a22c91SHong Zhang    Concepts: matrices^equality between
2286a22c91SHong Zhang @*/
2386a22c91SHong Zhang PetscErrorCode MatMultEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
2486a22c91SHong Zhang {
2586a22c91SHong Zhang   PetscErrorCode ierr;
2686a22c91SHong Zhang   Vec            x,s1,s2;
2786a22c91SHong Zhang   PetscRandom    rctx;
2886a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
2986a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
3086a22c91SHong Zhang 
3186a22c91SHong Zhang   PetscFunctionBegin;
3286a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
3386a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
3486a22c91SHong 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);
35cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
36cb5d8e9eSHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
37cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
3886a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
3986a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
40cb5d8e9eSHong Zhang 
41cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
42cb5d8e9eSHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
43cb5d8e9eSHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
44cb5d8e9eSHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
4586a22c91SHong Zhang 
464eb6d288SHong Zhang   *flg = PETSC_TRUE;
4786a22c91SHong Zhang   for (k=0; k<n; k++) {
4886a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
4986a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
5086a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
5186a22c91SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
5286a22c91SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
5386a22c91SHong Zhang     r1 -= r2;
5486a22c91SHong Zhang     if (r1<-tol || r1>tol) {
5586a22c91SHong Zhang       *flg = PETSC_FALSE;
56*63879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMult() %g\n",k,r1);
574eb6d288SHong Zhang       break;
5886a22c91SHong Zhang     }
5986a22c91SHong Zhang   }
6086a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
6186a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
6286a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
6386a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
6486a22c91SHong Zhang   PetscFunctionReturn(0);
6586a22c91SHong Zhang }
6686a22c91SHong Zhang 
6786a22c91SHong Zhang #undef __FUNCT__
6886a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual"
6986a22c91SHong Zhang /*@
7086a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
7186a22c91SHong Zhang 
7286a22c91SHong Zhang    Collective on Mat
7386a22c91SHong Zhang 
7486a22c91SHong Zhang    Input Parameters:
7586a22c91SHong Zhang +  A - the first matrix
7686a22c91SHong Zhang -  B - the second matrix
7786a22c91SHong Zhang -  n - number of random vectors to be tested
7886a22c91SHong Zhang 
7986a22c91SHong Zhang    Output Parameter:
8086a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
8186a22c91SHong Zhang 
8286a22c91SHong Zhang    Level: intermediate
8386a22c91SHong Zhang 
8486a22c91SHong Zhang    Concepts: matrices^equality between
8586a22c91SHong Zhang @*/
8686a22c91SHong Zhang PetscErrorCode MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
8786a22c91SHong Zhang {
8886a22c91SHong Zhang   PetscErrorCode ierr;
8986a22c91SHong Zhang   Vec            x,y,s1,s2;
9086a22c91SHong Zhang   PetscRandom    rctx;
9186a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
9286a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
9386a22c91SHong Zhang 
9486a22c91SHong Zhang   PetscFunctionBegin;
9586a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
9686a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
9786a22c91SHong 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);
98cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
99cb5d8e9eSHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
100cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
10186a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
10286a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
103*63879571SHong Zhang 
104*63879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
105*63879571SHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
106*63879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
107*63879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
108*63879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
10986a22c91SHong Zhang 
1104eb6d288SHong Zhang   *flg = PETSC_TRUE;
11186a22c91SHong Zhang   for (k=0; k<n; k++) {
11286a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
11386a22c91SHong Zhang     ierr = VecSetRandom(rctx,y);CHKERRQ(ierr);
11486a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
11586a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
11686a22c91SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
11786a22c91SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
11886a22c91SHong Zhang     r1 -= r2;
11986a22c91SHong Zhang     if (r1<-tol || r1>tol) {
12086a22c91SHong Zhang       *flg = PETSC_FALSE;
121*63879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMultAdd() %g\n",k,r1);
122*63879571SHong Zhang       break;
123*63879571SHong Zhang     }
124*63879571SHong Zhang   }
125*63879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
126*63879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
127*63879571SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
128*63879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
129*63879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
130*63879571SHong Zhang   PetscFunctionReturn(0);
131*63879571SHong Zhang }
132*63879571SHong Zhang 
133*63879571SHong Zhang #undef __FUNCT__
134*63879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual"
135*63879571SHong Zhang /*@
136*63879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
137*63879571SHong Zhang 
138*63879571SHong Zhang    Collective on Mat
139*63879571SHong Zhang 
140*63879571SHong Zhang    Input Parameters:
141*63879571SHong Zhang +  A - the first matrix
142*63879571SHong Zhang -  B - the second matrix
143*63879571SHong Zhang -  n - number of random vectors to be tested
144*63879571SHong Zhang 
145*63879571SHong Zhang    Output Parameter:
146*63879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
147*63879571SHong Zhang 
148*63879571SHong Zhang    Level: intermediate
149*63879571SHong Zhang 
150*63879571SHong Zhang    Concepts: matrices^equality between
151*63879571SHong Zhang @*/
152*63879571SHong Zhang PetscErrorCode MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
153*63879571SHong Zhang {
154*63879571SHong Zhang   PetscErrorCode ierr;
155*63879571SHong Zhang   Vec            x,s1,s2;
156*63879571SHong Zhang   PetscRandom    rctx;
157*63879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
158*63879571SHong Zhang   PetscInt       am,an,bm,bn,k;
159*63879571SHong Zhang 
160*63879571SHong Zhang   PetscFunctionBegin;
161*63879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
162*63879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
163*63879571SHong 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);
164*63879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
165*63879571SHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
166*63879571SHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
167*63879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
168*63879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
169*63879571SHong Zhang 
170*63879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
171*63879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
172*63879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
173*63879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
174*63879571SHong Zhang 
175*63879571SHong Zhang   *flg = PETSC_TRUE;
176*63879571SHong Zhang   for (k=0; k<n; k++) {
177*63879571SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
178*63879571SHong Zhang     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
179*63879571SHong Zhang     ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr);
180*63879571SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
181*63879571SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
182*63879571SHong Zhang     r1 -= r2;
183*63879571SHong Zhang     if (r1<-tol || r1>tol) {
184*63879571SHong Zhang       *flg = PETSC_FALSE;
185*63879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMultTranspose() %g\n",k,r1);
186*63879571SHong Zhang       break;
187*63879571SHong Zhang     }
188*63879571SHong Zhang   }
189*63879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
190*63879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
191*63879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
192*63879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
193*63879571SHong Zhang   PetscFunctionReturn(0);
194*63879571SHong Zhang }
195*63879571SHong Zhang 
196*63879571SHong Zhang #undef __FUNCT__
197*63879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual"
198*63879571SHong Zhang /*@
199*63879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
200*63879571SHong Zhang 
201*63879571SHong Zhang    Collective on Mat
202*63879571SHong Zhang 
203*63879571SHong Zhang    Input Parameters:
204*63879571SHong Zhang +  A - the first matrix
205*63879571SHong Zhang -  B - the second matrix
206*63879571SHong Zhang -  n - number of random vectors to be tested
207*63879571SHong Zhang 
208*63879571SHong Zhang    Output Parameter:
209*63879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
210*63879571SHong Zhang 
211*63879571SHong Zhang    Level: intermediate
212*63879571SHong Zhang 
213*63879571SHong Zhang    Concepts: matrices^equality between
214*63879571SHong Zhang @*/
215*63879571SHong Zhang PetscErrorCode MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
216*63879571SHong Zhang {
217*63879571SHong Zhang   PetscErrorCode ierr;
218*63879571SHong Zhang   Vec            x,y,s1,s2;
219*63879571SHong Zhang   PetscRandom    rctx;
220*63879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
221*63879571SHong Zhang   PetscInt       am,an,bm,bn,k;
222*63879571SHong Zhang 
223*63879571SHong Zhang   PetscFunctionBegin;
224*63879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
225*63879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
226*63879571SHong 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);
227*63879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
228*63879571SHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
229*63879571SHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
230*63879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
231*63879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
232*63879571SHong Zhang 
233*63879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
234*63879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
235*63879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
236*63879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
237*63879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
238*63879571SHong Zhang 
239*63879571SHong Zhang   *flg = PETSC_TRUE;
240*63879571SHong Zhang   for (k=0; k<n; k++) {
241*63879571SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
242*63879571SHong Zhang     ierr = VecSetRandom(rctx,y);CHKERRQ(ierr);
243*63879571SHong Zhang     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
244*63879571SHong Zhang     ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr);
245*63879571SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
246*63879571SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
247*63879571SHong Zhang     r1 -= r2;
248*63879571SHong Zhang     if (r1<-tol || r1>tol) {
249*63879571SHong Zhang       *flg = PETSC_FALSE;
250*63879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMultTransposeAdd() %g\n",k,r1);
2514eb6d288SHong Zhang       break;
25286a22c91SHong Zhang     }
25386a22c91SHong Zhang   }
25486a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
25586a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
25686a22c91SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
25786a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
25886a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
25986a22c91SHong Zhang   PetscFunctionReturn(0);
26086a22c91SHong Zhang }
261