static char help[] ="Tests MatPtAP() for MPIMAIJ and MPIAIJ \n "; #include int main(int argc,char **argv) { PetscErrorCode ierr; 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; ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; comm = PETSC_COMM_WORLD; ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); M = 10; N = 10; Z = 10; dof = 10; ierr = PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL);CHKERRQ(ierr); /* Set up distributed array for fine grid */ if (!Test_3D) { ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm);CHKERRQ(ierr); } else { ierr = 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);CHKERRQ(ierr); } ierr = DMSetFromOptions(coarsedm);CHKERRQ(ierr); ierr = DMSetUp(coarsedm);CHKERRQ(ierr); /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ ierr = DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm);CHKERRQ(ierr); /*------------------------------------------------------------*/ ierr = DMSetMatType(finedm,MATAIJ);CHKERRQ(ierr); ierr = DMCreateMatrix(finedm,&A);CHKERRQ(ierr); /* set val=one to A */ if (size == 1) { ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); if (flg) { ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr); for (i=0; i