1 2 static char help[] = "Tests MatGetColumnNorms()/Sums()/Means() for matrix read from file."; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat A; 9 PetscErrorCode ierr; 10 PetscReal *reductions_real; 11 PetscScalar *reductions_scalar; 12 char file[PETSC_MAX_PATH_LEN]; 13 PetscBool flg; 14 PetscViewer fd; 15 PetscInt n; 16 PetscMPIInt rank; 17 18 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 20 ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr); 21 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 22 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 23 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 24 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 25 ierr = MatLoad(A,fd);CHKERRQ(ierr); 26 ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 27 28 ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 29 ierr = PetscMalloc1(n,&reductions_real);CHKERRQ(ierr); 30 ierr = PetscMalloc1(n,&reductions_scalar);CHKERRQ(ierr); 31 32 ierr = MatGetColumnNorms(A,NORM_2,reductions_real);CHKERRQ(ierr); 33 if (rank == 0) { 34 ierr = PetscPrintf(PETSC_COMM_SELF,"NORM_2:\n");CHKERRQ(ierr); 35 ierr = PetscRealView(n,reductions_real,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 36 } 37 38 ierr = MatGetColumnNorms(A,NORM_1,reductions_real);CHKERRQ(ierr); 39 if (rank == 0) { 40 ierr = PetscPrintf(PETSC_COMM_SELF,"NORM_1:\n");CHKERRQ(ierr); 41 ierr = PetscRealView(n,reductions_real,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 42 } 43 44 ierr = MatGetColumnNorms(A,NORM_INFINITY,reductions_real);CHKERRQ(ierr); 45 if (rank == 0) { 46 ierr = PetscPrintf(PETSC_COMM_SELF,"NORM_INFINITY:\n");CHKERRQ(ierr); 47 ierr = PetscRealView(n,reductions_real,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 48 } 49 50 ierr = MatGetColumnSums(A,reductions_scalar);CHKERRQ(ierr); 51 if (!rank) { 52 ierr = PetscPrintf(PETSC_COMM_SELF,"REDUCTION_SUM:\n");CHKERRQ(ierr); 53 ierr = PetscScalarView(n,reductions_scalar,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 54 } 55 56 ierr = MatGetColumnMeans(A,reductions_scalar);CHKERRQ(ierr); 57 if (!rank) { 58 ierr = PetscPrintf(PETSC_COMM_SELF,"REDUCTION_MEAN:\n");CHKERRQ(ierr); 59 ierr = PetscScalarView(n,reductions_scalar,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 60 } 61 62 ierr = PetscFree(reductions_real);CHKERRQ(ierr); 63 ierr = PetscFree(reductions_scalar);CHKERRQ(ierr); 64 ierr = MatDestroy(&A);CHKERRQ(ierr); 65 ierr = PetscFinalize(); 66 return ierr; 67 } 68 69 /*TEST 70 71 test: 72 suffix: 1 73 nsize: 2 74 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 75 args: -f ${DATAFILESPATH}/matrices/small -mat_type aij 76 output_file: output/ex138.out 77 78 test: 79 suffix: 2 80 nsize: {{1 2}} 81 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 82 args: -f ${DATAFILESPATH}/matrices/small -mat_type baij -matload_block_size {{2 3}} 83 output_file: output/ex138.out 84 85 test: 86 suffix: complex 87 nsize: 2 88 requires: datafilespath complex double !defined(PETSC_USE_64BIT_INDICES) 89 args: -f ${DATAFILESPATH}/matrices/nimrod/small_112905 -mat_type aij 90 output_file: output/ex138_complex.out 91 filter: grep -E "\ 0:|1340:|1344:" 92 93 TEST*/ 94