xref: /petsc/src/mat/tests/ex27.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static char help[]= "Test MatSetRandom on MATMPIAIJ matrices\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    Adapted from an example Contributed-by: Jakub Kruzik <jakub.kruzik@vsb.cz>
5c4762a1bSJed Brown */
6c4762a1bSJed Brown #include <petscmat.h>
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            A[2];
10c4762a1bSJed Brown   PetscReal      nrm,tol=10*PETSC_SMALL;
11c4762a1bSJed Brown   PetscRandom    rctx;
12c4762a1bSJed Brown 
13*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
14*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   /* Call MatSetRandom on unassembled matrices */
17*9566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,20,20,3,NULL,3,NULL,&A[0]));
18*9566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,20,20,3,NULL,3,NULL,&A[1]));
19*9566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(A[0],rctx));
20*9566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(A[1],rctx));
21c4762a1bSJed Brown 
22*9566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A[0],1.0,A[1],DIFFERENT_NONZERO_PATTERN));
23*9566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A[0],-1.0,A[0],SAME_NONZERO_PATTERN));
24*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A[0],NORM_1,&nrm));
25*9566063dSJacob Faibussowitsch   if (nrm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm(), norm1=: %g\n",(double)nrm));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /* Call MatSetRandom on assembled matrices */
28*9566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(A[0],rctx));
29*9566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(A[1],rctx));
30c4762a1bSJed Brown 
31*9566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A[0],1.0,A[1],DIFFERENT_NONZERO_PATTERN));
32*9566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A[0],-1.0,A[0],SAME_NONZERO_PATTERN));
33*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A[0],NORM_1,&nrm));
34*9566063dSJacob Faibussowitsch   if (nrm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm(), norm1=: %g\n",(double)nrm));
35c4762a1bSJed Brown 
36*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A[0]));
37*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A[1]));
38*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
39*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
40b122ec5aSJacob Faibussowitsch   return 0;
41c4762a1bSJed Brown }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown /*TEST
44c4762a1bSJed Brown    test:
45c4762a1bSJed Brown       nsize: 3
46c4762a1bSJed Brown TEST*/
47