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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* Call MatSetRandom on unassembled matrices */ 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,20,20,3,NULL,3,NULL,&A[0])); 185f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,20,20,3,NULL,3,NULL,&A[1])); 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A[0],rctx)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A[1],rctx)); 21c4762a1bSJed Brown 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A[0],1.0,A[1],DIFFERENT_NONZERO_PATTERN)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A[0],-1.0,A[0],SAME_NONZERO_PATTERN)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A[0],NORM_1,&nrm)); 255f80ce2aSJacob Faibussowitsch if (nrm > tol) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm(), norm1=: %g\n",(double)nrm)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* Call MatSetRandom on assembled matrices */ 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A[0],rctx)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A[1],rctx)); 30c4762a1bSJed Brown 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A[0],1.0,A[1],DIFFERENT_NONZERO_PATTERN)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A[0],-1.0,A[0],SAME_NONZERO_PATTERN)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A[0],NORM_1,&nrm)); 345f80ce2aSJacob Faibussowitsch if (nrm > tol) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm(), norm1=: %g\n",(double)nrm)); 35c4762a1bSJed Brown 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A[0])); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A[1])); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 39*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 40*b122ec5aSJacob Faibussowitsch return 0; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown /*TEST 44c4762a1bSJed Brown test: 45c4762a1bSJed Brown nsize: 3 46c4762a1bSJed Brown TEST*/ 47