static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n"; /* Example: mpiexec -n ./ex104 -mat_type elemental */ #include int main(int argc,char **argv) { Mat A,B,C,D; PetscInt i,M=10,N=5,j,nrows,ncols,am,an,rstart,rend; PetscRandom r; PetscBool equal,Aiselemental; PetscReal fill = 1.0; IS isrows,iscols; const PetscInt *rows,*cols; PetscScalar *v,rval; #if defined(PETSC_HAVE_ELEMENTAL) PetscBool Test_MatMatMult=PETSC_TRUE; #else PetscBool Test_MatMatMult=PETSC_FALSE; #endif PetscMPIInt size; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); PetscCall(MatSetType(A,MATDENSE)); PetscCall(MatSetFromOptions(A)); PetscCall(MatSetUp(A)); PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r)); PetscCall(PetscRandomSetFromOptions(r)); /* Set local matrix entries */ PetscCall(MatGetOwnershipIS(A,&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