static char help[] = "Tests ScaLAPACK interface.\n\n"; #include int main(int argc,char **args) { Mat Cdense,Caij,B,C,Ct,Asub; Vec d; PetscInt i,j,M = 5,N,mb = 2,nb,nrows,ncols,mloc,nloc; const PetscInt *rows,*cols; IS isrows,iscols; PetscErrorCode ierr; PetscScalar *v; PetscMPIInt rank,color; PetscReal Cnorm; PetscBool flg,mats_view=PETSC_FALSE; MPI_Comm subcomm; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); N = M; CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mb",&mb,NULL)); nb = mb; CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nb",&nb,NULL)); CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); CHKERRQ(MatSetType(C,MATSCALAPACK)); mloc = PETSC_DECIDE; CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M)); nloc = PETSC_DECIDE; CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N)); CHKERRQ(MatSetSizes(C,mloc,nloc,M,N)); CHKERRQ(MatScaLAPACKSetBlockSizes(C,mb,nb)); CHKERRQ(MatSetFromOptions(C)); CHKERRQ(MatSetUp(C)); /*CHKERRQ(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&C)); */ CHKERRQ(MatGetOwnershipIS(C,&isrows,&iscols)); CHKERRQ(ISGetLocalSize(isrows,&nrows)); CHKERRQ(ISGetIndices(isrows,&rows)); CHKERRQ(ISGetLocalSize(iscols,&ncols)); CHKERRQ(ISGetIndices(iscols,&cols)); CHKERRQ(PetscMalloc1(nrows*ncols,&v)); for (i=0;iN) CHKERRQ(MatCreateVecs(C,&d,NULL)); else CHKERRQ(MatCreateVecs(C,NULL,&d)); CHKERRQ(MatGetDiagonal(C,d)); if (mats_view) { CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); CHKERRQ(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); } if (M>N) { CHKERRQ(MatDiagonalScale(C,NULL,d)); } else { CHKERRQ(MatDiagonalScale(C,d,NULL)); } if (mats_view) { CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); } /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); CHKERRQ(MatSetType(B,MATSCALAPACK)); CHKERRQ(MatSetSizes(B,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE)); CHKERRQ(MatScaLAPACKSetBlockSizes(B,mb,nb)); CHKERRQ(MatSetFromOptions(B)); CHKERRQ(MatSetUp(B)); /* CHKERRQ(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&B)); */ CHKERRQ(MatGetOwnershipIS(B,&isrows,&iscols)); CHKERRQ(ISGetLocalSize(isrows,&nrows)); CHKERRQ(ISGetIndices(isrows,&rows)); CHKERRQ(ISGetLocalSize(iscols,&ncols)); CHKERRQ(ISGetIndices(iscols,&cols)); CHKERRQ(PetscMalloc1(nrows*ncols,&v)); for (i=0;i