static void do_test(MAP *rmap, MAP *cmap) { MAT *A; VEC *x, *y; INT i, j, m0, n0; INT *cols = phgAlloc(cmap->nglobal * sizeof(*cols)); FLOAT *data = phgAlloc(cmap->nglobal * sizeof(*data)); A = phgMapCreateMat(rmap, cmap); m0 = A->rmap->partition[A->rmap->rank]; n0 = A->cmap->partition[A->cmap->rank]; /* Matrix entries: A(I,J) = 1 + (I-1) + (J-1) */ for (i = m0; i < m0 + A->rmap->nlocal; i++) { #if 0 /* Test MatAddEntry */ for (j = 0; j < A->cmap->nglobal; j++) phgMatAddGlobalEntry(A, i, j, 1.0 + i + j); #else /* Test MatAddEntries */ for (j = 0; j < A->cmap->nglobal; j++) { cols[j] = j; data[j] = 1.0 + i + j; } phgMatAddGlobalEntries(A, 1, &i, A->cmap->nglobal, cols, data); #endif } phgFree(cols); phgFree(data); phgMatAssemble(A); phgInfo(-1, "y = A * x\n"); x = phgMapCreateVec(A->cmap, 1); for (i = 0; i < x->map->nlocal; i++) x->data[i] = 1.0 + i + n0; phgVecAssemble(x); y = phgMatVec(MAT_OP_N, 1.0, A, x, 0.0, NULL); for (i = 0; i < y->map->nlocal; i++) phgInfo(-1, " y->data[%d] = %lg\n", i + m0, (double)y->data[i]); phgVecDestroy(&x); phgVecDestroy(&y); phgInfo(-1, "y = A' * x\n"); x = phgMapCreateVec(A->rmap, 1); for (i = 0; i < x->map->nlocal; i++) x->data[i] = 1.0 + i + m0; phgVecAssemble(x); y = phgMatVec(MAT_OP_T, 1.0, A, x, 0.0, NULL); for (i = 0; i < y->map->nlocal; i++) phgInfo(-1, " y->data[%d] = %lg\n", i + n0, (double)y->data[i]); phgVecDestroy(&x); phgVecDestroy(&y); phgMatDestroy(&A); }
/* callback function which computes A+UU^t */ static int funcB(MAT_OP op, MAT *B, VEC *x, VEC *y) /* computes y = op(A) * x */ { MAT *A = ((void **)B->mv_data)[0]; VEC *U = ((void **)B->mv_data)[1]; INT i, j, k, n = U->map->nlocal, K = U->nvec; /* compute y = op(A)*x */ phgMatVec(op, 1.0, A, x, 0.0, &y); if (K == 0) return 0; if (op != MAT_OP_D) { /* compute y += U*U^t*x */ FLOAT *z = phgCalloc(K, sizeof(*z)); /* z = U^t*x */ for (j = 0; j < K; j++) for (i = 0; i < n; i++) z[j] += U->data[j * n + i] * x->data[i]; #if USE_MPI if (U->map->nprocs > 1) { FLOAT *tmp = phgCalloc(K, sizeof(*tmp)); MPI_Allreduce(z, tmp, K, PHG_MPI_FLOAT, MPI_SUM, U->map->comm); phgFree(z); z = tmp; } #endif /* USE_MPI */ /* y += U*z */ for (i = 0; i < n; i++) for (j = 0; j < K; j++) y->data[i] += U->data[j * n + i] * z[j]; phgFree(z); return 0; } for (i = 0; i < n; i++) { FLOAT d = 0.0; for (j = 0; j < K; j++) for (k = 0; k < K; k++) d += U->data[j * n + i] * U->data[k * n + i]; y->data[i] += d * x->data[i]; } return 0; }
/* * Jacobi smoother. * * Implement using phgMatvec * ref: phgSolverJacobi, solver-phg.c * */ void mg_Jacobi(MAT *A, VEC *x, VEC *b, int nsmooth, void *ctx) { VEC *r; INT i, k; FLOAT dx, omega = _p->smooth_damp; if (A->diag == NULL) phgMatSetupDiagonal(A); r = phgMapCreateVec(b->map, 1); for (k = 0; k < nsmooth; k++) { phgVecCopy(b, &r); phgMatVec(MAT_OP_N, -1.0, A, x, 1.0, &r); for (i = 0; i < A->cmap->nlocal; i++) { assert(fabs(A->diag[i]) > 1e-12); dx = r->data[i] / A->diag[i]; x->data[i] += omega * dx; } } phgVecDestroy(&r); return; }
int main(int argc, char *argv[]) { MAT *A, *B; VEC *U = NULL, *x; FLOAT *C; SOLVER *solver, *solver1 = NULL; INT i, j, k, n, *pvt, N = 1000, K = 2; char *main_opts = NULL, *sub_opts = NULL; phgOptionsRegisterInt("-n", "N value", &N); phgOptionsRegisterInt("-k", "K value", &K); phgOptionsRegisterString("-main_solver_opts", "Options for the main solver", &main_opts); phgOptionsRegisterString("-sub_solver_opts", "Options for the subsolver", &sub_opts); /* a direct solver is preferable for the sparse matrix */ phgOptionsPreset("-solver mumps"); phgInit(&argc, &argv); phgPrintf( "----------------------------------------------------------------------------\n" "This code solves (A+UU^t)x=b using the Sherman-Morrison-Woodbury formula.\n" "Note: may use the following to disable use of the Sherman-Morrison-Woodbury\n" "algorithm and change to the default solver instead:\n" " -preonly_pc_type solver -preonly_pc_opts \"-solver_maxit 2000\"\n" "----------------------------------------------------------------------------\n" ); phgPrintf("Generating the linear system: N = %"dFMT", K = %"dFMT"\n", N, K); /* A is a distributed NxN SPD tridiagonal matrix (A = [-1, 2, -1]) */ n = N / phgNProcs + (phgRank < (N % phgNProcs) ? 1 : 0); A = phgMatCreate(phgComm, n, N); phgPrintf(" Generating matrix A.\n"); for (i = 0; i < n; i++) { /* diagonal */ phgMatAddEntry(A, i, i, 2.0); /* diagonal - 1 */ if (i > 0) phgMatAddEntry(A, i, i - 1, -1.0); else if (phgRank > 0) phgMatAddLGEntry(A, i, A->rmap->partition[phgRank] - 1, -1.0); /* diagonal + 1 */ if (i < n - 1) phgMatAddEntry(A, i, i + 1, -1.0); else if (phgRank < phgNProcs - 1) phgMatAddLGEntry(A, i, A->rmap->partition[phgRank] + n, -1.0); } phgMatAssemble(A); /* U is a K-component vector */ U = phgMapCreateVec(A->rmap, K); phgVecRandomize(U, 123); /* solver1 is the solver for A */ phgOptionsPush(); phgOptionsSetOptions(sub_opts); solver1 = phgMat2Solver(SOLVER_DEFAULT, A); phgOptionsPop(); /* x is a scratch vector */ x = phgMapCreateVec(A->rmap, 1); /* C is a KxK dense matrix, pvt is an integer array, they store the LU * factorization of (I + U^t*inv(A)*U) */ phgPrintf(" Generating the dense matrix I+U^t*inv(A)*U.\n"); C = phgCalloc(K * K, sizeof(*C)); pvt = phgAlloc(K * sizeof(*pvt)); for (i = 0; i < K; i++) { for (j = 0; j < n; j++) { solver1->rhs->data[j] = U->data[i * n + j]; x->data[j] = 0.0; } solver1->rhs->assembled = TRUE; phgSolverVecSolve(solver1, FALSE, x); for (j = 0; j < K; j++) for (k = 0; k < n; k++) C[i * K + j] += U->data[j * n + k] * x->data[k]; } #if USE_MPI if (U->map->nprocs > 1) { FLOAT *tmp = phgAlloc(K * K * sizeof(*tmp)); MPI_Allreduce(C, tmp, K * K, PHG_MPI_FLOAT, MPI_SUM, U->map->comm); phgFree(C); C = tmp; } #endif /* USE_MPI */ for (i = 0; i < K; i++) C[i * K + i] += 1.0; phgPrintf(" Factorizing the dense matrix I+U^t*inv(A)*U.\n"); phgSolverDenseLU(K, C, pvt); /* B is a matrix-free matrix representing A + U*U^t, B->mv_data is used * to pass A, U, solver1, C and pvt to callback functions */ B = phgMapCreateMatrixFreeMat(A->rmap, A->cmap, funcB, /* arguments carried over to CB functions */ A, U, solver1, C, pvt, NULL); /* solver is a PreOnly solver for B whose pc_proc is set to sherman(). * * Note: can also use pcg, gmres, or petsc for this solver, in this case * the solution obtained with the Sherman-Morisson formula is iteratively * refined. */ phgOptionsPush(); phgOptionsSetOptions("-solver preonly"); phgOptionsSetOptions(main_opts); solver = phgMat2Solver(SOLVER_DEFAULT, B); phgSolverSetPC(solver, solver, sherman); phgOptionsPop(); for (i = 0; i < n; i++) x->data[i] = 1.0; phgMatVec(MAT_OP_N, 1.0, B, x, 0.0, &solver->rhs); phgPrintf("Solving the linear system.\n"); /* reset initial solution to zero */ memset(x->data, 0, n * sizeof(*x->data)); phgSolverVecSolve(solver, TRUE, x); for (i = 0; i < n; i++) solver->rhs->data[i] = 1.0; phgVecAXPBY(-1.0, solver->rhs, 1.0, &x); phgPrintf("Checking the result: |x - x_exact| / |x_exact| = %lg\n", (double)phgVecNorm2(x, 0, NULL) / sqrt((double)N)); phgSolverDestroy(&solver); phgSolverDestroy(&solver1); phgMatDestroy(&A); phgMatDestroy(&B); phgVecDestroy(&U); phgVecDestroy(&x); phgFree(C); phgFree(pvt); phgFinalize(); return 0; }
/* * Preconditioning procedure: * p = Qp^-1 Fp Ap^-1 * q * u = F^-1 ( bu - Bt * p ) * */ static void pc_proc(void *ctx, VEC *b0, VEC **x0) { SOLVER *pc_solver = (SOLVER *) ctx; GRID *g = matF->rmap->dofs[0]->g; VEC *xu, *xp, *xu2, *xp2; int Nu = matF->rmap->nlocal; int Np = solver_Ap->mat->rmap->nlocal; INT verb = phgVerbosity; FLOAT *rhsF, *rhsAp, *rhsQp; double t; if (phgVerbosity > 0) phgVerbosity--; /* save old rhs */ rhsF = solver_F->rhs->data; rhsAp = solver_Ap->rhs->data; rhsQp = solver_Qp->rhs->data; xu = phgMapCreateVec(matF->rmap, 1); xu2 = phgMapCreateVec(matF->rmap, 1); xp = phgMapCreateVec(solver_Ap->mat->rmap, 1); xp2 = phgMapCreateVec(solver_Ap->mat->rmap, 1); /* Ap^-1 * q */ t = phgGetTime(NULL); solver_Ap->rhs->data = xp->data; solver_Ap->rhs->assembled = TRUE; solver_Ap->rhs_updated = TRUE; memcpy(xp->data, b0->data + Nu, sizeof(*xp->data) * Np); bzero(xp2->data, sizeof(*xp->data) * Np); phgSolverVecSolve(solver_Ap, FALSE, xp2); memcpy(xp->data, xp2->data, sizeof(*xp->data) * Np); if (verb > 0) phgPrintf("\t Ap: nits = %d, residual = %0.4le [%0.4lgMB %0.4lfs]\n", solver_Ap->nits, (double)solver_Ap->residual, phgMemoryUsage(g, NULL) / (1024.0 * 1024.0), phgGetTime(NULL) - t); /* Fp Ap^-1 * q */ phgMatVec(MAT_OP_N, 1.0, matFp, xp, 0., &xp2); memcpy(xp->data, xp2->data, sizeof(*xp->data) * Np); /* xp = Qp^-1 Fp Ap^-1 * q */ t = phgGetTime(NULL); solver_Qp->rhs->data = xp->data; solver_Qp->rhs->assembled = TRUE; solver_Qp->rhs_updated = TRUE; bzero(xp2->data, sizeof(*xp->data) * Np); phgSolverVecSolve(solver_Qp, FALSE, xp2); if (verb > 0) phgPrintf("\t Qp: nits = %d, residual = %0.4le [%0.4lgMB %0.4lfs]\n", solver_Qp->nits, (double)solver_Qp->residual, phgMemoryUsage(g, NULL) / (1024.0 * 1024.0), phgGetTime(NULL) - t); /* xu = F^-1 ( bu - Bt * xp ) */ t = phgGetTime(NULL); memcpy(xu->data, b0->data, sizeof(*xu->data) * Nu); phgMatVec(MAT_OP_N, -1.0, matBt, xp2, 1., &xu); if (solver_F->mat->cmap->nlocal == Nu) { /* use F in the preconditioning matrix */ solver_F->rhs->data = xu->data; solver_F->rhs->assembled = TRUE; solver_F->rhs_updated = TRUE; bzero(xu2->data, sizeof(*xu->data) * Nu); phgSolverVecSolve(solver_F, FALSE, xu2); } else { pc_proc1(solver_F, xu, &xu2); } if (verb > 0) phgPrintf("\t F: nits = %d, residual = %0.4le [%0.4lgMB %0.4lfs]\n", solver_F->nits, (double)solver_F->residual, phgMemoryUsage(g, NULL) / (1024.0 * 1024.0), phgGetTime(NULL) - t); /* Copy xu, xp to x */ memcpy((*x0)->data, xu2->data, sizeof(*xu->data) * Nu); memcpy((*x0)->data + Nu, xp2->data, sizeof(*xp->data) * Np); /* restore rhs */ solver_F->rhs->data = rhsF; solver_Ap->rhs->data = rhsAp; solver_Qp->rhs->data = rhsQp; phgVecDestroy(&xu); phgVecDestroy(&xu2); phgVecDestroy(&xp); phgVecDestroy(&xp2); phgVerbosity = verb; return; }
int main(int argc, char *argv[]) { GRID *g; DOF *u_h; MAT *A, *A0, *B; MAP *map; INT i; size_t nnz, mem, mem_peak; VEC *x, *y0, *y1, *y2; double t0, t1, dnz, dnz1, mflops, mop; char *fn = "../test/cube.dat"; FLOAT mem_max = 300; INT refine = 0; phgOptionsRegisterFilename("-mesh_file", "Mesh file", (char **)&fn); phgOptionsRegisterInt("-loop_count", "Loop count", &loop_count); phgOptionsRegisterInt("-refine", "Refinement level", &refine); phgOptionsRegisterFloat("-mem_max", "Maximum memory", &mem_max); phgInit(&argc, &argv); g = phgNewGrid(-1); if (!phgImport(g, fn, FALSE)) phgError(1, "can't read file \"%s\".\n", fn); phgRefineAllElements(g, refine); u_h = phgDofNew(g, DOF_DEFAULT, 1, "u_h", DofNoAction); while (TRUE) { phgPrintf("\n"); if (phgBalanceGrid(g, 1.2, 1, NULL, 0.)) phgPrintf("Repartition mesh, %d submeshes, load imbalance: %lg\n", g->nprocs, (double)g->lif); map = phgMapCreate(u_h, NULL); A = phgMapCreateMat(map, map); A->handle_bdry_eqns = TRUE; build_matrix(A, u_h); phgMatAssemble(A); /* Note: A is unsymmetric (A' != A) if boundary entries not removed */ phgMatRemoveBoundaryEntries(A); #if 0 /* test block matrix operation */ A0 = phgMatCreateBlockMatrix(g->comm, 1, 1, &A, NULL); #else A0 = A; #endif phgPrintf("%d DOF, %d elems, %d submeshes, matrix size: %d, LIF: %lg\n", DofGetDataCountGlobal(u_h), g->nleaf_global, g->nprocs, A->rmap->nglobal, (double)g->lif); /* test PHG mat-vec multiply */ x = phgMapCreateVec(A->cmap, 1); y1 = phgMapCreateVec(A->rmap, 1); phgVecRandomize(x, 123); phgMatVec(MAT_OP_N, 1.0, A0, x, 0.0, &y1); phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); for (i = 0; i < loop_count; i++) { phgMatVec(MAT_OP_N, 1.0, A0, x, 0.0, &y1); } t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); y0 = phgVecCopy(y1, NULL); nnz = A->nnz_d + A->nnz_o; #if USE_MPI dnz1 = nnz; MPI_Reduce(&dnz1, &dnz, 1, MPI_DOUBLE, MPI_SUM, 0, g->comm); #else dnz = nnz; #endif mop = loop_count * (dnz + dnz - A->rmap->nlocal) * 1e-6; phgPrintf("\n"); t1 -= t0; phgPrintf(" PHG: time %0.4lf, nnz %0.16lg, %0.2lfMF (%0.2lfMF)\n", t1, dnz, mop / (t1 == 0 ? 1. : t1), mflops); /* test trans(A)*x */ phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); for (i = 0; i < loop_count; i++) { phgMatVec(MAT_OP_T, 1.0, A0, x, 0.0, &y1); } t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); t1 -= t0; phgPrintf(" A'*x: time %0.4lf, nnz %0.16lg, %0.2lfMF (%0.2lfMF), " "err: %le\n", t1, dnz, mop / (t1 == 0 ? 1. : t1), mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y0, 1.0, &y1), 0, NULL)); /* time A * trans(A) */ phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); B = phgMatMat(MAT_OP_N, MAT_OP_N, 1.0, A, A, 0.0, NULL); t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); nnz = B->nnz_d + B->nnz_o; #if USE_MPI dnz1 = nnz; MPI_Reduce(&dnz1, &dnz, 1, MPI_DOUBLE, MPI_SUM, 0, g->comm); #else dnz = nnz; #endif /* compare B*x <--> A*A*x */ y2 = phgMatVec(MAT_OP_N, 1.0, B, x, 0.0, NULL); phgMatVec(MAT_OP_N, 1.0, A0, y0, 0.0, &y1); phgMatDestroy(&B); t1 -= t0; phgPrintf(" A*A: time %0.4lf, nnz %0.16lg, %0.2lfMF, err: %le\n", t1, dnz, mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y1, 1.0, &y2), 0, NULL)); #if USE_PETSC { Mat ma, mb; MatInfo info; Vec va, vb, vc; PetscScalar *vec; ma = phgPetscCreateMatAIJ(A); MatGetVecs(ma, PETSC_NULL, &va); VecDuplicate(va, &vb); VecGetArray(va, &vec); memcpy(vec, x->data, x->map->nlocal * sizeof(*vec)); VecRestoreArray(va, &vec); MatMult(ma, va, vb); phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); for (i = 0; i < loop_count; i++) { MatMult(ma, va, vb); } t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); VecGetArray(vb, &vec); memcpy(y1->data, vec, x->map->nlocal * sizeof(*vec)); VecRestoreArray(vb, &vec); MatGetInfo(ma, MAT_GLOBAL_SUM, &info); /*phgPrintf(" --------------------------------------------" "-------------------------\n");*/ phgPrintf("\n"); t1 -= t0; dnz = info.nz_used; phgPrintf(" PETSc: time %0.4lf, nnz %0.16lg, %0.2lfMF (%0.2lfMF), " "err: %le\n", t1, dnz, mop / (t1==0 ? 1.:t1), mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y0, 1.0, &y1), 0, NULL)); phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); for (i = 0; i < loop_count; i++) { MatMultTranspose(ma, va, vb); } t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); VecGetArray(vb, &vec); memcpy(y1->data, vec, x->map->nlocal * sizeof(*vec)); VecRestoreArray(vb, &vec); t1 -= t0; phgPrintf(" A'*x: time %0.4lf, nnz %0.16lg, %0.2lfMF (%0.2lfMF), " "err: %le\n", t1, dnz, mop / (t1==0 ? 1.:t1), mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y0, 1.0, &y1), 0, NULL)); phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); MatMatMult(ma, ma, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &mb); t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); t1 -= t0; MatGetInfo(mb, MAT_GLOBAL_SUM, &info); dnz = info.nz_used; VecDuplicate(va, &vc); /* compare B*x <--> A*A*x */ MatMult(ma, vb, vc); MatMult(mb, va, vb); VecGetArray(vb, &vec); memcpy(y1->data, vec, x->map->nlocal * sizeof(*vec)); VecRestoreArray(vb, &vec); VecGetArray(vc, &vec); memcpy(y2->data, vec, x->map->nlocal * sizeof(*vec)); VecRestoreArray(vc, &vec); phgPrintf(" A*A: time %0.4lf, nnz %0.16lg, %0.2lfMF, err: %le\n", t1, dnz, mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y1, 1.0, &y2), 0, NULL)); phgPetscMatDestroy(&mb); phgPetscMatDestroy(&ma); phgPetscVecDestroy(&va); phgPetscVecDestroy(&vb); phgPetscVecDestroy(&vc); } #endif /* USE_PETSC */ #if USE_HYPRE { HYPRE_IJMatrix ma; HYPRE_IJVector va, vb, vc; HYPRE_ParCSRMatrix par_ma; hypre_ParCSRMatrix *par_mb; HYPRE_ParVector par_va, par_vb, par_vc; HYPRE_Int offset, *ni, start, end; assert(sizeof(INT)==sizeof(int) && sizeof(FLOAT)==sizeof(double)); setup_hypre_mat(A, &ma); ni = phgAlloc(2 * A->rmap->nlocal * sizeof(*ni)); offset = A->cmap->partition[A->cmap->rank]; for (i = 0; i < A->rmap->nlocal; i++) ni[i] = i + offset; HYPRE_IJVectorCreate(g->comm, offset, offset + A->rmap->nlocal - 1, &va); HYPRE_IJVectorCreate(g->comm, offset, offset + A->rmap->nlocal - 1, &vb); HYPRE_IJVectorCreate(g->comm, offset, offset + A->rmap->nlocal - 1, &vc); HYPRE_IJVectorSetObjectType(va, HYPRE_PARCSR); HYPRE_IJVectorSetObjectType(vb, HYPRE_PARCSR); HYPRE_IJVectorSetObjectType(vc, HYPRE_PARCSR); HYPRE_IJVectorSetMaxOffProcElmts(va, 0); HYPRE_IJVectorSetMaxOffProcElmts(vb, 0); HYPRE_IJVectorSetMaxOffProcElmts(vc, 0); HYPRE_IJVectorInitialize(va); HYPRE_IJVectorInitialize(vb); HYPRE_IJVectorInitialize(vc); HYPRE_IJMatrixGetObject(ma, (void **)(void *)&par_ma); HYPRE_IJVectorGetObject(va, (void **)(void *)&par_va); HYPRE_IJVectorGetObject(vb, (void **)(void *)&par_vb); HYPRE_IJVectorGetObject(vc, (void **)(void *)&par_vc); HYPRE_IJVectorSetValues(va, A->cmap->nlocal, ni, (double *)x->data); HYPRE_IJVectorAssemble(va); HYPRE_IJVectorAssemble(vb); HYPRE_IJVectorAssemble(vc); HYPRE_IJMatrixGetRowCounts(ma, A->cmap->nlocal, ni, ni + A->rmap->nlocal); for (i = 0, nnz = 0; i < A->rmap->nlocal; i++) nnz += ni[A->rmap->nlocal + i]; #if USE_MPI dnz1 = nnz; MPI_Reduce(&dnz1, &dnz, 1, MPI_DOUBLE, MPI_SUM, 0, g->comm); #else dnz = nnz; #endif HYPRE_ParCSRMatrixMatvec(1.0, par_ma, par_va, 0.0, par_vb); phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); for (i = 0; i < loop_count; i++) { HYPRE_ParCSRMatrixMatvec(1.0, par_ma, par_va, 0.0, par_vb); } t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); HYPRE_IJVectorGetValues(vb, A->rmap->nlocal, ni, (double*)y1->data); /*phgPrintf(" --------------------------------------------" "-------------------------\n");*/ phgPrintf("\n"); t1 -= t0; phgPrintf(" HYPRE: time %0.4lf, nnz %0.16lg, %0.2lfMF (%0.2lfMF), " "err: %le\n", t1, dnz, mop / (t1==0 ? 1.:t1), mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y0, 1.0, &y1), 0, NULL)); phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); for (i = 0; i < loop_count; i++) { HYPRE_ParCSRMatrixMatvecT(1.0, par_ma, par_va, 0.0, par_vb); } t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); HYPRE_IJVectorGetValues(vb, A->rmap->nlocal, ni, (double*)y1->data); t1 -= t0; phgPrintf(" A'*x: time %0.4lf, nnz %0.16lg, %0.2lfMF (%0.2lfMF), " "err: %le\n", t1, dnz, mop / (t1==0 ? 1.:t1), mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y0, 1.0, &y1), 0, NULL)); phgPerfGetMflops(g, NULL, NULL); /* reset flops counter */ t0 = phgGetTime(NULL); /* Note: 'HYPRE_ParCSRMatrix' is currently typedef'ed to * 'hypre_ParCSRMatrix *' */ par_mb = hypre_ParMatmul((hypre_ParCSRMatrix *)par_ma, (hypre_ParCSRMatrix *)par_ma); t1 = phgGetTime(NULL); mflops = phgPerfGetMflops(g, NULL, NULL); start = hypre_ParCSRMatrixFirstRowIndex(par_mb); end = hypre_ParCSRMatrixLastRowIndex(par_mb) + 1; for (i = start, nnz = 0; i < end; i++) { HYPRE_Int ncols; hypre_ParCSRMatrixGetRow(par_mb, i, &ncols, NULL, NULL); hypre_ParCSRMatrixRestoreRow(par_mb, i, &ncols, NULL, NULL); nnz += ncols; } #if USE_MPI dnz1 = nnz; MPI_Reduce(&dnz1, &dnz, 1, MPI_DOUBLE, MPI_SUM, 0, g->comm); #else dnz = nnz; #endif /* compare B*x <--> A*A*x */ HYPRE_ParCSRMatrixMatvec(1.0, par_ma, par_vb, 0.0, par_vc); HYPRE_ParCSRMatrixMatvec(1.0, (void *)par_mb, par_va, 0.0, par_vb); HYPRE_IJVectorGetValues(vb, A->rmap->nlocal, ni, (double*)y1->data); HYPRE_IJVectorGetValues(vc, A->rmap->nlocal, ni, (double*)y2->data); hypre_ParCSRMatrixDestroy((par_mb)); t1 -= t0; phgPrintf(" A*A: time %0.4lf, nnz %0.16lg, %0.2lfMF, err: %le\n", t1, dnz, mflops, (double)phgVecNorm2(phgVecAXPBY(-1.0, y1, 1.0, &y2), 0, NULL)); phgFree(ni); HYPRE_IJMatrixDestroy(ma); HYPRE_IJVectorDestroy(va); HYPRE_IJVectorDestroy(vb); HYPRE_IJVectorDestroy(vc); } #endif /* USE_HYPRE */ if (A0 != A) phgMatDestroy(&A0); #if 0 if (A->rmap->nglobal > 1000) { VEC *v = phgMapCreateVec(A->rmap, 3); for (i = 0; i < v->map->nlocal; i++) { v->data[i + 0 * v->map->nlocal] = 1 * (i + v->map->partition[g->rank]); v->data[i + 1 * v->map->nlocal] = 2 * (i + v->map->partition[g->rank]); v->data[i + 2 * v->map->nlocal] = 3 * (i + v->map->partition[g->rank]); } phgMatDumpMATLAB(A, "A", "A.m"); phgVecDumpMATLAB(v, "v", "v.m"); phgFinalize(); exit(0); } #endif phgMatDestroy(&A); phgVecDestroy(&x); phgVecDestroy(&y0); phgVecDestroy(&y1); phgVecDestroy(&y2); phgMapDestroy(&map); mem = phgMemoryUsage(g, &mem_peak); dnz = mem / (1024.0 * 1024.0); dnz1 = mem_peak / (1024.0 * 1024.0); /*phgPrintf(" --------------------------------------------" "-------------------------\n");*/ phgPrintf("\n"); phgPrintf(" Memory: current %0.4lgMB, peak %0.4lgMB\n", dnz, dnz1); #if 0 { static int loop_count = 0; if (++loop_count == 4) break; } #endif if (mem_peak > 1024 * (size_t)1024 * mem_max) break; phgRefineAllElements(g, 1); } phgDofFree(&u_h); phgFreeGrid(&g); phgFinalize(); return 0; }