Actual source code: ex174.cxx

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: static char help[] = "Tests MatConvert(), MatLoad() for MATELEMENTAL interface.\n\n";
  3: /*
  4:  Example:
  5:    mpiexec -n <np> ./ex173 -fA <A_data> -fB <B_data> -orig_mat_type <type> -orig_mat_type <mat_type>
  6: */

  8:  #include <petscmat.h>
  9:  #include <petscmatelemental.h>

 11: int main(int argc,char **args)
 12: {
 13:   Mat            A,Ae,B,Be;
 15:   PetscViewer    view;
 16:   char           file[2][PETSC_MAX_PATH_LEN];
 17:   PetscBool      flg,flgB,isElemental,isDense,isAij,isSbaij;
 18:   PetscScalar    one = 1.0;
 19:   PetscMPIInt    rank,size;
 20:   PetscInt       M,N;

 22:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 23:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 24:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 26:   /* Load PETSc matrices */
 27:   PetscOptionsGetString(NULL,NULL,"-fA",file[0],PETSC_MAX_PATH_LEN,NULL);
 28:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&view);
 29:   MatCreate(PETSC_COMM_WORLD,&A);
 30:   MatSetOptionsPrefix(A,"orig_");
 31:   MatSetType(A,MATAIJ); 
 32:   MatSetFromOptions(A);
 33:   MatLoad(A,view);
 34:   PetscViewerDestroy(&view);

 36:   PetscOptionsGetString(NULL,NULL,"-fB",file[1],PETSC_MAX_PATH_LEN,&flgB);
 37:   if (flgB) {
 38:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&view);
 39:     MatCreate(PETSC_COMM_WORLD,&B);
 40:     MatSetOptionsPrefix(B,"orig_");
 41:     MatSetType(B,MATAIJ);
 42:     MatSetFromOptions(B);
 43:     MatLoad(B,view);
 44:     PetscViewerDestroy(&view);
 45:   } else {
 46:     /* Create matrix B = I */
 47:     PetscInt rstart,rend,i;
 48:     MatGetSize(A,&M,&N);
 49:     MatGetOwnershipRange(A,&rstart,&rend);

 51:     MatCreate(PETSC_COMM_WORLD,&B);
 52:     MatSetOptionsPrefix(B,"orig_");
 53:     MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
 54:     MatSetType(B,MATAIJ);
 55:     MatSetFromOptions(B);
 56:     MatSetUp(B);
 57:     for (i=rstart; i<rend; i++) {
 58:       MatSetValues(B,1,&i,1,&i,&one,ADD_VALUES);
 59:     }
 60:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 61:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 62:   }

 64:   PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&isElemental);
 65:   if (isElemental) {
 66:     Ae = A;
 67:     Be = B;
 68:     isDense = isAij = isSbaij = PETSC_FALSE;
 69:   } else { /* Convert AIJ/DENSE/SBAIJ matrices into Elemental matrices */
 70:     if (size == 1) {
 71:       PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense);
 72:       PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAij);
 73:       PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSbaij);
 74:     } else {
 75:       PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense);
 76:       PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAij);
 77:        PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&isSbaij);
 78:     }

 80:     if (!rank) {
 81:       if (isDense) {
 82:         printf(" Convert DENSE matrices A and B into Elemental matrix... \n");
 83:       } else if (isAij) {
 84:         printf(" Convert AIJ matrices A and B into Elemental matrix... \n");
 85:       } else if (isSbaij) {
 86:         printf(" Convert SBAIJ matrices A and B into Elemental matrix... \n");
 87:       } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not supported yet");
 88:     }
 89:     MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae);
 90:     MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &Be);

 92:     /* Test accuracy */
 93:     MatMultEqual(A,Ae,5,&flg);
 94:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != A_elemental.");
 95:     MatMultEqual(B,Be,5,&flg);
 96:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"B != B_elemental.");
 97:   }

 99:   if (!isElemental) {
100:     MatDestroy(&Ae);
101:     MatDestroy(&Be);

103:     /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
104:     MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A);
105:     //MatView(A,PETSC_VIEWER_STDOUT_WORLD);
106:   }

108:   MatDestroy(&A);
109:   MatDestroy(&B);
110:   PetscFinalize();
111:   return ierr;
112: }


115: /*TEST

117:    build:
118:       requires: elemental

120:    test:
121:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
122:       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
123:       output_file: output/ex174.out

125:    test:
126:       suffix: 2
127:       nsize: 8
128:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
129:       args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij
130:       output_file: output/ex174.out

132:    test:
133:       suffix: 2_dense
134:       nsize: 8
135:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
136:       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
137:       output_file: output/ex174_dense.out

139:    test:
140:       suffix: 2_elemental
141:       nsize: 8
142:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
143:       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
144:       output_file: output/ex174_elemental.out

146:    test:
147:       suffix: 2_sbaij
148:       nsize: 8
149:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
150:       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
151:       output_file: output/ex174_sbaij.out

153:    test:
154:       suffix: complex
155:       requires: complex double datafilespath !define(PETSC_USE_64BIT_INDICES)
156:       args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
157:       output_file: output/ex174.out

159:    test:
160:       suffix: complex_2
161:       nsize: 4
162:       requires: complex double datafilespath !define(PETSC_USE_64BIT_INDICES)
163:       args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905
164:       output_file: output/ex174.out

166:    test:
167:       suffix: dense
168:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
169:       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

171:    test:
172:       suffix: elemental
173:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
174:       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

176:    test:
177:       suffix: sbaij
178:       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
179:       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

181: TEST*/