static char help[] = "Test interface of Elemental. \n\n"; #include int main(int argc,char **args) { Mat C,Caij; PetscInt i,j,m = 5,n,nrows,ncols; const PetscInt *rows,*cols; IS isrows,iscols; PetscErrorCode ierr; PetscBool flg,Test_MatMatMult=PETSC_FALSE,mats_view=PETSC_FALSE; PetscScalar *v; PetscMPIInt rank,size; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); ierr = PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);CHKERRQ(ierr); /* Get local block or element size*/ ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); n = m; ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); ierr = MatSetType(C,MATELEMENTAL);CHKERRQ(ierr); ierr = MatSetFromOptions(C);CHKERRQ(ierr); ierr = MatSetUp(C);CHKERRQ(ierr); ierr = PetscOptionsHasName(NULL,NULL,"-row_oriented",&flg);CHKERRQ(ierr); if (flg) {ierr = MatSetOption(C,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);} ierr = MatGetOwnershipIS(C,&isrows,&iscols);CHKERRQ(ierr); ierr = PetscOptionsHasName(NULL,NULL,"-Cexp_view_ownership",&flg);CHKERRQ(ierr); if (flg) { /* View ownership of explicit C */ IS tmp; ierr = PetscPrintf(PETSC_COMM_WORLD,"Ownership of explicit C:\n");CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"Row index set:\n");CHKERRQ(ierr); ierr = ISOnComm(isrows,PETSC_COMM_WORLD,PETSC_USE_POINTER,&tmp);CHKERRQ(ierr); ierr = ISView(tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = ISDestroy(&tmp);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"Column index set:\n");CHKERRQ(ierr); ierr = ISOnComm(iscols,PETSC_COMM_WORLD,PETSC_USE_POINTER,&tmp);CHKERRQ(ierr); ierr = ISView(tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = ISDestroy(&tmp);CHKERRQ(ierr); } /* Set local matrix entries */ 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