1 2 static char help[] = "Tests various routines in MatMPIBAIJ format.\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat A; 9 PetscInt m=2,bs=1,M,row,col,start,end,i,j,k; 10 PetscErrorCode ierr; 11 PetscMPIInt rank,size; 12 PetscScalar data=100; 13 PetscBool flg; 14 15 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 16 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 17 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 18 19 /* Test MatSetValues() and MatGetValues() */ 20 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 21 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL)); 22 23 M = m*bs*size; 24 CHKERRQ(MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,M,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A)); 25 26 CHKERRQ(MatGetOwnershipRange(A,&start,&end)); 27 CHKERRQ(PetscOptionsHasName(NULL,NULL,"-column_oriented",&flg)); 28 if (flg) { 29 CHKERRQ(MatSetOption(A,MAT_ROW_ORIENTED,PETSC_FALSE)); 30 } 31 32 /* inproc assembly */ 33 for (row=start; row<end; row++) { 34 for (col=start; col<end; col++,data+=1) { 35 CHKERRQ(MatSetValues(A,1,&row,1,&col,&data,INSERT_VALUES)); 36 } 37 } 38 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 39 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 40 41 /* offproc assembly */ 42 data = 5.0; 43 row = (M+start-1)%M; 44 for (col=0; col<M; col++) { 45 CHKERRQ(MatSetValues(A,1,&row,1,&col,&data,ADD_VALUES)); 46 } 47 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 48 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 49 50 /* Test MatSetValuesBlocked() */ 51 CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_setvaluesblocked",&flg)); 52 if (flg) { 53 PetscScalar *bval; 54 row /= bs; 55 col = start/bs; 56 CHKERRQ(PetscMalloc1(bs*bs,&bval)); 57 k = 1; 58 /* row oriented - default */ 59 for (i=0; i<bs; i++) { 60 for (j=0; j<bs; j++) { 61 bval[i*bs+j] = (PetscScalar)k; k++; 62 } 63 } 64 CHKERRQ(MatSetValuesBlocked(A,1,&row,1,&col,bval,INSERT_VALUES)); 65 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 66 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 67 CHKERRQ(PetscFree(bval)); 68 } 69 70 CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 71 CHKERRQ(MatDestroy(&A)); 72 ierr = PetscFinalize(); 73 return ierr; 74 } 75 76 /*TEST 77 78 test: 79 suffix: 1 80 nsize: 3 81 args: -mat_block_size 2 -test_setvaluesblocked 82 83 test: 84 suffix: 2 85 nsize: 3 86 args: -mat_block_size 2 -test_setvaluesblocked -column_oriented 87 88 test: 89 suffix: 3 90 nsize: 3 91 args: -mat_block_size 1 -test_setvaluesblocked 92 93 test: 94 suffix: 4 95 nsize: 3 96 args: -mat_block_size 1 -test_setvaluesblocked -column_oriented 97 98 TEST*/ 99