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