xref: /petsc/src/mat/tests/ex174.cxx (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1 static char help[] = "Tests MatConvert(), MatLoad() for MATELEMENTAL interface.\n\n";
2 /*
3  Example:
4    mpiexec -n <np> ./ex173 -fA <A_data> -fB <B_data> -orig_mat_type <type> -orig_mat_type <mat_type>
5 */
6 
7 #include <petscmat.h>
8 #include <petscmatelemental.h>
9 
main(int argc,char ** args)10 int main(int argc, char **args)
11 {
12   Mat         A, Ae, B, Be;
13   PetscViewer view;
14   char        file[2][PETSC_MAX_PATH_LEN];
15   PetscBool   flg, flgB, isElemental, isDense, isAij, isSbaij;
16   PetscScalar one = 1.0;
17   PetscMPIInt rank, size;
18   PetscInt    M, N;
19 
20   PetscFunctionBeginUser;
21   PetscCall(PetscInitialize(&argc, &args, nullptr, help));
22   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24 
25   /* Load PETSc matrices */
26   PetscCall(PetscOptionsGetString(NULL, NULL, "-fA", file[0], sizeof(file[0]), NULL));
27   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[0], FILE_MODE_READ, &view));
28   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
29   PetscCall(MatSetOptionsPrefix(A, "orig_"));
30   PetscCall(MatSetType(A, MATAIJ));
31   PetscCall(MatSetFromOptions(A));
32   PetscCall(MatLoad(A, view));
33   PetscCall(PetscViewerDestroy(&view));
34 
35   PetscCall(PetscOptionsGetString(NULL, NULL, "-fB", file[1], sizeof(file[1]), &flgB));
36   if (flgB) {
37     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[1], FILE_MODE_READ, &view));
38     PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
39     PetscCall(MatSetOptionsPrefix(B, "orig_"));
40     PetscCall(MatSetType(B, MATAIJ));
41     PetscCall(MatSetFromOptions(B));
42     PetscCall(MatLoad(B, view));
43     PetscCall(PetscViewerDestroy(&view));
44   } else {
45     /* Create matrix B = I */
46     PetscInt rstart, rend, i;
47     PetscCall(MatGetSize(A, &M, &N));
48     PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
49 
50     PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
51     PetscCall(MatSetOptionsPrefix(B, "orig_"));
52     PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N));
53     PetscCall(MatSetType(B, MATAIJ));
54     PetscCall(MatSetFromOptions(B));
55     PetscCall(MatSetUp(B));
56     for (i = rstart; i < rend; i++) PetscCall(MatSetValues(B, 1, &i, 1, &i, &one, ADD_VALUES));
57     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
58     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
59   }
60 
61   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATELEMENTAL, &isElemental));
62   if (isElemental) {
63     Ae      = A;
64     Be      = B;
65     isDense = isAij = isSbaij = PETSC_FALSE;
66   } else { /* Convert AIJ/DENSE/SBAIJ matrices into Elemental matrices */
67     if (size == 1) {
68       PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQDENSE, &isDense));
69       PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isAij));
70       PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSbaij));
71     } else {
72       PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIDENSE, &isDense));
73       PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isAij));
74       PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &isSbaij));
75     }
76 
77     if (rank == 0) {
78       if (isDense) {
79         printf(" Convert DENSE matrices A and B into Elemental matrix... \n");
80       } else if (isAij) {
81         printf(" Convert AIJ matrices A and B into Elemental matrix... \n");
82       } else if (isSbaij) {
83         printf(" Convert SBAIJ matrices A and B into Elemental matrix... \n");
84       } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not supported yet");
85     }
86     PetscCall(MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae));
87     PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &Be));
88 
89     /* Test accuracy */
90     PetscCall(MatMultEqual(A, Ae, 5, &flg));
91     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "A != A_elemental.");
92     PetscCall(MatMultEqual(B, Be, 5, &flg));
93     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "B != B_elemental.");
94   }
95 
96   if (!isElemental) {
97     PetscCall(MatDestroy(&Ae));
98     PetscCall(MatDestroy(&Be));
99 
100     /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
101     PetscCall(MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A));
102     //PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
103   }
104 
105   PetscCall(MatDestroy(&A));
106   PetscCall(MatDestroy(&B));
107   PetscCall(PetscFinalize());
108   return 0;
109 }
110 
111 /*TEST
112 
113    build:
114       requires: elemental
115 
116    test:
117       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
118       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
119       output_file: output/ex174.out
120 
121    test:
122       suffix: 2
123       nsize: 8
124       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
125       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
126       output_file: output/ex174.out
127 
128    test:
129       suffix: 2_dense
130       nsize: 8
131       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
132       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type dense
133       output_file: output/ex174_dense.out
134 
135    test:
136       suffix: 2_elemental
137       nsize: 8
138       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
139       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type elemental
140       output_file: output/empty.out
141 
142    test:
143       suffix: 2_sbaij
144       nsize: 8
145       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
146       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B -orig_mat_type sbaij
147       output_file: output/ex174_sbaij.out
148 
149    test:
150       suffix: complex
151       requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
152       args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
153       output_file: output/ex174.out
154 
155    test:
156       suffix: complex_2
157       nsize: 4
158       requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
159       args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
160       output_file: output/ex174.out
161 
162    test:
163       suffix: dense
164       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
165       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type dense
166 
167    test:
168       suffix: elemental
169       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
170       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type elemental
171       output_file: output/empty.out
172 
173    test:
174       suffix: sbaij
175       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
176       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B -orig_mat_type sbaij
177 
178 TEST*/
179