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