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 PetscErrorCode ierr; 13c4762a1bSJed Brown Vec X, Y; 14c4762a1bSJed Brown Mat A,B,Af; 15c4762a1bSJed Brown PetscBool flg; 16e5f3c4beSJose E. Roman PetscReal xnorm,ynorm,anorm; 17c4762a1bSJed Brown 18c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19c4762a1bSJed Brown 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,20,20,3.0,&A)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&X,&Y)); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 23c4762a1bSJed Brown 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(X,NULL)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(X,NORM_2,&xnorm)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,X,Y)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Y,NORM_2,&ynorm)); 282c71b3e2SJacob 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); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(A,5.0)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,.5)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_FROBENIUS,&anorm)); 332c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(anorm - 4.0) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm 4.0 actual norm %g",(double)anorm); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* Convert to AIJ (exercises MatGetRow/MatRestoreRow) */ 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,B,10,&flg)); 38*5f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMult\n")); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,B,10,&flg)); 40*5f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultAdd\n")); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeEqual(A,B,10,&flg)); 42*5f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTranspose\n")); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(A,B,10,&flg)); 44*5f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTransposeAdd\n")); 45c4762a1bSJed Brown 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,Y)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&Af)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(Af,A,NULL,NULL,NULL)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(Af,A,NULL)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(Af,X,Y)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Y,NORM_2,&ynorm)); 522c71b3e2SJacob 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); 53c4762a1bSJed Brown 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Af)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown ierr = PetscFinalize(); 61c4762a1bSJed Brown return ierr; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /*TEST 65c4762a1bSJed Brown 66c4762a1bSJed Brown test: 67c4762a1bSJed Brown nsize: 2 68c4762a1bSJed Brown 69c4762a1bSJed Brown TEST*/ 70