186a22c91SHong Zhang 2af0996ceSBarry Smith #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; 28*2a0b5429SBarry Smith PetscReal r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON; 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 40ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 41c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 422a7a6963SBarry Smith ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr); 43cb5d8e9eSHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 4486a22c91SHong Zhang 454eb6d288SHong Zhang *flg = PETSC_TRUE; 4686a22c91SHong Zhang for (k=0; k<n; k++) { 47abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 4886a22c91SHong Zhang ierr = MatMult(A,x,s1);CHKERRQ(ierr); 4986a22c91SHong Zhang ierr = MatMult(B,x,s2);CHKERRQ(ierr); 50e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 51e6e6b1bdSHong Zhang if (r2 < tol) { 52e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 53e6e6b1bdSHong Zhang } else { 542dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 55e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 56e6e6b1bdSHong Zhang r1 /= r2; 57e6e6b1bdSHong Zhang } 58e6e6b1bdSHong Zhang if (r1 > tol) { 5986a22c91SHong Zhang *flg = PETSC_FALSE; 6057622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %D-th MatMult() %g\n",k,(double)r1);CHKERRQ(ierr); 614eb6d288SHong Zhang break; 6286a22c91SHong Zhang } 6386a22c91SHong Zhang } 646bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 656bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 666bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 676bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 6886a22c91SHong Zhang PetscFunctionReturn(0); 6986a22c91SHong Zhang } 7086a22c91SHong Zhang 7186a22c91SHong Zhang #undef __FUNCT__ 7286a22c91SHong Zhang #define __FUNCT__ "MatMultAddEqual" 7386a22c91SHong Zhang /*@ 7486a22c91SHong Zhang MatMultAddEqual - Compares matrix-vector products of two matrices. 7586a22c91SHong Zhang 7686a22c91SHong Zhang Collective on Mat 7786a22c91SHong Zhang 7886a22c91SHong Zhang Input Parameters: 7986a22c91SHong Zhang + A - the first matrix 8086a22c91SHong Zhang - B - the second matrix 8186a22c91SHong Zhang - n - number of random vectors to be tested 8286a22c91SHong Zhang 8386a22c91SHong Zhang Output Parameter: 8486a22c91SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 8586a22c91SHong Zhang 8686a22c91SHong Zhang Level: intermediate 8786a22c91SHong Zhang 8886a22c91SHong Zhang Concepts: matrices^equality between 8986a22c91SHong Zhang @*/ 907087cfbeSBarry Smith PetscErrorCode MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 9186a22c91SHong Zhang { 9286a22c91SHong Zhang PetscErrorCode ierr; 9386a22c91SHong Zhang Vec x,y,s1,s2; 9486a22c91SHong Zhang PetscRandom rctx; 95*2a0b5429SBarry Smith PetscReal r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON; 9686a22c91SHong Zhang PetscInt am,an,bm,bn,k; 97e6e6b1bdSHong Zhang PetscScalar none = -1.0; 9886a22c91SHong Zhang 9986a22c91SHong Zhang PetscFunctionBegin; 10086a22c91SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 10186a22c91SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 102e32f2f54SBarry 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); 103cb5d8e9eSHong Zhang PetscCheckSameComm(A,1,B,2); 104ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 105c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1062a7a6963SBarry Smith ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr); 10763879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 10863879571SHong Zhang ierr = VecDuplicate(s1,&y);CHKERRQ(ierr); 10986a22c91SHong Zhang 1104eb6d288SHong Zhang *flg = PETSC_TRUE; 11186a22c91SHong Zhang for (k=0; k<n; k++) { 112abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 113abb0e124SMatthew Knepley ierr = VecSetRandom(y,rctx);CHKERRQ(ierr); 11486a22c91SHong Zhang ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr); 11586a22c91SHong Zhang ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr); 116e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 117e6e6b1bdSHong Zhang if (r2 < tol) { 118e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 119e6e6b1bdSHong Zhang } else { 1202dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 121e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 122e6e6b1bdSHong Zhang r1 /= r2; 123e6e6b1bdSHong Zhang } 124e6e6b1bdSHong Zhang if (r1 > tol) { 12586a22c91SHong Zhang *flg = PETSC_FALSE; 12657622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultAdd() %g\n",k,(double)r1);CHKERRQ(ierr); 12763879571SHong Zhang break; 12863879571SHong Zhang } 12963879571SHong Zhang } 1306bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 1316bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1326bf464f9SBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 1336bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 1346bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 13563879571SHong Zhang PetscFunctionReturn(0); 13663879571SHong Zhang } 13763879571SHong Zhang 13863879571SHong Zhang #undef __FUNCT__ 13963879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual" 14063879571SHong Zhang /*@ 14163879571SHong Zhang MatMultTransposeEqual - Compares matrix-vector products of two matrices. 14263879571SHong Zhang 14363879571SHong Zhang Collective on Mat 14463879571SHong Zhang 14563879571SHong Zhang Input Parameters: 14663879571SHong Zhang + A - the first matrix 14763879571SHong Zhang - B - the second matrix 14863879571SHong Zhang - n - number of random vectors to be tested 14963879571SHong Zhang 15063879571SHong Zhang Output Parameter: 15163879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 15263879571SHong Zhang 15363879571SHong Zhang Level: intermediate 15463879571SHong Zhang 15563879571SHong Zhang Concepts: matrices^equality between 15663879571SHong Zhang @*/ 1577087cfbeSBarry Smith PetscErrorCode MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 15863879571SHong Zhang { 15963879571SHong Zhang PetscErrorCode ierr; 16063879571SHong Zhang Vec x,s1,s2; 16163879571SHong Zhang PetscRandom rctx; 162*2a0b5429SBarry Smith PetscReal r1,r2,tol= PETSC_SQRT_MACHINE_EPSILON; 16363879571SHong Zhang PetscInt am,an,bm,bn,k; 164e6e6b1bdSHong Zhang PetscScalar none = -1.0; 16563879571SHong Zhang 16663879571SHong Zhang PetscFunctionBegin; 16763879571SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 16863879571SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 169e32f2f54SBarry 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); 17063879571SHong Zhang PetscCheckSameComm(A,1,B,2); 171ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 172c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1732a7a6963SBarry Smith ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr); 17463879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 17563879571SHong Zhang 17663879571SHong Zhang *flg = PETSC_TRUE; 17763879571SHong Zhang for (k=0; k<n; k++) { 178abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 17963879571SHong Zhang ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr); 18063879571SHong Zhang ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr); 181e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 182e6e6b1bdSHong Zhang if (r2 < tol) { 183e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 184e6e6b1bdSHong Zhang } else { 1852dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 186e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 187e6e6b1bdSHong Zhang r1 /= r2; 188e6e6b1bdSHong Zhang } 189e6e6b1bdSHong Zhang if (r1 > tol) { 19063879571SHong Zhang *flg = PETSC_FALSE; 19157622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultTranspose() %g\n",k,(double)r1);CHKERRQ(ierr); 19263879571SHong Zhang break; 19363879571SHong Zhang } 19463879571SHong Zhang } 1956bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 1966bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1976bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 1986bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 19963879571SHong Zhang PetscFunctionReturn(0); 20063879571SHong Zhang } 20163879571SHong Zhang 20263879571SHong Zhang #undef __FUNCT__ 20363879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual" 20463879571SHong Zhang /*@ 20563879571SHong Zhang MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 20663879571SHong Zhang 20763879571SHong Zhang Collective on Mat 20863879571SHong Zhang 20963879571SHong Zhang Input Parameters: 21063879571SHong Zhang + A - the first matrix 21163879571SHong Zhang - B - the second matrix 21263879571SHong Zhang - n - number of random vectors to be tested 21363879571SHong Zhang 21463879571SHong Zhang Output Parameter: 21563879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 21663879571SHong Zhang 21763879571SHong Zhang Level: intermediate 21863879571SHong Zhang 21963879571SHong Zhang Concepts: matrices^equality between 22063879571SHong Zhang @*/ 2217087cfbeSBarry Smith PetscErrorCode MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 22263879571SHong Zhang { 22363879571SHong Zhang PetscErrorCode ierr; 22463879571SHong Zhang Vec x,y,s1,s2; 22563879571SHong Zhang PetscRandom rctx; 226*2a0b5429SBarry Smith PetscReal r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON; 22763879571SHong Zhang PetscInt am,an,bm,bn,k; 228e6e6b1bdSHong Zhang PetscScalar none = -1.0; 22963879571SHong Zhang 23063879571SHong Zhang PetscFunctionBegin; 23163879571SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 23263879571SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 233e32f2f54SBarry 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); 23463879571SHong Zhang PetscCheckSameComm(A,1,B,2); 235ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 236c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 2372a7a6963SBarry Smith ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr); 23863879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 23963879571SHong Zhang ierr = VecDuplicate(s1,&y);CHKERRQ(ierr); 24063879571SHong Zhang 24163879571SHong Zhang *flg = PETSC_TRUE; 24263879571SHong Zhang for (k=0; k<n; k++) { 243abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 244abb0e124SMatthew Knepley ierr = VecSetRandom(y,rctx);CHKERRQ(ierr); 24563879571SHong Zhang ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr); 24663879571SHong Zhang ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr); 247e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 248e6e6b1bdSHong Zhang if (r2 < tol) { 249e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 250e6e6b1bdSHong Zhang } else { 2512dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 252e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 253e6e6b1bdSHong Zhang r1 /= r2; 254e6e6b1bdSHong Zhang } 255e6e6b1bdSHong Zhang if (r1 > tol) { 25663879571SHong Zhang *flg = PETSC_FALSE; 25757622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %g\n",k,(double)r1);CHKERRQ(ierr); 2584eb6d288SHong Zhang break; 25986a22c91SHong Zhang } 26086a22c91SHong Zhang } 2616bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 2626bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 2636bf464f9SBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 2646bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 2656bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 26686a22c91SHong Zhang PetscFunctionReturn(0); 26786a22c91SHong Zhang } 268