xref: /petsc/src/mat/utils/multequal.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
286a22c91SHong Zhang 
37c4f633dSBarry 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;
34*0700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
35*0700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,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);
404fcfbe5fSBarry Smith 
414fcfbe5fSBarry Smith #if defined(PETSC_USE_SCALAR_SINGLE)
424fcfbe5fSBarry Smith   tol  = 1.e-5;
434fcfbe5fSBarry Smith #endif
447adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
45c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
467adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
4786a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
4886a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
49cb5d8e9eSHong Zhang 
507adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
51cb5d8e9eSHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
52cb5d8e9eSHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
53cb5d8e9eSHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
5486a22c91SHong Zhang 
554eb6d288SHong Zhang   *flg = PETSC_TRUE;
5686a22c91SHong Zhang   for (k=0; k<n; k++) {
57abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
5886a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
5986a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
60e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
61e6e6b1bdSHong Zhang     if (r2 < tol){
62e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
63e6e6b1bdSHong Zhang     } else {
642dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
65e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
66e6e6b1bdSHong Zhang       r1 /= r2;
67e6e6b1bdSHong Zhang     }
68e6e6b1bdSHong Zhang     if (r1 > tol) {
6986a22c91SHong Zhang       *flg = PETSC_FALSE;
701e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %D-th MatMult() %G\n",k,r1);CHKERRQ(ierr);
714eb6d288SHong Zhang       break;
7286a22c91SHong Zhang     }
7386a22c91SHong Zhang   }
7486a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
7586a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
7686a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
7786a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
7886a22c91SHong Zhang   PetscFunctionReturn(0);
7986a22c91SHong Zhang }
8086a22c91SHong Zhang 
8186a22c91SHong Zhang #undef __FUNCT__
8286a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual"
8386a22c91SHong Zhang /*@
8486a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
8586a22c91SHong Zhang 
8686a22c91SHong Zhang    Collective on Mat
8786a22c91SHong Zhang 
8886a22c91SHong Zhang    Input Parameters:
8986a22c91SHong Zhang +  A - the first matrix
9086a22c91SHong Zhang -  B - the second matrix
9186a22c91SHong Zhang -  n - number of random vectors to be tested
9286a22c91SHong Zhang 
9386a22c91SHong Zhang    Output Parameter:
9486a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
9586a22c91SHong Zhang 
9686a22c91SHong Zhang    Level: intermediate
9786a22c91SHong Zhang 
9886a22c91SHong Zhang    Concepts: matrices^equality between
9986a22c91SHong Zhang @*/
100be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
10186a22c91SHong Zhang {
10286a22c91SHong Zhang   PetscErrorCode ierr;
10386a22c91SHong Zhang   Vec            x,y,s1,s2;
10486a22c91SHong Zhang   PetscRandom    rctx;
10586a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
10686a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
107e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
10886a22c91SHong Zhang 
10986a22c91SHong Zhang   PetscFunctionBegin;
11086a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
11186a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
11286a22c91SHong 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);
113cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
1147adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
115c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1167adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
11786a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
11886a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
11963879571SHong Zhang 
1207adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
12163879571SHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
12263879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
12363879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
12463879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
12586a22c91SHong Zhang 
1264eb6d288SHong Zhang   *flg = PETSC_TRUE;
12786a22c91SHong Zhang   for (k=0; k<n; k++) {
128abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
129abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
13086a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
13186a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
132e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
133e6e6b1bdSHong Zhang     if (r2 < tol){
134e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
135e6e6b1bdSHong Zhang     } else {
1362dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
137e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
138e6e6b1bdSHong Zhang       r1 /= r2;
139e6e6b1bdSHong Zhang     }
140e6e6b1bdSHong Zhang     if (r1 > tol) {
14186a22c91SHong Zhang       *flg = PETSC_FALSE;
1421e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultAdd() %G\n",k,r1);CHKERRQ(ierr);
14363879571SHong Zhang       break;
14463879571SHong Zhang     }
14563879571SHong Zhang   }
14663879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
14763879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
14863879571SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
14963879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
15063879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
15163879571SHong Zhang   PetscFunctionReturn(0);
15263879571SHong Zhang }
15363879571SHong Zhang 
15463879571SHong Zhang #undef __FUNCT__
15563879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual"
15663879571SHong Zhang /*@
15763879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
15863879571SHong Zhang 
15963879571SHong Zhang    Collective on Mat
16063879571SHong Zhang 
16163879571SHong Zhang    Input Parameters:
16263879571SHong Zhang +  A - the first matrix
16363879571SHong Zhang -  B - the second matrix
16463879571SHong Zhang -  n - number of random vectors to be tested
16563879571SHong Zhang 
16663879571SHong Zhang    Output Parameter:
16763879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
16863879571SHong Zhang 
16963879571SHong Zhang    Level: intermediate
17063879571SHong Zhang 
17163879571SHong Zhang    Concepts: matrices^equality between
17263879571SHong Zhang @*/
173be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
17463879571SHong Zhang {
17563879571SHong Zhang   PetscErrorCode ierr;
17663879571SHong Zhang   Vec            x,s1,s2;
17763879571SHong Zhang   PetscRandom    rctx;
17863879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
17963879571SHong Zhang   PetscInt       am,an,bm,bn,k;
180e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
18163879571SHong Zhang 
18263879571SHong Zhang   PetscFunctionBegin;
18363879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
18463879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
18563879571SHong 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);
18663879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
1877adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
188c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1897adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
19063879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
19163879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
19263879571SHong Zhang 
1937adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
19463879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
19563879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
19663879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
19763879571SHong Zhang 
19863879571SHong Zhang   *flg = PETSC_TRUE;
19963879571SHong Zhang   for (k=0; k<n; k++) {
200abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
20163879571SHong Zhang     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
20263879571SHong Zhang     ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr);
203e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
204e6e6b1bdSHong Zhang     if (r2 < tol){
205e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
206e6e6b1bdSHong Zhang     } else {
2072dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
208e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
209e6e6b1bdSHong Zhang       r1 /= r2;
210e6e6b1bdSHong Zhang     }
211e6e6b1bdSHong Zhang     if (r1 > tol) {
21263879571SHong Zhang       *flg = PETSC_FALSE;
2131e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTranspose() %G\n",k,r1);CHKERRQ(ierr);
21463879571SHong Zhang       break;
21563879571SHong Zhang     }
21663879571SHong Zhang   }
21763879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
21863879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
21963879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
22063879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
22163879571SHong Zhang   PetscFunctionReturn(0);
22263879571SHong Zhang }
22363879571SHong Zhang 
22463879571SHong Zhang #undef __FUNCT__
22563879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual"
22663879571SHong Zhang /*@
22763879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
22863879571SHong Zhang 
22963879571SHong Zhang    Collective on Mat
23063879571SHong Zhang 
23163879571SHong Zhang    Input Parameters:
23263879571SHong Zhang +  A - the first matrix
23363879571SHong Zhang -  B - the second matrix
23463879571SHong Zhang -  n - number of random vectors to be tested
23563879571SHong Zhang 
23663879571SHong Zhang    Output Parameter:
23763879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
23863879571SHong Zhang 
23963879571SHong Zhang    Level: intermediate
24063879571SHong Zhang 
24163879571SHong Zhang    Concepts: matrices^equality between
24263879571SHong Zhang @*/
243be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
24463879571SHong Zhang {
24563879571SHong Zhang   PetscErrorCode ierr;
24663879571SHong Zhang   Vec            x,y,s1,s2;
24763879571SHong Zhang   PetscRandom    rctx;
24863879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
24963879571SHong Zhang   PetscInt       am,an,bm,bn,k;
250e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
25163879571SHong Zhang 
25263879571SHong Zhang   PetscFunctionBegin;
25363879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
25463879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
25563879571SHong 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);
25663879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
2577adad957SLisandro Dalcin   ierr = PetscRandomCreate(((PetscObject)A)->comm,&rctx);CHKERRQ(ierr);
258c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
2597adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&x);CHKERRQ(ierr);
26063879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
26163879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
26263879571SHong Zhang 
2637adad957SLisandro Dalcin   ierr = VecCreate(((PetscObject)A)->comm,&s1);CHKERRQ(ierr);
26463879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
26563879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
26663879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
26763879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
26863879571SHong Zhang 
26963879571SHong Zhang   *flg = PETSC_TRUE;
27063879571SHong Zhang   for (k=0; k<n; k++) {
271abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
272abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
27363879571SHong Zhang     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
27463879571SHong Zhang     ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr);
275e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
276e6e6b1bdSHong Zhang     if (r2 < tol){
277e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
278e6e6b1bdSHong Zhang     } else {
2792dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
280e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
281e6e6b1bdSHong Zhang       r1 /= r2;
282e6e6b1bdSHong Zhang     }
283e6e6b1bdSHong Zhang     if (r1 > tol) {
28463879571SHong Zhang       *flg = PETSC_FALSE;
2851e2582c4SBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %G\n",k,r1);CHKERRQ(ierr);
2864eb6d288SHong Zhang       break;
28786a22c91SHong Zhang     }
28886a22c91SHong Zhang   }
28986a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
29086a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
29186a22c91SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
29286a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
29386a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
29486a22c91SHong Zhang   PetscFunctionReturn(0);
29586a22c91SHong Zhang }
296