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