1c4762a1bSJed Brown static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown Mat C, A; 8c4762a1bSJed Brown PetscInt i, j, m = 5, n = 5, Ii, J; 9c4762a1bSJed Brown PetscScalar v; 10c4762a1bSJed Brown IS perm, iperm; 11c4762a1bSJed Brown 12327415f7SBarry Smith PetscFunctionBeginUser; 13*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 149566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C)); 15c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 16c4762a1bSJed Brown for (i = 0; i < m; i++) { 17c4762a1bSJed Brown for (j = 0; j < n; j++) { 189371c9d4SSatish Balay v = -1.0; 199371c9d4SSatish Balay Ii = j + n * i; 209371c9d4SSatish Balay if (i > 0) { 219371c9d4SSatish Balay J = Ii - n; 229371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 239371c9d4SSatish Balay } 249371c9d4SSatish Balay if (i < m - 1) { 259371c9d4SSatish Balay J = Ii + n; 269371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 279371c9d4SSatish Balay } 289371c9d4SSatish Balay if (j > 0) { 299371c9d4SSatish Balay J = Ii - 1; 309371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 319371c9d4SSatish Balay } 329371c9d4SSatish Balay if (j < n - 1) { 339371c9d4SSatish Balay J = Ii + 1; 349371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 359371c9d4SSatish Balay } 369371c9d4SSatish Balay v = 4.0; 379371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown } 409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 42c4762a1bSJed Brown 439566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSAME, MAT_INITIAL_MATRIX, &A)); 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm)); 469566063dSJacob Faibussowitsch PetscCall(ISView(perm, PETSC_VIEWER_STDOUT_SELF)); 479566063dSJacob Faibussowitsch PetscCall(ISView(iperm, PETSC_VIEWER_STDOUT_SELF)); 489566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 549566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 55b122ec5aSJacob Faibussowitsch return 0; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 58c4762a1bSJed Brown /*TEST 59c4762a1bSJed Brown 60c4762a1bSJed Brown test: 618cc725e6SPierre Jolivet filter: grep -v " MPI process" 62c4762a1bSJed Brown 63c4762a1bSJed Brown TEST*/ 64