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