1d24d4204SJose E. Roman static char help[] = "Tests MatConvert(), MatLoad() for MATSCALAPACK interface.\n\n";
2d24d4204SJose E. Roman /*
3d24d4204SJose E. Roman Example:
4d24d4204SJose E. Roman mpiexec -n <np> ./ex244 -fA <A_data> -fB <B_data> -orig_mat_type <type> -orig_mat_type <mat_type>
5d24d4204SJose E. Roman */
6d24d4204SJose E. Roman
7d24d4204SJose E. Roman #include <petscmat.h>
8d24d4204SJose E. Roman
main(int argc,char ** args)9d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
10d71ae5a4SJacob Faibussowitsch {
11d24d4204SJose E. Roman Mat A, Ae, B, Be;
12d24d4204SJose E. Roman PetscViewer view;
13d24d4204SJose E. Roman char file[2][PETSC_MAX_PATH_LEN];
14d24d4204SJose E. Roman PetscBool flg, flgB, isScaLAPACK, isDense, isAij, isSbaij;
15d24d4204SJose E. Roman PetscScalar one = 1.0;
16d24d4204SJose E. Roman PetscMPIInt rank, size;
17d24d4204SJose E. Roman PetscInt M, N;
18d24d4204SJose E. Roman
19327415f7SBarry Smith PetscFunctionBeginUser;
20c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, nullptr, help));
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23d24d4204SJose E. Roman
24d24d4204SJose E. Roman /* Load PETSc matrices */
259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-fA", file[0], PETSC_MAX_PATH_LEN, NULL));
269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[0], FILE_MODE_READ, &view));
279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
289566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A, "orig_"));
299566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ));
309566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
319566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view));
329566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
33d24d4204SJose E. Roman
343ba16761SJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-fB", file[1], PETSC_MAX_PATH_LEN, &flgB));
35d24d4204SJose E. Roman if (flgB) {
369566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[1], FILE_MODE_READ, &view));
379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
389566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "orig_"));
399566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ));
409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B));
419566063dSJacob Faibussowitsch PetscCall(MatLoad(B, view));
429566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
43d24d4204SJose E. Roman } else {
44d24d4204SJose E. Roman /* Create matrix B = I */
45d24d4204SJose E. Roman PetscInt rstart, rend, i;
469566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
479566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
48d24d4204SJose E. Roman
499566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
509566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "orig_"));
519566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N));
529566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ));
539566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B));
549566063dSJacob Faibussowitsch PetscCall(MatSetUp(B));
5548a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(B, 1, &i, 1, &i, &one, ADD_VALUES));
569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
58d24d4204SJose E. Roman }
59d24d4204SJose E. Roman
609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSCALAPACK, &isScaLAPACK));
61d24d4204SJose E. Roman if (isScaLAPACK) {
62d24d4204SJose E. Roman Ae = A;
63d24d4204SJose E. Roman Be = B;
64d24d4204SJose E. Roman isDense = isAij = isSbaij = PETSC_FALSE;
65d24d4204SJose E. Roman } else { /* Convert AIJ/DENSE/SBAIJ matrices into ScaLAPACK matrices */
66d24d4204SJose E. Roman if (size == 1) {
679566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQDENSE, &isDense));
689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isAij));
699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSbaij));
70d24d4204SJose E. Roman } else {
719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIDENSE, &isDense));
729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isAij));
739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &isSbaij));
74d24d4204SJose E. Roman }
75d24d4204SJose E. Roman
76dd400576SPatrick Sanan if (rank == 0) {
77d24d4204SJose E. Roman if (isDense) {
78d24d4204SJose E. Roman printf(" Convert DENSE matrices A and B into ScaLAPACK matrix... \n");
79d24d4204SJose E. Roman } else if (isAij) {
80d24d4204SJose E. Roman printf(" Convert AIJ matrices A and B into ScaLAPACK matrix... \n");
81d24d4204SJose E. Roman } else if (isSbaij) {
82d24d4204SJose E. Roman printf(" Convert SBAIJ matrices A and B into ScaLAPACK matrix... \n");
83d24d4204SJose E. Roman } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not supported yet");
84d24d4204SJose E. Roman }
859566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSCALAPACK, MAT_INITIAL_MATRIX, &Ae));
869566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATSCALAPACK, MAT_INITIAL_MATRIX, &Be));
87d24d4204SJose E. Roman
88d24d4204SJose E. Roman /* Test accuracy */
899566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, Ae, 5, &flg));
9028b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "A != A_scalapack.");
919566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B, Be, 5, &flg));
9228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "B != B_scalapack.");
93d24d4204SJose E. Roman }
94d24d4204SJose E. Roman
95d24d4204SJose E. Roman if (!isScaLAPACK) {
969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ae));
979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Be));
98d24d4204SJose E. Roman
99d24d4204SJose E. Roman /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
1009566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSCALAPACK, MAT_INPLACE_MATRIX, &A));
1019566063dSJacob Faibussowitsch //PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
102d24d4204SJose E. Roman }
103d24d4204SJose E. Roman
1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
1069566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
107b122ec5aSJacob Faibussowitsch return 0;
108d24d4204SJose E. Roman }
109d24d4204SJose E. Roman
110d24d4204SJose E. Roman /*TEST
111d24d4204SJose E. Roman
112d24d4204SJose E. Roman build:
113d24d4204SJose E. Roman requires: scalapack
114d24d4204SJose E. Roman
115d24d4204SJose E. Roman test:
116dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
117d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
118d24d4204SJose E. Roman output_file: output/ex244.out
119d24d4204SJose E. Roman
120d24d4204SJose E. Roman test:
121d24d4204SJose E. Roman suffix: 2
122d24d4204SJose E. Roman nsize: 8
123dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
124d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
125d24d4204SJose E. Roman output_file: output/ex244.out
126d24d4204SJose E. Roman
127d24d4204SJose E. Roman test:
128d24d4204SJose E. Roman suffix: 2_dense
129d24d4204SJose E. Roman nsize: 8
1308a33ae86SSatish Balay requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) !defined(PETSCTEST_VALGRIND)
131d24d4204SJose E. Roman 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
132d24d4204SJose E. Roman output_file: output/ex244_dense.out
133d24d4204SJose E. Roman
134d24d4204SJose E. Roman test:
135d24d4204SJose E. Roman suffix: 2_scalapack
136d24d4204SJose E. Roman nsize: 8
1378a33ae86SSatish Balay requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) !defined(PETSCTEST_VALGRIND)
138d24d4204SJose E. Roman 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 scalapack
139*3886731fSPierre Jolivet output_file: output/empty.out
140d24d4204SJose E. Roman
141d24d4204SJose E. Roman test:
142d24d4204SJose E. Roman suffix: 2_sbaij
143d24d4204SJose E. Roman nsize: 8
144dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
145d24d4204SJose E. Roman 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
146d24d4204SJose E. Roman output_file: output/ex244_sbaij.out
147d24d4204SJose E. Roman
148d24d4204SJose E. Roman test:
149d24d4204SJose E. Roman suffix: complex
150dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
151d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
152d24d4204SJose E. Roman output_file: output/ex244.out
153d24d4204SJose E. Roman
154d24d4204SJose E. Roman test:
155d24d4204SJose E. Roman suffix: complex_2
156d24d4204SJose E. Roman nsize: 4
157dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
158d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
159d24d4204SJose E. Roman output_file: output/ex244.out
160d24d4204SJose E. Roman
161d24d4204SJose E. Roman test:
162d24d4204SJose E. Roman suffix: dense
163dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
164d24d4204SJose E. Roman 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
165d24d4204SJose E. Roman
166d24d4204SJose E. Roman test:
167d24d4204SJose E. Roman suffix: scalapack
168dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
169d24d4204SJose E. Roman 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 scalapack
170*3886731fSPierre Jolivet output_file: output/empty.out
171d24d4204SJose E. Roman
172d24d4204SJose E. Roman test:
173d24d4204SJose E. Roman suffix: sbaij
174dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
175d24d4204SJose E. Roman 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
176d24d4204SJose E. Roman
177d24d4204SJose E. Roman TEST*/
178