xref: /petsc/src/mat/utils/multequal.c (revision 57622a8ec272e2944ba90348a8875f087fbb2bc6)
186a22c91SHong Zhang 
2b45d2f2cSJed Brown #include <petsc-private/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 @*/
237087cfbeSBarry Smith PetscErrorCode  MatMultEqual(Mat A,Mat B,PetscInt n,PetscBool  *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;
30e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
3186a22c91SHong Zhang 
3286a22c91SHong Zhang   PetscFunctionBegin;
330700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
340700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,2);
3586a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
3686a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
37e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
38cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
394fcfbe5fSBarry Smith 
40ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
414fcfbe5fSBarry Smith   tol = 1.e-5;
424fcfbe5fSBarry Smith #endif
43ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
44c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
45c0dedaeaSBarry Smith   ierr = MatGetVecs(A,&x,&s1);CHKERRQ(ierr);
46cb5d8e9eSHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
4786a22c91SHong Zhang 
484eb6d288SHong Zhang   *flg = PETSC_TRUE;
4986a22c91SHong Zhang   for (k=0; k<n; k++) {
50abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
5186a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
5286a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
53e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
54e6e6b1bdSHong Zhang     if (r2 < tol) {
55e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
56e6e6b1bdSHong Zhang     } else {
572dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
58e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
59e6e6b1bdSHong Zhang       r1  /= r2;
60e6e6b1bdSHong Zhang     }
61e6e6b1bdSHong Zhang     if (r1 > tol) {
6286a22c91SHong Zhang       *flg = PETSC_FALSE;
63*57622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %D-th MatMult() %g\n",k,(double)r1);CHKERRQ(ierr);
644eb6d288SHong Zhang       break;
6586a22c91SHong Zhang     }
6686a22c91SHong Zhang   }
676bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
686bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
696bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
706bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
7186a22c91SHong Zhang   PetscFunctionReturn(0);
7286a22c91SHong Zhang }
7386a22c91SHong Zhang 
7486a22c91SHong Zhang #undef __FUNCT__
7586a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual"
7686a22c91SHong Zhang /*@
7786a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
7886a22c91SHong Zhang 
7986a22c91SHong Zhang    Collective on Mat
8086a22c91SHong Zhang 
8186a22c91SHong Zhang    Input Parameters:
8286a22c91SHong Zhang +  A - the first matrix
8386a22c91SHong Zhang -  B - the second matrix
8486a22c91SHong Zhang -  n - number of random vectors to be tested
8586a22c91SHong Zhang 
8686a22c91SHong Zhang    Output Parameter:
8786a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
8886a22c91SHong Zhang 
8986a22c91SHong Zhang    Level: intermediate
9086a22c91SHong Zhang 
9186a22c91SHong Zhang    Concepts: matrices^equality between
9286a22c91SHong Zhang @*/
937087cfbeSBarry Smith PetscErrorCode  MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscBool  *flg)
9486a22c91SHong Zhang {
9586a22c91SHong Zhang   PetscErrorCode ierr;
9686a22c91SHong Zhang   Vec            x,y,s1,s2;
9786a22c91SHong Zhang   PetscRandom    rctx;
9886a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
9986a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
100e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
10186a22c91SHong Zhang 
10286a22c91SHong Zhang   PetscFunctionBegin;
10386a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
10486a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
105e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
106cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
107ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
108c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
109c0dedaeaSBarry Smith   ierr = MatGetVecs(A,&x,&s1);CHKERRQ(ierr);
11063879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
11163879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
11286a22c91SHong Zhang 
1134eb6d288SHong Zhang   *flg = PETSC_TRUE;
11486a22c91SHong Zhang   for (k=0; k<n; k++) {
115abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
116abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
11786a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
11886a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
119e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
120e6e6b1bdSHong Zhang     if (r2 < tol) {
121e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
122e6e6b1bdSHong Zhang     } else {
1232dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
124e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
125e6e6b1bdSHong Zhang       r1  /= r2;
126e6e6b1bdSHong Zhang     }
127e6e6b1bdSHong Zhang     if (r1 > tol) {
12886a22c91SHong Zhang       *flg = PETSC_FALSE;
129*57622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultAdd() %g\n",k,(double)r1);CHKERRQ(ierr);
13063879571SHong Zhang       break;
13163879571SHong Zhang     }
13263879571SHong Zhang   }
1336bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
1346bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1356bf464f9SBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
1366bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
1376bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
13863879571SHong Zhang   PetscFunctionReturn(0);
13963879571SHong Zhang }
14063879571SHong Zhang 
14163879571SHong Zhang #undef __FUNCT__
14263879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual"
14363879571SHong Zhang /*@
14463879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
14563879571SHong Zhang 
14663879571SHong Zhang    Collective on Mat
14763879571SHong Zhang 
14863879571SHong Zhang    Input Parameters:
14963879571SHong Zhang +  A - the first matrix
15063879571SHong Zhang -  B - the second matrix
15163879571SHong Zhang -  n - number of random vectors to be tested
15263879571SHong Zhang 
15363879571SHong Zhang    Output Parameter:
15463879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
15563879571SHong Zhang 
15663879571SHong Zhang    Level: intermediate
15763879571SHong Zhang 
15863879571SHong Zhang    Concepts: matrices^equality between
15963879571SHong Zhang @*/
1607087cfbeSBarry Smith PetscErrorCode  MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscBool  *flg)
16163879571SHong Zhang {
16263879571SHong Zhang   PetscErrorCode ierr;
16363879571SHong Zhang   Vec            x,s1,s2;
16463879571SHong Zhang   PetscRandom    rctx;
16563879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
16663879571SHong Zhang   PetscInt       am,an,bm,bn,k;
167e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
16863879571SHong Zhang 
16963879571SHong Zhang   PetscFunctionBegin;
17063879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
17163879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
172e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
17363879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
174ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
175c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
176c0dedaeaSBarry Smith   ierr = MatGetVecs(A,&s1,&x);CHKERRQ(ierr);
17763879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
17863879571SHong Zhang 
17963879571SHong Zhang   *flg = PETSC_TRUE;
18063879571SHong Zhang   for (k=0; k<n; k++) {
181abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
18263879571SHong Zhang     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
18363879571SHong Zhang     ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr);
184e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
185e6e6b1bdSHong Zhang     if (r2 < tol) {
186e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
187e6e6b1bdSHong Zhang     } else {
1882dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
189e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
190e6e6b1bdSHong Zhang       r1  /= r2;
191e6e6b1bdSHong Zhang     }
192e6e6b1bdSHong Zhang     if (r1 > tol) {
19363879571SHong Zhang       *flg = PETSC_FALSE;
194*57622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTranspose() %g\n",k,(double)r1);CHKERRQ(ierr);
19563879571SHong Zhang       break;
19663879571SHong Zhang     }
19763879571SHong Zhang   }
1986bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
1996bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
2006bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
2016bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
20263879571SHong Zhang   PetscFunctionReturn(0);
20363879571SHong Zhang }
20463879571SHong Zhang 
20563879571SHong Zhang #undef __FUNCT__
20663879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual"
20763879571SHong Zhang /*@
20863879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
20963879571SHong Zhang 
21063879571SHong Zhang    Collective on Mat
21163879571SHong Zhang 
21263879571SHong Zhang    Input Parameters:
21363879571SHong Zhang +  A - the first matrix
21463879571SHong Zhang -  B - the second matrix
21563879571SHong Zhang -  n - number of random vectors to be tested
21663879571SHong Zhang 
21763879571SHong Zhang    Output Parameter:
21863879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
21963879571SHong Zhang 
22063879571SHong Zhang    Level: intermediate
22163879571SHong Zhang 
22263879571SHong Zhang    Concepts: matrices^equality between
22363879571SHong Zhang @*/
2247087cfbeSBarry Smith PetscErrorCode  MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscBool  *flg)
22563879571SHong Zhang {
22663879571SHong Zhang   PetscErrorCode ierr;
22763879571SHong Zhang   Vec            x,y,s1,s2;
22863879571SHong Zhang   PetscRandom    rctx;
22963879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
23063879571SHong Zhang   PetscInt       am,an,bm,bn,k;
231e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
23263879571SHong Zhang 
23363879571SHong Zhang   PetscFunctionBegin;
23463879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
23563879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
236e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
23763879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
238ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
239c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
240c0dedaeaSBarry Smith   ierr = MatGetVecs(A,&s1,&x);CHKERRQ(ierr);
24163879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
24263879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
24363879571SHong Zhang 
24463879571SHong Zhang   *flg = PETSC_TRUE;
24563879571SHong Zhang   for (k=0; k<n; k++) {
246abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
247abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
24863879571SHong Zhang     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
24963879571SHong Zhang     ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr);
250e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
251e6e6b1bdSHong Zhang     if (r2 < tol) {
252e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
253e6e6b1bdSHong Zhang     } else {
2542dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
255e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
256e6e6b1bdSHong Zhang       r1  /= r2;
257e6e6b1bdSHong Zhang     }
258e6e6b1bdSHong Zhang     if (r1 > tol) {
25963879571SHong Zhang       *flg = PETSC_FALSE;
260*57622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %g\n",k,(double)r1);CHKERRQ(ierr);
2614eb6d288SHong Zhang       break;
26286a22c91SHong Zhang     }
26386a22c91SHong Zhang   }
2646bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
2656bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
2666bf464f9SBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
2676bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
2686bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
26986a22c91SHong Zhang   PetscFunctionReturn(0);
27086a22c91SHong Zhang }
271