static char help[] ="Tests MatPtAP() for MPIMAIJ and MPIAIJ \n "; #include int main(int argc,char **argv) { DM coarsedm,finedm; PetscMPIInt size,rank; PetscInt M,N,Z,i,nrows; PetscScalar one = 1.0; PetscReal fill=2.0; Mat A,P,C; PetscScalar *array,alpha; PetscBool Test_3D=PETSC_FALSE,flg; const PetscInt *ia,*ja; PetscInt dof; MPI_Comm comm; PetscFunctionBeginUser; PetscCall(PetscInitialize(&argc,&argv,NULL,help)); comm = PETSC_COMM_WORLD; PetscCallMPI(MPI_Comm_rank(comm,&rank)); PetscCallMPI(MPI_Comm_size(comm,&size)); M = 10; N = 10; Z = 10; dof = 10; PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL)); /* Set up distributed array for fine grid */ if (!Test_3D) { PetscCall(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm)); } else { PetscCall(DMDACreate3d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,Z,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,NULL,&coarsedm)); } PetscCall(DMSetFromOptions(coarsedm)); PetscCall(DMSetUp(coarsedm)); /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ PetscCall(DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm)); /*------------------------------------------------------------*/ PetscCall(DMSetMatType(finedm,MATAIJ)); PetscCall(DMCreateMatrix(finedm,&A)); /* set val=one to A */ if (size == 1) { PetscCall(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); if (flg) { PetscCall(MatSeqAIJGetArray(A,&array)); for (i=0; i