xref: /petsc/src/mat/tests/ex64.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Saves 4by4 block matrix.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat         A;
8c4762a1bSJed Brown   PetscInt    i, j;
9c4762a1bSJed Brown   PetscMPIInt size;
10c4762a1bSJed Brown   PetscViewer fd;
11c4762a1bSJed Brown   PetscScalar values[16], one = 1.0;
12c4762a1bSJed Brown   Vec         x;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
15*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, m "Can only run on one processor");
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   /*
20c4762a1bSJed Brown      Open binary file.  Note that we use FILE_MODE_WRITE to indicate
21c4762a1bSJed Brown      writing to this file.
22c4762a1bSJed Brown   */
239566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "4by4", FILE_MODE_WRITE, &fd));
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD, 4, 12, 12, 0, 0, &A));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   for (i = 0; i < 16; i++) values[i] = i;
28c4762a1bSJed Brown   for (i = 0; i < 4; i++) values[4 * i + i] += 5;
299371c9d4SSatish Balay   i = 0;
309371c9d4SSatish Balay   j = 0;
319566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 1, &i, 1, &j, values, INSERT_VALUES));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   for (i = 0; i < 16; i++) values[i] = i;
349371c9d4SSatish Balay   i = 0;
359371c9d4SSatish Balay   j = 2;
369566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 1, &i, 1, &j, values, INSERT_VALUES));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   for (i = 0; i < 16; i++) values[i] = i;
399371c9d4SSatish Balay   i = 1;
409371c9d4SSatish Balay   j = 0;
419566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 1, &i, 1, &j, values, INSERT_VALUES));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   for (i = 0; i < 16; i++) values[i] = i;
449371c9d4SSatish Balay   for (i = 0; i < 4; i++) values[4 * i + i] += 6;
459371c9d4SSatish Balay   i = 1;
469371c9d4SSatish Balay   j = 1;
479566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 1, &i, 1, &j, values, INSERT_VALUES));
48c4762a1bSJed Brown 
499371c9d4SSatish Balay   for (i = 0; i < 16; i++) values[i] = i;
509371c9d4SSatish Balay   i = 2;
519371c9d4SSatish Balay   j = 0;
529371c9d4SSatish Balay   PetscCall(MatSetValuesBlocked(A, 1, &i, 1, &j, values, INSERT_VALUES));
539371c9d4SSatish Balay 
549371c9d4SSatish Balay   for (i = 0; i < 16; i++) values[i] = i;
559371c9d4SSatish Balay   for (i = 0; i < 4; i++) values[4 * i + i] += 7;
569371c9d4SSatish Balay   i = 2;
579371c9d4SSatish Balay   j = 2;
589566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 1, &i, 1, &j, values, INSERT_VALUES));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
629566063dSJacob Faibussowitsch   PetscCall(MatView(A, fd));
639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 12, &x));
669566063dSJacob Faibussowitsch   PetscCall(VecSet(x, one));
679566063dSJacob Faibussowitsch   PetscCall(VecView(x, fd));
689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
72b122ec5aSJacob Faibussowitsch   return 0;
73c4762a1bSJed Brown }
74