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