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 PetscMPIInt NP; 12c4762a1bSJed Brown MPI_Comm comm; 13c4762a1bSJed Brown PetscInt *rows; 14c4762a1bSJed Brown 15*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 16c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 175f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&NP)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 192c71b3e2SJacob Faibussowitsch PetscCheckFalse(M < 6,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Matrix has to have more than 6 columns"); 20c4762a1bSJed Brown /* Hypre matrix */ 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&B)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,M)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATHYPRE)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatHYPRESetPreallocation(B,9,NULL,9,NULL)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* PETSc AIJ matrix */ 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&A)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,M)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,9,NULL)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,9,NULL,9,NULL)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /*Set Values */ 34c4762a1bSJed Brown for (i=0; i<M; i++) { 35c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,5}; 36c4762a1bSJed Brown PetscScalar vals[6] = {0}; 37c4762a1bSJed Brown PetscScalar value[] = {100}; 38c4762a1bSJed Brown for (j=0; j<6; j++) 39c4762a1bSJed Brown vals[j] = ((PetscReal)j)/NP; 40c4762a1bSJed Brown 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,6,cols,vals,ADD_VALUES)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,1,&i,value,ADD_VALUES)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,6,cols,vals,ADD_VALUES)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&i,value,ADD_VALUES)); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* MAT_FLUSH_ASSEMBLY currently not supported */ 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Compare A and B */ 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_INFINITY,&err)); 572c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues %g",err); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* MatZeroRows */ 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(M, &rows)); 62c4762a1bSJed Brown for (i=0; i<M; i++) rows[i] = i; 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRows(B, M, rows, 10.0, NULL, NULL)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRows(A, M, rows, 10.0,NULL, NULL)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_INFINITY,&err)); 692c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatZeroRows %g",err); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rows)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Test MatZeroEntries */ 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(B)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_INFINITY,&err)); 772c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatZeroEntries %g",err); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Insert Values */ 81c4762a1bSJed Brown for (i=0; i<M; i++) { 82c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,5}; 83c4762a1bSJed Brown PetscScalar vals[6] = {0}; 84c4762a1bSJed Brown PetscScalar value[] = {100}; 85c4762a1bSJed Brown 86c4762a1bSJed Brown for (j=0; j<6; j++) 87c4762a1bSJed Brown vals[j] = ((PetscReal)j)/NP; 88c4762a1bSJed Brown 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,6,cols,vals,INSERT_VALUES)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,1,&i,value,INSERT_VALUES)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,6,cols,vals,INSERT_VALUES)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&i,value,INSERT_VALUES)); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* MAT_FLUSH_ASSEMBLY currently not supported */ 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* Rows are not sorted with HYPRE so we need an intermediate sort 102c4762a1bSJed Brown They use a temporary buffer, so we can sort inplace the const memory */ 103c4762a1bSJed Brown { 104c4762a1bSJed Brown const PetscInt *idxA,*idxB; 105c4762a1bSJed Brown const PetscScalar *vA, *vB; 106c4762a1bSJed Brown PetscInt rstart, rend, nzA, nzB; 107c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,-5}; 108c4762a1bSJed Brown PetscInt *rows; 109c4762a1bSJed Brown PetscScalar *valuesA, *valuesB; 110c4762a1bSJed Brown PetscBool flg; 111c4762a1bSJed Brown 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 113c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(A,i,&nzA,&idxA,&vA)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(B,i,&nzB,&idxB,&vB)); 1162c71b3e2SJacob Faibussowitsch PetscCheckFalse(nzA!=nzB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT, nzA-nzB); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSortIntWithScalarArray(nzB,(PetscInt*)idxB,(PetscScalar*)vB)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycmp(idxA,idxB,nzA,&flg)); 11928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT " (indices)",i); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycmp(vA,vB,nzA,&flg)); 12128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT " (values)",i); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(A,i,&nzA,&idxA,&vA)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(B,i,&nzB,&idxB,&vB)); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc3((rend-rstart)*6,&valuesA,(rend-rstart)*6,&valuesB,rend-rstart,&rows)); 128c4762a1bSJed Brown for (i=rstart; i<rend; i++) rows[i-rstart] =i; 129c4762a1bSJed Brown 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetValues(A,rend-rstart,rows,6,cols,valuesA)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetValues(B,rend-rstart,rows,6,cols,valuesB)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown for (i=0; i<(rend-rstart); i++) { 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycmp(valuesA + 6*i,valuesB + 6*i,6,&flg)); 13528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetValues %" PetscInt_FMT,i + rstart); 136c4762a1bSJed Brown } 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(valuesA,valuesB,rows)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Compare A and B */ 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_INFINITY,&err)); 1442c71b3e2SJacob Faibussowitsch PetscCheckFalse(err > PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues with INSERT_VALUES %g",err); 145c4762a1bSJed Brown 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 149c4762a1bSJed Brown 150*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 151*b122ec5aSJacob Faibussowitsch return 0; 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154c4762a1bSJed Brown /*TEST 155c4762a1bSJed Brown 156c4762a1bSJed Brown build: 157c4762a1bSJed Brown requires: hypre 158c4762a1bSJed Brown 159c4762a1bSJed Brown test: 160c4762a1bSJed Brown suffix: 1 161263f2b91SStefano Zampini requires: !defined(PETSC_HAVE_HYPRE_DEVICE) 162c4762a1bSJed Brown 163c4762a1bSJed Brown test: 164263f2b91SStefano Zampini suffix: 2 165263f2b91SStefano Zampini requires: !defined(PETSC_HAVE_HYPRE_DEVICE) 166c4762a1bSJed Brown output_file: output/ex225_1.out 167c4762a1bSJed Brown nsize: 2 168c4762a1bSJed Brown 169c4762a1bSJed Brown TEST*/ 170