1c4762a1bSJed Brown static char help[] = "Tests MatMPI{AIJ,BAIJ,SBAIJ}SetPreallocationCSR\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **argv) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown PetscErrorCode ierr; 8c4762a1bSJed Brown PetscInt ia[3]={0,2,4}; 9c4762a1bSJed Brown PetscInt ja[4]={0,1,0,1}; 10c4762a1bSJed Brown PetscScalar c[16]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15}; 11c4762a1bSJed Brown PetscInt ia2[5]={0,4,8,12,16}; 12c4762a1bSJed Brown PetscInt ja2[16]={0,1,2,3,0,1,2,3,0,1,2,3,0,1,2,3}; 13c4762a1bSJed Brown PetscScalar c2[16]={0,1,4,5,2,3,6,7,8,9,12,13,10,11,14,15}; 14c4762a1bSJed Brown PetscMPIInt size,rank; 15c4762a1bSJed Brown Mat ssbaij; 16c4762a1bSJed Brown PetscBool rect = PETSC_FALSE; 17c4762a1bSJed Brown 18c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 19ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 20ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 21*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size < 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is an example with more then one processors"); 22c4762a1bSJed Brown if (rank) { 23c4762a1bSJed Brown PetscInt i; 24c4762a1bSJed Brown for (i = 0; i < 3; i++) ia[i] = 0; 25c4762a1bSJed Brown for (i = 0; i < 5; i++) ia2[i] = 0; 26c4762a1bSJed Brown } 27c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-rect",&rect,NULL);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&ssbaij);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSetBlockSize(ssbaij,2);CHKERRQ(ierr); 30c4762a1bSJed Brown if (rect) { 31c4762a1bSJed Brown ierr = MatSetType(ssbaij,MATMPIBAIJ);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatSetSizes(ssbaij,4,6,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 33c4762a1bSJed Brown } else { 34c4762a1bSJed Brown ierr = MatSetType(ssbaij,MATMPISBAIJ);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MatSetSizes(ssbaij,4,4,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown ierr = MatSetFromOptions(ssbaij);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocationCSR(ssbaij,ia2,ja2,c2);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatMPIBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatMPISBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatViewFromOptions(ssbaij,NULL,"-view");CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = MatDestroy(&ssbaij);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscFinalize(); 44c4762a1bSJed Brown return ierr; 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown /*TEST 48c4762a1bSJed Brown 49c4762a1bSJed Brown test: 50c4762a1bSJed Brown filter: grep -v type | sed -e "s/\.//g" 51c4762a1bSJed Brown suffix: aijbaij_csr 52c4762a1bSJed Brown nsize: 2 53c4762a1bSJed Brown args: -mat_type {{aij baij}} -view -rect {{0 1}} 54c4762a1bSJed Brown 55c4762a1bSJed Brown test: 56c4762a1bSJed Brown filter: sed -e "s/\.//g" 57c4762a1bSJed Brown suffix: sbaij_csr 58c4762a1bSJed Brown nsize: 2 59c4762a1bSJed Brown args: -mat_type sbaij -view -rect {{0 1}} 60c4762a1bSJed Brown 61c4762a1bSJed Brown TEST*/ 62