xref: /petsc/src/mat/tests/ex20.c (revision c6a7a37075f8bf8d34d92c4910d42445b7a3482d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests converting a matrix to another format with MatConvert().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         C, A;
9c4762a1bSJed Brown   PetscInt    i, j, m = 5, n = 4, Ii, J;
10c4762a1bSJed Brown   PetscMPIInt rank, size;
11c4762a1bSJed Brown   PetscScalar v;
12c4762a1bSJed Brown   char        mtype[256];
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   /* This example does not work correctly for np > 2 */
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
20be096a46SBarry Smith   PetscCheck(size <= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Use np <= 2");
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
239566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
27c4762a1bSJed Brown   for (i = 0; i < m; i++) {
28c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
299371c9d4SSatish Balay       v  = -1.0;
309371c9d4SSatish Balay       Ii = j + n * i;
319371c9d4SSatish Balay       if (i > 0) {
329371c9d4SSatish Balay         J = Ii - n;
339371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
349371c9d4SSatish Balay       }
359371c9d4SSatish Balay       if (i < m - 1) {
369371c9d4SSatish Balay         J = Ii + n;
379371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
389371c9d4SSatish Balay       }
399371c9d4SSatish Balay       if (j > 0) {
409371c9d4SSatish Balay         J = Ii - 1;
419371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
429371c9d4SSatish Balay       }
439371c9d4SSatish Balay       if (j < n - 1) {
449371c9d4SSatish Balay         J = Ii + 1;
459371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
469371c9d4SSatish Balay       }
479371c9d4SSatish Balay       v = 4.0;
489371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
49c4762a1bSJed Brown     }
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
549566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO));
559566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
569566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
579566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
58c4762a1bSJed Brown 
59*c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(mtype, MATSAME, sizeof(mtype)));
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-conv_mat_type", mtype, sizeof(mtype), NULL));
619566063dSJacob Faibussowitsch   PetscCall(MatConvert(C, mtype, MAT_INITIAL_MATRIX, &A));
629566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO));
639566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_IMPL));
669566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
679566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* Free data structures */
709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
74b122ec5aSJacob Faibussowitsch   return 0;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*TEST
78c4762a1bSJed Brown 
79c4762a1bSJed Brown    test:
80c4762a1bSJed Brown       args: -conv_mat_type seqaij
81c4762a1bSJed Brown 
82c4762a1bSJed Brown TEST*/
83