static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.\n\n";

#include <petscmat.h>

int main(int argc, char **args)
{
  Mat         C, A;
  PetscInt    i, j, m = 5, n = 5, Ii, J;
  PetscScalar v;
  IS          perm, iperm;

  PetscFunctionBeginUser;
  PetscCall(PetscInitialize(&argc, &args, NULL, help));
  PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C));
  /* create the matrix for the five point stencil, YET AGAIN*/
  for (i = 0; i < m; i++) {
    for (j = 0; j < n; j++) {
      v  = -1.0;
      Ii = j + n * i;
      if (i > 0) {
        J = Ii - n;
        PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
      }
      if (i < m - 1) {
        J = Ii + n;
        PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
      }
      if (j > 0) {
        J = Ii - 1;
        PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
      }
      if (j < n - 1) {
        J = Ii + 1;
        PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
      }
      v = 4.0;
      PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
    }
  }
  PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
  PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));

  PetscCall(MatConvert(C, MATSAME, MAT_INITIAL_MATRIX, &A));

  PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
  PetscCall(ISView(perm, PETSC_VIEWER_STDOUT_SELF));
  PetscCall(ISView(iperm, PETSC_VIEWER_STDOUT_SELF));
  PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF));

  PetscCall(ISDestroy(&perm));
  PetscCall(ISDestroy(&iperm));
  PetscCall(MatDestroy(&C));
  PetscCall(MatDestroy(&A));
  PetscCall(PetscFinalize());
  return 0;
}

/*TEST

   test:
      filter: grep -v " MPI process"

TEST*/
