1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Reads a PETSc matrix and vector from a file and reorders it.\n\ 3c4762a1bSJed Brown -f0 <input_file> : first file to load (small system)\n\ 4c4762a1bSJed Brown -f1 <input_file> : second file to load (larger system)\n\n"; 5c4762a1bSJed Brown 6c4762a1bSJed Brown /*T 7c4762a1bSJed Brown Concepts: Mat^ordering a matrix - loading a binary matrix and vector; 8c4762a1bSJed Brown Concepts: Mat^loading a binary matrix and vector; 9c4762a1bSJed Brown Concepts: Vectors^loading a binary vector; 10c4762a1bSJed Brown Concepts: PetscLog^preloading executable 11c4762a1bSJed Brown Processors: 1 12c4762a1bSJed Brown T*/ 13c4762a1bSJed Brown 14c4762a1bSJed Brown /* 15c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. 16c4762a1bSJed Brown automatically includes: 17c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 18c4762a1bSJed Brown petscmat.h - matrices 19c4762a1bSJed Brown petscis.h - index sets petscviewer.h - viewers 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown #include <petscmat.h> 22c4762a1bSJed Brown 23c4762a1bSJed Brown int main(int argc,char **args) 24c4762a1bSJed Brown { 25c4762a1bSJed Brown Mat A; /* matrix */ 26c4762a1bSJed Brown PetscViewer fd; /* viewer */ 27c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 28c4762a1bSJed Brown IS isrow,iscol; /* row and column permutations */ 29c4762a1bSJed Brown PetscErrorCode ierr; 30c4762a1bSJed Brown MatOrderingType rtype = MATORDERINGRCM; 31c4762a1bSJed Brown PetscBool flg; 32c4762a1bSJed Brown 33c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown Determine files from which we read the two linear systems 36c4762a1bSJed Brown (matrix and right-hand-side vector). 37c4762a1bSJed Brown */ 38589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr); 39c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* ----------------------------------------------------------- 42c4762a1bSJed Brown Beginning of loop 43c4762a1bSJed Brown ----------------------------------------------------------- */ 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown Loop through the reordering 2 times. 46c4762a1bSJed Brown - The intention here is to preload and solve a small system; 47c4762a1bSJed Brown then load another (larger) system and solve it as well. 48c4762a1bSJed Brown This process preloads the instructions with the smaller 49c4762a1bSJed Brown system so that more accurate performance monitoring (via 50c4762a1bSJed Brown -log_view) can be done with the larger one (that actually 51c4762a1bSJed Brown is the system of interest). 52c4762a1bSJed Brown */ 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* - - - - - - - - - - - New Stage - - - - - - - - - - - - - 55c4762a1bSJed Brown Load system i 56c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* 59c4762a1bSJed Brown Open binary file. Note that we use FILE_MODE_READ to indicate 60c4762a1bSJed Brown reading from this file. 61c4762a1bSJed Brown */ 62c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown Load the matrix; then destroy the viewer. 66c4762a1bSJed Brown */ 67c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 71c4762a1bSJed Brown 72c4762a1bSJed Brown ierr = MatGetOrdering(A,rtype,&isrow,&iscol);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = ISView(isrow,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* 76c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 77c4762a1bSJed Brown are no longer needed. 78c4762a1bSJed Brown */ 79c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = ISDestroy(&iscol);CHKERRQ(ierr); 82c4762a1bSJed Brown 83c4762a1bSJed Brown ierr = PetscFinalize(); 84c4762a1bSJed Brown return ierr; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87c4762a1bSJed Brown /*TEST 88c4762a1bSJed Brown 89c4762a1bSJed Brown test: 90*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 91c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/medium -viewer_binary_skip_info 92c4762a1bSJed Brown 93c4762a1bSJed Brown TEST*/ 94