static char help[] = "Test Hypre matrix APIs\n"; #include int main(int argc,char **args) { Mat A, B, C; PetscReal err; PetscInt i,j,M = 20; PetscMPIInt NP; MPI_Comm comm; PetscInt *rows; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); comm = PETSC_COMM_WORLD; PetscCallMPI(MPI_Comm_size(comm,&NP)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); PetscCheckFalse(M < 6,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Matrix has to have more than 6 columns"); /* Hypre matrix */ PetscCall(MatCreate(comm,&B)); PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,M)); PetscCall(MatSetType(B,MATHYPRE)); PetscCall(MatHYPRESetPreallocation(B,9,NULL,9,NULL)); /* PETSc AIJ matrix */ PetscCall(MatCreate(comm,&A)); PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,M)); PetscCall(MatSetType(A,MATAIJ)); PetscCall(MatSeqAIJSetPreallocation(A,9,NULL)); PetscCall(MatMPIAIJSetPreallocation(A,9,NULL,9,NULL)); /*Set Values */ for (i=0; i PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues %g",err); PetscCall(MatDestroy(&C)); /* MatZeroRows */ PetscCall(PetscMalloc1(M, &rows)); for (i=0; i PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatZeroRows %g",err); PetscCall(MatDestroy(&C)); PetscCall(PetscFree(rows)); /* Test MatZeroEntries */ PetscCall(MatZeroEntries(B)); PetscCall(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); PetscCall(MatNorm(C,NORM_INFINITY,&err)); PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatZeroEntries %g",err); PetscCall(MatDestroy(&C)); /* Insert Values */ for (i=0; i PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues with INSERT_VALUES %g",err); PetscCall(MatDestroy(&A)); PetscCall(MatDestroy(&B)); PetscCall(MatDestroy(&C)); PetscCall(PetscFinalize()); return 0; } /*TEST build: requires: hypre test: suffix: 1 requires: !defined(PETSC_HAVE_HYPRE_DEVICE) test: suffix: 2 requires: !defined(PETSC_HAVE_HYPRE_DEVICE) output_file: output/ex225_1.out nsize: 2 TEST*/