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