1c4762a1bSJed Brown static char help[] = "Tests MatConvert(), MatLoad() for MATELEMENTAL interface.\n\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown Example:
4c4762a1bSJed Brown mpiexec -n <np> ./ex173 -fA <A_data> -fB <B_data> -orig_mat_type <type> -orig_mat_type <mat_type>
5c4762a1bSJed Brown */
6c4762a1bSJed Brown
7c4762a1bSJed Brown #include <petscmat.h>
8c4762a1bSJed Brown #include <petscmatelemental.h>
9c4762a1bSJed Brown
main(int argc,char ** args)10d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
11d71ae5a4SJacob Faibussowitsch {
12c4762a1bSJed Brown Mat A, Ae, B, Be;
13c4762a1bSJed Brown PetscViewer view;
14c4762a1bSJed Brown char file[2][PETSC_MAX_PATH_LEN];
15c4762a1bSJed Brown PetscBool flg, flgB, isElemental, isDense, isAij, isSbaij;
16c4762a1bSJed Brown PetscScalar one = 1.0;
17c4762a1bSJed Brown PetscMPIInt rank, size;
18c4762a1bSJed Brown PetscInt M, N;
19c4762a1bSJed Brown
20327415f7SBarry Smith PetscFunctionBeginUser;
2107e90273SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, nullptr, help));
229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24c4762a1bSJed Brown
25c4762a1bSJed Brown /* Load PETSc matrices */
269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-fA", file[0], sizeof(file[0]), NULL));
279566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[0], FILE_MODE_READ, &view));
289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
299566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A, "orig_"));
309566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ));
319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
329566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view));
339566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
34c4762a1bSJed Brown
353ba16761SJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-fB", file[1], sizeof(file[1]), &flgB));
36c4762a1bSJed Brown if (flgB) {
379566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[1], FILE_MODE_READ, &view));
389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
399566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "orig_"));
409566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ));
419566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B));
429566063dSJacob Faibussowitsch PetscCall(MatLoad(B, view));
439566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
44c4762a1bSJed Brown } else {
45c4762a1bSJed Brown /* Create matrix B = I */
46c4762a1bSJed Brown PetscInt rstart, rend, i;
479566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
489566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
49c4762a1bSJed Brown
509566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
519566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "orig_"));
529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N));
539566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ));
549566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B));
559566063dSJacob Faibussowitsch PetscCall(MatSetUp(B));
5648a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(B, 1, &i, 1, &i, &one, ADD_VALUES));
579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
59c4762a1bSJed Brown }
60c4762a1bSJed Brown
619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATELEMENTAL, &isElemental));
62c4762a1bSJed Brown if (isElemental) {
63c4762a1bSJed Brown Ae = A;
64c4762a1bSJed Brown Be = B;
65c4762a1bSJed Brown isDense = isAij = isSbaij = PETSC_FALSE;
66c4762a1bSJed Brown } else { /* Convert AIJ/DENSE/SBAIJ matrices into Elemental matrices */
67c4762a1bSJed Brown if (size == 1) {
689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQDENSE, &isDense));
699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isAij));
709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSbaij));
71c4762a1bSJed Brown } else {
729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIDENSE, &isDense));
739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isAij));
749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &isSbaij));
75c4762a1bSJed Brown }
76c4762a1bSJed Brown
77dd400576SPatrick Sanan if (rank == 0) {
78c4762a1bSJed Brown if (isDense) {
79c4762a1bSJed Brown printf(" Convert DENSE matrices A and B into Elemental matrix... \n");
80c4762a1bSJed Brown } else if (isAij) {
81c4762a1bSJed Brown printf(" Convert AIJ matrices A and B into Elemental matrix... \n");
82c4762a1bSJed Brown } else if (isSbaij) {
83c4762a1bSJed Brown printf(" Convert SBAIJ matrices A and B into Elemental matrix... \n");
84c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not supported yet");
85c4762a1bSJed Brown }
869566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae));
879566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &Be));
88c4762a1bSJed Brown
89c4762a1bSJed Brown /* Test accuracy */
909566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, Ae, 5, &flg));
9128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "A != A_elemental.");
929566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B, Be, 5, &flg));
9328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "B != B_elemental.");
94c4762a1bSJed Brown }
95c4762a1bSJed Brown
96c4762a1bSJed Brown if (!isElemental) {
979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ae));
989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Be));
99c4762a1bSJed Brown
100c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
1019566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A));
1029566063dSJacob Faibussowitsch //PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
103c4762a1bSJed Brown }
104c4762a1bSJed Brown
1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
1079566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
108b122ec5aSJacob Faibussowitsch return 0;
109c4762a1bSJed Brown }
110c4762a1bSJed Brown
111c4762a1bSJed Brown /*TEST
112c4762a1bSJed Brown
113c4762a1bSJed Brown build:
114c4762a1bSJed Brown requires: elemental
115c4762a1bSJed Brown
116c4762a1bSJed Brown test:
117dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
118c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
119c4762a1bSJed Brown output_file: output/ex174.out
120c4762a1bSJed Brown
121c4762a1bSJed Brown test:
122c4762a1bSJed Brown suffix: 2
123c4762a1bSJed Brown nsize: 8
124dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
125c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
126c4762a1bSJed Brown output_file: output/ex174.out
127c4762a1bSJed Brown
128c4762a1bSJed Brown test:
129c4762a1bSJed Brown suffix: 2_dense
130c4762a1bSJed Brown nsize: 8
131dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
132c4762a1bSJed Brown 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
133c4762a1bSJed Brown output_file: output/ex174_dense.out
134c4762a1bSJed Brown
135c4762a1bSJed Brown test:
136c4762a1bSJed Brown suffix: 2_elemental
137c4762a1bSJed Brown nsize: 8
138dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
139c4762a1bSJed Brown 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*3886731fSPierre Jolivet output_file: output/empty.out
141c4762a1bSJed Brown
142c4762a1bSJed Brown test:
143c4762a1bSJed Brown suffix: 2_sbaij
144c4762a1bSJed Brown nsize: 8
145dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
146c4762a1bSJed Brown 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
147c4762a1bSJed Brown output_file: output/ex174_sbaij.out
148c4762a1bSJed Brown
149c4762a1bSJed Brown test:
150c4762a1bSJed Brown suffix: complex
151dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
152c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
153c4762a1bSJed Brown output_file: output/ex174.out
154c4762a1bSJed Brown
155c4762a1bSJed Brown test:
156c4762a1bSJed Brown suffix: complex_2
157c4762a1bSJed Brown nsize: 4
158dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
159c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
160c4762a1bSJed Brown output_file: output/ex174.out
161c4762a1bSJed Brown
162c4762a1bSJed Brown test:
163c4762a1bSJed Brown suffix: dense
164dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
165c4762a1bSJed Brown 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
166c4762a1bSJed Brown
167c4762a1bSJed Brown test:
168c4762a1bSJed Brown suffix: elemental
169dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
170c4762a1bSJed Brown 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*3886731fSPierre Jolivet output_file: output/empty.out
172c4762a1bSJed Brown
173c4762a1bSJed Brown test:
174c4762a1bSJed Brown suffix: sbaij
175dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
176c4762a1bSJed Brown 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
177c4762a1bSJed Brown
178c4762a1bSJed Brown TEST*/
179