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; 282a0b5429SBarry 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); 39ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 40c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 412a7a6963SBarry Smith ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr); 42cb5d8e9eSHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 4386a22c91SHong Zhang 444eb6d288SHong Zhang *flg = PETSC_TRUE; 4586a22c91SHong Zhang for (k=0; k<n; k++) { 46abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 4786a22c91SHong Zhang ierr = MatMult(A,x,s1);CHKERRQ(ierr); 4886a22c91SHong Zhang ierr = MatMult(B,x,s2);CHKERRQ(ierr); 49e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 50e6e6b1bdSHong Zhang if (r2 < tol) { 51e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 52e6e6b1bdSHong Zhang } else { 532dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 54e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 55e6e6b1bdSHong Zhang r1 /= r2; 56e6e6b1bdSHong Zhang } 57e6e6b1bdSHong Zhang if (r1 > tol) { 5886a22c91SHong Zhang *flg = PETSC_FALSE; 5957622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %D-th MatMult() %g\n",k,(double)r1);CHKERRQ(ierr); 604eb6d288SHong Zhang break; 6186a22c91SHong Zhang } 6286a22c91SHong Zhang } 636bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 646bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 656bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 666bf464f9SBarry Smith 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 @*/ 897087cfbeSBarry Smith PetscErrorCode MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 9086a22c91SHong Zhang { 9186a22c91SHong Zhang PetscErrorCode ierr; 9286a22c91SHong Zhang Vec x,y,s1,s2; 9386a22c91SHong Zhang PetscRandom rctx; 942a0b5429SBarry Smith PetscReal r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON; 9586a22c91SHong Zhang PetscInt am,an,bm,bn,k; 96e6e6b1bdSHong Zhang PetscScalar none = -1.0; 9786a22c91SHong Zhang 9886a22c91SHong Zhang PetscFunctionBegin; 9986a22c91SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 10086a22c91SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 101e32f2f54SBarry 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); 102cb5d8e9eSHong Zhang PetscCheckSameComm(A,1,B,2); 103ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 104c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1052a7a6963SBarry Smith ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr); 10663879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 10763879571SHong Zhang ierr = VecDuplicate(s1,&y);CHKERRQ(ierr); 10886a22c91SHong Zhang 1094eb6d288SHong Zhang *flg = PETSC_TRUE; 11086a22c91SHong Zhang for (k=0; k<n; k++) { 111abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 112abb0e124SMatthew Knepley ierr = VecSetRandom(y,rctx);CHKERRQ(ierr); 11386a22c91SHong Zhang ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr); 11486a22c91SHong Zhang ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr); 115e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 116e6e6b1bdSHong Zhang if (r2 < tol) { 117e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 118e6e6b1bdSHong Zhang } else { 1192dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 120e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 121e6e6b1bdSHong Zhang r1 /= r2; 122e6e6b1bdSHong Zhang } 123e6e6b1bdSHong Zhang if (r1 > tol) { 12486a22c91SHong Zhang *flg = PETSC_FALSE; 12557622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultAdd() %g\n",k,(double)r1);CHKERRQ(ierr); 12663879571SHong Zhang break; 12763879571SHong Zhang } 12863879571SHong Zhang } 1296bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 1306bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1316bf464f9SBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 1326bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 1336bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 13463879571SHong Zhang PetscFunctionReturn(0); 13563879571SHong Zhang } 13663879571SHong Zhang 13763879571SHong Zhang #undef __FUNCT__ 13863879571SHong Zhang #define __FUNCT__ "MatMultTransposeEqual" 13963879571SHong Zhang /*@ 14063879571SHong Zhang MatMultTransposeEqual - Compares matrix-vector products of two matrices. 14163879571SHong Zhang 14263879571SHong Zhang Collective on Mat 14363879571SHong Zhang 14463879571SHong Zhang Input Parameters: 14563879571SHong Zhang + A - the first matrix 14663879571SHong Zhang - B - the second matrix 14763879571SHong Zhang - n - number of random vectors to be tested 14863879571SHong Zhang 14963879571SHong Zhang Output Parameter: 15063879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 15163879571SHong Zhang 15263879571SHong Zhang Level: intermediate 15363879571SHong Zhang 15463879571SHong Zhang Concepts: matrices^equality between 15563879571SHong Zhang @*/ 1567087cfbeSBarry Smith PetscErrorCode MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 15763879571SHong Zhang { 15863879571SHong Zhang PetscErrorCode ierr; 15963879571SHong Zhang Vec x,s1,s2; 16063879571SHong Zhang PetscRandom rctx; 1612a0b5429SBarry Smith PetscReal r1,r2,tol= PETSC_SQRT_MACHINE_EPSILON; 16263879571SHong Zhang PetscInt am,an,bm,bn,k; 163e6e6b1bdSHong Zhang PetscScalar none = -1.0; 16463879571SHong Zhang 16563879571SHong Zhang PetscFunctionBegin; 16663879571SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 16763879571SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 168e32f2f54SBarry 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); 16963879571SHong Zhang PetscCheckSameComm(A,1,B,2); 170ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 171c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1722a7a6963SBarry Smith ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr); 17363879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 17463879571SHong Zhang 17563879571SHong Zhang *flg = PETSC_TRUE; 17663879571SHong Zhang for (k=0; k<n; k++) { 177abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 17863879571SHong Zhang ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr); 17963879571SHong Zhang ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr); 180e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 181e6e6b1bdSHong Zhang if (r2 < tol) { 182e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 183e6e6b1bdSHong Zhang } else { 1842dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 185e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 186e6e6b1bdSHong Zhang r1 /= r2; 187e6e6b1bdSHong Zhang } 188e6e6b1bdSHong Zhang if (r1 > tol) { 18963879571SHong Zhang *flg = PETSC_FALSE; 19057622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultTranspose() %g\n",k,(double)r1);CHKERRQ(ierr); 19163879571SHong Zhang break; 19263879571SHong Zhang } 19363879571SHong Zhang } 1946bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 1956bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1966bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 1976bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 19863879571SHong Zhang PetscFunctionReturn(0); 19963879571SHong Zhang } 20063879571SHong Zhang 20163879571SHong Zhang #undef __FUNCT__ 20263879571SHong Zhang #define __FUNCT__ "MatMultTransposeAddEqual" 20363879571SHong Zhang /*@ 20463879571SHong Zhang MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 20563879571SHong Zhang 20663879571SHong Zhang Collective on Mat 20763879571SHong Zhang 20863879571SHong Zhang Input Parameters: 20963879571SHong Zhang + A - the first matrix 21063879571SHong Zhang - B - the second matrix 21163879571SHong Zhang - n - number of random vectors to be tested 21263879571SHong Zhang 21363879571SHong Zhang Output Parameter: 21463879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 21563879571SHong Zhang 21663879571SHong Zhang Level: intermediate 21763879571SHong Zhang 21863879571SHong Zhang Concepts: matrices^equality between 21963879571SHong Zhang @*/ 2207087cfbeSBarry Smith PetscErrorCode MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 22163879571SHong Zhang { 22263879571SHong Zhang PetscErrorCode ierr; 22363879571SHong Zhang Vec x,y,s1,s2; 22463879571SHong Zhang PetscRandom rctx; 2252a0b5429SBarry Smith PetscReal r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON; 22663879571SHong Zhang PetscInt am,an,bm,bn,k; 227e6e6b1bdSHong Zhang PetscScalar none = -1.0; 22863879571SHong Zhang 22963879571SHong Zhang PetscFunctionBegin; 23063879571SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 23163879571SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 232e32f2f54SBarry 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); 23363879571SHong Zhang PetscCheckSameComm(A,1,B,2); 234ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 235c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 2362a7a6963SBarry Smith ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr); 23763879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 23863879571SHong Zhang ierr = VecDuplicate(s1,&y);CHKERRQ(ierr); 23963879571SHong Zhang 24063879571SHong Zhang *flg = PETSC_TRUE; 24163879571SHong Zhang for (k=0; k<n; k++) { 242abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 243abb0e124SMatthew Knepley ierr = VecSetRandom(y,rctx);CHKERRQ(ierr); 24463879571SHong Zhang ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr); 24563879571SHong Zhang ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr); 246e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 247e6e6b1bdSHong Zhang if (r2 < tol) { 248e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 249e6e6b1bdSHong Zhang } else { 2502dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 251e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 252e6e6b1bdSHong Zhang r1 /= r2; 253e6e6b1bdSHong Zhang } 254e6e6b1bdSHong Zhang if (r1 > tol) { 25563879571SHong Zhang *flg = PETSC_FALSE; 25657622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %g\n",k,(double)r1);CHKERRQ(ierr); 2574eb6d288SHong Zhang break; 25886a22c91SHong Zhang } 25986a22c91SHong Zhang } 2606bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 2616bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 2626bf464f9SBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 2636bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 2646bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 26586a22c91SHong Zhang PetscFunctionReturn(0); 26686a22c91SHong Zhang } 267*a52676ecSHong Zhang 268*a52676ecSHong Zhang #undef __FUNCT__ 269*a52676ecSHong Zhang #define __FUNCT__ "MatMatMultEqual" 270*a52676ecSHong Zhang /*@ 271*a52676ecSHong Zhang MatMatMultEqual - Test A*B*x = C*x for n random vector x 272*a52676ecSHong Zhang 273*a52676ecSHong Zhang Collective on Mat 274*a52676ecSHong Zhang 275*a52676ecSHong Zhang Input Parameters: 276*a52676ecSHong Zhang + A - the first matrix 277*a52676ecSHong Zhang - B - the second matrix 278*a52676ecSHong Zhang - C - the third matrix 279*a52676ecSHong Zhang - n - number of random vectors to be tested 280*a52676ecSHong Zhang 281*a52676ecSHong Zhang Output Parameter: 282*a52676ecSHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 283*a52676ecSHong Zhang 284*a52676ecSHong Zhang Level: intermediate 285*a52676ecSHong Zhang 286*a52676ecSHong Zhang Concepts: matrices^equality between 287*a52676ecSHong Zhang @*/ 288*a52676ecSHong Zhang PetscErrorCode MatMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg) 289*a52676ecSHong Zhang { 290*a52676ecSHong Zhang PetscErrorCode ierr; 291*a52676ecSHong Zhang Vec x,s1,s2,s3; 292*a52676ecSHong Zhang PetscRandom rctx; 293*a52676ecSHong Zhang PetscReal r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON; 294*a52676ecSHong Zhang PetscInt am,an,bm,bn,cm,cn,k; 295*a52676ecSHong Zhang PetscScalar none = -1.0; 296*a52676ecSHong Zhang 297*a52676ecSHong Zhang PetscFunctionBegin; 298*a52676ecSHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 299*a52676ecSHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 300*a52676ecSHong Zhang ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr); 301*a52676ecSHong Zhang if (an != bm || am != cm || bn != cn) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %D %D %D %D",am,an,bm,bn,cm, cn); 302*a52676ecSHong Zhang 303*a52676ecSHong Zhang ierr = PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx);CHKERRQ(ierr); 304*a52676ecSHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 305*a52676ecSHong Zhang ierr = MatCreateVecs(B,&x,&s1);CHKERRQ(ierr); 306*a52676ecSHong Zhang ierr = MatCreateVecs(C,NULL,&s2);CHKERRQ(ierr); 307*a52676ecSHong Zhang ierr = VecDuplicate(s2,&s3);CHKERRQ(ierr); 308*a52676ecSHong Zhang 309*a52676ecSHong Zhang *flg = PETSC_TRUE; 310*a52676ecSHong Zhang for (k=0; k<n; k++) { 311*a52676ecSHong Zhang ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 312*a52676ecSHong Zhang ierr = MatMult(B,x,s1);CHKERRQ(ierr); 313*a52676ecSHong Zhang ierr = MatMult(A,s1,s2);CHKERRQ(ierr); 314*a52676ecSHong Zhang ierr = MatMult(C,x,s3);CHKERRQ(ierr); 315*a52676ecSHong Zhang 316*a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 317*a52676ecSHong Zhang if (r2 < tol) { 318*a52676ecSHong Zhang ierr = VecNorm(s3,NORM_INFINITY,&r1);CHKERRQ(ierr); 319*a52676ecSHong Zhang } else { 320*a52676ecSHong Zhang ierr = VecAXPY(s2,none,s3);CHKERRQ(ierr); 321*a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 322*a52676ecSHong Zhang r1 /= r2; 323*a52676ecSHong Zhang } 324*a52676ecSHong Zhang if (r1 > tol) { 325*a52676ecSHong Zhang *flg = PETSC_FALSE; 326*a52676ecSHong Zhang ierr = PetscInfo2(A,"Error: %D-th MatMatMult() %g\n",k,(double)r1);CHKERRQ(ierr); 327*a52676ecSHong Zhang break; 328*a52676ecSHong Zhang } 329*a52676ecSHong Zhang } 330*a52676ecSHong Zhang ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 331*a52676ecSHong Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 332*a52676ecSHong Zhang ierr = VecDestroy(&s1);CHKERRQ(ierr); 333*a52676ecSHong Zhang ierr = VecDestroy(&s2);CHKERRQ(ierr); 334*a52676ecSHong Zhang ierr = VecDestroy(&s3);CHKERRQ(ierr); 335*a52676ecSHong Zhang PetscFunctionReturn(0); 336*a52676ecSHong Zhang } 337*a52676ecSHong Zhang 338*a52676ecSHong Zhang #undef __FUNCT__ 339*a52676ecSHong Zhang #define __FUNCT__ "MatTransposeMatMultEqual" 340*a52676ecSHong Zhang /*@ 341*a52676ecSHong Zhang MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x 342*a52676ecSHong Zhang 343*a52676ecSHong Zhang Collective on Mat 344*a52676ecSHong Zhang 345*a52676ecSHong Zhang Input Parameters: 346*a52676ecSHong Zhang + A - the first matrix 347*a52676ecSHong Zhang - B - the second matrix 348*a52676ecSHong Zhang - C - the third matrix 349*a52676ecSHong Zhang - n - number of random vectors to be tested 350*a52676ecSHong Zhang 351*a52676ecSHong Zhang Output Parameter: 352*a52676ecSHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 353*a52676ecSHong Zhang 354*a52676ecSHong Zhang Level: intermediate 355*a52676ecSHong Zhang 356*a52676ecSHong Zhang Concepts: matrices^equality between 357*a52676ecSHong Zhang @*/ 358*a52676ecSHong Zhang PetscErrorCode MatTransposeMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg) 359*a52676ecSHong Zhang { 360*a52676ecSHong Zhang PetscErrorCode ierr; 361*a52676ecSHong Zhang Vec x,s1,s2,s3; 362*a52676ecSHong Zhang PetscRandom rctx; 363*a52676ecSHong Zhang PetscReal r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON; 364*a52676ecSHong Zhang PetscInt am,an,bm,bn,cm,cn,k; 365*a52676ecSHong Zhang PetscScalar none = -1.0; 366*a52676ecSHong Zhang 367*a52676ecSHong Zhang PetscFunctionBegin; 368*a52676ecSHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 369*a52676ecSHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 370*a52676ecSHong Zhang ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr); 371*a52676ecSHong Zhang if (am != bm || an != cm || bn != cn) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %D %D %D %D",am,an,bm,bn,cm, cn); 372*a52676ecSHong Zhang 373*a52676ecSHong Zhang ierr = PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx);CHKERRQ(ierr); 374*a52676ecSHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 375*a52676ecSHong Zhang ierr = MatCreateVecs(B,&x,&s1);CHKERRQ(ierr); 376*a52676ecSHong Zhang ierr = MatCreateVecs(C,NULL,&s2);CHKERRQ(ierr); 377*a52676ecSHong Zhang ierr = VecDuplicate(s2,&s3);CHKERRQ(ierr); 378*a52676ecSHong Zhang 379*a52676ecSHong Zhang *flg = PETSC_TRUE; 380*a52676ecSHong Zhang for (k=0; k<n; k++) { 381*a52676ecSHong Zhang ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 382*a52676ecSHong Zhang ierr = MatMult(B,x,s1);CHKERRQ(ierr); 383*a52676ecSHong Zhang ierr = MatMultTranspose(A,s1,s2);CHKERRQ(ierr); 384*a52676ecSHong Zhang ierr = MatMult(C,x,s3);CHKERRQ(ierr); 385*a52676ecSHong Zhang 386*a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 387*a52676ecSHong Zhang if (r2 < tol) { 388*a52676ecSHong Zhang ierr = VecNorm(s3,NORM_INFINITY,&r1);CHKERRQ(ierr); 389*a52676ecSHong Zhang } else { 390*a52676ecSHong Zhang ierr = VecAXPY(s2,none,s3);CHKERRQ(ierr); 391*a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 392*a52676ecSHong Zhang r1 /= r2; 393*a52676ecSHong Zhang } 394*a52676ecSHong Zhang if (r1 > tol) { 395*a52676ecSHong Zhang *flg = PETSC_FALSE; 396*a52676ecSHong Zhang ierr = PetscInfo2(A,"Error: %D-th MatTransposeMatMult() %g\n",k,(double)r1);CHKERRQ(ierr); 397*a52676ecSHong Zhang break; 398*a52676ecSHong Zhang } 399*a52676ecSHong Zhang } 400*a52676ecSHong Zhang ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 401*a52676ecSHong Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 402*a52676ecSHong Zhang ierr = VecDestroy(&s1);CHKERRQ(ierr); 403*a52676ecSHong Zhang ierr = VecDestroy(&s2);CHKERRQ(ierr); 404*a52676ecSHong Zhang ierr = VecDestroy(&s3);CHKERRQ(ierr); 405*a52676ecSHong Zhang PetscFunctionReturn(0); 406*a52676ecSHong Zhang } 407