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; PetscScalar *v; PetscMPIInt rank,color; PetscReal Cnorm; PetscBool flg,mats_view=PETSC_FALSE; MPI_Comm subcomm; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); N = M; PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-mb",&mb,NULL)); nb = mb; PetscCall(PetscOptionsGetInt(NULL,NULL,"-nb",&nb,NULL)); PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); PetscCall(MatSetType(C,MATSCALAPACK)); mloc = PETSC_DECIDE; PetscCall(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M)); nloc = PETSC_DECIDE; PetscCall(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N)); PetscCall(MatSetSizes(C,mloc,nloc,M,N)); PetscCall(MatScaLAPACKSetBlockSizes(C,mb,nb)); PetscCall(MatSetFromOptions(C)); PetscCall(MatSetUp(C)); /*PetscCall(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&C)); */ PetscCall(MatGetOwnershipIS(C,&isrows,&iscols)); PetscCall(ISGetLocalSize(isrows,&nrows)); PetscCall(ISGetIndices(isrows,&rows)); PetscCall(ISGetLocalSize(iscols,&ncols)); PetscCall(ISGetIndices(iscols,&cols)); PetscCall(PetscMalloc1(nrows*ncols,&v)); for (i=0;iN) PetscCall(MatCreateVecs(C,&d,NULL)); else PetscCall(MatCreateVecs(C,NULL,&d)); PetscCall(MatGetDiagonal(C,d)); if (mats_view) { PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); PetscCall(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); } if (M>N) { PetscCall(MatDiagonalScale(C,NULL,d)); } else { PetscCall(MatDiagonalScale(C,d,NULL)); } if (mats_view) { PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); } /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); PetscCall(MatSetType(B,MATSCALAPACK)); PetscCall(MatSetSizes(B,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE)); PetscCall(MatScaLAPACKSetBlockSizes(B,mb,nb)); PetscCall(MatSetFromOptions(B)); PetscCall(MatSetUp(B)); /* PetscCall(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&B)); */ PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); PetscCall(ISGetLocalSize(isrows,&nrows)); PetscCall(ISGetIndices(isrows,&rows)); PetscCall(ISGetLocalSize(iscols,&ncols)); PetscCall(ISGetIndices(iscols,&cols)); PetscCall(PetscMalloc1(nrows*ncols,&v)); for (i=0;i