1c4762a1bSJed Brown static char help[] = "Tests MatLoad() for dense matrix with uneven dimensions set in program\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **args) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown Mat A; 8c4762a1bSJed Brown PetscViewer fd; 9c4762a1bSJed Brown PetscErrorCode ierr; 10c4762a1bSJed Brown PetscMPIInt rank; 11c4762a1bSJed Brown PetscScalar *Av; 12c4762a1bSJed Brown PetscInt i; 13c4762a1bSJed Brown 14c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 15*ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 16c4762a1bSJed Brown 17c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,6,6,12,12,NULL,&A);CHKERRQ(ierr); 18c4762a1bSJed Brown ierr = MatDenseGetArray(A,&Av);CHKERRQ(ierr); 19c4762a1bSJed Brown for (i=0; i<6*12; i++) Av[i] = (PetscScalar) i; 20c4762a1bSJed Brown ierr = MatDenseRestoreArray(A,&Av);CHKERRQ(ierr); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Load matrices */ 23c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_WRITE,&fd);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = PetscViewerPushFormat(fd,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = MatView(A,fd);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = PetscViewerPopFormat(fd);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 29c4762a1bSJed Brown 30c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 32c4762a1bSJed Brown if (!rank) { 33c4762a1bSJed Brown ierr = MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 34c4762a1bSJed Brown } else { 35c4762a1bSJed Brown ierr = MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_READ,&fd);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = PetscFinalize(); 46c4762a1bSJed Brown return ierr; 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown 50c4762a1bSJed Brown /*TEST 51c4762a1bSJed Brown 52c4762a1bSJed Brown test: 53c4762a1bSJed Brown nsize: 2 54c4762a1bSJed Brown filter: grep -v alloced 55c4762a1bSJed Brown 56c4762a1bSJed Brown TEST*/ 57