static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n"; #include int main(int argc,char **argv) { Mat A,B,C,D,AT; PetscInt i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an; PetscScalar v; PetscErrorCode ierr; PetscRandom r; PetscBool equal=PETSC_FALSE,flg; PetscReal fill = 1.0,norm; PetscMPIInt size; ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL)); CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r)); CHKERRQ(PetscRandomSetFromOptions(r)); /* Create a aij matrix A */ M = N = m*n; CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); CHKERRQ(MatSetType(A,MATAIJ)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL)); CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend)); am = Iend - Istart; for (Ii=Istart; Ii0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (i0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (j