xref: /petsc/src/mat/utils/multequal.c (revision 2dcb1b2a648e6798dfef446393ab7396b79ac556)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
286a22c91SHong Zhang 
386a22c91SHong Zhang #include "src/mat/matimpl.h"  /*I   "petscmat.h"  I*/
486a22c91SHong Zhang 
586a22c91SHong Zhang #undef __FUNCT__
686a22c91SHong Zhang #define __FUNCT__ "MatMultEqual"
786a22c91SHong Zhang /*@
886a22c91SHong Zhang    MatMultEqual - Compares matrix-vector products of two matrices.
986a22c91SHong Zhang 
1086a22c91SHong Zhang    Collective on Mat
1186a22c91SHong Zhang 
1286a22c91SHong Zhang    Input Parameters:
1386a22c91SHong Zhang +  A - the first matrix
1486a22c91SHong Zhang -  B - the second matrix
1586a22c91SHong Zhang -  n - number of random vectors to be tested
1686a22c91SHong Zhang 
1786a22c91SHong Zhang    Output Parameter:
1886a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
1986a22c91SHong Zhang 
2086a22c91SHong Zhang    Level: intermediate
2186a22c91SHong Zhang 
2286a22c91SHong Zhang    Concepts: matrices^equality between
2386a22c91SHong Zhang @*/
24be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
2586a22c91SHong Zhang {
2686a22c91SHong Zhang   PetscErrorCode ierr;
2786a22c91SHong Zhang   Vec            x,s1,s2;
2886a22c91SHong Zhang   PetscRandom    rctx;
2986a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
3086a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
31e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
3286a22c91SHong Zhang 
3386a22c91SHong Zhang   PetscFunctionBegin;
344b497936SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
354b497936SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
3686a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
3786a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
3886a22c91SHong 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);
39cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
40cb5d8e9eSHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
41cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
4286a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
4386a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
44cb5d8e9eSHong Zhang 
45cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
46cb5d8e9eSHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
47cb5d8e9eSHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
48cb5d8e9eSHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
4986a22c91SHong Zhang 
504eb6d288SHong Zhang   *flg = PETSC_TRUE;
5186a22c91SHong Zhang   for (k=0; k<n; k++) {
5286a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
5386a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
5486a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
55e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
56e6e6b1bdSHong Zhang     if (r2 < tol){
57e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
58e6e6b1bdSHong Zhang     } else {
59*2dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
60e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
61e6e6b1bdSHong Zhang       r1 /= r2;
62e6e6b1bdSHong Zhang     }
63e6e6b1bdSHong Zhang     if (r1 > tol) {
6486a22c91SHong Zhang       *flg = PETSC_FALSE;
65e6e6b1bdSHong Zhang       ierr = PetscLogInfo((0,"Error: %D-th MatMult() %g\n",k,r1));CHKERRQ(ierr);
664eb6d288SHong Zhang       break;
6786a22c91SHong Zhang     }
6886a22c91SHong Zhang   }
6986a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
7086a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
7186a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
7286a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
7386a22c91SHong Zhang   PetscFunctionReturn(0);
7486a22c91SHong Zhang }
7586a22c91SHong Zhang 
7686a22c91SHong Zhang #undef __FUNCT__
7786a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual"
7886a22c91SHong Zhang /*@
7986a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
8086a22c91SHong Zhang 
8186a22c91SHong Zhang    Collective on Mat
8286a22c91SHong Zhang 
8386a22c91SHong Zhang    Input Parameters:
8486a22c91SHong Zhang +  A - the first matrix
8586a22c91SHong Zhang -  B - the second matrix
8686a22c91SHong Zhang -  n - number of random vectors to be tested
8786a22c91SHong Zhang 
8886a22c91SHong Zhang    Output Parameter:
8986a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
9086a22c91SHong Zhang 
9186a22c91SHong Zhang    Level: intermediate
9286a22c91SHong Zhang 
9386a22c91SHong Zhang    Concepts: matrices^equality between
9486a22c91SHong Zhang @*/
95be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
9686a22c91SHong Zhang {
9786a22c91SHong Zhang   PetscErrorCode ierr;
9886a22c91SHong Zhang   Vec            x,y,s1,s2;
9986a22c91SHong Zhang   PetscRandom    rctx;
10086a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
10186a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
102e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
10386a22c91SHong Zhang 
10486a22c91SHong Zhang   PetscFunctionBegin;
10586a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
10686a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
10786a22c91SHong 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);
108cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
109cb5d8e9eSHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
110cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
11186a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
11286a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
11363879571SHong Zhang 
11463879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
11563879571SHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
11663879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
11763879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
11863879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
11986a22c91SHong Zhang 
1204eb6d288SHong Zhang   *flg = PETSC_TRUE;
12186a22c91SHong Zhang   for (k=0; k<n; k++) {
12286a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
12386a22c91SHong Zhang     ierr = VecSetRandom(rctx,y);CHKERRQ(ierr);
12486a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
12586a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
126e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
127e6e6b1bdSHong Zhang     if (r2 < tol){
128e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
129e6e6b1bdSHong Zhang     } else {
130*2dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
131e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
132e6e6b1bdSHong Zhang       r1 /= r2;
133e6e6b1bdSHong Zhang     }
134e6e6b1bdSHong Zhang     if (r1 > tol) {
13586a22c91SHong Zhang       *flg = PETSC_FALSE;
136e6e6b1bdSHong Zhang       ierr = PetscLogInfo((0,"Error: %d-th MatMultAdd() %g\n",k,r1));CHKERRQ(ierr);
13763879571SHong Zhang       break;
13863879571SHong Zhang     }
13963879571SHong Zhang   }
14063879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
14163879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
14263879571SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
14363879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
14463879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
14563879571SHong Zhang   PetscFunctionReturn(0);
14663879571SHong Zhang }
14763879571SHong Zhang 
14863879571SHong Zhang #undef __FUNCT__
14963879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual"
15063879571SHong Zhang /*@
15163879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
15263879571SHong Zhang 
15363879571SHong Zhang    Collective on Mat
15463879571SHong Zhang 
15563879571SHong Zhang    Input Parameters:
15663879571SHong Zhang +  A - the first matrix
15763879571SHong Zhang -  B - the second matrix
15863879571SHong Zhang -  n - number of random vectors to be tested
15963879571SHong Zhang 
16063879571SHong Zhang    Output Parameter:
16163879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
16263879571SHong Zhang 
16363879571SHong Zhang    Level: intermediate
16463879571SHong Zhang 
16563879571SHong Zhang    Concepts: matrices^equality between
16663879571SHong Zhang @*/
167be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
16863879571SHong Zhang {
16963879571SHong Zhang   PetscErrorCode ierr;
17063879571SHong Zhang   Vec            x,s1,s2;
17163879571SHong Zhang   PetscRandom    rctx;
17263879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
17363879571SHong Zhang   PetscInt       am,an,bm,bn,k;
174e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
17563879571SHong Zhang 
17663879571SHong Zhang   PetscFunctionBegin;
17763879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
17863879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
17963879571SHong 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);
18063879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
18163879571SHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
18263879571SHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
18363879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
18463879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
18563879571SHong Zhang 
18663879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
18763879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
18863879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
18963879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
19063879571SHong Zhang 
19163879571SHong Zhang   *flg = PETSC_TRUE;
19263879571SHong Zhang   for (k=0; k<n; k++) {
19363879571SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
19463879571SHong Zhang     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
19563879571SHong Zhang     ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr);
196e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
197e6e6b1bdSHong Zhang     if (r2 < tol){
198e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
199e6e6b1bdSHong Zhang     } else {
200*2dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
201e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
202e6e6b1bdSHong Zhang       r1 /= r2;
203e6e6b1bdSHong Zhang     }
204e6e6b1bdSHong Zhang     if (r1 > tol) {
20563879571SHong Zhang       *flg = PETSC_FALSE;
206e6e6b1bdSHong Zhang       ierr = PetscLogInfo((0,"Error: %d-th MatMultTranspose() %g\n",k,r1));CHKERRQ(ierr);
20763879571SHong Zhang       break;
20863879571SHong Zhang     }
20963879571SHong Zhang   }
21063879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
21163879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
21263879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
21363879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
21463879571SHong Zhang   PetscFunctionReturn(0);
21563879571SHong Zhang }
21663879571SHong Zhang 
21763879571SHong Zhang #undef __FUNCT__
21863879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual"
21963879571SHong Zhang /*@
22063879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
22163879571SHong Zhang 
22263879571SHong Zhang    Collective on Mat
22363879571SHong Zhang 
22463879571SHong Zhang    Input Parameters:
22563879571SHong Zhang +  A - the first matrix
22663879571SHong Zhang -  B - the second matrix
22763879571SHong Zhang -  n - number of random vectors to be tested
22863879571SHong Zhang 
22963879571SHong Zhang    Output Parameter:
23063879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
23163879571SHong Zhang 
23263879571SHong Zhang    Level: intermediate
23363879571SHong Zhang 
23463879571SHong Zhang    Concepts: matrices^equality between
23563879571SHong Zhang @*/
236be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
23763879571SHong Zhang {
23863879571SHong Zhang   PetscErrorCode ierr;
23963879571SHong Zhang   Vec            x,y,s1,s2;
24063879571SHong Zhang   PetscRandom    rctx;
24163879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
24263879571SHong Zhang   PetscInt       am,an,bm,bn,k;
243e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
24463879571SHong Zhang 
24563879571SHong Zhang   PetscFunctionBegin;
24663879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
24763879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
24863879571SHong 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);
24963879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
25063879571SHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
25163879571SHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
25263879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
25363879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
25463879571SHong Zhang 
25563879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
25663879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
25763879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
25863879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
25963879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
26063879571SHong Zhang 
26163879571SHong Zhang   *flg = PETSC_TRUE;
26263879571SHong Zhang   for (k=0; k<n; k++) {
26363879571SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
26463879571SHong Zhang     ierr = VecSetRandom(rctx,y);CHKERRQ(ierr);
26563879571SHong Zhang     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
26663879571SHong Zhang     ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr);
267e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
268e6e6b1bdSHong Zhang     if (r2 < tol){
269e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
270e6e6b1bdSHong Zhang     } else {
271*2dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
272e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
273e6e6b1bdSHong Zhang       r1 /= r2;
274e6e6b1bdSHong Zhang     }
275e6e6b1bdSHong Zhang     if (r1 > tol) {
27663879571SHong Zhang       *flg = PETSC_FALSE;
277e6e6b1bdSHong Zhang       ierr = PetscLogInfo((0,"Error: %d-th MatMultTransposeAdd() %g\n",k,r1));CHKERRQ(ierr);
2784eb6d288SHong Zhang       break;
27986a22c91SHong Zhang     }
28086a22c91SHong Zhang   }
28186a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
28286a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
28386a22c91SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
28486a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
28586a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
28686a22c91SHong Zhang   PetscFunctionReturn(0);
28786a22c91SHong Zhang }
288