1 static char help[] = "Tests MATHTOOL\n\n"; 2 3 #include <petscmat.h> 4 5 static PetscScalar GenEntry(PetscInt sdim,PetscInt i,PetscInt j,void *ctx) 6 { 7 PetscInt d; 8 PetscReal diff = 0.0,*coords = (PetscReal*)(ctx); 9 10 PetscFunctionBeginUser; 11 for (d = 0; d < sdim; d++) { diff += (coords[i*sdim+d] - coords[j*sdim+d]) * (coords[i*sdim+d] - coords[j*sdim+d]); } 12 PetscFunctionReturn(1.0/(1.0e-2 + PetscSqrtReal(diff))); 13 } 14 15 static PetscScalar GenEntryRectangular(PetscInt sdim,PetscInt i,PetscInt j,void *ctx) 16 { 17 PetscInt d; 18 PetscReal diff = 0.0,**coords = (PetscReal**)(ctx); 19 20 PetscFunctionBeginUser; 21 for (d = 0; d < sdim; d++) { diff += (coords[0][i*sdim+d] - coords[1][j*sdim+d]) * (coords[0][i*sdim+d] - coords[1][j*sdim+d]); } 22 PetscFunctionReturn(1.0/(1.0e-2 + PetscSqrtReal(diff))); 23 } 24 25 int main(int argc,char **argv) 26 { 27 Mat A,AT,D,B,P,R,RT; 28 PetscInt m = 100,dim = 3,M,K = 10,begin,n = 0,N; 29 PetscMPIInt size; 30 PetscScalar *ptr; 31 PetscReal *coords,*gcoords,*scoords,*gscoords,*(ctx[2]),norm,epsilon; 32 MatHtoolKernel kernel = GenEntry; 33 PetscBool flg,sym = PETSC_FALSE; 34 PetscRandom rdm; 35 PetscErrorCode ierr; 36 37 ierr = PetscInitialize(&argc,&argv,(char*)NULL,help);if (ierr) return ierr; 38 ierr = PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL);CHKERRQ(ierr); 39 ierr = PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL);CHKERRQ(ierr); 40 ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr); 41 ierr = PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL);CHKERRQ(ierr); 42 ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL);CHKERRQ(ierr); 43 ierr = PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL);CHKERRQ(ierr); 44 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 45 M = size*m; 46 ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 47 ierr = PetscMalloc1(m*dim,&coords);CHKERRQ(ierr); 48 ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); 49 ierr = PetscRandomGetValuesReal(rdm,m*dim,coords);CHKERRQ(ierr); 50 ierr = PetscCalloc1(M*dim,&gcoords);CHKERRQ(ierr); 51 ierr = MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr); 52 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 53 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54 ierr = MatSetRandom(B,rdm);CHKERRQ(ierr); 55 ierr = MatGetOwnershipRange(B,&begin,NULL);CHKERRQ(ierr); 56 ierr = PetscArraycpy(gcoords+begin*dim,coords,m*dim);CHKERRQ(ierr); 57 ierr = MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 58 ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A);CHKERRQ(ierr); 59 ierr = MatSetOption(A,MAT_SYMMETRIC,sym);CHKERRQ(ierr); 60 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 61 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 62 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63 ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 64 if (PetscAbsReal(epsilon) >= PETSC_SMALL) { /* when there is compression, it is more difficult to check against MATDENSE, so just compare symmetric and nonsymmetric assemblies */ 65 PetscReal relative; 66 ierr = MatDestroy(&B);CHKERRQ(ierr); 67 ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B);CHKERRQ(ierr); 68 ierr = MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym);CHKERRQ(ierr); 69 ierr = MatSetFromOptions(B);CHKERRQ(ierr); 70 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 72 ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr); 73 ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P);CHKERRQ(ierr); 74 ierr = MatNorm(P,NORM_FROBENIUS,&relative);CHKERRQ(ierr); 75 ierr = MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr); 76 ierr = MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 77 ierr = MatNorm(R,NORM_INFINITY,&norm);CHKERRQ(ierr); 78 if (PetscAbsReal(norm/relative) > epsilon) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon); 79 ierr = MatDestroy(&B);CHKERRQ(ierr); 80 ierr = MatDestroy(&R);CHKERRQ(ierr); 81 ierr = MatDestroy(&P);CHKERRQ(ierr); 82 } else { 83 ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr); 84 ierr = MatViewFromOptions(D,NULL,"-D_view");CHKERRQ(ierr); 85 ierr = MatMultEqual(A,D,10,&flg);CHKERRQ(ierr); 86 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 87 ierr = MatMultAddEqual(A,D,10,&flg);CHKERRQ(ierr); 88 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx"); 89 ierr = MatGetOwnershipRange(B,&begin,NULL);CHKERRQ(ierr); 90 ierr = MatDenseGetArrayWrite(D,&ptr);CHKERRQ(ierr); 91 for (PetscInt i = 0; i < m; ++i) 92 for (PetscInt j = 0; j < M; ++j) ptr[i + j*m] = GenEntry(dim,begin+i,j,gcoords); 93 ierr = MatDenseRestoreArrayWrite(D,&ptr);CHKERRQ(ierr); 94 ierr = MatMultEqual(A,D,10,&flg);CHKERRQ(ierr); 95 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 96 ierr = MatTranspose(D,MAT_INPLACE_MATRIX,&D);CHKERRQ(ierr); 97 ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr); 98 ierr = MatMultEqual(AT,D,10,&flg);CHKERRQ(ierr); 99 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 100 ierr = MatTranspose(A,MAT_REUSE_MATRIX,&AT);CHKERRQ(ierr); 101 ierr = MatMultEqual(AT,D,10,&flg);CHKERRQ(ierr); 102 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 103 ierr = MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 104 ierr = MatNorm(D,NORM_INFINITY,&norm);CHKERRQ(ierr); 105 if (PetscAbsReal(norm) > PETSC_SMALL) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A-D|| = %g (> %g)",(double)norm,(double)PETSC_SMALL); 106 ierr = MatDestroy(&AT);CHKERRQ(ierr); 107 ierr = MatDestroy(&D);CHKERRQ(ierr); 108 ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P);CHKERRQ(ierr); 109 ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 110 ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 111 ierr = MatMatMultEqual(A,B,P,10,&flg);CHKERRQ(ierr); 112 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px"); 113 ierr = MatDestroy(&B);CHKERRQ(ierr); 114 ierr = MatDestroy(&P);CHKERRQ(ierr); 115 if (n) { 116 ierr = PetscMalloc1(n*dim,&scoords);CHKERRQ(ierr); 117 ierr = PetscRandomGetValuesReal(rdm,n*dim,scoords);CHKERRQ(ierr); 118 N = n; 119 ierr = MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 120 ierr = PetscCalloc1(N*dim,&gscoords);CHKERRQ(ierr); 121 ierr = MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 122 ierr = PetscArraycpy(gscoords+begin*dim,scoords,n*dim);CHKERRQ(ierr); 123 ierr = MPIU_Allreduce(MPI_IN_PLACE,gscoords,N*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 124 kernel = GenEntryRectangular; 125 ctx[0] = gcoords; 126 ctx[1] = gscoords; 127 ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R);CHKERRQ(ierr); 128 ierr = MatSetFromOptions(R);CHKERRQ(ierr); 129 ierr = MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 130 ierr = MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 131 ierr = MatViewFromOptions(R,NULL,"-R_view");CHKERRQ(ierr); 132 ierr = MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr); 133 ierr = MatViewFromOptions(D,NULL,"-D_view");CHKERRQ(ierr); 134 ierr = MatMultEqual(R,D,10,&flg);CHKERRQ(ierr); 135 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx"); 136 ierr = MatTranspose(D,MAT_INPLACE_MATRIX,&D);CHKERRQ(ierr); 137 ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&RT);CHKERRQ(ierr); 138 ierr = MatMultEqual(RT,D,10,&flg);CHKERRQ(ierr); 139 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 140 ierr = MatTranspose(R,MAT_REUSE_MATRIX,&RT);CHKERRQ(ierr); 141 ierr = MatMultEqual(RT,D,10,&flg);CHKERRQ(ierr); 142 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 143 ierr = MatDestroy(&RT);CHKERRQ(ierr); 144 ierr = MatDestroy(&D);CHKERRQ(ierr); 145 ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B);CHKERRQ(ierr); 146 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 148 ierr = MatSetRandom(B,rdm);CHKERRQ(ierr); 149 ierr = MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P);CHKERRQ(ierr); 150 ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151 ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 152 ierr = MatMatMultEqual(R,B,P,10,&flg);CHKERRQ(ierr); 153 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px"); 154 ierr = MatDestroy(&B);CHKERRQ(ierr); 155 ierr = MatDestroy(&P);CHKERRQ(ierr); 156 ierr = MatDestroy(&R);CHKERRQ(ierr); 157 ierr = PetscFree(gscoords);CHKERRQ(ierr); 158 ierr = PetscFree(scoords);CHKERRQ(ierr); 159 } 160 } 161 ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 162 ierr = MatDestroy(&A);CHKERRQ(ierr); 163 ierr = PetscFree(gcoords);CHKERRQ(ierr); 164 ierr = PetscFree(coords);CHKERRQ(ierr); 165 ierr = PetscFinalize(); 166 return ierr; 167 } 168 169 /*TEST 170 171 build: 172 requires: htool 173 174 test: 175 requires: htool 176 suffix: 1 177 nsize: 4 178 args: -m_local 80 -n_local 25 -mat_htool_epsilon 1.0e-11 -symmetric {{false true}shared output} 179 output_file: output/ex101.out 180 181 test: 182 requires: htool 183 suffix: 2 184 nsize: 4 185 args: -m_local 120 -mat_htool_epsilon 1.0e-2 -mat_htool_compressor {{sympartialACA fullACA SVD}shared output} 186 output_file: output/ex101.out 187 188 TEST*/ 189