186a22c91SHong Zhang 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 386a22c91SHong Zhang 486a22c91SHong Zhang /*@ 586a22c91SHong Zhang MatMultEqual - Compares matrix-vector products of two matrices. 686a22c91SHong Zhang 786a22c91SHong Zhang Collective on Mat 886a22c91SHong Zhang 986a22c91SHong Zhang Input Parameters: 1086a22c91SHong Zhang + A - the first matrix 1186a22c91SHong Zhang - B - the second matrix 1286a22c91SHong Zhang - n - number of random vectors to be tested 1386a22c91SHong Zhang 1486a22c91SHong Zhang Output Parameter: 1586a22c91SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 1686a22c91SHong Zhang 1786a22c91SHong Zhang Level: intermediate 1886a22c91SHong Zhang 1986a22c91SHong Zhang Concepts: matrices^equality between 2086a22c91SHong Zhang @*/ 217087cfbeSBarry Smith PetscErrorCode MatMultEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 2286a22c91SHong Zhang { 2386a22c91SHong Zhang PetscErrorCode ierr; 2486a22c91SHong Zhang Vec x,s1,s2; 2586a22c91SHong Zhang PetscRandom rctx; 262a0b5429SBarry Smith PetscReal r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON; 2786a22c91SHong Zhang PetscInt am,an,bm,bn,k; 28e6e6b1bdSHong Zhang PetscScalar none = -1.0; 2986a22c91SHong Zhang 3086a22c91SHong Zhang PetscFunctionBegin; 310700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 320700a824SBarry Smith PetscValidHeaderSpecific(B,MAT_CLASSID,2); 3386a22c91SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 3486a22c91SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 35e32f2f54SBarry 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); 36cb5d8e9eSHong Zhang PetscCheckSameComm(A,1,B,2); 37ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 38c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 392a7a6963SBarry Smith ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr); 40cb5d8e9eSHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 4186a22c91SHong Zhang 424eb6d288SHong Zhang *flg = PETSC_TRUE; 4386a22c91SHong Zhang for (k=0; k<n; k++) { 44abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 4586a22c91SHong Zhang ierr = MatMult(A,x,s1);CHKERRQ(ierr); 4686a22c91SHong Zhang ierr = MatMult(B,x,s2);CHKERRQ(ierr); 47e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 48e6e6b1bdSHong Zhang if (r2 < tol) { 49e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 50e6e6b1bdSHong Zhang } else { 512dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 52e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 53e6e6b1bdSHong Zhang r1 /= r2; 54e6e6b1bdSHong Zhang } 55e6e6b1bdSHong Zhang if (r1 > tol) { 5686a22c91SHong Zhang *flg = PETSC_FALSE; 5757622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %D-th MatMult() %g\n",k,(double)r1);CHKERRQ(ierr); 584eb6d288SHong Zhang break; 5986a22c91SHong Zhang } 6086a22c91SHong Zhang } 616bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 626bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 636bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 646bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 6586a22c91SHong Zhang PetscFunctionReturn(0); 6686a22c91SHong Zhang } 6786a22c91SHong Zhang 6886a22c91SHong Zhang /*@ 6986a22c91SHong Zhang MatMultAddEqual - Compares matrix-vector products of two matrices. 7086a22c91SHong Zhang 7186a22c91SHong Zhang Collective on Mat 7286a22c91SHong Zhang 7386a22c91SHong Zhang Input Parameters: 7486a22c91SHong Zhang + A - the first matrix 7586a22c91SHong Zhang - B - the second matrix 7686a22c91SHong Zhang - n - number of random vectors to be tested 7786a22c91SHong Zhang 7886a22c91SHong Zhang Output Parameter: 7986a22c91SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 8086a22c91SHong Zhang 8186a22c91SHong Zhang Level: intermediate 8286a22c91SHong Zhang 8386a22c91SHong Zhang Concepts: matrices^equality between 8486a22c91SHong Zhang @*/ 857087cfbeSBarry Smith PetscErrorCode MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 8686a22c91SHong Zhang { 8786a22c91SHong Zhang PetscErrorCode ierr; 8886a22c91SHong Zhang Vec x,y,s1,s2; 8986a22c91SHong Zhang PetscRandom rctx; 902a0b5429SBarry Smith PetscReal r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON; 9186a22c91SHong Zhang PetscInt am,an,bm,bn,k; 92e6e6b1bdSHong Zhang PetscScalar none = -1.0; 9386a22c91SHong Zhang 9486a22c91SHong Zhang PetscFunctionBegin; 9586a22c91SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 9686a22c91SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 97e32f2f54SBarry 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); 98cb5d8e9eSHong Zhang PetscCheckSameComm(A,1,B,2); 99ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 100c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1012a7a6963SBarry Smith ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr); 10263879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 10363879571SHong Zhang ierr = VecDuplicate(s1,&y);CHKERRQ(ierr); 10486a22c91SHong Zhang 1054eb6d288SHong Zhang *flg = PETSC_TRUE; 10686a22c91SHong Zhang for (k=0; k<n; k++) { 107abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 108abb0e124SMatthew Knepley ierr = VecSetRandom(y,rctx);CHKERRQ(ierr); 10986a22c91SHong Zhang ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr); 11086a22c91SHong Zhang ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr); 111e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 112e6e6b1bdSHong Zhang if (r2 < tol) { 113e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 114e6e6b1bdSHong Zhang } else { 1152dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 116e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 117e6e6b1bdSHong Zhang r1 /= r2; 118e6e6b1bdSHong Zhang } 119e6e6b1bdSHong Zhang if (r1 > tol) { 12086a22c91SHong Zhang *flg = PETSC_FALSE; 12157622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultAdd() %g\n",k,(double)r1);CHKERRQ(ierr); 12263879571SHong Zhang break; 12363879571SHong Zhang } 12463879571SHong Zhang } 1256bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 1266bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1276bf464f9SBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 1286bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 1296bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 13063879571SHong Zhang PetscFunctionReturn(0); 13163879571SHong Zhang } 13263879571SHong Zhang 13363879571SHong Zhang /*@ 13463879571SHong Zhang MatMultTransposeEqual - Compares matrix-vector products of two matrices. 13563879571SHong Zhang 13663879571SHong Zhang Collective on Mat 13763879571SHong Zhang 13863879571SHong Zhang Input Parameters: 13963879571SHong Zhang + A - the first matrix 14063879571SHong Zhang - B - the second matrix 14163879571SHong Zhang - n - number of random vectors to be tested 14263879571SHong Zhang 14363879571SHong Zhang Output Parameter: 14463879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 14563879571SHong Zhang 14663879571SHong Zhang Level: intermediate 14763879571SHong Zhang 14863879571SHong Zhang Concepts: matrices^equality between 14963879571SHong Zhang @*/ 1507087cfbeSBarry Smith PetscErrorCode MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 15163879571SHong Zhang { 15263879571SHong Zhang PetscErrorCode ierr; 15363879571SHong Zhang Vec x,s1,s2; 15463879571SHong Zhang PetscRandom rctx; 1552a0b5429SBarry Smith PetscReal r1,r2,tol= PETSC_SQRT_MACHINE_EPSILON; 15663879571SHong Zhang PetscInt am,an,bm,bn,k; 157e6e6b1bdSHong Zhang PetscScalar none = -1.0; 15863879571SHong Zhang 15963879571SHong Zhang PetscFunctionBegin; 16063879571SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 16163879571SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 162e32f2f54SBarry 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); 16363879571SHong Zhang PetscCheckSameComm(A,1,B,2); 164ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 165c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1662a7a6963SBarry Smith ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr); 16763879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 16863879571SHong Zhang 16963879571SHong Zhang *flg = PETSC_TRUE; 17063879571SHong Zhang for (k=0; k<n; k++) { 171abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 17263879571SHong Zhang ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr); 17363879571SHong Zhang ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr); 174e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 175e6e6b1bdSHong Zhang if (r2 < tol) { 176e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 177e6e6b1bdSHong Zhang } else { 1782dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 179e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 180e6e6b1bdSHong Zhang r1 /= r2; 181e6e6b1bdSHong Zhang } 182e6e6b1bdSHong Zhang if (r1 > tol) { 18363879571SHong Zhang *flg = PETSC_FALSE; 18457622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultTranspose() %g\n",k,(double)r1);CHKERRQ(ierr); 18563879571SHong Zhang break; 18663879571SHong Zhang } 18763879571SHong Zhang } 1886bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 1896bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1906bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 1916bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 19263879571SHong Zhang PetscFunctionReturn(0); 19363879571SHong Zhang } 19463879571SHong Zhang 19563879571SHong Zhang /*@ 19663879571SHong Zhang MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 19763879571SHong Zhang 19863879571SHong Zhang Collective on Mat 19963879571SHong Zhang 20063879571SHong Zhang Input Parameters: 20163879571SHong Zhang + A - the first matrix 20263879571SHong Zhang - B - the second matrix 20363879571SHong Zhang - n - number of random vectors to be tested 20463879571SHong Zhang 20563879571SHong Zhang Output Parameter: 20663879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 20763879571SHong Zhang 20863879571SHong Zhang Level: intermediate 20963879571SHong Zhang 21063879571SHong Zhang Concepts: matrices^equality between 21163879571SHong Zhang @*/ 2127087cfbeSBarry Smith PetscErrorCode MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg) 21363879571SHong Zhang { 21463879571SHong Zhang PetscErrorCode ierr; 21563879571SHong Zhang Vec x,y,s1,s2; 21663879571SHong Zhang PetscRandom rctx; 2172a0b5429SBarry Smith PetscReal r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON; 21863879571SHong Zhang PetscInt am,an,bm,bn,k; 219e6e6b1bdSHong Zhang PetscScalar none = -1.0; 22063879571SHong Zhang 22163879571SHong Zhang PetscFunctionBegin; 22263879571SHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 22363879571SHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 224e32f2f54SBarry 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); 22563879571SHong Zhang PetscCheckSameComm(A,1,B,2); 226ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr); 227c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 2282a7a6963SBarry Smith ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr); 22963879571SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 23063879571SHong Zhang ierr = VecDuplicate(s1,&y);CHKERRQ(ierr); 23163879571SHong Zhang 23263879571SHong Zhang *flg = PETSC_TRUE; 23363879571SHong Zhang for (k=0; k<n; k++) { 234abb0e124SMatthew Knepley ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 235abb0e124SMatthew Knepley ierr = VecSetRandom(y,rctx);CHKERRQ(ierr); 23663879571SHong Zhang ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr); 23763879571SHong Zhang ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr); 238e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 239e6e6b1bdSHong Zhang if (r2 < tol) { 240e6e6b1bdSHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr); 241e6e6b1bdSHong Zhang } else { 2422dcb1b2aSMatthew Knepley ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr); 243e6e6b1bdSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 244e6e6b1bdSHong Zhang r1 /= r2; 245e6e6b1bdSHong Zhang } 246e6e6b1bdSHong Zhang if (r1 > tol) { 24763879571SHong Zhang *flg = PETSC_FALSE; 24857622a8eSBarry Smith ierr = PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %g\n",k,(double)r1);CHKERRQ(ierr); 2494eb6d288SHong Zhang break; 25086a22c91SHong Zhang } 25186a22c91SHong Zhang } 2526bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 2536bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 2546bf464f9SBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 2556bf464f9SBarry Smith ierr = VecDestroy(&s1);CHKERRQ(ierr); 2566bf464f9SBarry Smith ierr = VecDestroy(&s2);CHKERRQ(ierr); 25786a22c91SHong Zhang PetscFunctionReturn(0); 25886a22c91SHong Zhang } 259a52676ecSHong Zhang 260a52676ecSHong Zhang /*@ 261a52676ecSHong Zhang MatMatMultEqual - Test A*B*x = C*x for n random vector x 262a52676ecSHong Zhang 263a52676ecSHong Zhang Collective on Mat 264a52676ecSHong Zhang 265a52676ecSHong Zhang Input Parameters: 266a52676ecSHong Zhang + A - the first matrix 267a52676ecSHong Zhang - B - the second matrix 268a52676ecSHong Zhang - C - the third matrix 269a52676ecSHong Zhang - n - number of random vectors to be tested 270a52676ecSHong Zhang 271a52676ecSHong Zhang Output Parameter: 272a52676ecSHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 273a52676ecSHong Zhang 274a52676ecSHong Zhang Level: intermediate 275a52676ecSHong Zhang 276a52676ecSHong Zhang Concepts: matrices^equality between 277a52676ecSHong Zhang @*/ 278a52676ecSHong Zhang PetscErrorCode MatMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg) 279a52676ecSHong Zhang { 280a52676ecSHong Zhang PetscErrorCode ierr; 281a52676ecSHong Zhang Vec x,s1,s2,s3; 282a52676ecSHong Zhang PetscRandom rctx; 283a52676ecSHong Zhang PetscReal r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON; 284a52676ecSHong Zhang PetscInt am,an,bm,bn,cm,cn,k; 285a52676ecSHong Zhang PetscScalar none = -1.0; 286a52676ecSHong Zhang 287a52676ecSHong Zhang PetscFunctionBegin; 288a52676ecSHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 289a52676ecSHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 290a52676ecSHong Zhang ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr); 291a52676ecSHong 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); 292a52676ecSHong Zhang 293a52676ecSHong Zhang ierr = PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx);CHKERRQ(ierr); 294a52676ecSHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 295a52676ecSHong Zhang ierr = MatCreateVecs(B,&x,&s1);CHKERRQ(ierr); 296a52676ecSHong Zhang ierr = MatCreateVecs(C,NULL,&s2);CHKERRQ(ierr); 297a52676ecSHong Zhang ierr = VecDuplicate(s2,&s3);CHKERRQ(ierr); 298a52676ecSHong Zhang 299a52676ecSHong Zhang *flg = PETSC_TRUE; 300a52676ecSHong Zhang for (k=0; k<n; k++) { 301a52676ecSHong Zhang ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 302a52676ecSHong Zhang ierr = MatMult(B,x,s1);CHKERRQ(ierr); 303a52676ecSHong Zhang ierr = MatMult(A,s1,s2);CHKERRQ(ierr); 304a52676ecSHong Zhang ierr = MatMult(C,x,s3);CHKERRQ(ierr); 305a52676ecSHong Zhang 306a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 307a52676ecSHong Zhang if (r2 < tol) { 308a52676ecSHong Zhang ierr = VecNorm(s3,NORM_INFINITY,&r1);CHKERRQ(ierr); 309a52676ecSHong Zhang } else { 310a52676ecSHong Zhang ierr = VecAXPY(s2,none,s3);CHKERRQ(ierr); 311a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 312a52676ecSHong Zhang r1 /= r2; 313a52676ecSHong Zhang } 314a52676ecSHong Zhang if (r1 > tol) { 315a52676ecSHong Zhang *flg = PETSC_FALSE; 316a52676ecSHong Zhang ierr = PetscInfo2(A,"Error: %D-th MatMatMult() %g\n",k,(double)r1);CHKERRQ(ierr); 317a52676ecSHong Zhang break; 318a52676ecSHong Zhang } 319a52676ecSHong Zhang } 320a52676ecSHong Zhang ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 321a52676ecSHong Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 322a52676ecSHong Zhang ierr = VecDestroy(&s1);CHKERRQ(ierr); 323a52676ecSHong Zhang ierr = VecDestroy(&s2);CHKERRQ(ierr); 324a52676ecSHong Zhang ierr = VecDestroy(&s3);CHKERRQ(ierr); 325a52676ecSHong Zhang PetscFunctionReturn(0); 326a52676ecSHong Zhang } 327a52676ecSHong Zhang 328a52676ecSHong Zhang /*@ 329a52676ecSHong Zhang MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x 330a52676ecSHong Zhang 331a52676ecSHong Zhang Collective on Mat 332a52676ecSHong Zhang 333a52676ecSHong Zhang Input Parameters: 334a52676ecSHong Zhang + A - the first matrix 335a52676ecSHong Zhang - B - the second matrix 336a52676ecSHong Zhang - C - the third matrix 337a52676ecSHong Zhang - n - number of random vectors to be tested 338a52676ecSHong Zhang 339a52676ecSHong Zhang Output Parameter: 340a52676ecSHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 341a52676ecSHong Zhang 342a52676ecSHong Zhang Level: intermediate 343a52676ecSHong Zhang 344a52676ecSHong Zhang Concepts: matrices^equality between 345a52676ecSHong Zhang @*/ 346a52676ecSHong Zhang PetscErrorCode MatTransposeMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg) 347a52676ecSHong Zhang { 348a52676ecSHong Zhang PetscErrorCode ierr; 349a52676ecSHong Zhang Vec x,s1,s2,s3; 350a52676ecSHong Zhang PetscRandom rctx; 351a52676ecSHong Zhang PetscReal r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON; 352a52676ecSHong Zhang PetscInt am,an,bm,bn,cm,cn,k; 353a52676ecSHong Zhang PetscScalar none = -1.0; 354a52676ecSHong Zhang 355a52676ecSHong Zhang PetscFunctionBegin; 356a52676ecSHong Zhang ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 357a52676ecSHong Zhang ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr); 358a52676ecSHong Zhang ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr); 359a52676ecSHong 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); 360a52676ecSHong Zhang 361a52676ecSHong Zhang ierr = PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx);CHKERRQ(ierr); 362a52676ecSHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 363a52676ecSHong Zhang ierr = MatCreateVecs(B,&x,&s1);CHKERRQ(ierr); 364a52676ecSHong Zhang ierr = MatCreateVecs(C,NULL,&s2);CHKERRQ(ierr); 365a52676ecSHong Zhang ierr = VecDuplicate(s2,&s3);CHKERRQ(ierr); 366a52676ecSHong Zhang 367a52676ecSHong Zhang *flg = PETSC_TRUE; 368a52676ecSHong Zhang for (k=0; k<n; k++) { 369a52676ecSHong Zhang ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 370a52676ecSHong Zhang ierr = MatMult(B,x,s1);CHKERRQ(ierr); 371a52676ecSHong Zhang ierr = MatMultTranspose(A,s1,s2);CHKERRQ(ierr); 372a52676ecSHong Zhang ierr = MatMult(C,x,s3);CHKERRQ(ierr); 373a52676ecSHong Zhang 374a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr); 375a52676ecSHong Zhang if (r2 < tol) { 376a52676ecSHong Zhang ierr = VecNorm(s3,NORM_INFINITY,&r1);CHKERRQ(ierr); 377a52676ecSHong Zhang } else { 378a52676ecSHong Zhang ierr = VecAXPY(s2,none,s3);CHKERRQ(ierr); 379a52676ecSHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr); 380a52676ecSHong Zhang r1 /= r2; 381a52676ecSHong Zhang } 382a52676ecSHong Zhang if (r1 > tol) { 383a52676ecSHong Zhang *flg = PETSC_FALSE; 384a52676ecSHong Zhang ierr = PetscInfo2(A,"Error: %D-th MatTransposeMatMult() %g\n",k,(double)r1);CHKERRQ(ierr); 385a52676ecSHong Zhang break; 386a52676ecSHong Zhang } 387a52676ecSHong Zhang } 388a52676ecSHong Zhang ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 389a52676ecSHong Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 390a52676ecSHong Zhang ierr = VecDestroy(&s1);CHKERRQ(ierr); 391a52676ecSHong Zhang ierr = VecDestroy(&s2);CHKERRQ(ierr); 392a52676ecSHong Zhang ierr = VecDestroy(&s3);CHKERRQ(ierr); 393a52676ecSHong Zhang PetscFunctionReturn(0); 394a52676ecSHong Zhang } 395*86919fd6SHong Zhang 396*86919fd6SHong Zhang /*@ 397*86919fd6SHong Zhang MatIsLinear - Check if a shell matrix A is a linear operator. 398*86919fd6SHong Zhang 399*86919fd6SHong Zhang Collective on Mat 400*86919fd6SHong Zhang 401*86919fd6SHong Zhang Input Parameters: 402*86919fd6SHong Zhang + A - the shell matrix 403*86919fd6SHong Zhang - n - number of random vectors to be tested 404*86919fd6SHong Zhang 405*86919fd6SHong Zhang Output Parameter: 406*86919fd6SHong Zhang . flg - PETSC_TRUE if the shell matrix is linear; PETSC_FALSE otherwise. 407*86919fd6SHong Zhang 408*86919fd6SHong Zhang Level: intermediate 409*86919fd6SHong Zhang @*/ 410*86919fd6SHong Zhang PetscErrorCode MatIsLinear(Mat A,PetscInt n,PetscBool *flg) 411*86919fd6SHong Zhang { 412*86919fd6SHong Zhang PetscErrorCode ierr; 413*86919fd6SHong Zhang Vec x,y,s1,s2; 414*86919fd6SHong Zhang PetscRandom rctx; 415*86919fd6SHong Zhang PetscScalar a; 416*86919fd6SHong Zhang PetscInt k; 417*86919fd6SHong Zhang PetscReal norm,normA; 418*86919fd6SHong Zhang MPI_Comm comm; 419*86919fd6SHong Zhang PetscMPIInt rank; 420*86919fd6SHong Zhang 421*86919fd6SHong Zhang PetscFunctionBegin; 422*86919fd6SHong Zhang PetscValidHeaderSpecific(A,MAT_CLASSID,1); 423*86919fd6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 424*86919fd6SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 425*86919fd6SHong Zhang 426*86919fd6SHong Zhang ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr); 427*86919fd6SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 428*86919fd6SHong Zhang ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr); 429*86919fd6SHong Zhang ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 430*86919fd6SHong Zhang ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr); 431*86919fd6SHong Zhang 432*86919fd6SHong Zhang *flg = PETSC_TRUE; 433*86919fd6SHong Zhang for (k=0; k<n; k++) { 434*86919fd6SHong Zhang ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 435*86919fd6SHong Zhang if (!rank) { 436*86919fd6SHong Zhang ierr = PetscRandomGetValue(rctx,&a);CHKERRQ(ierr); 437*86919fd6SHong Zhang } 438*86919fd6SHong Zhang ierr = MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm);CHKERRQ(ierr); 439*86919fd6SHong Zhang 440*86919fd6SHong Zhang /* s2 = a*A*x + A*y */ 441*86919fd6SHong Zhang ierr = MatMult(A,y,s2);CHKERRQ(ierr); /* s2 = A*y */ 442*86919fd6SHong Zhang ierr = MatMult(A,x,s1);CHKERRQ(ierr); /* s1 = A*x */ 443*86919fd6SHong Zhang ierr = VecAXPY(s2,a,s1);CHKERRQ(ierr); /* s2 = a s1 + s2 */ 444*86919fd6SHong Zhang 445*86919fd6SHong Zhang /* s1 = A * (a x + y) */ 446*86919fd6SHong Zhang ierr = VecAXPY(y,a,x);CHKERRQ(ierr); /* y = a x + y */ 447*86919fd6SHong Zhang ierr = MatMult(A,y,s1);CHKERRQ(ierr); 448*86919fd6SHong Zhang ierr = VecNorm(s1,NORM_INFINITY,&normA);CHKERRQ(ierr); 449*86919fd6SHong Zhang 450*86919fd6SHong Zhang ierr = VecAXPY(s2,-1.0,s1);CHKERRQ(ierr); /* s2 = - s1 + s2 */ 451*86919fd6SHong Zhang ierr = VecNorm(s2,NORM_INFINITY,&norm);CHKERRQ(ierr); 452*86919fd6SHong Zhang if (norm > 10*normA*PETSC_MACHINE_EPSILON) { 453*86919fd6SHong Zhang *flg = PETSC_FALSE; 454*86919fd6SHong Zhang ierr = PetscInfo2(A,"Error: %D-th || A*(a x + y) - (a*A*x + A*y)|| = %g\n",k,norm);CHKERRQ(ierr); 455*86919fd6SHong Zhang break; 456*86919fd6SHong Zhang } 457*86919fd6SHong Zhang } 458*86919fd6SHong Zhang ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 459*86919fd6SHong Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 460*86919fd6SHong Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 461*86919fd6SHong Zhang ierr = VecDestroy(&s1);CHKERRQ(ierr); 462*86919fd6SHong Zhang ierr = VecDestroy(&s2);CHKERRQ(ierr); 463*86919fd6SHong Zhang PetscFunctionReturn(0); 464*86919fd6SHong Zhang } 465