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