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