xref: /petsc/src/mat/utils/multequal.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
286a22c91SHong Zhang 
3*7c4f633dSBarry Smith #include "private/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);
407adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
41c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
427adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
4386a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
4486a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
45cb5d8e9eSHong Zhang 
467adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
47cb5d8e9eSHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
48cb5d8e9eSHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
49cb5d8e9eSHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
5086a22c91SHong Zhang 
514eb6d288SHong Zhang   *flg = PETSC_TRUE;
5286a22c91SHong Zhang   for (k=0; k<n; k++) {
53abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
5486a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
5586a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
56e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
57e6e6b1bdSHong Zhang     if (r2 < tol){
58e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
59e6e6b1bdSHong Zhang     } else {
602dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
61e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
62e6e6b1bdSHong Zhang       r1 /= r2;
63e6e6b1bdSHong Zhang     }
64e6e6b1bdSHong Zhang     if (r1 > tol) {
6586a22c91SHong Zhang       *flg = PETSC_FALSE;
661e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %D-th MatMult() %G\n",k,r1);CHKERRQ(ierr);
674eb6d288SHong Zhang       break;
6886a22c91SHong Zhang     }
6986a22c91SHong Zhang   }
7086a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
7186a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
7286a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
7386a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
7486a22c91SHong Zhang   PetscFunctionReturn(0);
7586a22c91SHong Zhang }
7686a22c91SHong Zhang 
7786a22c91SHong Zhang #undef __FUNCT__
7886a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual"
7986a22c91SHong Zhang /*@
8086a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
8186a22c91SHong Zhang 
8286a22c91SHong Zhang    Collective on Mat
8386a22c91SHong Zhang 
8486a22c91SHong Zhang    Input Parameters:
8586a22c91SHong Zhang +  A - the first matrix
8686a22c91SHong Zhang -  B - the second matrix
8786a22c91SHong Zhang -  n - number of random vectors to be tested
8886a22c91SHong Zhang 
8986a22c91SHong Zhang    Output Parameter:
9086a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
9186a22c91SHong Zhang 
9286a22c91SHong Zhang    Level: intermediate
9386a22c91SHong Zhang 
9486a22c91SHong Zhang    Concepts: matrices^equality between
9586a22c91SHong Zhang @*/
96be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
9786a22c91SHong Zhang {
9886a22c91SHong Zhang   PetscErrorCode ierr;
9986a22c91SHong Zhang   Vec            x,y,s1,s2;
10086a22c91SHong Zhang   PetscRandom    rctx;
10186a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
10286a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
103e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
10486a22c91SHong Zhang 
10586a22c91SHong Zhang   PetscFunctionBegin;
10686a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
10786a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
10886a22c91SHong 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);
109cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
1107adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
111c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1127adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
11386a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
11486a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
11563879571SHong Zhang 
1167adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
11763879571SHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
11863879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
11963879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
12063879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
12186a22c91SHong Zhang 
1224eb6d288SHong Zhang   *flg = PETSC_TRUE;
12386a22c91SHong Zhang   for (k=0; k<n; k++) {
124abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
125abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
12686a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
12786a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
128e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
129e6e6b1bdSHong Zhang     if (r2 < tol){
130e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
131e6e6b1bdSHong Zhang     } else {
1322dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
133e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
134e6e6b1bdSHong Zhang       r1 /= r2;
135e6e6b1bdSHong Zhang     }
136e6e6b1bdSHong Zhang     if (r1 > tol) {
13786a22c91SHong Zhang       *flg = PETSC_FALSE;
1381e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultAdd() %G\n",k,r1);CHKERRQ(ierr);
13963879571SHong Zhang       break;
14063879571SHong Zhang     }
14163879571SHong Zhang   }
14263879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
14363879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
14463879571SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
14563879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
14663879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
14763879571SHong Zhang   PetscFunctionReturn(0);
14863879571SHong Zhang }
14963879571SHong Zhang 
15063879571SHong Zhang #undef __FUNCT__
15163879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual"
15263879571SHong Zhang /*@
15363879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
15463879571SHong Zhang 
15563879571SHong Zhang    Collective on Mat
15663879571SHong Zhang 
15763879571SHong Zhang    Input Parameters:
15863879571SHong Zhang +  A - the first matrix
15963879571SHong Zhang -  B - the second matrix
16063879571SHong Zhang -  n - number of random vectors to be tested
16163879571SHong Zhang 
16263879571SHong Zhang    Output Parameter:
16363879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
16463879571SHong Zhang 
16563879571SHong Zhang    Level: intermediate
16663879571SHong Zhang 
16763879571SHong Zhang    Concepts: matrices^equality between
16863879571SHong Zhang @*/
169be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
17063879571SHong Zhang {
17163879571SHong Zhang   PetscErrorCode ierr;
17263879571SHong Zhang   Vec            x,s1,s2;
17363879571SHong Zhang   PetscRandom    rctx;
17463879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
17563879571SHong Zhang   PetscInt       am,an,bm,bn,k;
176e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
17763879571SHong Zhang 
17863879571SHong Zhang   PetscFunctionBegin;
17963879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
18063879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
18163879571SHong 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);
18263879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
1837adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
184c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1857adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
18663879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
18763879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
18863879571SHong Zhang 
1897adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
19063879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
19163879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
19263879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
19363879571SHong Zhang 
19463879571SHong Zhang   *flg = PETSC_TRUE;
19563879571SHong Zhang   for (k=0; k<n; k++) {
196abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
19763879571SHong Zhang     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
19863879571SHong Zhang     ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr);
199e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
200e6e6b1bdSHong Zhang     if (r2 < tol){
201e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
202e6e6b1bdSHong Zhang     } else {
2032dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
204e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
205e6e6b1bdSHong Zhang       r1 /= r2;
206e6e6b1bdSHong Zhang     }
207e6e6b1bdSHong Zhang     if (r1 > tol) {
20863879571SHong Zhang       *flg = PETSC_FALSE;
2091e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTranspose() %G\n",k,r1);CHKERRQ(ierr);
21063879571SHong Zhang       break;
21163879571SHong Zhang     }
21263879571SHong Zhang   }
21363879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
21463879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
21563879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
21663879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
21763879571SHong Zhang   PetscFunctionReturn(0);
21863879571SHong Zhang }
21963879571SHong Zhang 
22063879571SHong Zhang #undef __FUNCT__
22163879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual"
22263879571SHong Zhang /*@
22363879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
22463879571SHong Zhang 
22563879571SHong Zhang    Collective on Mat
22663879571SHong Zhang 
22763879571SHong Zhang    Input Parameters:
22863879571SHong Zhang +  A - the first matrix
22963879571SHong Zhang -  B - the second matrix
23063879571SHong Zhang -  n - number of random vectors to be tested
23163879571SHong Zhang 
23263879571SHong Zhang    Output Parameter:
23363879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
23463879571SHong Zhang 
23563879571SHong Zhang    Level: intermediate
23663879571SHong Zhang 
23763879571SHong Zhang    Concepts: matrices^equality between
23863879571SHong Zhang @*/
239be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
24063879571SHong Zhang {
24163879571SHong Zhang   PetscErrorCode ierr;
24263879571SHong Zhang   Vec            x,y,s1,s2;
24363879571SHong Zhang   PetscRandom    rctx;
24463879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
24563879571SHong Zhang   PetscInt       am,an,bm,bn,k;
246e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
24763879571SHong Zhang 
24863879571SHong Zhang   PetscFunctionBegin;
24963879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
25063879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
25163879571SHong 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);
25263879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
2537adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
254c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
2557adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
25663879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
25763879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
25863879571SHong Zhang 
2597adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
26063879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
26163879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
26263879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
26363879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
26463879571SHong Zhang 
26563879571SHong Zhang   *flg = PETSC_TRUE;
26663879571SHong Zhang   for (k=0; k<n; k++) {
267abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
268abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
26963879571SHong Zhang     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
27063879571SHong Zhang     ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr);
271e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
272e6e6b1bdSHong Zhang     if (r2 < tol){
273e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
274e6e6b1bdSHong Zhang     } else {
2752dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
276e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
277e6e6b1bdSHong Zhang       r1 /= r2;
278e6e6b1bdSHong Zhang     }
279e6e6b1bdSHong Zhang     if (r1 > tol) {
28063879571SHong Zhang       *flg = PETSC_FALSE;
2811e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %G\n",k,r1);CHKERRQ(ierr);
2824eb6d288SHong Zhang       break;
28386a22c91SHong Zhang     }
28486a22c91SHong Zhang   }
28586a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
28686a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
28786a22c91SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
28886a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
28986a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
29086a22c91SHong Zhang   PetscFunctionReturn(0);
29186a22c91SHong Zhang }
292