1c4762a1bSJed Brown static char help[] = "Tests MatCreateConstantDiagonal().\n" 2c4762a1bSJed Brown "\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown /*T 7c4762a1bSJed Brown Concepts: Mat 8c4762a1bSJed Brown T*/ 9c4762a1bSJed Brown 10c4762a1bSJed Brown int main(int argc, char **args) 11c4762a1bSJed Brown { 12c4762a1bSJed Brown Vec X, Y; 13c4762a1bSJed Brown Mat A,B,Af; 14c4762a1bSJed Brown PetscBool flg; 15e5f3c4beSJose E. Roman PetscReal xnorm,ynorm,anorm; 16c4762a1bSJed Brown 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 18c4762a1bSJed Brown 199566063dSJacob Faibussowitsch PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,20,20,3.0,&A)); 209566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&X,&Y)); 219566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 22c4762a1bSJed Brown 239566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X,NULL)); 249566063dSJacob Faibussowitsch PetscCall(VecNorm(X,NORM_2,&xnorm)); 259566063dSJacob Faibussowitsch PetscCall(MatMult(A,X,Y)); 269566063dSJacob Faibussowitsch PetscCall(VecNorm(Y,NORM_2,&ynorm)); 27*08401ef6SPierre Jolivet PetscCheck(PetscAbsReal(ynorm - 3*xnorm) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm %g actual norm %g",(double)(3*xnorm),(double)ynorm); 289566063dSJacob Faibussowitsch PetscCall(MatShift(A,5.0)); 299566063dSJacob Faibussowitsch PetscCall(MatScale(A,.5)); 309566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 319566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_FROBENIUS,&anorm)); 32*08401ef6SPierre Jolivet PetscCheck(PetscAbsReal(anorm - 4.0) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm 4.0 actual norm %g",(double)anorm); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Convert to AIJ (exercises MatGetRow/MatRestoreRow) */ 359566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 369566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,B,10,&flg)); 379566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMult\n")); 389566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A,B,10,&flg)); 399566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultAdd\n")); 409566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(A,B,10,&flg)); 419566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTranspose\n")); 429566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAddEqual(A,B,10,&flg)); 439566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTransposeAdd\n")); 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,Y)); 469566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&Af)); 479566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(Af,A,NULL,NULL,NULL)); 489566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(Af,A,NULL)); 499566063dSJacob Faibussowitsch PetscCall(MatSolve(Af,X,Y)); 509566063dSJacob Faibussowitsch PetscCall(VecNorm(Y,NORM_2,&ynorm)); 51*08401ef6SPierre Jolivet PetscCheck(PetscAbsReal(ynorm - xnorm/4) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm %g actual norm %g",(double)(.25*xnorm),(double)ynorm); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Af)); 569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 60b122ec5aSJacob Faibussowitsch return 0; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /*TEST 64c4762a1bSJed Brown 65c4762a1bSJed Brown test: 66c4762a1bSJed Brown nsize: 2 67c4762a1bSJed Brown 68c4762a1bSJed Brown TEST*/ 69