1 static char help[] = "Tests MatCreateConstantDiagonal().\n" 2 "\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc, char **args) 7 { 8 Vec X, Y; 9 Mat A,B,Af; 10 PetscBool flg; 11 PetscReal xnorm,ynorm,anorm; 12 13 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 14 15 PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,20,20,3.0,&A)); 16 PetscCall(MatCreateVecs(A,&X,&Y)); 17 PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 18 19 PetscCall(VecSetRandom(X,NULL)); 20 PetscCall(VecNorm(X,NORM_2,&xnorm)); 21 PetscCall(MatMult(A,X,Y)); 22 PetscCall(VecNorm(Y,NORM_2,&ynorm)); 23 PetscCheckFalse(PetscAbsReal(ynorm - 3*xnorm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm %g actual norm %g",(double)(3*xnorm),(double)ynorm); 24 PetscCall(MatShift(A,5.0)); 25 PetscCall(MatScale(A,.5)); 26 PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 27 PetscCall(MatNorm(A,NORM_FROBENIUS,&anorm)); 28 PetscCheckFalse(PetscAbsReal(anorm - 4.0) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm 4.0 actual norm %g",(double)anorm); 29 30 /* Convert to AIJ (exercises MatGetRow/MatRestoreRow) */ 31 PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 32 PetscCall(MatMultEqual(A,B,10,&flg)); 33 if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMult\n")); 34 PetscCall(MatMultAddEqual(A,B,10,&flg)); 35 if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultAdd\n")); 36 PetscCall(MatMultTransposeEqual(A,B,10,&flg)); 37 if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTranspose\n")); 38 PetscCall(MatMultTransposeAddEqual(A,B,10,&flg)); 39 if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTransposeAdd\n")); 40 41 PetscCall(MatGetDiagonal(A,Y)); 42 PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&Af)); 43 PetscCall(MatLUFactorSymbolic(Af,A,NULL,NULL,NULL)); 44 PetscCall(MatLUFactorNumeric(Af,A,NULL)); 45 PetscCall(MatSolve(Af,X,Y)); 46 PetscCall(VecNorm(Y,NORM_2,&ynorm)); 47 PetscCheckFalse(PetscAbsReal(ynorm - xnorm/4) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm %g actual norm %g",(double)(.25*xnorm),(double)ynorm); 48 49 PetscCall(MatDestroy(&A)); 50 PetscCall(MatDestroy(&B)); 51 PetscCall(MatDestroy(&Af)); 52 PetscCall(VecDestroy(&X)); 53 PetscCall(VecDestroy(&Y)); 54 55 PetscCall(PetscFinalize()); 56 return 0; 57 } 58 59 /*TEST 60 61 test: 62 nsize: 2 63 64 TEST*/ 65