1 2 static char help[] = "Tests converting a SBAIJ matrix to BAIJ format with MatConvert. Modified from ex55.c\n\n"; 3 4 #include <petscmat.h> 5 /* Usage: ./ex141 -mat_view */ 6 7 int main(int argc,char **args) 8 { 9 Mat C,B; 10 PetscErrorCode ierr; 11 PetscInt i,bs=2,mbs,m,block,d_nz=6,col[3]; 12 PetscMPIInt size; 13 char file[PETSC_MAX_PATH_LEN]; 14 PetscViewer fd; 15 PetscBool equal,loadmat; 16 PetscScalar value[4]; 17 PetscInt ridx[2],cidx[2]; 18 19 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20 ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&loadmat);CHKERRQ(ierr); 21 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 22 if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 23 24 /* input matrix C */ 25 if (loadmat) { 26 /* Open binary file. Load a sbaij matrix, then destroy the viewer. */ 27 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 28 ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 29 ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 30 ierr = MatLoad(C,fd);CHKERRQ(ierr); 31 ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 32 } else { /* Create a sbaij mat with bs>1 */ 33 mbs =8; 34 ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr); 35 m = mbs*bs; 36 ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 37 ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr); 38 ierr = MatSetType(C,MATSBAIJ);CHKERRQ(ierr); 39 ierr = MatSetFromOptions(C);CHKERRQ(ierr); 40 ierr = MatSeqSBAIJSetPreallocation(C,bs,d_nz,NULL);CHKERRQ(ierr); 41 ierr = MatSetUp(C);CHKERRQ(ierr); 42 ierr = MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 43 44 for (block=0; block<mbs; block++) { 45 /* diagonal blocks */ 46 value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 47 for (i=1+block*bs; i<bs-1+block*bs; i++) { 48 col[0] = i-1; col[1] = i; col[2] = i+1; 49 ierr = MatSetValues(C,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 50 } 51 i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 52 53 value[0]=-1.0; value[1]=4.0; 54 ierr = MatSetValues(C,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 55 56 i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 57 58 value[0]=4.0; value[1] = -1.0; 59 ierr = MatSetValues(C,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 60 } 61 /* off-diagonal blocks */ 62 value[0]=-1.0; value[1] = -0.1; value[2] = 0.0; value[3] = -1.0; /* row-oriented */ 63 for (block=0; block<mbs-1; block++) { 64 for (i=0; i<bs; i++) { 65 ridx[i] = block*bs+i; cidx[i] = (block+1)*bs+i; 66 } 67 ierr = MatSetValues(C,bs,ridx,bs,cidx,value,INSERT_VALUES);CHKERRQ(ierr); 68 } 69 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71 } 72 73 /* convert C to BAIJ format */ 74 ierr = MatConvert(C,MATSEQBAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 75 ierr = MatMultEqual(B,C,10,&equal);CHKERRQ(ierr); 76 if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert fails!"); 77 78 ierr = MatDestroy(&B);CHKERRQ(ierr); 79 ierr = MatDestroy(&C);CHKERRQ(ierr); 80 ierr = PetscFinalize(); 81 return ierr; 82 } 83 84 /*TEST 85 86 test: 87 output_file: output/ex141.out 88 89 TEST*/ 90