1c4762a1bSJed Brown static char help[] = "Tests converting a SBAIJ matrix to BAIJ format with MatConvert. Modified from ex55.c\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown /* Usage: ./ex141 -mat_view */
5c4762a1bSJed Brown
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown Mat C, B;
9c4762a1bSJed Brown PetscInt i, bs = 2, mbs, m, block, d_nz = 6, col[3];
10c4762a1bSJed Brown PetscMPIInt size;
11c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN];
12c4762a1bSJed Brown PetscViewer fd;
13c4762a1bSJed Brown PetscBool equal, loadmat;
14c4762a1bSJed Brown PetscScalar value[4];
15c4762a1bSJed Brown PetscInt ridx[2], cidx[2];
16c4762a1bSJed Brown
17327415f7SBarry Smith PetscFunctionBeginUser;
18c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &loadmat));
209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
21be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
22c4762a1bSJed Brown
23c4762a1bSJed Brown /* input matrix C */
24c4762a1bSJed Brown if (loadmat) {
25c4762a1bSJed Brown /* Open binary file. Load a sbaij matrix, then destroy the viewer. */
269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
289566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATSEQSBAIJ));
299566063dSJacob Faibussowitsch PetscCall(MatLoad(C, fd));
309566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd));
31c4762a1bSJed Brown } else { /* Create a sbaij mat with bs>1 */
32c4762a1bSJed Brown mbs = 8;
339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL));
34c4762a1bSJed Brown m = mbs * bs;
359566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m, m));
379566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATSBAIJ));
389566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C));
399566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(C, bs, d_nz, NULL));
409566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
419566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
42c4762a1bSJed Brown
43c4762a1bSJed Brown for (block = 0; block < mbs; block++) {
44c4762a1bSJed Brown /* diagonal blocks */
459371c9d4SSatish Balay value[0] = -1.0;
469371c9d4SSatish Balay value[1] = 4.0;
479371c9d4SSatish Balay value[2] = -1.0;
48c4762a1bSJed Brown for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) {
499371c9d4SSatish Balay col[0] = i - 1;
509371c9d4SSatish Balay col[1] = i;
519371c9d4SSatish Balay col[2] = i + 1;
529566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 3, col, value, INSERT_VALUES));
53c4762a1bSJed Brown }
549371c9d4SSatish Balay i = bs - 1 + block * bs;
559371c9d4SSatish Balay col[0] = bs - 2 + block * bs;
569371c9d4SSatish Balay col[1] = bs - 1 + block * bs;
57c4762a1bSJed Brown
589371c9d4SSatish Balay value[0] = -1.0;
599371c9d4SSatish Balay value[1] = 4.0;
609566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 2, col, value, INSERT_VALUES));
61c4762a1bSJed Brown
629371c9d4SSatish Balay i = 0 + block * bs;
639371c9d4SSatish Balay col[0] = 0 + block * bs;
649371c9d4SSatish Balay col[1] = 1 + block * bs;
65c4762a1bSJed Brown
669371c9d4SSatish Balay value[0] = 4.0;
679371c9d4SSatish Balay value[1] = -1.0;
689566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 2, col, value, INSERT_VALUES));
69c4762a1bSJed Brown }
70c4762a1bSJed Brown /* off-diagonal blocks */
719371c9d4SSatish Balay value[0] = -1.0;
729371c9d4SSatish Balay value[1] = -0.1;
739371c9d4SSatish Balay value[2] = 0.0;
749371c9d4SSatish Balay value[3] = -1.0; /* row-oriented */
75c4762a1bSJed Brown for (block = 0; block < mbs - 1; block++) {
76c4762a1bSJed Brown for (i = 0; i < bs; i++) {
779371c9d4SSatish Balay ridx[i] = block * bs + i;
789371c9d4SSatish Balay cidx[i] = (block + 1) * bs + i;
79c4762a1bSJed Brown }
809566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, bs, ridx, bs, cidx, value, INSERT_VALUES));
81c4762a1bSJed Brown }
829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
84c4762a1bSJed Brown }
85c4762a1bSJed Brown
86c4762a1bSJed Brown /* convert C to BAIJ format */
879566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSEQBAIJ, MAT_INITIAL_MATRIX, &B));
889566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B, C, 10, &equal));
8928b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatConvert fails!");
90c4762a1bSJed Brown
919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
939566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
94b122ec5aSJacob Faibussowitsch return 0;
95c4762a1bSJed Brown }
96c4762a1bSJed Brown
97c4762a1bSJed Brown /*TEST
98c4762a1bSJed Brown
99c4762a1bSJed Brown test:
100*3886731fSPierre Jolivet output_file: output/empty.out
101c4762a1bSJed Brown
102c4762a1bSJed Brown TEST*/
103