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