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> 7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 8d71ae5a4SJacob Faibussowitsch { 9c4762a1bSJed Brown Mat A[2]; 10c4762a1bSJed Brown PetscReal nrm, tol = 10 * PETSC_SMALL; 11c4762a1bSJed Brown PetscRandom rctx; 12c4762a1bSJed Brown 13327415f7SBarry Smith PetscFunctionBeginUser; 14c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 159566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* Call MatSetRandom on unassembled matrices */ 189566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 20, 20, 3, NULL, 3, NULL, &A[0])); 199566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 20, 20, 3, NULL, 3, NULL, &A[1])); 209566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A[0], rctx)); 219566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A[1], rctx)); 22c4762a1bSJed Brown 239566063dSJacob Faibussowitsch PetscCall(MatAXPY(A[0], 1.0, A[1], DIFFERENT_NONZERO_PATTERN)); 249566063dSJacob Faibussowitsch PetscCall(MatAXPY(A[0], -1.0, A[0], SAME_NONZERO_PATTERN)); 259566063dSJacob Faibussowitsch PetscCall(MatNorm(A[0], NORM_1, &nrm)); 269566063dSJacob Faibussowitsch if (nrm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: MatNorm(), norm1=: %g\n", (double)nrm)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Call MatSetRandom on assembled matrices */ 299566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A[0], rctx)); 309566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A[1], rctx)); 31c4762a1bSJed Brown 329566063dSJacob Faibussowitsch PetscCall(MatAXPY(A[0], 1.0, A[1], DIFFERENT_NONZERO_PATTERN)); 339566063dSJacob Faibussowitsch PetscCall(MatAXPY(A[0], -1.0, A[0], SAME_NONZERO_PATTERN)); 349566063dSJacob Faibussowitsch PetscCall(MatNorm(A[0], NORM_1, &nrm)); 359566063dSJacob Faibussowitsch if (nrm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: MatNorm(), norm1=: %g\n", (double)nrm)); 36c4762a1bSJed Brown 379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A[0])); 389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A[1])); 399566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 41b122ec5aSJacob Faibussowitsch return 0; 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown /*TEST 45c4762a1bSJed Brown test: 46c4762a1bSJed Brown nsize: 3 47*3886731fSPierre Jolivet output_file: output/empty.out 48c4762a1bSJed Brown TEST*/ 49