static char help[] = "Tests Elemental interface.\n\n"; #include int main(int argc,char **args) { Mat Cdense,B,C,Ct; Vec d; PetscInt i,j,m = 5,n,nrows,ncols; const PetscInt *rows,*cols; IS isrows,iscols; PetscScalar *v; PetscMPIInt rank,size; PetscReal Cnorm; PetscBool flg,mats_view=PETSC_FALSE; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); n = m; PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); PetscCall(MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE)); PetscCall(MatSetType(C,MATELEMENTAL)); PetscCall(MatSetFromOptions(C)); PetscCall(MatSetUp(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)); #if defined(PETSC_USE_COMPLEX) PetscRandom rand; PetscScalar rval; PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); PetscCall(PetscRandomSetFromOptions(rand)); for (i=0; in ? n : m,PETSC_DECIDE)); PetscCall(VecSetFromOptions(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(MatSetSizes(B,m,n,PETSC_DECIDE,PETSC_DECIDE)); PetscCall(MatSetType(B,MATELEMENTAL)); PetscCall(MatSetFromOptions(B)); PetscCall(MatSetUp(B)); PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); PetscCall(ISGetLocalSize(isrows,&nrows)); PetscCall(ISGetIndices(isrows,&rows)); PetscCall(ISGetLocalSize(iscols,&ncols)); PetscCall(ISGetIndices(iscols,&cols)); for (i=0; i