Actual source code: ex52.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Tests the vatious routines in MatMPIBAIJ format.\n";


  5: #include <petscmat.h>

  9: int main(int argc,char **args)
 10: {
 11:   Mat            A;
 12:   PetscInt       m=2,bs=1,M,row,col,start,end,i,j,k;
 14:   PetscMPIInt    rank,size;
 15:   PetscScalar    data=100;
 16:   PetscBool      flg;

 18:   PetscInitialize(&argc,&args,(char*)0,help);
 19:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 20:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 22:   /* Test MatSetValues() and MatGetValues() */
 23:   PetscOptionsGetInt(NULL,"-mat_block_size",&bs,NULL);
 24:   PetscOptionsGetInt(NULL,"-mat_size",&m,NULL);

 26:   M    = m*bs*size;
 27:   MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,M,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A);

 29:   MatGetOwnershipRange(A,&start,&end);
 30:   PetscOptionsHasName(NULL,"-column_oriented",&flg);
 31:   if (flg) {
 32:     MatSetOption(A,MAT_ROW_ORIENTED,PETSC_FALSE);
 33:   }

 35:   /* inproc assembly */
 36:   for (row=start; row<end; row++) {
 37:     for (col=start; col<end; col++,data+=1) {
 38:       MatSetValues(A,1,&row,1,&col,&data,INSERT_VALUES);
 39:     }
 40:   }
 41:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 42:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 44:   /* offproc assembly */
 45:   data = 5.0;
 46:   row  = (M+start-1)%M;
 47:   for (col=0; col<M; col++) {
 48:     MatSetValues(A,1,&row,1,&col,&data,ADD_VALUES);
 49:   }
 50:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 51:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 53:   /* Test MatSetValuesBlocked() */
 54:   PetscOptionsHasName(NULL,"-test_setvaluesblocked",&flg);
 55:   if (flg) {
 56:     PetscScalar *bval;
 57:     row /= bs;
 58:     col  = start/bs;
 59:     PetscMalloc1(bs*bs,&bval);
 60:     /* PetscPrintf(PETSC_COMM_SELF,"[%d] Set offproc blockvalues, blockrow %d, blockcol: %d\n",rank,row,col); */
 61:     k = 1;
 62:     /* row oriented - defalt */
 63:     for (i=0; i<bs; i++) {
 64:       for (j=0; j<bs; j++) {
 65:         bval[i*bs+j] = (PetscScalar)k; k++;
 66:       }
 67:     }
 68:     MatSetValuesBlocked(A,1,&row,1,&col,bval,INSERT_VALUES);
 69:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 70:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 71:     PetscFree(bval);
 72:   }

 74:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 75:   MatDestroy(&A);
 76:   PetscFinalize();
 77:   return 0;
 78: }