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 17*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 18c4762a1bSJed Brown 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,20,20,3.0,&A)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&X,&Y)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 22c4762a1bSJed Brown 235f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(X,NULL)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(X,NORM_2,&xnorm)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,X,Y)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Y,NORM_2,&ynorm)); 272c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(ynorm - 3*xnorm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm %g actual norm %g",(double)(3*xnorm),(double)ynorm); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(A,5.0)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,.5)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_FROBENIUS,&anorm)); 322c71b3e2SJacob Faibussowitsch PetscCheckFalse(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) */ 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,B,10,&flg)); 375f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMult\n")); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,B,10,&flg)); 395f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultAdd\n")); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeEqual(A,B,10,&flg)); 415f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTranspose\n")); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(A,B,10,&flg)); 435f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTransposeAdd\n")); 44c4762a1bSJed Brown 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,Y)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&Af)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(Af,A,NULL,NULL,NULL)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(Af,A,NULL)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(Af,X,Y)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Y,NORM_2,&ynorm)); 512c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Af)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 58c4762a1bSJed Brown 59*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 60*b122ec5aSJacob 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