xref: /petsc/src/mat/tests/ex185.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown static char help[] = "Tests MatCreateConstantDiagonal().\n"
2c4762a1bSJed Brown                      "\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Vec       X, Y;
9c4762a1bSJed Brown   Mat       A, B, Af;
10c4762a1bSJed Brown   PetscBool flg;
11e5f3c4beSJose E. Roman   PetscReal xnorm, ynorm, anorm;
12c4762a1bSJed Brown 
13327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
15c4762a1bSJed Brown 
169566063dSJacob Faibussowitsch   PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, 20, 20, 3.0, &A));
179566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &X, &Y));
189566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
19c4762a1bSJed Brown 
209566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(X, NULL));
219566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &xnorm));
229566063dSJacob Faibussowitsch   PetscCall(MatMult(A, X, Y));
239566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
2408401ef6SPierre Jolivet   PetscCheck(PetscAbsReal(ynorm - 3 * xnorm) <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Expected norm %g actual norm %g", (double)(3 * xnorm), (double)ynorm);
259566063dSJacob Faibussowitsch   PetscCall(MatShift(A, 5.0));
269566063dSJacob Faibussowitsch   PetscCall(MatScale(A, .5));
279566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
289566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_FROBENIUS, &anorm));
2908401ef6SPierre Jolivet   PetscCheck(PetscAbsReal(anorm - 4.0) <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Expected norm 4.0 actual norm %g", (double)anorm);
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* Convert to AIJ (exercises MatGetRow/MatRestoreRow) */
329566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
339566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A, B, 10, &flg));
349566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatMult\n"));
359566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(A, B, 10, &flg));
369566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatMultAdd\n"));
379566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeEqual(A, B, 10, &flg));
389566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatMultTranspose\n"));
399566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAddEqual(A, B, 10, &flg));
409566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatMultTransposeAdd\n"));
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A, Y));
439566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &Af));
449566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(Af, A, NULL, NULL, NULL));
459566063dSJacob Faibussowitsch   PetscCall(MatLUFactorNumeric(Af, A, NULL));
469566063dSJacob Faibussowitsch   PetscCall(MatSolve(Af, X, Y));
479566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
4808401ef6SPierre Jolivet   PetscCheck(PetscAbsReal(ynorm - xnorm / 4) <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Expected norm %g actual norm %g", (double)(.25 * xnorm), (double)ynorm);
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Af));
539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
57b122ec5aSJacob Faibussowitsch   return 0;
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown /*TEST
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   test:
63c4762a1bSJed Brown     nsize: 2
64c4762a1bSJed Brown 
65c4762a1bSJed Brown TEST*/
66