1 static char help[] = "Tests MATHARA\n\n"; 2 3 #include <petscmat.h> 4 #include <petscsf.h> 5 6 static PetscScalar GenEntry_Symm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx) 7 { 8 PetscInt d; 9 PetscReal clength = sdim == 3 ? 0.2 : 0.1; 10 PetscReal dist, diff = 0.0; 11 12 for (d = 0; d < sdim; d++) { diff += (x[d] - y[d]) * (x[d] - y[d]); } 13 dist = PetscSqrtReal(diff); 14 return PetscExpReal(-dist / clength); 15 } 16 17 static PetscScalar GenEntry_Unsymm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx) 18 { 19 PetscInt d; 20 PetscReal clength = sdim == 3 ? 0.2 : 0.1; 21 PetscReal dist, diff = 0.0, nx = 0.0, ny = 0.0; 22 23 for (d = 0; d < sdim; d++) { nx += x[d]*x[d]; } 24 for (d = 0; d < sdim; d++) { ny += y[d]*y[d]; } 25 for (d = 0; d < sdim; d++) { diff += (x[d] - y[d]) * (x[d] - y[d]); } 26 dist = PetscSqrtReal(diff); 27 return nx > ny ? PetscExpReal(-dist / clength) : PetscExpReal(-dist / clength) + 1.; 28 } 29 30 int main(int argc,char **argv) 31 { 32 Mat A,B,C,D; 33 Vec v,x,y,Ax,Ay,Bx,By; 34 PetscRandom r; 35 PetscLayout map; 36 PetscScalar *Adata = NULL, *Cdata = NULL, scale = 1.0; 37 PetscReal *coords,nA,nD,nB,err,nX,norms[3]; 38 PetscInt N, n = 64, dim = 1, i, j, nrhs = 11, lda = 0, ldc = 0, nt, ntrials = 2; 39 PetscMPIInt size,rank; 40 PetscBool testlayout = PETSC_FALSE,flg,symm = PETSC_FALSE, Asymm = PETSC_TRUE, kernel = PETSC_TRUE; 41 PetscBool checkexpl = PETSC_FALSE, agpu = PETSC_FALSE, bgpu = PETSC_FALSE, cgpu = PETSC_FALSE; 42 PetscBool testtrans, testnorm, randommat = PETSC_TRUE; 43 void (*approxnormfunc)(void); 44 void (*Anormfunc)(void); 45 PetscErrorCode ierr; 46 47 PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_MULTIPLE; 48 ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr; 49 ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 50 ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr); 51 ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr); 52 ierr = PetscOptionsGetInt(NULL,NULL,"-lda",&lda,NULL);CHKERRQ(ierr); 53 ierr = PetscOptionsGetInt(NULL,NULL,"-ldc",&ldc,NULL);CHKERRQ(ierr); 54 ierr = PetscOptionsGetInt(NULL,NULL,"-matmattrials",&ntrials,NULL);CHKERRQ(ierr); 55 ierr = PetscOptionsGetBool(NULL,NULL,"-randommat",&randommat,NULL);CHKERRQ(ierr); 56 ierr = PetscOptionsGetBool(NULL,NULL,"-testlayout",&testlayout,NULL);CHKERRQ(ierr); 57 ierr = PetscOptionsGetBool(NULL,NULL,"-Asymm",&Asymm,NULL);CHKERRQ(ierr); 58 ierr = PetscOptionsGetBool(NULL,NULL,"-symm",&symm,NULL);CHKERRQ(ierr); 59 ierr = PetscOptionsGetBool(NULL,NULL,"-kernel",&kernel,NULL);CHKERRQ(ierr); 60 ierr = PetscOptionsGetBool(NULL,NULL,"-checkexpl",&checkexpl,NULL);CHKERRQ(ierr); 61 ierr = PetscOptionsGetBool(NULL,NULL,"-agpu",&agpu,NULL);CHKERRQ(ierr); 62 ierr = PetscOptionsGetBool(NULL,NULL,"-bgpu",&bgpu,NULL);CHKERRQ(ierr); 63 ierr = PetscOptionsGetBool(NULL,NULL,"-cgpu",&cgpu,NULL);CHKERRQ(ierr); 64 ierr = PetscOptionsGetScalar(NULL,NULL,"-scale",&scale,NULL);CHKERRQ(ierr); 65 if (!Asymm) symm = PETSC_FALSE; 66 67 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 68 /* MatMultTranspose for nonsymmetric matrices not implemented */ 69 testtrans = (PetscBool)(size == 1 || symm); 70 testnorm = (PetscBool)(size == 1 || symm); 71 72 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 73 ierr = PetscLayoutCreate(PETSC_COMM_WORLD,&map);CHKERRQ(ierr); 74 if (testlayout) { 75 if (rank%2) n = PetscMax(2*n-5*rank,0); 76 else n = 2*n+rank; 77 } 78 ierr = PetscLayoutSetLocalSize(map,n);CHKERRQ(ierr); 79 ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 80 ierr = PetscLayoutGetSize(map,&N);CHKERRQ(ierr); 81 ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 82 83 ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr); 84 ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 85 if (lda) { 86 ierr = PetscMalloc1(N*(n+lda),&Adata);CHKERRQ(ierr); 87 } 88 ierr = MatCreateDense(PETSC_COMM_WORLD,n,n,N,N,Adata,&A);CHKERRQ(ierr); 89 ierr = MatDenseSetLDA(A,n+lda);CHKERRQ(ierr); 90 ierr = PetscMalloc1(n*dim,&coords);CHKERRQ(ierr); 91 for (j = 0; j < n; j++) { 92 for (i = 0; i < dim; i++) { 93 PetscScalar a; 94 95 ierr = PetscRandomGetValue(r,&a);CHKERRQ(ierr); 96 coords[j*dim + i] = PetscRealPart(a); 97 } 98 } 99 if (kernel || !randommat) { 100 MatHaraKernel k = Asymm ? GenEntry_Symm : GenEntry_Unsymm; 101 PetscInt ist,ien; 102 PetscReal *gcoords; 103 104 if (size > 1) { /* replicated coords so that we can populate the dense matrix */ 105 PetscSF sf; 106 MPI_Datatype dtype; 107 108 ierr = MPI_Type_contiguous(dim,MPIU_REAL,&dtype);CHKERRQ(ierr); 109 ierr = MPI_Type_commit(&dtype);CHKERRQ(ierr); 110 111 ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); 112 ierr = MatGetLayouts(A,&map,NULL);CHKERRQ(ierr); 113 ierr = PetscSFSetGraphWithPattern(sf,map,PETSCSF_PATTERN_ALLGATHER);CHKERRQ(ierr); 114 ierr = PetscMalloc1(dim*N,&gcoords);CHKERRQ(ierr); 115 ierr = PetscSFBcastBegin(sf,dtype,coords,gcoords);CHKERRQ(ierr); 116 ierr = PetscSFBcastEnd(sf,dtype,coords,gcoords);CHKERRQ(ierr); 117 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 118 ierr = MPI_Type_free(&dtype);CHKERRQ(ierr); 119 } else gcoords = (PetscReal*)coords; 120 121 ierr = MatGetOwnershipRange(A,&ist,&ien);CHKERRQ(ierr); 122 for (i = ist; i < ien; i++) { 123 for (j = 0; j < N; j++) { 124 ierr = MatSetValue(A,i,j,(*k)(dim,gcoords + i*dim,gcoords + j*dim,NULL),INSERT_VALUES);CHKERRQ(ierr); 125 } 126 } 127 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 128 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129 if (kernel) { 130 ierr = MatCreateHaraFromKernel(PETSC_COMM_WORLD,n,n,N,N,dim,coords,k,NULL,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B);CHKERRQ(ierr); 131 } else { 132 ierr = MatCreateHaraFromMat(A,dim,coords,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B);CHKERRQ(ierr); 133 } 134 if (gcoords != coords) { ierr = PetscFree(gcoords);CHKERRQ(ierr); } 135 } else { 136 ierr = MatSetRandom(A,r);CHKERRQ(ierr); 137 if (Asymm) { 138 ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 139 ierr = MatAXPY(A,1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 140 ierr = MatDestroy(&B);CHKERRQ(ierr); 141 ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 142 } 143 ierr = MatCreateHaraFromMat(A,dim,coords,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B);CHKERRQ(ierr); 144 } 145 if (agpu) { 146 ierr = MatConvert(A,MATDENSECUDA,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 147 } 148 ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 149 150 ierr = MatSetOption(B,MAT_SYMMETRIC,symm);CHKERRQ(ierr); 151 ierr = PetscFree(coords);CHKERRQ(ierr); 152 153 /* assemble the H-matrix */ 154 ierr = MatBindToCPU(B,(PetscBool)!bgpu);CHKERRQ(ierr); 155 ierr = MatSetFromOptions(B);CHKERRQ(ierr); 156 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158 ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr); 159 160 /* Test MatScale */ 161 ierr = MatScale(A,scale);CHKERRQ(ierr); 162 ierr = MatScale(B,scale);CHKERRQ(ierr); 163 164 /* Test MatMult */ 165 ierr = MatCreateVecs(A,&Ax,&Ay);CHKERRQ(ierr); 166 ierr = MatCreateVecs(B,&Bx,&By);CHKERRQ(ierr); 167 ierr = VecSetRandom(Ax,r);CHKERRQ(ierr); 168 ierr = VecCopy(Ax,Bx);CHKERRQ(ierr); 169 ierr = MatMult(A,Ax,Ay);CHKERRQ(ierr); 170 ierr = MatMult(B,Bx,By);CHKERRQ(ierr); 171 ierr = VecViewFromOptions(Ay,NULL,"-mult_vec_view");CHKERRQ(ierr); 172 ierr = VecViewFromOptions(By,NULL,"-mult_vec_view");CHKERRQ(ierr); 173 ierr = VecNorm(Ay,NORM_INFINITY,&nX);CHKERRQ(ierr); 174 ierr = VecAXPY(Ay,-1.0,By);CHKERRQ(ierr); 175 ierr = VecViewFromOptions(Ay,NULL,"-mult_vec_view");CHKERRQ(ierr); 176 ierr = VecNorm(Ay,NORM_INFINITY,&err);CHKERRQ(ierr); 177 ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMult err %g\n",err/nX);CHKERRQ(ierr); 178 ierr = VecScale(By,-1.0);CHKERRQ(ierr); 179 ierr = MatMultAdd(B,Bx,By,By);CHKERRQ(ierr); 180 ierr = VecNorm(By,NORM_INFINITY,&err);CHKERRQ(ierr); 181 ierr = VecViewFromOptions(By,NULL,"-mult_vec_view");CHKERRQ(ierr); 182 if (err > PETSC_SMALL) { 183 ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultAdd err %g\n",err);CHKERRQ(ierr); 184 } 185 186 /* Test MatNorm */ 187 ierr = MatNorm(A,NORM_INFINITY,&norms[0]);CHKERRQ(ierr); 188 ierr = MatNorm(A,NORM_1,&norms[1]);CHKERRQ(ierr); 189 norms[2] = -1.; /* NORM_2 not supported */ 190 ierr = PetscPrintf(PETSC_COMM_WORLD,"A Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]);CHKERRQ(ierr); 191 ierr = MatGetOperation(A,MATOP_NORM,&Anormfunc);CHKERRQ(ierr); 192 ierr = MatGetOperation(B,MATOP_NORM,&approxnormfunc);CHKERRQ(ierr); 193 ierr = MatSetOperation(A,MATOP_NORM,approxnormfunc);CHKERRQ(ierr); 194 ierr = MatNorm(A,NORM_INFINITY,&norms[0]);CHKERRQ(ierr); 195 ierr = MatNorm(A,NORM_1,&norms[1]);CHKERRQ(ierr); 196 ierr = MatNorm(A,NORM_2,&norms[2]);CHKERRQ(ierr); 197 ierr = PetscPrintf(PETSC_COMM_WORLD,"A Approx Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]);CHKERRQ(ierr); 198 if (testnorm) { 199 ierr = MatNorm(B,NORM_INFINITY,&norms[0]);CHKERRQ(ierr); 200 ierr = MatNorm(B,NORM_1,&norms[1]);CHKERRQ(ierr); 201 ierr = MatNorm(B,NORM_2,&norms[2]);CHKERRQ(ierr); 202 } else { 203 norms[0] = -1.; 204 norms[1] = -1.; 205 norms[2] = -1.; 206 } 207 ierr = PetscPrintf(PETSC_COMM_WORLD,"B Approx Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]);CHKERRQ(ierr); 208 ierr = MatSetOperation(A,MATOP_NORM,Anormfunc);CHKERRQ(ierr); 209 210 /* Test MatDuplicate */ 211 ierr = MatDuplicate(B,MAT_COPY_VALUES,&D);CHKERRQ(ierr); 212 ierr = MatMultEqual(B,D,10,&flg);CHKERRQ(ierr); 213 if (!flg) { 214 ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMult error after MatDuplicate\n");CHKERRQ(ierr); 215 } 216 if (testtrans) { 217 ierr = MatMultTransposeEqual(B,D,10,&flg);CHKERRQ(ierr); 218 if (!flg) { 219 ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultTranpose error after MatDuplicate\n");CHKERRQ(ierr); 220 } 221 } 222 ierr = MatDestroy(&D);CHKERRQ(ierr); 223 224 if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 225 ierr = VecSetRandom(Ay,r);CHKERRQ(ierr); 226 ierr = VecCopy(Ay,By);CHKERRQ(ierr); 227 ierr = MatMultTranspose(A,Ay,Ax);CHKERRQ(ierr); 228 ierr = MatMultTranspose(B,By,Bx);CHKERRQ(ierr); 229 ierr = VecViewFromOptions(Ax,NULL,"-multtrans_vec_view");CHKERRQ(ierr); 230 ierr = VecViewFromOptions(Bx,NULL,"-multtrans_vec_view");CHKERRQ(ierr); 231 ierr = VecNorm(Ax,NORM_INFINITY,&nX);CHKERRQ(ierr); 232 ierr = VecAXPY(Ax,-1.0,Bx);CHKERRQ(ierr); 233 ierr = VecViewFromOptions(Ax,NULL,"-multtrans_vec_view");CHKERRQ(ierr); 234 ierr = VecNorm(Ax,NORM_INFINITY,&err);CHKERRQ(ierr); 235 ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultTranspose err %g\n",err/nX);CHKERRQ(ierr); 236 ierr = VecScale(Bx,-1.0);CHKERRQ(ierr); 237 ierr = MatMultTransposeAdd(B,By,Bx,Bx);CHKERRQ(ierr); 238 ierr = VecNorm(Bx,NORM_INFINITY,&err);CHKERRQ(ierr); 239 if (err > PETSC_SMALL) { 240 ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultTransposeAdd err %g\n",err);CHKERRQ(ierr); 241 } 242 } 243 ierr = VecDestroy(&Ax);CHKERRQ(ierr); 244 ierr = VecDestroy(&Ay);CHKERRQ(ierr); 245 ierr = VecDestroy(&Bx);CHKERRQ(ierr); 246 ierr = VecDestroy(&By);CHKERRQ(ierr); 247 248 /* Test MatMatMult */ 249 if (ldc) { 250 ierr = PetscMalloc1(nrhs*(n+ldc),&Cdata);CHKERRQ(ierr); 251 } 252 ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,N,nrhs,Cdata,&C);CHKERRQ(ierr); 253 ierr = MatDenseSetLDA(C,n+ldc);CHKERRQ(ierr); 254 ierr = MatSetRandom(C,r);CHKERRQ(ierr); 255 if (cgpu) { 256 ierr = MatConvert(C,MATDENSECUDA,MAT_INPLACE_MATRIX,&C);CHKERRQ(ierr); 257 } 258 for (nt = 0; nt < ntrials; nt++) { 259 ierr = MatMatMult(B,C,nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D);CHKERRQ(ierr); 260 ierr = MatViewFromOptions(D,NULL,"-bc_view");CHKERRQ(ierr); 261 ierr = PetscObjectBaseTypeCompareAny((PetscObject)D,&flg,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 262 if (flg) { 263 ierr = MatCreateVecs(B,&x,&y);CHKERRQ(ierr); 264 ierr = MatCreateVecs(D,NULL,&v);CHKERRQ(ierr); 265 for (i = 0; i < nrhs; i++) { 266 ierr = MatGetColumnVector(D,v,i);CHKERRQ(ierr); 267 ierr = MatGetColumnVector(C,x,i);CHKERRQ(ierr); 268 ierr = MatMult(B,x,y);CHKERRQ(ierr); 269 ierr = VecAXPY(y,-1.0,v);CHKERRQ(ierr); 270 ierr = VecNorm(y,NORM_INFINITY,&err);CHKERRQ(ierr); 271 if (err > PETSC_SMALL) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMat err %D %g\n",i,err);CHKERRQ(ierr); } 272 } 273 ierr = VecDestroy(&y);CHKERRQ(ierr); 274 ierr = VecDestroy(&x);CHKERRQ(ierr); 275 ierr = VecDestroy(&v);CHKERRQ(ierr); 276 } 277 } 278 ierr = MatDestroy(&D);CHKERRQ(ierr); 279 280 /* Test MatTransposeMatMult */ 281 if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 282 for (nt = 0; nt < ntrials; nt++) { 283 ierr = MatTransposeMatMult(B,C,nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D);CHKERRQ(ierr); 284 ierr = MatViewFromOptions(D,NULL,"-btc_view");CHKERRQ(ierr); 285 ierr = PetscObjectBaseTypeCompareAny((PetscObject)D,&flg,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 286 if (flg) { 287 ierr = MatCreateVecs(B,&y,&x);CHKERRQ(ierr); 288 ierr = MatCreateVecs(D,NULL,&v);CHKERRQ(ierr); 289 for (i = 0; i < nrhs; i++) { 290 ierr = MatGetColumnVector(D,v,i);CHKERRQ(ierr); 291 ierr = MatGetColumnVector(C,x,i);CHKERRQ(ierr); 292 ierr = MatMultTranspose(B,x,y);CHKERRQ(ierr); 293 ierr = VecAXPY(y,-1.0,v);CHKERRQ(ierr); 294 ierr = VecNorm(y,NORM_INFINITY,&err);CHKERRQ(ierr); 295 if (err > PETSC_SMALL) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatTransMat err %D %g\n",i,err);CHKERRQ(ierr); } 296 } 297 ierr = VecDestroy(&y);CHKERRQ(ierr); 298 ierr = VecDestroy(&x);CHKERRQ(ierr); 299 ierr = VecDestroy(&v);CHKERRQ(ierr); 300 } 301 } 302 ierr = MatDestroy(&D);CHKERRQ(ierr); 303 } 304 305 /* check explicit operator */ 306 if (checkexpl) { 307 Mat Be, Bet; 308 309 ierr = MatComputeOperator(B,MATDENSE,&D);CHKERRQ(ierr); 310 ierr = MatDuplicate(D,MAT_COPY_VALUES,&Be);CHKERRQ(ierr); 311 ierr = MatNorm(D,NORM_FROBENIUS,&nB);CHKERRQ(ierr); 312 ierr = MatViewFromOptions(D,NULL,"-expl_view");CHKERRQ(ierr); 313 ierr = MatAXPY(D,-1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 314 ierr = MatViewFromOptions(D,NULL,"-diff_view");CHKERRQ(ierr); 315 ierr = MatNorm(D,NORM_FROBENIUS,&nD);CHKERRQ(ierr); 316 ierr = MatNorm(A,NORM_FROBENIUS,&nA);CHKERRQ(ierr); 317 ierr = PetscPrintf(PETSC_COMM_WORLD,"Approximation error %g (%g / %g, %g)\n",nD/nA,nD,nA,nB);CHKERRQ(ierr); 318 ierr = MatDestroy(&D);CHKERRQ(ierr); 319 320 if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 321 ierr = MatTranspose(A,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 322 ierr = MatComputeOperatorTranspose(B,MATDENSE,&D);CHKERRQ(ierr); 323 ierr = MatDuplicate(D,MAT_COPY_VALUES,&Bet);CHKERRQ(ierr); 324 ierr = MatNorm(D,NORM_FROBENIUS,&nB);CHKERRQ(ierr); 325 ierr = MatViewFromOptions(D,NULL,"-expl_trans_view");CHKERRQ(ierr); 326 ierr = MatAXPY(D,-1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 327 ierr = MatViewFromOptions(D,NULL,"-diff_trans_view");CHKERRQ(ierr); 328 ierr = MatNorm(D,NORM_FROBENIUS,&nD);CHKERRQ(ierr); 329 ierr = MatNorm(A,NORM_FROBENIUS,&nA);CHKERRQ(ierr); 330 ierr = PetscPrintf(PETSC_COMM_WORLD,"Approximation error transpose %g (%g / %g, %g)\n",nD/nA,nD,nA,nB);CHKERRQ(ierr); 331 ierr = MatDestroy(&D);CHKERRQ(ierr); 332 333 ierr = MatTranspose(Bet,MAT_INPLACE_MATRIX,&Bet);CHKERRQ(ierr); 334 ierr = MatAXPY(Be,-1.0,Bet,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 335 ierr = MatViewFromOptions(Be,NULL,"-diff_expl_view");CHKERRQ(ierr); 336 ierr = MatNorm(Be,NORM_FROBENIUS,&nB);CHKERRQ(ierr); 337 ierr = PetscPrintf(PETSC_COMM_WORLD,"Approximation error B - (B^T)^T %g\n",nB);CHKERRQ(ierr); 338 ierr = MatDestroy(&Be);CHKERRQ(ierr); 339 ierr = MatDestroy(&Bet);CHKERRQ(ierr); 340 } 341 } 342 ierr = MatDestroy(&A);CHKERRQ(ierr); 343 ierr = MatDestroy(&B);CHKERRQ(ierr); 344 ierr = MatDestroy(&C);CHKERRQ(ierr); 345 ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 346 ierr = PetscFree(Cdata);CHKERRQ(ierr); 347 ierr = PetscFree(Adata);CHKERRQ(ierr); 348 ierr = PetscFinalize(); 349 return ierr; 350 } 351 352 /*TEST 353 354 build: 355 requires: hara 356 357 #tests from kernel 358 test: 359 requires: hara 360 nsize: 1 361 suffix: 1 362 args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 0 363 364 test: 365 requires: hara 366 nsize: 1 367 suffix: 1_ld 368 output_file: output/ex66_1.out 369 args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -symm 0 -checkexpl -bgpu 0 370 371 test: 372 requires: hara cuda 373 nsize: 1 374 suffix: 1_cuda 375 output_file: output/ex66_1.out 376 args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 1 377 378 test: 379 requires: hara cuda 380 nsize: 1 381 suffix: 1_cuda_ld 382 output_file: output/ex66_1.out 383 args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -symm 0 -checkexpl -bgpu 1 384 385 test: 386 requires: hara define(PETSC_HAVE_MPI_INIT_THREAD) 387 nsize: 2 388 suffix: 1_par 389 args: -n 32 -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu 0 -cgpu 0 390 391 test: 392 requires: hara cuda define(PETSC_HAVE_MPI_INIT_THREAD) 393 nsize: 2 394 suffix: 1_par_cuda 395 args: -n 32 -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu {{0 1}} -cgpu {{0 1}} 396 output_file: output/ex66_1_par.out 397 398 #tests from matrix sampling (parallel or unsymmetric not supported) 399 test: 400 requires: hara 401 nsize: 1 402 suffix: 2 403 args: -n {{17 33}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu 0 404 405 test: 406 requires: hara cuda 407 nsize: 1 408 suffix: 2_cuda 409 output_file: output/ex66_2.out 410 args: -n {{17 33}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu {{0 1}} -agpu {{0 1}} 411 412 TEST*/ 413