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 PetscMPIInt rank; 10c4762a1bSJed Brown PetscScalar *Av; 11c4762a1bSJed Brown PetscInt i; 12c4762a1bSJed Brown 13*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 145f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 15c4762a1bSJed Brown 165f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,6,6,12,12,NULL,&A)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(A,&Av)); 18c4762a1bSJed Brown for (i=0; i<6*12; i++) Av[i] = (PetscScalar) i; 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(A,&Av)); 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* Load matrices */ 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_WRITE,&fd)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(fd,PETSC_VIEWER_NATIVE)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,fd)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(fd)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 28c4762a1bSJed Brown 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATDENSE)); 31dd400576SPatrick Sanan if (rank == 0) { 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE)); 33c4762a1bSJed Brown } else { 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE)); 35c4762a1bSJed Brown } 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_READ,&fd)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 44*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 45*b122ec5aSJacob Faibussowitsch return 0; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown /*TEST 49c4762a1bSJed Brown 50c4762a1bSJed Brown test: 51c4762a1bSJed Brown nsize: 2 52c4762a1bSJed Brown filter: grep -v alloced 53c4762a1bSJed Brown 54c4762a1bSJed Brown TEST*/ 55