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; PetscRandom r; PetscBool equal=PETSC_FALSE,flg; PetscReal fill = 1.0,norm; PetscMPIInt size; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL)); PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r)); PetscCall(PetscRandomSetFromOptions(r)); /* Create a aij matrix A */ M = N = m*n; PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); PetscCall(MatSetType(A,MATAIJ)); PetscCall(MatSetFromOptions(A)); PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); PetscCall(MatSeqAIJSetPreallocation(A,5,NULL)); PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); am = Iend - Istart; for (Ii=Istart; Ii0) {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (i0) {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (j