1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Test the use of MatSetValuesBlocked(), MatZeroRows() for rectangular MatBAIJ matrix, test MatSetValuesBlocked() for MatSBAIJ matrix (-test_mat_sbaij)."; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A; 9c4762a1bSJed Brown PetscInt bs=3,m=4,n=6,i,j,val = 10,row[2],col[3],eval,rstart; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown PetscMPIInt size,rank; 12c4762a1bSJed Brown PetscScalar x[6][9],y[3][3],one=1.0; 13c4762a1bSJed Brown PetscBool flg,testsbaij=PETSC_FALSE; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 16*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 17*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 18c4762a1bSJed Brown 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_mat_sbaij",&testsbaij)); 20c4762a1bSJed Brown 21c4762a1bSJed Brown if (testsbaij) { 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A)); 23c4762a1bSJed Brown } else { 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A)); 25c4762a1bSJed Brown } 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 27c4762a1bSJed Brown eval = 9; 28c4762a1bSJed Brown 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-ass_extern",&flg)); 30c4762a1bSJed Brown if (flg && (size != 1)) rstart = m*((rank+1)%size); 31c4762a1bSJed Brown else rstart = m*(rank); 32c4762a1bSJed Brown 33c4762a1bSJed Brown row[0] =rstart+0; row[1] =rstart+2; 34c4762a1bSJed Brown col[0] =rstart+0; col[1] =rstart+1; col[2] =rstart+3; 35c4762a1bSJed Brown for (i=0; i<6; i++) { 36c4762a1bSJed Brown for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++; 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES)); 40c4762a1bSJed Brown 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown This option does not work for rectangular matrices 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 47c4762a1bSJed Brown */ 48c4762a1bSJed Brown 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Do another MatSetValues to test the case when only one local block is specified */ 52c4762a1bSJed Brown for (i=0; i<3; i++) { 53c4762a1bSJed Brown for (j =0; j<3; j++) y[i][j] = (PetscScalar)(10 + i*eval + j); 54c4762a1bSJed Brown } 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,1,row,1,col,&y[0][0],INSERT_VALUES)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 58c4762a1bSJed Brown 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-zero_rows",&flg)); 60c4762a1bSJed Brown if (flg) { 61c4762a1bSJed Brown col[0] = rstart*bs+0; 62c4762a1bSJed Brown col[1] = rstart*bs+1; 63c4762a1bSJed Brown col[2] = rstart*bs+2; 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRows(A,3,col,one,0,0)); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 68c4762a1bSJed Brown 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 70c4762a1bSJed Brown ierr = PetscFinalize(); 71c4762a1bSJed Brown return ierr; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /*TEST 75c4762a1bSJed Brown 76c4762a1bSJed Brown test: 77c4762a1bSJed Brown filter: grep -v "MPI processes" 78c4762a1bSJed Brown 79c4762a1bSJed Brown test: 80c4762a1bSJed Brown suffix: 4 81c4762a1bSJed Brown nsize: 3 82c4762a1bSJed Brown args: -ass_extern 83c4762a1bSJed Brown filter: grep -v "MPI processes" 84c4762a1bSJed Brown 85c4762a1bSJed Brown test: 86c4762a1bSJed Brown suffix: 5 87c4762a1bSJed Brown nsize: 3 88c4762a1bSJed Brown args: -ass_extern -zero_rows 89c4762a1bSJed Brown filter: grep -v "MPI processes" 90c4762a1bSJed Brown 91c4762a1bSJed Brown TEST*/ 92