xref: /petsc/src/mat/tests/ex174.cxx (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscViewer view;
15c4762a1bSJed Brown   char        file[2][PETSC_MAX_PATH_LEN];
16c4762a1bSJed Brown   PetscBool   flg,flgB,isElemental,isDense,isAij,isSbaij;
17c4762a1bSJed Brown   PetscScalar one = 1.0;
18c4762a1bSJed Brown   PetscMPIInt rank,size;
19c4762a1bSJed Brown   PetscInt    M,N;
20c4762a1bSJed Brown 
21*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* Load PETSc matrices */
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fA",file[0],sizeof(file[0]),NULL));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&view));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(A,"orig_"));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,view));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
34c4762a1bSJed Brown 
35589a23caSBarry Smith   PetscOptionsGetString(NULL,NULL,"-fB",file[1],sizeof(file[1]),&flgB);
36c4762a1bSJed Brown   if (flgB) {
375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&view));
385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOptionsPrefix(B,"orig_"));
405f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(B,MATAIJ));
415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(B));
425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLoad(B,view));
435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&view));
44c4762a1bSJed Brown   } else {
45c4762a1bSJed Brown     /* Create matrix B = I */
46c4762a1bSJed Brown     PetscInt rstart,rend,i;
475f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(A,&M,&N));
485f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
49c4762a1bSJed Brown 
505f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
515f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOptionsPrefix(B,"orig_"));
525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N));
535f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(B,MATAIJ));
545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(B));
555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetUp(B));
56c4762a1bSJed Brown     for (i=rstart; i<rend; i++) {
575f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(B,1,&i,1,&i,&one,ADD_VALUES));
58c4762a1bSJed Brown     }
595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
605f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown 
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&isElemental));
64c4762a1bSJed Brown   if (isElemental) {
65c4762a1bSJed Brown     Ae = A;
66c4762a1bSJed Brown     Be = B;
67c4762a1bSJed Brown     isDense = isAij = isSbaij = PETSC_FALSE;
68c4762a1bSJed Brown   } else { /* Convert AIJ/DENSE/SBAIJ matrices into Elemental matrices */
69c4762a1bSJed Brown     if (size == 1) {
705f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense));
715f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAij));
725f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSbaij));
73c4762a1bSJed Brown     } else {
745f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense));
755f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAij));
765f80ce2aSJacob Faibussowitsch        CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&isSbaij));
77c4762a1bSJed Brown     }
78c4762a1bSJed Brown 
79dd400576SPatrick Sanan     if (rank == 0) {
80c4762a1bSJed Brown       if (isDense) {
81c4762a1bSJed Brown         printf(" Convert DENSE matrices A and B into Elemental matrix... \n");
82c4762a1bSJed Brown       } else if (isAij) {
83c4762a1bSJed Brown         printf(" Convert AIJ matrices A and B into Elemental matrix... \n");
84c4762a1bSJed Brown       } else if (isSbaij) {
85c4762a1bSJed Brown         printf(" Convert SBAIJ matrices A and B into Elemental matrix... \n");
86c4762a1bSJed Brown       } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not supported yet");
87c4762a1bSJed Brown     }
885f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae));
895f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &Be));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown     /* Test accuracy */
925f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,Ae,5,&flg));
9328b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != A_elemental.");
945f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(B,Be,5,&flg));
9528b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"B != B_elemental.");
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   if (!isElemental) {
995f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Ae));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Be));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
1035f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A));
1045f80ce2aSJacob Faibussowitsch     //CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
109*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
110*b122ec5aSJacob Faibussowitsch   return 0;
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown /*TEST
114c4762a1bSJed Brown 
115c4762a1bSJed Brown    build:
116c4762a1bSJed Brown       requires: elemental
117c4762a1bSJed Brown 
118c4762a1bSJed Brown    test:
119dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
120c4762a1bSJed Brown       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
121c4762a1bSJed Brown       output_file: output/ex174.out
122c4762a1bSJed Brown 
123c4762a1bSJed Brown    test:
124c4762a1bSJed Brown       suffix: 2
125c4762a1bSJed Brown       nsize: 8
126dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
127c4762a1bSJed Brown       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
128c4762a1bSJed Brown       output_file: output/ex174.out
129c4762a1bSJed Brown 
130c4762a1bSJed Brown    test:
131c4762a1bSJed Brown       suffix: 2_dense
132c4762a1bSJed Brown       nsize: 8
133dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
134c4762a1bSJed 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
135c4762a1bSJed Brown       output_file: output/ex174_dense.out
136c4762a1bSJed Brown 
137c4762a1bSJed Brown    test:
138c4762a1bSJed Brown       suffix: 2_elemental
139c4762a1bSJed Brown       nsize: 8
140dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
141c4762a1bSJed 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
142c4762a1bSJed Brown       output_file: output/ex174_elemental.out
143c4762a1bSJed Brown 
144c4762a1bSJed Brown    test:
145c4762a1bSJed Brown       suffix: 2_sbaij
146c4762a1bSJed Brown       nsize: 8
147dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
148c4762a1bSJed 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
149c4762a1bSJed Brown       output_file: output/ex174_sbaij.out
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    test:
152c4762a1bSJed Brown       suffix: complex
153dfd57a17SPierre Jolivet       requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
154c4762a1bSJed Brown       args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
155c4762a1bSJed Brown       output_file: output/ex174.out
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    test:
158c4762a1bSJed Brown       suffix: complex_2
159c4762a1bSJed Brown       nsize: 4
160dfd57a17SPierre Jolivet       requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
161c4762a1bSJed Brown       args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
162c4762a1bSJed Brown       output_file: output/ex174.out
163c4762a1bSJed Brown 
164c4762a1bSJed Brown    test:
165c4762a1bSJed Brown       suffix: dense
166dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
167c4762a1bSJed 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
168c4762a1bSJed Brown 
169c4762a1bSJed Brown    test:
170c4762a1bSJed Brown       suffix: elemental
171dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
172c4762a1bSJed 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
173c4762a1bSJed Brown 
174c4762a1bSJed Brown    test:
175c4762a1bSJed Brown       suffix: sbaij
176dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
177c4762a1bSJed 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
178c4762a1bSJed Brown 
179c4762a1bSJed Brown TEST*/
180