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