xref: /petsc/src/mat/utils/multequal.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris 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 @*/
24*be1d678aSKris 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;
3186a22c91SHong Zhang 
3286a22c91SHong Zhang   PetscFunctionBegin;
334b497936SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
344b497936SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
3586a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
3686a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
3786a22c91SHong 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);
38cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
39cb5d8e9eSHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
40cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
4186a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
4286a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
43cb5d8e9eSHong Zhang 
44cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
45cb5d8e9eSHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
46cb5d8e9eSHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
47cb5d8e9eSHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
4886a22c91SHong Zhang 
494eb6d288SHong Zhang   *flg = PETSC_TRUE;
5086a22c91SHong Zhang   for (k=0; k<n; k++) {
5186a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
5286a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
5386a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
5486a22c91SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
5586a22c91SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
5686a22c91SHong Zhang     r1 -= r2;
5786a22c91SHong Zhang     if (r1<-tol || r1>tol) {
5886a22c91SHong Zhang       *flg = PETSC_FALSE;
5963879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMult() %g\n",k,r1);
604eb6d288SHong Zhang       break;
6186a22c91SHong Zhang     }
6286a22c91SHong Zhang   }
6386a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
6486a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
6586a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
6686a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
6786a22c91SHong Zhang   PetscFunctionReturn(0);
6886a22c91SHong Zhang }
6986a22c91SHong Zhang 
7086a22c91SHong Zhang #undef __FUNCT__
7186a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual"
7286a22c91SHong Zhang /*@
7386a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
7486a22c91SHong Zhang 
7586a22c91SHong Zhang    Collective on Mat
7686a22c91SHong Zhang 
7786a22c91SHong Zhang    Input Parameters:
7886a22c91SHong Zhang +  A - the first matrix
7986a22c91SHong Zhang -  B - the second matrix
8086a22c91SHong Zhang -  n - number of random vectors to be tested
8186a22c91SHong Zhang 
8286a22c91SHong Zhang    Output Parameter:
8386a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
8486a22c91SHong Zhang 
8586a22c91SHong Zhang    Level: intermediate
8686a22c91SHong Zhang 
8786a22c91SHong Zhang    Concepts: matrices^equality between
8886a22c91SHong Zhang @*/
89*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
9086a22c91SHong Zhang {
9186a22c91SHong Zhang   PetscErrorCode ierr;
9286a22c91SHong Zhang   Vec            x,y,s1,s2;
9386a22c91SHong Zhang   PetscRandom    rctx;
9486a22c91SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
9586a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
9686a22c91SHong Zhang 
9786a22c91SHong Zhang   PetscFunctionBegin;
9886a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
9986a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
10086a22c91SHong 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);
101cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
102cb5d8e9eSHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
103cb5d8e9eSHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
10486a22c91SHong Zhang   ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr);
10586a22c91SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
10663879571SHong Zhang 
10763879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
10863879571SHong Zhang   ierr = VecSetSizes(s1,am,PETSC_DECIDE);CHKERRQ(ierr);
10963879571SHong Zhang   ierr = VecSetFromOptions(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++) {
11586a22c91SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
11686a22c91SHong Zhang     ierr = VecSetRandom(rctx,y);CHKERRQ(ierr);
11786a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
11886a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
11986a22c91SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
12086a22c91SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
12186a22c91SHong Zhang     r1 -= r2;
12286a22c91SHong Zhang     if (r1<-tol || r1>tol) {
12386a22c91SHong Zhang       *flg = PETSC_FALSE;
12463879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMultAdd() %g\n",k,r1);
12563879571SHong Zhang       break;
12663879571SHong Zhang     }
12763879571SHong Zhang   }
12863879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
12963879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
13063879571SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
13163879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
13263879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
13363879571SHong Zhang   PetscFunctionReturn(0);
13463879571SHong Zhang }
13563879571SHong Zhang 
13663879571SHong Zhang #undef __FUNCT__
13763879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual"
13863879571SHong Zhang /*@
13963879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
14063879571SHong Zhang 
14163879571SHong Zhang    Collective on Mat
14263879571SHong Zhang 
14363879571SHong Zhang    Input Parameters:
14463879571SHong Zhang +  A - the first matrix
14563879571SHong Zhang -  B - the second matrix
14663879571SHong Zhang -  n - number of random vectors to be tested
14763879571SHong Zhang 
14863879571SHong Zhang    Output Parameter:
14963879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
15063879571SHong Zhang 
15163879571SHong Zhang    Level: intermediate
15263879571SHong Zhang 
15363879571SHong Zhang    Concepts: matrices^equality between
15463879571SHong Zhang @*/
155*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
15663879571SHong Zhang {
15763879571SHong Zhang   PetscErrorCode ierr;
15863879571SHong Zhang   Vec            x,s1,s2;
15963879571SHong Zhang   PetscRandom    rctx;
16063879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
16163879571SHong Zhang   PetscInt       am,an,bm,bn,k;
16263879571SHong Zhang 
16363879571SHong Zhang   PetscFunctionBegin;
16463879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
16563879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
16663879571SHong 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);
16763879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
16863879571SHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
16963879571SHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
17063879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
17163879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
17263879571SHong Zhang 
17363879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
17463879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
17563879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
17663879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
17763879571SHong Zhang 
17863879571SHong Zhang   *flg = PETSC_TRUE;
17963879571SHong Zhang   for (k=0; k<n; k++) {
18063879571SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
18163879571SHong Zhang     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
18263879571SHong Zhang     ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr);
18363879571SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
18463879571SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
18563879571SHong Zhang     r1 -= r2;
18663879571SHong Zhang     if (r1<-tol || r1>tol) {
18763879571SHong Zhang       *flg = PETSC_FALSE;
18863879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMultTranspose() %g\n",k,r1);
18963879571SHong Zhang       break;
19063879571SHong Zhang     }
19163879571SHong Zhang   }
19263879571SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
19363879571SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
19463879571SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
19563879571SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
19663879571SHong Zhang   PetscFunctionReturn(0);
19763879571SHong Zhang }
19863879571SHong Zhang 
19963879571SHong Zhang #undef __FUNCT__
20063879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual"
20163879571SHong Zhang /*@
20263879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
20363879571SHong Zhang 
20463879571SHong Zhang    Collective on Mat
20563879571SHong Zhang 
20663879571SHong Zhang    Input Parameters:
20763879571SHong Zhang +  A - the first matrix
20863879571SHong Zhang -  B - the second matrix
20963879571SHong Zhang -  n - number of random vectors to be tested
21063879571SHong Zhang 
21163879571SHong Zhang    Output Parameter:
21263879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
21363879571SHong Zhang 
21463879571SHong Zhang    Level: intermediate
21563879571SHong Zhang 
21663879571SHong Zhang    Concepts: matrices^equality between
21763879571SHong Zhang @*/
218*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscTruth *flg)
21963879571SHong Zhang {
22063879571SHong Zhang   PetscErrorCode ierr;
22163879571SHong Zhang   Vec            x,y,s1,s2;
22263879571SHong Zhang   PetscRandom    rctx;
22363879571SHong Zhang   PetscReal      r1,r2,tol=1.e-10;
22463879571SHong Zhang   PetscInt       am,an,bm,bn,k;
22563879571SHong Zhang 
22663879571SHong Zhang   PetscFunctionBegin;
22763879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
22863879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
22963879571SHong 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);
23063879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
23163879571SHong Zhang   ierr = PetscRandomCreate(A->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
23263879571SHong Zhang   ierr = VecCreate(A->comm,&x);CHKERRQ(ierr);
23363879571SHong Zhang   ierr = VecSetSizes(x,am,PETSC_DECIDE);CHKERRQ(ierr);
23463879571SHong Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
23563879571SHong Zhang 
23663879571SHong Zhang   ierr = VecCreate(A->comm,&s1);CHKERRQ(ierr);
23763879571SHong Zhang   ierr = VecSetSizes(s1,an,PETSC_DECIDE);CHKERRQ(ierr);
23863879571SHong Zhang   ierr = VecSetFromOptions(s1);CHKERRQ(ierr);
23963879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
24063879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
24163879571SHong Zhang 
24263879571SHong Zhang   *flg = PETSC_TRUE;
24363879571SHong Zhang   for (k=0; k<n; k++) {
24463879571SHong Zhang     ierr = VecSetRandom(rctx,x);CHKERRQ(ierr);
24563879571SHong Zhang     ierr = VecSetRandom(rctx,y);CHKERRQ(ierr);
24663879571SHong Zhang     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
24763879571SHong Zhang     ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr);
24863879571SHong Zhang     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
24963879571SHong Zhang     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
25063879571SHong Zhang     r1 -= r2;
25163879571SHong Zhang     if (r1<-tol || r1>tol) {
25263879571SHong Zhang       *flg = PETSC_FALSE;
25363879571SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: %d-th MatMultTransposeAdd() %g\n",k,r1);
2544eb6d288SHong Zhang       break;
25586a22c91SHong Zhang     }
25686a22c91SHong Zhang   }
25786a22c91SHong Zhang   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
25886a22c91SHong Zhang   ierr = VecDestroy(x);CHKERRQ(ierr);
25986a22c91SHong Zhang   ierr = VecDestroy(y);CHKERRQ(ierr);
26086a22c91SHong Zhang   ierr = VecDestroy(s1);CHKERRQ(ierr);
26186a22c91SHong Zhang   ierr = VecDestroy(s2);CHKERRQ(ierr);
26286a22c91SHong Zhang   PetscFunctionReturn(0);
26386a22c91SHong Zhang }
264