xref: /petsc/src/mat/tests/ex27.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscReal      nrm,tol=10*PETSC_SMALL;
12c4762a1bSJed Brown   PetscRandom    rctx;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   /* Call MatSetRandom on unassembled matrices */
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,20,20,3,NULL,3,NULL,&A[0]));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,20,20,3,NULL,3,NULL,&A[1]));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(A[0],rctx));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(A[1],rctx));
22c4762a1bSJed Brown 
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A[0],1.0,A[1],DIFFERENT_NONZERO_PATTERN));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A[0],-1.0,A[0],SAME_NONZERO_PATTERN));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A[0],NORM_1,&nrm));
26*5f80ce2aSJacob Faibussowitsch   if (nrm > tol) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm(), norm1=: %g\n",(double)nrm));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Call MatSetRandom on assembled matrices */
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(A[0],rctx));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(A[1],rctx));
31c4762a1bSJed Brown 
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A[0],1.0,A[1],DIFFERENT_NONZERO_PATTERN));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A[0],-1.0,A[0],SAME_NONZERO_PATTERN));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A[0],NORM_1,&nrm));
35*5f80ce2aSJacob Faibussowitsch   if (nrm > tol) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm(), norm1=: %g\n",(double)nrm));
36c4762a1bSJed Brown 
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A[0]));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A[1]));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rctx));
40c4762a1bSJed Brown   ierr = PetscFinalize();
41c4762a1bSJed Brown   return ierr;
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*TEST
45c4762a1bSJed Brown    test:
46c4762a1bSJed Brown       nsize: 3
47c4762a1bSJed Brown TEST*/
48