153022affSStefano Zampini static char help[] = "Tests MATH2OPUS\n\n"; 2c4fb69feSStefano Zampini 3c4fb69feSStefano Zampini #include <petscmat.h> 4c4fb69feSStefano Zampini #include <petscsf.h> 5c4fb69feSStefano Zampini 6d71ae5a4SJacob Faibussowitsch static PetscScalar GenEntry_Symm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx) 7d71ae5a4SJacob Faibussowitsch { 8c4fb69feSStefano Zampini PetscInt d; 9c4fb69feSStefano Zampini PetscReal clength = sdim == 3 ? 0.2 : 0.1; 10c4fb69feSStefano Zampini PetscReal dist, diff = 0.0; 11c4fb69feSStefano Zampini 12ad540459SPierre Jolivet for (d = 0; d < sdim; d++) diff += (x[d] - y[d]) * (x[d] - y[d]); 13c4fb69feSStefano Zampini dist = PetscSqrtReal(diff); 14c4fb69feSStefano Zampini return PetscExpReal(-dist / clength); 15c4fb69feSStefano Zampini } 16c4fb69feSStefano Zampini 17d71ae5a4SJacob Faibussowitsch static PetscScalar GenEntry_Unsymm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx) 18d71ae5a4SJacob Faibussowitsch { 19c4fb69feSStefano Zampini PetscInt d; 20c4fb69feSStefano Zampini PetscReal clength = sdim == 3 ? 0.2 : 0.1; 21c4fb69feSStefano Zampini PetscReal dist, diff = 0.0, nx = 0.0, ny = 0.0; 22c4fb69feSStefano Zampini 23ad540459SPierre Jolivet for (d = 0; d < sdim; d++) nx += x[d] * x[d]; 24ad540459SPierre Jolivet for (d = 0; d < sdim; d++) ny += y[d] * y[d]; 25ad540459SPierre Jolivet for (d = 0; d < sdim; d++) diff += (x[d] - y[d]) * (x[d] - y[d]); 26c4fb69feSStefano Zampini dist = PetscSqrtReal(diff); 27c4fb69feSStefano Zampini return nx > ny ? PetscExpReal(-dist / clength) : PetscExpReal(-dist / clength) + 1.; 28c4fb69feSStefano Zampini } 29c4fb69feSStefano Zampini 30d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 31d71ae5a4SJacob Faibussowitsch { 32c4fb69feSStefano Zampini Mat A, B, C, D; 33c4fb69feSStefano Zampini Vec v, x, y, Ax, Ay, Bx, By; 34c4fb69feSStefano Zampini PetscRandom r; 35c4fb69feSStefano Zampini PetscLayout map; 36c4fb69feSStefano Zampini PetscScalar *Adata = NULL, *Cdata = NULL, scale = 1.0; 37c4fb69feSStefano Zampini PetscReal *coords, nA, nD, nB, err, nX, norms[3]; 38300d917bSStefano Zampini PetscInt N, n = 64, dim = 1, i, j, nrhs = 11, lda = 0, ldc = 0, ldu = 0, nlr = 7, nt, ntrials = 2; 39c4fb69feSStefano Zampini PetscMPIInt size, rank; 40c4fb69feSStefano Zampini PetscBool testlayout = PETSC_FALSE, flg, symm = PETSC_FALSE, Asymm = PETSC_TRUE, kernel = PETSC_TRUE; 4153022affSStefano Zampini PetscBool checkexpl = PETSC_FALSE, agpu = PETSC_FALSE, bgpu = PETSC_FALSE, cgpu = PETSC_FALSE, flgglob; 42300d917bSStefano Zampini PetscBool testtrans, testnorm, randommat = PETSC_TRUE, testorthog, testcompress, testhlru; 43c4fb69feSStefano Zampini void (*approxnormfunc)(void); 44c4fb69feSStefano Zampini void (*Anormfunc)(void); 45c4fb69feSStefano Zampini 4653022affSStefano Zampini #if defined(PETSC_HAVE_MPI_INIT_THREAD) 47c4fb69feSStefano Zampini PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_MULTIPLE; 4853022affSStefano Zampini #endif 49327415f7SBarry Smith PetscFunctionBeginUser; 509566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ng", &N, &flgglob)); 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL)); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-lda", &lda, NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ldc", &ldc, NULL)); 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nlr", &nlr, NULL)); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ldu", &ldu, NULL)); 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-matmattrials", &ntrials, NULL)); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-randommat", &randommat, NULL)); 619566063dSJacob Faibussowitsch if (!flgglob) PetscCall(PetscOptionsGetBool(NULL, NULL, "-testlayout", &testlayout, NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-Asymm", &Asymm, NULL)); 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-symm", &symm, NULL)); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-kernel", &kernel, NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-checkexpl", &checkexpl, NULL)); 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-agpu", &agpu, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-bgpu", &bgpu, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-cgpu", &cgpu, NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-scale", &scale, NULL)); 70c4fb69feSStefano Zampini if (!Asymm) symm = PETSC_FALSE; 71c4fb69feSStefano Zampini 729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 73300d917bSStefano Zampini 74300d917bSStefano Zampini /* Disable tests for unimplemented variants */ 75c4fb69feSStefano Zampini testtrans = (PetscBool)(size == 1 || symm); 76c4fb69feSStefano Zampini testnorm = (PetscBool)(size == 1 || symm); 7753022affSStefano Zampini testorthog = (PetscBool)(size == 1 || symm); 7853022affSStefano Zampini testcompress = (PetscBool)(size == 1 || symm); 79300d917bSStefano Zampini testhlru = (PetscBool)(size == 1); 80c4fb69feSStefano Zampini 819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 829566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(PETSC_COMM_WORLD, &map)); 83c4fb69feSStefano Zampini if (testlayout) { 84c4fb69feSStefano Zampini if (rank % 2) n = PetscMax(2 * n - 5 * rank, 0); 85c4fb69feSStefano Zampini else n = 2 * n + rank; 86c4fb69feSStefano Zampini } 8753022affSStefano Zampini if (!flgglob) { 889566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(map, n)); 899566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map)); 909566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 9153022affSStefano Zampini } else { 929566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(map, N)); 939566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map)); 949566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(map, &n)); 9553022affSStefano Zampini } 969566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&map)); 97c4fb69feSStefano Zampini 9848a46eb9SPierre Jolivet if (lda) PetscCall(PetscMalloc1(N * (n + lda), &Adata)); 999566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, n, n, N, N, Adata, &A)); 1009566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(A, n + lda)); 101c4fb69feSStefano Zampini 10253022affSStefano Zampini /* Create random points; these are replicated in order to populate a dense matrix and to compare sequential and dense runs 10353022affSStefano Zampini The constructor for MATH2OPUS can take as input the distributed coordinates and replicates them internally in case 10453022affSStefano Zampini a kernel construction is requested */ 1059566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r)); 1069566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 1079566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(r, 123456)); 1089566063dSJacob Faibussowitsch PetscCall(PetscRandomSeed(r)); 1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * dim, &coords)); 1109566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValuesReal(r, N * dim, coords)); 1119566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 11253022affSStefano Zampini 113c4fb69feSStefano Zampini if (kernel || !randommat) { 11453022affSStefano Zampini MatH2OpusKernel k = Asymm ? GenEntry_Symm : GenEntry_Unsymm; 115c4fb69feSStefano Zampini PetscInt ist, ien; 116c4fb69feSStefano Zampini 1179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &ist, &ien)); 118c4fb69feSStefano Zampini for (i = ist; i < ien; i++) { 11948a46eb9SPierre Jolivet for (j = 0; j < N; j++) PetscCall(MatSetValue(A, i, j, (*k)(dim, coords + i * dim, coords + j * dim, NULL), INSERT_VALUES)); 120c4fb69feSStefano Zampini } 1219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 123c4fb69feSStefano Zampini if (kernel) { 1249566063dSJacob Faibussowitsch PetscCall(MatCreateH2OpusFromKernel(PETSC_COMM_WORLD, n, n, N, N, dim, coords + ist * dim, PETSC_TRUE, k, NULL, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, &B)); 125c4fb69feSStefano Zampini } else { 1269566063dSJacob Faibussowitsch PetscCall(MatCreateH2OpusFromMat(A, dim, coords + ist * dim, PETSC_TRUE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, &B)); 127c4fb69feSStefano Zampini } 128c4fb69feSStefano Zampini } else { 12953022affSStefano Zampini PetscInt ist; 13053022affSStefano Zampini 1319566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &ist, NULL)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A, NULL)); 133c4fb69feSStefano Zampini if (Asymm) { 1349566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B)); 1359566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, 1.0, B, SAME_NONZERO_PATTERN)); 1369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1379566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 138c4fb69feSStefano Zampini } 1399566063dSJacob Faibussowitsch PetscCall(MatCreateH2OpusFromMat(A, dim, coords + ist * dim, PETSC_TRUE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, &B)); 140c4fb69feSStefano Zampini } 1419566063dSJacob Faibussowitsch PetscCall(PetscFree(coords)); 14248a46eb9SPierre Jolivet if (agpu) PetscCall(MatConvert(A, MATDENSECUDA, MAT_INPLACE_MATRIX, &A)); 1439566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-A_view")); 144c4fb69feSStefano Zampini 1459566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_SYMMETRIC, symm)); 146c4fb69feSStefano Zampini 147c4fb69feSStefano Zampini /* assemble the H-matrix */ 1489566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(B, (PetscBool)!bgpu)); 1499566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 1529566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(B, NULL, "-B_view")); 153c4fb69feSStefano Zampini 154c4fb69feSStefano Zampini /* Test MatScale */ 1559566063dSJacob Faibussowitsch PetscCall(MatScale(A, scale)); 1569566063dSJacob Faibussowitsch PetscCall(MatScale(B, scale)); 157c4fb69feSStefano Zampini 158c4fb69feSStefano Zampini /* Test MatMult */ 1599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Ax, &Ay)); 1609566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &Bx, &By)); 1619566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ax, NULL)); 1629566063dSJacob Faibussowitsch PetscCall(VecCopy(Ax, Bx)); 1639566063dSJacob Faibussowitsch PetscCall(MatMult(A, Ax, Ay)); 1649566063dSJacob Faibussowitsch PetscCall(MatMult(B, Bx, By)); 1659566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(Ay, NULL, "-mult_vec_view")); 1669566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(By, NULL, "-mult_vec_view")); 1679566063dSJacob Faibussowitsch PetscCall(VecNorm(Ay, NORM_INFINITY, &nX)); 1689566063dSJacob Faibussowitsch PetscCall(VecAXPY(Ay, -1.0, By)); 1699566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(Ay, NULL, "-mult_vec_view")); 1709566063dSJacob Faibussowitsch PetscCall(VecNorm(Ay, NORM_INFINITY, &err)); 1719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMult err %g\n", err / nX)); 1729566063dSJacob Faibussowitsch PetscCall(VecScale(By, -1.0)); 1739566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B, Bx, By, By)); 1749566063dSJacob Faibussowitsch PetscCall(VecNorm(By, NORM_INFINITY, &err)); 1759566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(By, NULL, "-mult_vec_view")); 17648a46eb9SPierre Jolivet if (err > 10. * PETSC_SMALL) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMultAdd err %g\n", err)); 177c4fb69feSStefano Zampini 178c4fb69feSStefano Zampini /* Test MatNorm */ 1799566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_INFINITY, &norms[0])); 1809566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_1, &norms[1])); 181c4fb69feSStefano Zampini norms[2] = -1.; /* NORM_2 not supported */ 1829566063dSJacob Faibussowitsch PetscCall(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])); 1839566063dSJacob Faibussowitsch PetscCall(MatGetOperation(A, MATOP_NORM, &Anormfunc)); 1849566063dSJacob Faibussowitsch PetscCall(MatGetOperation(B, MATOP_NORM, &approxnormfunc)); 1859566063dSJacob Faibussowitsch PetscCall(MatSetOperation(A, MATOP_NORM, approxnormfunc)); 1869566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_INFINITY, &norms[0])); 1879566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_1, &norms[1])); 1889566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_2, &norms[2])); 1899566063dSJacob Faibussowitsch PetscCall(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])); 190c4fb69feSStefano Zampini if (testnorm) { 1919566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_INFINITY, &norms[0])); 1929566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_1, &norms[1])); 1939566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_2, &norms[2])); 194c4fb69feSStefano Zampini } else { 195c4fb69feSStefano Zampini norms[0] = -1.; 196c4fb69feSStefano Zampini norms[1] = -1.; 197c4fb69feSStefano Zampini norms[2] = -1.; 198c4fb69feSStefano Zampini } 1999566063dSJacob Faibussowitsch PetscCall(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])); 2009566063dSJacob Faibussowitsch PetscCall(MatSetOperation(A, MATOP_NORM, Anormfunc)); 201c4fb69feSStefano Zampini 202c4fb69feSStefano Zampini /* Test MatDuplicate */ 2039566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &D)); 2049566063dSJacob Faibussowitsch PetscCall(MatSetOption(D, MAT_SYMMETRIC, symm)); 2059566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B, D, 10, &flg)); 20648a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMult error after MatDuplicate\n")); 207c4fb69feSStefano Zampini if (testtrans) { 2089566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(B, D, 10, &flg)); 20948a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMultTranspose error after MatDuplicate\n")); 210c4fb69feSStefano Zampini } 2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 212c4fb69feSStefano Zampini 213c4fb69feSStefano Zampini if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 2149566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ay, NULL)); 2159566063dSJacob Faibussowitsch PetscCall(VecCopy(Ay, By)); 2169566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Ay, Ax)); 2179566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, By, Bx)); 2189566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(Ax, NULL, "-multtrans_vec_view")); 2199566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(Bx, NULL, "-multtrans_vec_view")); 2209566063dSJacob Faibussowitsch PetscCall(VecNorm(Ax, NORM_INFINITY, &nX)); 2219566063dSJacob Faibussowitsch PetscCall(VecAXPY(Ax, -1.0, Bx)); 2229566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(Ax, NULL, "-multtrans_vec_view")); 2239566063dSJacob Faibussowitsch PetscCall(VecNorm(Ax, NORM_INFINITY, &err)); 2249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMultTranspose err %g\n", err / nX)); 2259566063dSJacob Faibussowitsch PetscCall(VecScale(Bx, -1.0)); 2269566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B, By, Bx, Bx)); 2279566063dSJacob Faibussowitsch PetscCall(VecNorm(Bx, NORM_INFINITY, &err)); 22848a46eb9SPierre Jolivet if (err > 10. * PETSC_SMALL) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMultTransposeAdd err %g\n", err)); 229c4fb69feSStefano Zampini } 2309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 2319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ay)); 2329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 2339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&By)); 234c4fb69feSStefano Zampini 235c4fb69feSStefano Zampini /* Test MatMatMult */ 23648a46eb9SPierre Jolivet if (ldc) PetscCall(PetscMalloc1(nrhs * (n + ldc), &Cdata)); 2379566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, n, PETSC_DECIDE, N, nrhs, Cdata, &C)); 2389566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(C, n + ldc)); 2399566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, NULL)); 24048a46eb9SPierre Jolivet if (cgpu) PetscCall(MatConvert(C, MATDENSECUDA, MAT_INPLACE_MATRIX, &C)); 241c4fb69feSStefano Zampini for (nt = 0; nt < ntrials; nt++) { 2429566063dSJacob Faibussowitsch PetscCall(MatMatMult(B, C, nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DEFAULT, &D)); 2439566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D, NULL, "-bc_view")); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)D, &flg, MATSEQDENSE, MATMPIDENSE, "")); 245c4fb69feSStefano Zampini if (flg) { 2469566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &x, &y)); 2479566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(D, NULL, &v)); 248c4fb69feSStefano Zampini for (i = 0; i < nrhs; i++) { 2499566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector(D, v, i)); 2509566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector(C, x, i)); 2519566063dSJacob Faibussowitsch PetscCall(MatMult(B, x, y)); 2529566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, v)); 2539566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_INFINITY, &err)); 2549566063dSJacob Faibussowitsch if (err > 10. * PETSC_SMALL) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMat err %" PetscInt_FMT " %g\n", i, err)); 255c4fb69feSStefano Zampini } 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 259c4fb69feSStefano Zampini } 260c4fb69feSStefano Zampini } 2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 262c4fb69feSStefano Zampini 263c4fb69feSStefano Zampini /* Test MatTransposeMatMult */ 264c4fb69feSStefano Zampini if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 265c4fb69feSStefano Zampini for (nt = 0; nt < ntrials; nt++) { 2669566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(B, C, nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DEFAULT, &D)); 2679566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D, NULL, "-btc_view")); 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)D, &flg, MATSEQDENSE, MATMPIDENSE, "")); 269c4fb69feSStefano Zampini if (flg) { 2709566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &y, &x)); 2719566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(D, NULL, &v)); 272c4fb69feSStefano Zampini for (i = 0; i < nrhs; i++) { 2739566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector(D, v, i)); 2749566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector(C, x, i)); 2759566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, x, y)); 2769566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, v)); 2779566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_INFINITY, &err)); 2789566063dSJacob Faibussowitsch if (err > 10. * PETSC_SMALL) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatTransMat err %" PetscInt_FMT " %g\n", i, err)); 279c4fb69feSStefano Zampini } 2809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 283c4fb69feSStefano Zampini } 284c4fb69feSStefano Zampini } 2859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 286c4fb69feSStefano Zampini } 287c4fb69feSStefano Zampini 288300d917bSStefano Zampini /* Test basis orthogonalization */ 28953022affSStefano Zampini if (testorthog) { 2909566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &D)); 2919566063dSJacob Faibussowitsch PetscCall(MatSetOption(D, MAT_SYMMETRIC, symm)); 2929566063dSJacob Faibussowitsch PetscCall(MatH2OpusOrthogonalize(D)); 2939566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B, D, 10, &flg)); 29448a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMult error after basis ortogonalization\n")); 2959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 29653022affSStefano Zampini } 29753022affSStefano Zampini 298300d917bSStefano Zampini /* Test matrix compression */ 29953022affSStefano Zampini if (testcompress) { 3009566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &D)); 3019566063dSJacob Faibussowitsch PetscCall(MatSetOption(D, MAT_SYMMETRIC, symm)); 3029566063dSJacob Faibussowitsch PetscCall(MatH2OpusCompress(D, PETSC_SMALL)); 3039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 30453022affSStefano Zampini } 30553022affSStefano Zampini 306300d917bSStefano Zampini /* Test low-rank update */ 307300d917bSStefano Zampini if (testhlru) { 308300d917bSStefano Zampini Mat U, V; 309300d917bSStefano Zampini PetscScalar *Udata = NULL, *Vdata = NULL; 310300d917bSStefano Zampini 311300d917bSStefano Zampini if (ldu) { 3129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlr * (n + ldu), &Udata)); 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlr * (n + ldu + 2), &Vdata)); 314300d917bSStefano Zampini } 3159566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &D)); 3169566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)D), n, PETSC_DECIDE, N, nlr, Udata, &U)); 3179566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(U, n + ldu)); 3189566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)D), n, PETSC_DECIDE, N, nlr, Vdata, &V)); 3199566063dSJacob Faibussowitsch if (ldu) PetscCall(MatDenseSetLDA(V, n + ldu + 2)); 3209566063dSJacob Faibussowitsch PetscCall(MatSetRandom(U, NULL)); 3219566063dSJacob Faibussowitsch PetscCall(MatSetRandom(V, NULL)); 3229566063dSJacob Faibussowitsch PetscCall(MatH2OpusLowRankUpdate(D, U, V, 0.5)); 3239566063dSJacob Faibussowitsch PetscCall(MatH2OpusLowRankUpdate(D, U, V, -0.5)); 3249566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B, D, 10, &flg)); 32548a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMult error after low-rank update\n")); 3269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 3279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&U)); 3289566063dSJacob Faibussowitsch PetscCall(PetscFree(Udata)); 3299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&V)); 3309566063dSJacob Faibussowitsch PetscCall(PetscFree(Vdata)); 331300d917bSStefano Zampini } 332300d917bSStefano Zampini 333c4fb69feSStefano Zampini /* check explicit operator */ 334c4fb69feSStefano Zampini if (checkexpl) { 335c4fb69feSStefano Zampini Mat Be, Bet; 336c4fb69feSStefano Zampini 3379566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(B, MATDENSE, &D)); 3389566063dSJacob Faibussowitsch PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &Be)); 3399566063dSJacob Faibussowitsch PetscCall(MatNorm(D, NORM_FROBENIUS, &nB)); 3409566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D, NULL, "-expl_view")); 3419566063dSJacob Faibussowitsch PetscCall(MatAXPY(D, -1.0, A, SAME_NONZERO_PATTERN)); 3429566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D, NULL, "-diff_view")); 3439566063dSJacob Faibussowitsch PetscCall(MatNorm(D, NORM_FROBENIUS, &nD)); 3449566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &nA)); 3459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Approximation error %g (%g / %g, %g)\n", nD / nA, nD, nA, nB)); 3469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 347c4fb69feSStefano Zampini 348c4fb69feSStefano Zampini if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 3499566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INPLACE_MATRIX, &A)); 3509566063dSJacob Faibussowitsch PetscCall(MatComputeOperatorTranspose(B, MATDENSE, &D)); 3519566063dSJacob Faibussowitsch PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &Bet)); 3529566063dSJacob Faibussowitsch PetscCall(MatNorm(D, NORM_FROBENIUS, &nB)); 3539566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D, NULL, "-expl_trans_view")); 3549566063dSJacob Faibussowitsch PetscCall(MatAXPY(D, -1.0, A, SAME_NONZERO_PATTERN)); 3559566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D, NULL, "-diff_trans_view")); 3569566063dSJacob Faibussowitsch PetscCall(MatNorm(D, NORM_FROBENIUS, &nD)); 3579566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &nA)); 3589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Approximation error transpose %g (%g / %g, %g)\n", nD / nA, nD, nA, nB)); 3599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 360c4fb69feSStefano Zampini 3619566063dSJacob Faibussowitsch PetscCall(MatTranspose(Bet, MAT_INPLACE_MATRIX, &Bet)); 3629566063dSJacob Faibussowitsch PetscCall(MatAXPY(Be, -1.0, Bet, SAME_NONZERO_PATTERN)); 3639566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Be, NULL, "-diff_expl_view")); 3649566063dSJacob Faibussowitsch PetscCall(MatNorm(Be, NORM_FROBENIUS, &nB)); 3659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Approximation error B - (B^T)^T %g\n", nB)); 3669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Be)); 3679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bet)); 368c4fb69feSStefano Zampini } 369c4fb69feSStefano Zampini } 3709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 3729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 3739566063dSJacob Faibussowitsch PetscCall(PetscFree(Cdata)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFree(Adata)); 3759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 376b122ec5aSJacob Faibussowitsch return 0; 377c4fb69feSStefano Zampini } 378c4fb69feSStefano Zampini 379c4fb69feSStefano Zampini /*TEST 380c4fb69feSStefano Zampini 381c4fb69feSStefano Zampini build: 38253022affSStefano Zampini requires: h2opus 383c4fb69feSStefano Zampini 384c4fb69feSStefano Zampini #tests from kernel 385c4fb69feSStefano Zampini test: 38653022affSStefano Zampini requires: h2opus 387c4fb69feSStefano Zampini nsize: 1 388c4fb69feSStefano Zampini suffix: 1 389c4fb69feSStefano Zampini args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 0 390c4fb69feSStefano Zampini 391c4fb69feSStefano Zampini test: 39253022affSStefano Zampini requires: h2opus 393c4fb69feSStefano Zampini nsize: 1 394c4fb69feSStefano Zampini suffix: 1_ld 395c4fb69feSStefano Zampini output_file: output/ex66_1.out 396300d917bSStefano Zampini args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -ldu 17 -symm 0 -checkexpl -bgpu 0 397c4fb69feSStefano Zampini 398c4fb69feSStefano Zampini test: 39953022affSStefano Zampini requires: h2opus cuda 400c4fb69feSStefano Zampini nsize: 1 401c4fb69feSStefano Zampini suffix: 1_cuda 402c4fb69feSStefano Zampini output_file: output/ex66_1.out 403c4fb69feSStefano Zampini args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 1 404c4fb69feSStefano Zampini 405c4fb69feSStefano Zampini test: 40653022affSStefano Zampini requires: h2opus cuda 407c4fb69feSStefano Zampini nsize: 1 408c4fb69feSStefano Zampini suffix: 1_cuda_ld 409c4fb69feSStefano Zampini output_file: output/ex66_1.out 410300d917bSStefano Zampini args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -ldu 17 -symm 0 -checkexpl -bgpu 1 411c4fb69feSStefano Zampini 412c4fb69feSStefano Zampini test: 41353022affSStefano Zampini requires: h2opus 41453022affSStefano Zampini nsize: {{2 3}} 415c4fb69feSStefano Zampini suffix: 1_par 41653022affSStefano Zampini args: -n 64 -symm -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu 0 -cgpu 0 417c4fb69feSStefano Zampini 418c4fb69feSStefano Zampini test: 41953022affSStefano Zampini requires: h2opus cuda 42053022affSStefano Zampini nsize: {{2 3}} 421c4fb69feSStefano Zampini suffix: 1_par_cuda 42253022affSStefano Zampini args: -n 64 -symm -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu {{0 1}} -cgpu {{0 1}} 423c4fb69feSStefano Zampini output_file: output/ex66_1_par.out 424c4fb69feSStefano Zampini 425c4fb69feSStefano Zampini #tests from matrix sampling (parallel or unsymmetric not supported) 426c4fb69feSStefano Zampini test: 42753022affSStefano Zampini requires: h2opus 428c4fb69feSStefano Zampini nsize: 1 429c4fb69feSStefano Zampini suffix: 2 430c4fb69feSStefano Zampini args: -n {{17 33}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu 0 431c4fb69feSStefano Zampini 432c4fb69feSStefano Zampini test: 43353022affSStefano Zampini requires: h2opus cuda 434c4fb69feSStefano Zampini nsize: 1 435c4fb69feSStefano Zampini suffix: 2_cuda 436c4fb69feSStefano Zampini output_file: output/ex66_2.out 437*283059ddSSatish Balay args: -n {{17 29}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu {{0 1}} -agpu {{0 1}} 438c4fb69feSStefano Zampini 43953022affSStefano Zampini #tests view operation 44053022affSStefano Zampini test: 44153022affSStefano Zampini requires: h2opus !cuda 4428cc725e6SPierre Jolivet filter: grep -v " MPI process" | grep -v "\[" | grep -v "\]" 44353022affSStefano Zampini nsize: {{1 2 3}} 44453022affSStefano Zampini suffix: view 44553022affSStefano Zampini args: -ng 64 -kernel 1 -dim 2 -symm 1 -checkexpl -B_view -mat_h2opus_leafsize 17 -mat_h2opus_normsamples 13 -mat_h2opus_indexmap_view ::ascii_matlab -mat_approximate_norm_samples 2 -mat_h2opus_normsamples 2 44653022affSStefano Zampini 44753022affSStefano Zampini test: 44853022affSStefano Zampini requires: h2opus cuda 4498cc725e6SPierre Jolivet filter: grep -v " MPI process" | grep -v "\[" | grep -v "\]" 45053022affSStefano Zampini nsize: {{1 2 3}} 45153022affSStefano Zampini suffix: view_cuda 45253022affSStefano Zampini args: -ng 64 -kernel 1 -dim 2 -symm 1 -checkexpl -bgpu -B_view -mat_h2opus_leafsize 17 -mat_h2opus_normsamples 13 -mat_h2opus_indexmap_view ::ascii_matlab -mat_approximate_norm_samples 2 -mat_h2opus_normsamples 2 45353022affSStefano Zampini 454c4fb69feSStefano Zampini TEST*/ 455