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*327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 17c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&NP)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 2008401ef6SPierre Jolivet PetscCheck(M >= 6,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Matrix has to have more than 6 columns"); 21c4762a1bSJed Brown /* Hypre matrix */ 229566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&B)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,M)); 249566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATHYPRE)); 259566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(B,9,NULL,9,NULL)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* PETSc AIJ matrix */ 289566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&A)); 299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,M)); 309566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJ)); 319566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,9,NULL)); 329566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,9,NULL,9,NULL)); 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 429566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,6,cols,vals,ADD_VALUES)); 439566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,1,&i,value,ADD_VALUES)); 449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,6,cols,vals,ADD_VALUES)); 459566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&i,value,ADD_VALUES)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* MAT_FLUSH_ASSEMBLY currently not supported */ 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Compare A and B */ 559566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 569566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN)); 579566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_INFINITY,&err)); 5808401ef6SPierre Jolivet PetscCheck(err <= PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues %g",err); 599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* MatZeroRows */ 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &rows)); 63c4762a1bSJed Brown for (i=0; i<M; i++) rows[i] = i; 649566063dSJacob Faibussowitsch PetscCall(MatZeroRows(B, M, rows, 10.0, NULL, NULL)); 659566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 669566063dSJacob Faibussowitsch PetscCall(MatZeroRows(A, M, rows, 10.0,NULL, NULL)); 679566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 689566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN)); 699566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_INFINITY,&err)); 7008401ef6SPierre Jolivet PetscCheck(err <= PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatZeroRows %g",err); 719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 729566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Test MatZeroEntries */ 759566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B)); 769566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 779566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_INFINITY,&err)); 7808401ef6SPierre Jolivet PetscCheck(err <= PETSC_SMALL,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatZeroEntries %g",err); 799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 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 909566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,6,cols,vals,INSERT_VALUES)); 919566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,1,&i,value,INSERT_VALUES)); 929566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,6,cols,vals,INSERT_VALUES)); 939566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&i,value,INSERT_VALUES)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* MAT_FLUSH_ASSEMBLY currently not supported */ 979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 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 1139566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 114c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 1159566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nzA,&idxA,&vA)); 1169566063dSJacob Faibussowitsch PetscCall(MatGetRow(B,i,&nzB,&idxB,&vB)); 11708401ef6SPierre Jolivet PetscCheck(nzA==nzB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT, nzA-nzB); 1189566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithScalarArray(nzB,(PetscInt*)idxB,(PetscScalar*)vB)); 1199566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(idxA,idxB,nzA,&flg)); 12028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT " (indices)",i); 1219566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(vA,vB,nzA,&flg)); 12228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetRow %" PetscInt_FMT " (values)",i); 1239566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nzA,&idxA,&vA)); 1249566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(B,i,&nzB,&idxB,&vB)); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 1279566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 1289566063dSJacob Faibussowitsch PetscCall(PetscCalloc3((rend-rstart)*6,&valuesA,(rend-rstart)*6,&valuesB,rend-rstart,&rows)); 129c4762a1bSJed Brown for (i=rstart; i<rend; i++) rows[i-rstart] =i; 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(MatGetValues(A,rend-rstart,rows,6,cols,valuesA)); 1329566063dSJacob Faibussowitsch PetscCall(MatGetValues(B,rend-rstart,rows,6,cols,valuesB)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown for (i=0; i<(rend-rstart); i++) { 1359566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(valuesA + 6*i,valuesB + 6*i,6,&flg)); 13628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error MatGetValues %" PetscInt_FMT,i + rstart); 137c4762a1bSJed Brown } 1389566063dSJacob Faibussowitsch PetscCall(PetscFree3(valuesA,valuesB,rows)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Compare A and B */ 1429566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C)); 1439566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN)); 1449566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_INFINITY,&err)); 14508401ef6SPierre Jolivet PetscCheck(err <= PETSC_SMALL,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatSetValues with INSERT_VALUES %g",err); 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 150c4762a1bSJed Brown 1519566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 152b122ec5aSJacob Faibussowitsch return 0; 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