static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n"; #include int main(int argc,char **argv) { Mat seqmat,mpimat; PetscMPIInt rank; PetscScalar value[3],*vals; PetscInt i,col[3],n=5,bs=1; PetscFunctionBeginUser; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); /* Create seqaij matrices of size (n+rank) by n */ PetscCall(MatCreate(PETSC_COMM_SELF,&seqmat)); PetscCall(MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs)); PetscCall(MatSetFromOptions(seqmat)); PetscCall(MatSeqAIJSetPreallocation(seqmat,3*bs,NULL)); PetscCall(MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL)); if (bs == 1) { value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; for (i=1; i