1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Test Hypre matrix APIs\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmathypre.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A, B, C; 9c4762a1bSJed Brown PetscReal err; 10c4762a1bSJed Brown PetscInt i,j,M = 20; 11c4762a1bSJed Brown PetscErrorCode ierr; 12c4762a1bSJed Brown PetscMPIInt NP; 13c4762a1bSJed Brown MPI_Comm comm; 14c4762a1bSJed Brown PetscInt *rows; 15c4762a1bSJed Brown 16c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 17c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 18ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&NP);CHKERRMPI(ierr); 19c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 20*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(M < 6,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Matrix has to have more than 6 columns"); 21c4762a1bSJed Brown /* Hypre matrix */ 22c4762a1bSJed Brown ierr = MatCreate(comm,&B);CHKERRQ(ierr); 23c4762a1bSJed Brown ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,M);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = MatSetType(B,MATHYPRE);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = MatHYPRESetPreallocation(B,9,NULL,9,NULL);CHKERRQ(ierr); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* PETSc AIJ matrix */ 28c4762a1bSJed Brown ierr = MatCreate(comm,&A);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,M);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,9,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,9,NULL,9,NULL);CHKERRQ(ierr); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /*Set Values */ 35c4762a1bSJed Brown for (i=0; i<M; i++) { 36c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,5}; 37c4762a1bSJed Brown PetscScalar vals[6] = {0}; 38c4762a1bSJed Brown PetscScalar value[] = {100}; 39c4762a1bSJed Brown for (j=0; j<6; j++) 40c4762a1bSJed Brown vals[j] = ((PetscReal)j)/NP; 41c4762a1bSJed Brown 42c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,6,cols,vals,ADD_VALUES);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,1,&i,value,ADD_VALUES);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,6,cols,vals,ADD_VALUES);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&i,value,ADD_VALUES);CHKERRQ(ierr); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* MAT_FLUSH_ASSEMBLY currently not supported */ 49c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Compare A and B */ 55c4762a1bSJed Brown ierr = MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 58*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues %g",err); 59c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* MatZeroRows */ 62c4762a1bSJed Brown ierr = PetscMalloc1(M, &rows);CHKERRQ(ierr); 63c4762a1bSJed Brown for (i=0; i<M; i++) rows[i] = i; 64c4762a1bSJed Brown ierr = MatZeroRows(B, M, rows, 10.0, NULL, NULL);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatZeroRows(A, M, rows, 10.0,NULL, NULL);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 70*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatZeroRows %g",err); 71c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = PetscFree(rows);CHKERRQ(ierr); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Test MatZeroEntries */ 75c4762a1bSJed Brown ierr = MatZeroEntries(B);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 78*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatZeroEntries %g",err); 79c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Insert Values */ 82c4762a1bSJed Brown for (i=0; i<M; i++) { 83c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,5}; 84c4762a1bSJed Brown PetscScalar vals[6] = {0}; 85c4762a1bSJed Brown PetscScalar value[] = {100}; 86c4762a1bSJed Brown 87c4762a1bSJed Brown for (j=0; j<6; j++) 88c4762a1bSJed Brown vals[j] = ((PetscReal)j)/NP; 89c4762a1bSJed Brown 90c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,6,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,1,&i,value,INSERT_VALUES);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,6,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&i,value,INSERT_VALUES);CHKERRQ(ierr); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* MAT_FLUSH_ASSEMBLY currently not supported */ 97c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Rows are not sorted with HYPRE so we need an intermediate sort 103c4762a1bSJed Brown They use a temporary buffer, so we can sort inplace the const memory */ 104c4762a1bSJed Brown { 105c4762a1bSJed Brown const PetscInt *idxA,*idxB; 106c4762a1bSJed Brown const PetscScalar *vA, *vB; 107c4762a1bSJed Brown PetscInt rstart, rend, nzA, nzB; 108c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,-5}; 109c4762a1bSJed Brown PetscInt *rows; 110c4762a1bSJed Brown PetscScalar *valuesA, *valuesB; 111c4762a1bSJed Brown PetscBool flg; 112c4762a1bSJed Brown 113c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 114c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 115c4762a1bSJed Brown ierr = MatGetRow(A,i,&nzA,&idxA,&vA);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = MatGetRow(B,i,&nzB,&idxB,&vB);CHKERRQ(ierr); 117*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nzA!=nzB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT, nzA-nzB); 118c4762a1bSJed Brown ierr = PetscSortIntWithScalarArray(nzB,(PetscInt*)idxB,(PetscScalar*)vB);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscArraycmp(idxA,idxB,nzA,&flg);CHKERRQ(ierr); 120*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT " (indices)",i); 121c4762a1bSJed Brown ierr = PetscArraycmp(vA,vB,nzA,&flg);CHKERRQ(ierr); 122*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT " (values)",i); 123c4762a1bSJed Brown ierr = MatRestoreRow(A,i,&nzA,&idxA,&vA);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = MatRestoreRow(B,i,&nzB,&idxB,&vB);CHKERRQ(ierr); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = PetscCalloc3((rend-rstart)*6,&valuesA,(rend-rstart)*6,&valuesB,rend-rstart,&rows);CHKERRQ(ierr); 129c4762a1bSJed Brown for (i=rstart; i<rend; i++) rows[i-rstart] =i; 130c4762a1bSJed Brown 131c4762a1bSJed Brown ierr = MatGetValues(A,rend-rstart,rows,6,cols,valuesA);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = MatGetValues(B,rend-rstart,rows,6,cols,valuesB);CHKERRQ(ierr); 133c4762a1bSJed Brown 134c4762a1bSJed Brown for (i=0; i<(rend-rstart); i++) { 135c4762a1bSJed Brown ierr = PetscArraycmp(valuesA + 6*i,valuesB + 6*i,6,&flg);CHKERRQ(ierr); 136*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetValues %" PetscInt_FMT,i + rstart); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown ierr = PetscFree3(valuesA,valuesB,rows);CHKERRQ(ierr); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Compare A and B */ 142c4762a1bSJed Brown ierr = MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 145*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues with INSERT_VALUES %g",err); 146c4762a1bSJed Brown 147c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown ierr = PetscFinalize(); 152c4762a1bSJed Brown return ierr; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /*TEST 156c4762a1bSJed Brown 157c4762a1bSJed Brown build: 158c4762a1bSJed Brown requires: hypre 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown suffix: 1 162263f2b91SStefano Zampini requires: !defined(PETSC_HAVE_HYPRE_DEVICE) 163c4762a1bSJed Brown 164c4762a1bSJed Brown test: 165263f2b91SStefano Zampini suffix: 2 166263f2b91SStefano Zampini requires: !defined(PETSC_HAVE_HYPRE_DEVICE) 167c4762a1bSJed Brown output_file: output/ex225_1.out 168c4762a1bSJed Brown nsize: 2 169c4762a1bSJed Brown 170c4762a1bSJed Brown TEST*/ 171