xref: /petsc/src/mat/tests/ex162.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Tests MatShift for SeqAIJ matrices with some missing diagonal entries\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat         A;
8c4762a1bSJed Brown   PetscInt    coli[4], row;
9c4762a1bSJed Brown   PetscScalar vali[4];
10c4762a1bSJed Brown   PetscMPIInt size;
11c4762a1bSJed Brown 
12327415f7SBarry Smith   PetscFunctionBeginUser;
13*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
15be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
16c4762a1bSJed Brown 
179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, 4, 4, 4, 4));
199566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATSEQAIJ));
209566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, 4, NULL));
21c4762a1bSJed Brown 
229371c9d4SSatish Balay   row     = 0;
239371c9d4SSatish Balay   coli[0] = 1;
249371c9d4SSatish Balay   coli[1] = 3;
259371c9d4SSatish Balay   vali[0] = 1.0;
269371c9d4SSatish Balay   vali[1] = 2.0;
279566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 2, coli, vali, ADD_VALUES));
28c4762a1bSJed Brown 
299371c9d4SSatish Balay   row     = 1;
309371c9d4SSatish Balay   coli[0] = 0;
319371c9d4SSatish Balay   coli[1] = 1;
329371c9d4SSatish Balay   coli[2] = 2;
339371c9d4SSatish Balay   coli[3] = 3;
349371c9d4SSatish Balay   vali[0] = 3.0;
359371c9d4SSatish Balay   vali[1] = 4.0;
369371c9d4SSatish Balay   vali[2] = 5.0;
379371c9d4SSatish Balay   vali[3] = 6.0;
389566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 4, coli, vali, ADD_VALUES));
39c4762a1bSJed Brown 
409371c9d4SSatish Balay   row     = 2;
419371c9d4SSatish Balay   coli[0] = 0;
429371c9d4SSatish Balay   coli[1] = 3;
439371c9d4SSatish Balay   vali[0] = 7.0;
449371c9d4SSatish Balay   vali[1] = 8.0;
459566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 2, coli, vali, ADD_VALUES));
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
499566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatShift(A, 0.0));
529566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
559566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
56b122ec5aSJacob Faibussowitsch   return 0;
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*TEST
60c4762a1bSJed Brown 
61c4762a1bSJed Brown    test:
62c4762a1bSJed Brown 
63c4762a1bSJed Brown TEST*/
64