1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests converting a SBAIJ matrix to BAIJ format with MatConvert. Modified from ex55.c\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown /* Usage: ./ex141 -mat_view */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown Mat C,B; 10c4762a1bSJed Brown PetscInt i,bs=2,mbs,m,block,d_nz=6,col[3]; 11c4762a1bSJed Brown PetscMPIInt size; 12c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 13c4762a1bSJed Brown PetscViewer fd; 14c4762a1bSJed Brown PetscBool equal,loadmat; 15c4762a1bSJed Brown PetscScalar value[4]; 16c4762a1bSJed Brown PetscInt ridx[2],cidx[2]; 17c4762a1bSJed Brown 18*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&loadmat)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 212c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* input matrix C */ 24c4762a1bSJed Brown if (loadmat) { 25c4762a1bSJed Brown /* Open binary file. Load a sbaij matrix, then destroy the viewer. */ 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATSEQSBAIJ)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(C,fd)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 31c4762a1bSJed Brown } else { /* Create a sbaij mat with bs>1 */ 32c4762a1bSJed Brown mbs =8; 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 34c4762a1bSJed Brown m = mbs*bs; 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATSBAIJ)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqSBAIJSetPreallocation(C,bs,d_nz,NULL)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown for (block=0; block<mbs; block++) { 44c4762a1bSJed Brown /* diagonal blocks */ 45c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 46c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 47c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,3,col,value,INSERT_VALUES)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 51c4762a1bSJed Brown 52c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 56c4762a1bSJed Brown 57c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES)); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown /* off-diagonal blocks */ 61c4762a1bSJed Brown value[0]=-1.0; value[1] = -0.1; value[2] = 0.0; value[3] = -1.0; /* row-oriented */ 62c4762a1bSJed Brown for (block=0; block<mbs-1; block++) { 63c4762a1bSJed Brown for (i=0; i<bs; i++) { 64c4762a1bSJed Brown ridx[i] = block*bs+i; cidx[i] = (block+1)*bs+i; 65c4762a1bSJed Brown } 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,bs,ridx,bs,cidx,value,INSERT_VALUES)); 67c4762a1bSJed Brown } 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* convert C to BAIJ format */ 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATSEQBAIJ,MAT_INITIAL_MATRIX,&B)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(B,C,10,&equal)); 7528b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert fails!"); 76c4762a1bSJed Brown 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 79*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 80*b122ec5aSJacob Faibussowitsch return 0; 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown /*TEST 84c4762a1bSJed Brown 85c4762a1bSJed Brown test: 86c4762a1bSJed Brown output_file: output/ex141.out 87c4762a1bSJed Brown 88c4762a1bSJed Brown TEST*/ 89