static void build_matrix(MAT *mat, DOF *u_h) { int N = u_h->type->nbas; /* number of basis functions in an element */ int i, j; GRID *g = u_h->g; ELEMENT *e; FLOAT A[N][N], buf[N]; INT I[N]; assert(u_h->dim == 1); ForAllElements(g, e) { /* compute \int \grad\phi_j \cdot \grad\phi_i making use of symmetry */ for (i = 0; i < N; i++) { I[i] = phgMapE2L(mat->cmap, 0, e, i); for (j = 0; j <= i; j++) A[j][i] = A[i][j] = phgQuadGradBasDotGradBas(e, u_h, j, u_h, i, QUAD_DEFAULT); } /* loop on basis functions */ for (i = 0; i < N; i++) { if (phgDofDirichletBC(u_h, e, i, NULL, buf, NULL, DOF_PROJ_NONE)) { phgMatAddEntries(mat, 1, I + i, N, I, buf); } else { /* interior node */ phgMatAddEntries(mat, 1, I + i, N, I, A[i]); } } } }
static void build_matrices(MAT *matA, MAT *matM, DOF *u_h, DOF *p_h) { int N = u_h->type->nbas; /* number of basis functions in an element */ int M = p_h->type->nbas; int i, j; GRID *g = u_h->g; ELEMENT *e; FLOAT A[N][N], B[N][N], C[M][N]; INT I[N], Ip[N]; INT k; assert(DofTypeDim(u_h) == Dim); assert(u_h->dim == 1); ForAllElements(g, e) { for (i = 0; i < N; i++) I[i] = phgMapE2L(matA->rmap, 0, e, i); for (i = 0; i < M; i++) Ip[i] = phgMapE2L(matA->rmap, 1, e, i); for (i = 0; i < N; i++) { for (k = 0; k < M; k++) { /* \int \grad\psi_k\cdot\phi_i */ C[k][i] = phgQuadGradBasDotBas(e, p_h, k, u_h, i, QUAD_DEFAULT); } for (j = 0; j <= i; j++) { /* \int \phi_j\cdot\phi_i */ B[j][i] = B[i][j] = phgQuadBasDotBas(e, u_h, j, u_h, i, QUAD_DEFAULT); /* \int \curl\phi_j\cdot\curl\phi_i */ A[j][i] = A[i][j] = phgQuadCurlBasDotCurlBas(e, u_h, j, u_h, i, QUAD_DEFAULT); } } /* loop on basis functions of u_h */ for (i = 0; i < N; i++) { if (phgDofDirichletBC(u_h, e, i, NULL, NULL, NULL, DOF_PROJ_CROSS)) continue; phgMatAddEntries(matA, 1, I + i, N, I, A[i]); phgMatAddEntries(matM, 1, I + i, N, I, B[i]); for (j = 0; j < M; j++) phgMatAddEntry(matA, I[i], Ip[j], C[j][i]); } /* loop on basis functions of p_h */ for (i = 0; i < M; i++) phgMatAddEntries(matA, 1, Ip + i, N, I, C[i]); } }
static void build_linear_system(SOLVER *solver, DOF *u_h, DOF *f_h) { int i, j; GRID *g = u_h->g; ELEMENT *e; assert(u_h->dim == 1); ForAllElements(g, e) { int N = DofGetNBas(u_h, e); /* number of bases in the element */ FLOAT A[N][N], rhs[N], buffer[N]; INT I[N]; /* compute \int \grad\phi_j \cdot \grad\phi_i making use of symmetry */ for (i = 0; i < N; i++) { I[i] = phgSolverMapE2L(solver, 0, e, i); for (j = 0; j <= i; j++) { A[j][i] = A[i][j] = /* stiffness */ phgQuadGradBasDotGradBas(e, u_h, j, u_h, i, QUAD_DEFAULT) + /* mass */ a * phgQuadBasDotBas(e, u_h, j, u_h, i, QUAD_DEFAULT); } } /* loop on basis functions */ for (i = 0; i < N; i++) { if (phgDofDirichletBC(u_h, e, i, func_u, buffer, rhs+i, DOF_PROJ_NONE)) { phgSolverAddMatrixEntries(solver, 1, I + i, N, I, buffer); } else { /* interior node */ /* right hand side = \int f * phi_i */ phgQuadDofTimesBas(e, f_h, u_h, i, QUAD_DEFAULT, rhs + i); phgSolverAddMatrixEntries(solver, 1, I + i, N, I, A[i]); } } phgSolverAddRHSEntries(solver, N, I, rhs); }
static void build_matrices(MAT *ma, MAT *mb, VEC *rhs_vec, DOF *u_h, DOF *f) /* build stiffness (ma) and mass (mb) matrices */ { int N = u_h->type->nbas; /* number of element basis functions */ int i, j, k, l; GRID *g = u_h->g; ELEMENT *e; FLOAT A[N][3][N][3], buffer[N], rhs[N][3], d; static FLOAT *B0 = NULL; FLOAT B[N][N]; INT I[N][3], J[3][N]; QUAD *quad = phgQuadGetQuad3D((u_h->type->order - 1) * 2); assert(u_h->type->dim == 1 && u_h->dim == 3); assert(mb == NULL || u_h->type->invariant); if (mb != NULL && B0 == NULL && g->nroot > 0) { /* (\int \phi_j\cdot\phi_i)/vol is independent of element */ FreeAtExit(B0); B0 = phgAlloc(N * N * sizeof(*B0)); e = g->roots; d = 1. / phgGeomGetVolume(g, e); for (i = 0; i < N; i++) for (j = 0; j <= i; j++) B0[i * N + j] = B0[i + j * N] = d * phgQuadBasDotBas(e, u_h, j, u_h, i, QUAD_DEFAULT); } ForAllElements(g, e) { d = phgGeomGetVolume(g, e) * rho; for (i = 0; i < N; i++) { for (k = 0; k < 3; k++) J[k][i] = I[i][k] = phgMapE2L(ma->rmap, 0, e, i * 3 + k); for (j = 0; j <= i; j++) { /* the stiffness matrix */ stiffness_matrix(e, u_h, i, j, E, nu, 3 * N, &A[i][0][j][0], quad); if (j != i) { for (k = 0; k < 3; k++) for (l = 0; l < 3; l++) A[j][k][i][l] = A[i][l][j][k]; } /* the mass matrix: \int \phi_i\cdot\phi_j */ if (mb != NULL) B[j][i] = B[i][j] = B0[i * N + j] * d; } } #if 0 /* print element stiffness matrix */ FLOAT *p = (void *)A; printf("vol = %0.16lg;\n", phgGeomGetVolume(g, e)); printf("a = [\n"); for (i = 0; i < 3 * N; i++) { for (j = 0; j < 3 * N; j++) printf("%0.16lg ", p[i * 3 * N + j]); printf("\n"); } printf("];\n"); printf("b = [\n"); for (i = 0; i < N; i++) for (k = 0; k < 3; k++) { for (j = 0; j < N; j++) for (l = 0; l < 3; l++) printf("%0.16lg ", l == k ? B[i][j] : 0.0); printf("\n"); } printf("];\n"); exit(0); #endif /* loop on basis functions */ for (i = 0; i < N; i++) { if (phgDofDirichletBC(u_h, e, i, func_g, buffer, rhs[i], DOF_PROJ_NONE)) { /* Dirichlet boundary */ for (k = 0; k < 3; k++) { phgMatAddEntries(ma, 1, I[i] + k, N, J[k], buffer); if (mb != NULL) phgMatAddEntries(mb, 1, I[i] + k, N, J[k], buffer); } } else { /* interior node */ phgMatAddEntries(ma, 3, I[i], 3 * N, I[0], &A[i][0][0][0]); if (mb == NULL) { phgQuadDofTimesBas(e, f, u_h, i, QUAD_DEFAULT, rhs[i]); for (k = 0; k < 3; k++) rhs[i][k] = -rhs[i][k]; } else { for (k = 0; k < 3; k++) phgMatAddEntries(mb, 1, I[i] + k, N, J[k], B[i]); } } } if (rhs_vec != NULL) phgVecAddEntries(rhs_vec, 0, N * 3, I[0], rhs[0]); }
/**************************************************************** * Build RHS which is the residual of the nonlinear system. ***************************************************************/ void phgNSBuildSolverURHS(NSSolver *ns) { GRID *g = ns->g; SIMPLEX *e; SOLVER *solver_u = ns->solver_u; int i, k, l, q, s; FLOAT *dt = ns->dt; BOOLEAN tstep_minus = (ns->u[-1] != NULL); VEC *vec_rhs = phgMapCreateVec(solver_u->rhs->map, 1); FLOAT Theta = _nsp->Theta, nu = _nsp->nu, Thet1; int viscosity_type = ns->viscosity_type; SURF_BAS *surf_bas = ns->surf_bas; DOF *surf_dof = surf_bas->dof; BOOLEAN *rotated = surf_bas->rotated; const FLOAT *Trans = DofData(surf_dof); #if STEADY_STATE assert(fabs(Theta - 1) < 1e-12); Thet1 = 0; Unused(Thet1); Unused(dt); #else Thet1 = 1 - Theta; Unused(dt); #endif /* STEADY_STATE */ phgPrintf(" DB_mask: ["); for (k = 0; k < Dim; k++) phgPrintf("%d ", ns->u[1]->DB_masks[k]); phgPrintf("] "); nu_max = -1e10; nu_min = +1e10; phgVecDisassemble(vec_rhs); ForAllElements(g, e) { int M = ns->u[1]->type->nbas; /* num of bases of Velocity */ int N = ns->p[1]->type->nbas; /* num of bases of Pressure */ int order = DofTypeOrder(ns->u[1], e) * 3 - 1; /* Note: * quad order is really high here, * highest order term (u \nabla u, phi) */ FLOAT bufu[M], bufp[N], rhsu[M][Dim], rhsp[N]; INT Iu[M][Dim], Ip[N]; QUAD *quad; FLOAT vol, area, det; const FLOAT *w, *p, *normal, **vu, *vu_queue[3], *vf[2], *gu[2], *vp[2], *vw, *vT; FLOAT *vf_cache[2]; vu = vu_queue + 1; quad = phgQuadGetQuad3D(order); vu[0] = phgQuadGetDofValues(e, ns->u[0], quad); /* u^{n} */ vp[0] = phgQuadGetDofValues(e, ns->p[0], quad); /* p^{n} */ gu[0] = phgQuadGetDofValues(e, ns->gradu[0], quad); /* grad u^{n} */ vw = phgQuadGetDofValues(e, ns->wind, quad); /* wind */ vT = phgQuadGetDofValues(e, ns->T[1], quad); /* T^{n} */ if (tstep_minus) { vu[-1] = phgQuadGetDofValues(e, ns->u[-1], quad); /* u^{n-1} */ } else { vu[-1] = vu[0]; } #if STEADY_STATE || TIME_DEP_NON vu[1] = phgQuadGetDofValues(e, ns->u[1], quad); /* u^{n+1} */ gu[1] = phgQuadGetDofValues(e, ns->gradu[1], quad); /* grad u^{n} */ vp[1] = phgQuadGetDofValues(e, ns->p[1], quad); /* p^{n+1} */ #else TIME_DEP_LINEAR_ENTRY; /* Unavailable */ #endif /* STEADY_STATE || TIME_DEP_NON */ Unused(l); Unused(vf); Unused(vf_cache); if (!_nsp->no_extern_source) { /* cache f values */ for (l = 0; l < 2; l++) { const FLOAT *cache; size_t cache_size; setFuncTime(ns->time[l]); /* set static time in ins-test.c */ /* cache f */ cache_size = Dim * quad->npoints * sizeof(FLOAT); cache = phgQuadGetFuncValues(g, e, Dim, func_f, quad); vf[l] = vf_cache[l] = phgAlloc(cache_size); memcpy(vf_cache[l], cache, cache_size); phgQuadGetFuncValues(NULL, NULL, 0, NULL, NULL); /* clear cache */ } } /* Global Matrix */ Bzero(rhsu); Bzero(rhsp); p = quad->points; w = quad->weights; for (q = 0; q < quad->npoints; q++) { phgGeomGetCurvedJacobianAtLambda(g, e, p, &det); vol = fabs(det / 6.); /* rhs u */ for (i = 0; i < M; i++) { /* interior node or Neumann */ const FLOAT *gi_u = phgQuadGetBasisValues(e, ns->u[1], i, quad) + q; /* phi_i */ const FLOAT *ggi_u = phgQuadGetBasisCurvedGradient(e, ns->u[1], i, quad, q); /* grad phi_i */ for (k = 0; k < Dim; k++) { #if ICE_BENCH_TEST nu = get_effective_viscosity(gu[1], 0, 0, viscosity_type); FLOAT eu[DDim]; MAT3_SYM(gu[1], eu); rhsu[i][k] += vol*(*w) * EQU_SCALING * (- nu * INNER_PRODUCT(eu+k*Dim, ggi_u) + (*vp[1]) * *(ggi_u+k) * LEN_SCALING * PRES_SCALING ); /* left */ if (k == Z_DIR) { const FLOAT rho = RHO_ICE; const FLOAT grav = GRAVITY; const FLOAT a = SEC_PER_YEAR; const FLOAT f = rho*grav * EQU_SCALING * LEN_SCALING2; Unused(a); rhsu[i][k] += vol*(*w) * (-f * (*gi_u) ); /* right */ } #elif ESIMINT_TEST || \ HEINO_TEST || \ TEST_CASE == ICE_GREEN_LAND nu = get_effective_viscosity(gu[1], *vT, 0, viscosity_type); FLOAT eu[DDim]; MAT3_SYM(gu[1], eu); rhsu[i][k] += vol*(*w) * EQU_SCALING * (- nu * INNER_PRODUCT(eu+k*Dim, ggi_u) + (*vp[1]) * *(ggi_u+k) * LEN_SCALING * PRES_SCALING ); /* left */ if (k == Z_DIR) { const FLOAT rho = RHO_ICE; const FLOAT grav = GRAVITY; const FLOAT a = SEC_PER_YEAR; const FLOAT f = rho*grav * EQU_SCALING * LEN_SCALING2; Unused(a); rhsu[i][k] += vol*(*w) * (-f * (*gi_u) ); /* right */ } #elif STEADY_STATE rhsu[i][k] += vol*(*w) * (- nu * INNER_PRODUCT(gu[1]+k*Dim, ggi_u) + (*vp[1]) * *(ggi_u+k) ); /* left */ if (!_nsp->no_extern_source) rhsu[i][k] += vol*(*w) * (*(vf[1]+k) * (*gi_u) ); /* right */ #elif TIME_DEP_NON rhsu[i][k] -= vol*(*w) * ((vu[1][k] - vu[0][k]) * (*gi_u) / dt[0] + Theta * (nu * INNER_PRODUCT(gu[1]+k*Dim, ggi_u) ) - (*vp[1]) * *(ggi_u+k) + Thet1 * (nu * INNER_PRODUCT(gu[0]+k*Dim, ggi_u) ) ); /* left */ if (!_nsp->no_extern_source) rhsu[i][k] += vol*(*w) * (Theta * *(vf[1]+k) * (*gi_u) + Thet1 * *(vf[0]+k) * (*gi_u) ); /* right */ #else TIME_DEP_LINEAR_ENTRY; /* Unavailable */ #endif /* STEADY_STATE */ } } /* rhs p */ for (i = 0; i < N; i++) { const FLOAT *gi_p = phgQuadGetBasisValues(e, ns->p[1], i, quad) + q; /* psi_i */ FLOAT divu1 = gu[1][0] + gu[1][4] + gu[1][8]; //FLOAT divu0 = gu[0][0] + gu[0][4] + gu[0][8]; rhsp[i] += vol*(*w) * (divu1 * (*gi_p) ); } if (tstep_minus) vu[-1] += Dim; #if STEADY_STATE || TIME_DEP_NON vu[1] += Dim; gu[1] += Dim*Dim; vp[1]++; #else TIME_DEP_LINEAR; /* Unavailable */ #endif /* STEADY_STATE || TIME_DEP_NON */ vu[0] += Dim; gu[0] += Dim * Dim; vp[0]++; vw += Dim; if (!_nsp->no_extern_source) { vf[0] += Dim; vf[1] += Dim; } vT++; w++; p += Dim + 1; } if (!_nsp->no_extern_source) { phgFree(vf_cache[0]); phgFree(vf_cache[1]); } normal = NULL; Unused(normal); area = 0; Unused(area); if (!_nsp->enclosed_flow) { /* slip boundary */ for (s = 0; s < NFace; s++) { if (e->bound_type[s] & INFLOW) { int v0, v1, v2; int nbas_face = NbasFace(ns->u[1]); SHORT bases[nbas_face]; FLOAT lambda[Dim + 1], x,y,z, beta; order = DofTypeOrder(ns->u[1], e) * 3 - 1; phgDofGetBasesOnFace(ns->u[1], e, s, bases); v0 = GetFaceVertex(s, 0); v1 = GetFaceVertex(s, 1); v2 = GetFaceVertex(s, 2); lambda[s] = 0.; area = phgGeomGetFaceArea(g, e, s); normal = phgGeomGetFaceOutNormal(g, e, s); quad = phgQuadGetQuad2D(order); p = quad->points; w = quad->weights; for (q = 0; q < quad->npoints; q++) { FLOAT vu[Dim]; lambda[v0] = *(p++); lambda[v1] = *(p++); lambda[v2] = *(p++); phgGeomLambda2XYZ(g, e, lambda, &x, &y, &z); func_beta(x, y, z, &beta); phgDofEval(ns->u[1], e, lambda, vu); for (i = 0; i < nbas_face; i++) { int ii = bases[i]; FLOAT gi_u = *ns->u[1]->type->BasFuncs(ns->u[1], e, ii, ii + 1, lambda); for (k = 0; k < Dim; k++) { #if STEADY_STATE rhus[ii][k] += 0.; #elif TIME_DEP_NON # if USE_SLIDING_BC abort(); rhsu[ii][k] += SIGN_FRICTION * area*(*w) * beta * vu[k] * (gi_u) * EQU_SCALING * LEN_SCALING; # else Unused(gi_u); # endif #else TIME_DEP_LINEAR_ENTRY; /* Unavailable */ #endif /* STEADY_STATE */ } } /* end of bas_i */ w++; } /* end of quad point */ } /* end of face outflow */ } /* end of all outflow face in element */ } /* end out flow boundary */ #if USE_SLIDING_BC /* Rotate bases */ for (i = 0; i < M; i++) { INT id = phgDofMapE2D(surf_dof, e, i * (Dim*Dim)) / (Dim*Dim); if (!rotated[id]) continue; const FLOAT *trans = Trans + id*(Dim*Dim); trans_left(&rhsu[i][0], 1, 1, trans); } #else Unused(Trans); Unused(rotated); #endif /* Map: Element -> system */ for (i = 0; i < M; i++) for (k = 0; k < Dim; k++) Iu[i][k] = phgMapE2L(solver_u->rhs->map, 0, e, i * Dim + k); for (i = 0; i < N; i++) Ip[i] = phgMapE2L(solver_u->rhs->map, 1, e, i); /* set velocity dirichlet bdry */ FLOAT tmp[Dim]; for (i = 0; i < M; i++) for (k = 0; k < Dim; k++) if (phgDofDirichletBC_(ns->u[1], e, i*Dim+k, NULL, bufu, tmp, DOF_PROJ_NONE)) { rhsu[i][k] = 0.; } #if STEADY_STATE || TIME_DEP_NON /* set pressure dirichlet bdry for pinned point */ for (i = 0; i < N; i++) if (phgDofDirichletBC(ns->p[1], e, i, NULL, bufp, &rhsp[i], DOF_PROJ_NONE)) { if (!_nsp->enclosed_flow) phgError(1, "No dirichlet bc for Unenclosed flow!\n"); if (_nsp->pin_node) { # if PIN_AT_ROOT if (g->rank != 0) phgError(1, "Pinned node only on rank 0!\n"); if (g, e->verts[i] != ns->pinned_node_id) phgError(1, "Build rhs: pinned node e:%d, bas:%d, [%d] and [%d] " "doesn't coincide when build RHS!\n", e->index, i, e->verts[i], ns->pinned_node_id); # else if (GlobalVertex(g, e->verts[i]) != ns->pinned_node_id) phgError(1, "Build rhs: pinned node e:%d, bas:%d, [%d] and [%d] " "doesn't coincide when build RHS!\n", e->index, i, e->verts[i], ns->pinned_node_id); # endif /* PIN_AT_ROOT */ } } #else TIME_DEP_LINEAR; /* Unavailable */ #endif /* STEADY_STATE || TIME_DEP_NON */ /* Global res */ phgVecAddEntries(vec_rhs, 0, M * Dim, Iu[0], &rhsu[0][0]); phgVecAddEntries(vec_rhs, 0, N, Ip, rhsp); } /* end element */
/*************************************************************************** * Build matrices *F, *Fu, *B, *Bt, *Fp, *Ap, and *Qp used * by the solvers. **************************************************************************/ static void build_matrices(SOLVER *solver, SOLVER *pc, DOF **dofs, MAT **mats, LTYPE type) { DOF *u = dofs[0], *p = dofs[1]; DOF *f, *pbc, *gn[3], *gradu, *divu; int M = u->type->nbas; /* num of bases of Velocity */ int N = p->type->nbas; /* num of bases of Pressure */ int i, j, k, l; GRID *g = u->g; ELEMENT *e; MAT *matB, *matC, *matAp, *matQp; FLOAT F[M][Dim][M][Dim], Fu[M][M], B[N][M][Dim], Bt[M][Dim][N], Ap[N][N], Fp[N][N], Qp[N][N], bufu[M], bufp[N], tmp[9]; INT Iu[M][Dim], Ju[Dim][M], Iu0[M], Ip[N]; /* Unpack Dofs */ unpackDof(dofs, 9, &u, &p, &gradu, &divu, &f, &pbc, &gn[0], &gn[1], &gn[2]); unpackMat(mats, 4, &matB, &matC, &matAp, &matQp); #if 0 #warning remove me! SOLVER *s = phgSolverCreate(SOLVER_PCG, p, NULL); #endif ForAllElements(g, e) { /* Map: Element -> system */ for (i = 0; i < M; i++) { for (k = 0; k < Dim; k++) Ju[k][i] = Iu[i][k] = phgMapE2L(matF->cmap, 0, e, i * Dim + k); Iu0[i] = phgMapE2L(matFu->cmap, 0, e, i); } for (i = 0; i < N; i++) { Ip[i] = phgMapE2L(matFp->cmap, 0, e, i); } /* A V W */ for (i = 0; i < M; i++) { for (j = 0; j < M; j++) { FLOAT m, a, v, w[9]; /* \phi_j \dot \phi_i */ m = phgQuadBasDotBas(e, u, j, u, i, QUAD_DEFAULT); /* (\grad \phi_j) \cdot (\grad \phi_i) */ a = nu * phgQuadGradBasDotGradBas(e, u, j, u, i, QUAD_DEFAULT); /* (u \cdot \grad) \phi_j \times \phi_i, 1 item */ v = phgQuadDofDotGradBasBas(e, u, u, j, i, QUAD_DEFAULT); Fu[i][j] = a + v; if (type == PICARD) { memset(w, 0, sizeof(w)); } else { /* \phi_j (\grad u) \times \phi_i, 9 items */ phgQuadGradDofBasDotBas(e, gradu, u, j, i, QUAD_DEFAULT, w); } for (k = 0; k < Dim; k++) { for (l = 0; l < Dim; l++) { F[i][k][j][l] = Theta * dt * *(w + k * Dim + l); if (k == l) F[i][k][j][k] += m + Theta * dt * (a + v); } } } } /* B Bt */ for (i = 0; i < M; i++) { for (j = 0; j < N; j++) { FLOAT bxyz[3]; phgQuadGradBasDotBas3(e, u, i, p, j, bxyz, QUAD_DEFAULT); for (k = 0; k < Dim; k++) Bt[i][k][j] = B[j][i][k] = -Theta * dt * bxyz[k]; } } /* Ap Qp; Fp(update) */ for (i = 0; i < N; i++) { for (j = 0; j < N; j++) { FLOAT ap, qp, cp; ap = phgQuadGradBasDotGradBas(e, p, j, p, i, QUAD_DEFAULT); Ap[i][j] = Theta * dt * ap; qp = phgQuadBasDotBas(e, p, j, p, i, QUAD_DEFAULT); Qp[i][j] = Theta * dt * qp; cp = phgQuadDofDotGradBasBas(e, u, p, j, i, QUAD_DEFAULT); Fp[i][j] = qp + Theta * dt * (nu * ap + cp); } } /* Global Matrix */ for (i = 0; i < M; i++) { /* Dirichle Boundary for velocity. */ if (phgDofDirichletBC(u, e, i, NULL, bufu, NULL, DOF_PROJ_NONE)) { for (k = 0; k < Dim; k++) phgMatAddEntries(matF, 1, Iu[i] + k, M, Ju[k], bufu); phgMatAddEntries(matFu, 1, Iu0 + i, M, Iu0, bufu); } else { /* interior node Or Neumann */ /* Matrix F */ phgMatAddEntries(matF, Dim, Iu[i], M * Dim, Iu[0], &(F[i][0][0][0])); /* Matrix Fu */ phgMatAddEntries(matFu, 1, Iu0 + i, M, Iu0, &(Fu[i][0])); /* Matrix Bt */ for (k = 0; k < Dim; k++) phgMatAddEntries(matBt, 1, &Iu[i][k], N, Ip, &Bt[i][k][0]); } } /* end of Block (1,1), (1,2) */ for (i = 0; i < N; i++) phgMatAddEntries(matB, 1, Ip + i, M * Dim, Iu[0], &B[i][0][0]); /* Matrix Ap Qp Fp */ for (i = 0; i < N; i++) { if (phgDofDirichletBC(pbc, e, i, NULL, bufp, tmp, DOF_PROJ_NONE)) { /* Dirichle Boundary for pressure PC. */ phgMatAddEntries(matAp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(matFp, 1, Ip + i, N, Ip, bufp); } else { /* interior node Or Neumann */ phgMatAddEntries(matAp, 1, Ip + i, N, Ip, Ap[i]); phgMatAddEntries(matFp, 1, Ip + i, N, Ip, Fp[i]); } #if 0 /* use only diagonal of the mass matrix in the preconditioner */ phgMatAddEntries(matQp, 1, Ip + i, 1, Ip + i, &Qp[i][i]); #else /* use full mass matrix in the preconditioner */ phgMatAddEntries(matQp, 1, Ip + i, N, Ip, Qp[i]); #endif #if 0 phgMatAddEntries(s->mat, 1, Ip + i, N, Ip, Qp[i]); #endif } } /* end element */ #if 0 VEC *v = phgMapCreateVec(s->rhs->map, 1); memcpy(s->rhs->data, solver->rhs->data + DofGetDataCount(u), s->rhs->map->nlocal * sizeof(FLOAT)); s->rhs_updated = TRUE; s->rhs->assembled = TRUE; s->rtol = 1e-10; s->maxit = 10000; phgSolverVecSolve(s, TRUE, v); phgPrintf("v = %e\n", phgVecNormInfty(v, 0, NULL)); phgVecDestroy(&v); phgSolverDestroy(&s); #endif return; }
/**************************************************************** * Build RHS which is the residual of the nonlinear system. ***************************************************************/ static void build_rhs(SOLVER *solver, SOLVER *pc, DOF **dofs, MAT **mats) { DOF *u = dofs[0], *p = dofs[1]; DOF *f, *pbc, *gn[3], *gradu, *divu, *lapu, *gradp, *f0; int M = u->type->nbas; /* num of bases of Velocity */ int i, k, s; GRID *g = u->g; ELEMENT *e; FLOAT bufu[M], resu[M][Dim], tmp[9]; INT Iu[M][Dim]; /* Unpack Dofs */ unpackDof(dofs, 9, &u, &p, &gradu, &divu, &f, &pbc, &gn[0], &gn[1], &gn[2]); lapu = phgDofDivergence(gradu, NULL, NULL, NULL); gradp = phgDofGradient(p, NULL, NULL, NULL); time -= dt; f0 = phgDofNew(g, DOF_HB6, 3, "p_n", func_f); time += dt; ForAllElements(g, e) { /* Map: Element -> system */ for (i = 0; i < M; i++) for (k = 0; k < Dim; k++) Iu[i][k] = phgMapE2L(solver->rhs->map, 0, e, i * Dim + k); /* Global Matrix */ bzero(resu, sizeof(resu)); for (i = 0; i < M; i++) { /* Dirichle Boundary for velocity. */ if (phgDofDirichletBC(u, e, i, func_u, bufu, &resu[i][0], DOF_PROJ_NONE)) { /* set velocity at Dirichlet bdry */ } else { /* interior node or Neumann */ /* (u(t_n), \phi) */ phgQuadDofTimesBas(e, u, u, i, QUAD_DEFAULT, tmp); for (k = 0; k < Dim; k++) resu[i][k] = tmp[k]; /* (f, \phi_i) */ phgQuadDofTimesBas(e, f0, u, i, QUAD_DEFAULT, tmp); for (k = 0; k < Dim; k++) resu[i][k] += dt * (1 - Theta) * tmp[k]; phgQuadDofTimesBas(e, f, u, i, QUAD_DEFAULT, tmp); for (k = 0; k < Dim; k++) resu[i][k] += dt * Theta * tmp[k]; /* -( ((u.\grad) u, \phi) */ phgQuadDofDotGradDofBas(e, u, gradu, i, QUAD_DEFAULT, tmp); for (k = 0; k < Dim; k++) resu[i][k] -= dt * (1 - Theta) * tmp[k]; /* +\nu ( lap(u(t_n)), \phi) */ phgQuadDofTimesBas(e, lapu, u, i, QUAD_DEFAULT, tmp); for (k = 0; k < Dim; k++) resu[i][k] += (1 - Theta) * dt * nu * tmp[k]; /* -(gradp(t_n), \phi) */ phgQuadDofTimesBas(e, gradp, u, i, QUAD_DEFAULT, tmp); for (k = 0; k < Dim; k++) resu[i][k] -= (1 - Theta) * dt * tmp[k]; } } /* end of Block (1,1), (1,2) */ /* Neumann Bdry */ for (s = 0; s < NFace; s++) { if (e->bound_type[s] & NEUMANN) { SHORT bases[NbasFace(u)]; phgDofGetBasesOnFace(u, e, s, bases); for (i = 0; i < NbasFace(u); i++) { if (phgDofGetElementBoundaryType(u, e, bases[i] * Dim) & DIRICHLET) { /* Dirichlet bas on Neumann face, do nothing */ } else if (phgDofGetElementBoundaryType(u, e, bases[i] * Dim) & NEUMANN) { for (k = 0; k < Dim; k++) resu[bases[i]][k] += dt * Theta * phgQuadFaceDofDotBas(e, s, gn[k], DOF_PROJ_DOT, u, bases[i], QUAD_DEFAULT); } else { fprintf(stderr, "Warning: unkown bdry!"); } } /* end of base on face */ } /* end of face neu */ } /* end of all neumann face in element */ /* Global res */ phgSolverAddRHSEntries(solver, M * Dim, Iu[0], &resu[0][0]); } /* end element */ solver->rhs_updated = FALSE; phgDofFree(&lapu); phgDofFree(&gradp); phgDofFree(&f0); return; }
static void build_matrices(MAT *matA, MAT *matM, MAT *matC, MAT *S, FLOAT s, DOF *u_h, DOF *p_h) /* S is used to store s*diag(M(p_h))^(-1) */ { int N = u_h->type->nbas; /* number of basis functions in an element */ int M = p_h->type->nbas; int i, j; GRID *g = u_h->g; ELEMENT *e; FLOAT A[N][N], B[N][N], C[N][M]; INT I[N], Ip[M]; INT k, n0; VEC *V = phgMapCreateVec(S->rmap, 1); phgVecDisassemble(V); /* for phgVecAddEntry */ ForAllElements(g, e) { for (i = 0; i < N; i++) { I[i] = phgMapE2L(matA->rmap, 0, e, i); for (k = 0; k < M; k++) { /* \int \grad\psi_k\cdot\phi_i */ C[i][k] = phgQuadGradBasDotBas(e, p_h, k, u_h, i, QUAD_DEFAULT); } for (j = 0; j <= i; j++) { /* \int \phi_i\cdot\phi_j */ B[j][i] = B[i][j] = phgQuadBasDotBas(e, u_h, j, u_h, i, QUAD_DEFAULT); /* \int \curl\phi_i\cdot\curl\phi_j */ A[j][i] = A[i][j] = phgQuadCurlBasDotCurlBas(e, u_h, j, u_h, i, QUAD_DEFAULT); } } for (i = 0; i < M; i++) { Ip[i] = phgMapE2L(matC->cmap, 0, e, i); if (Ip[i] < 0) /* boundary entry */ continue; phgVecAddEntry(V, 0, Ip[i], phgQuadBasDotBas(e, p_h, i, p_h, i, QUAD_DEFAULT)); } /* loop on basis functions */ for (i = 0; i < N; i++) { if (phgDofDirichletBC(u_h, e, i, NULL, NULL, NULL, DOF_PROJ_CROSS)) continue; phgMatAddEntries(matA, 1, I + i, N, I, A[i]); phgMatAddEntries(matM, 1, I + i, N, I, B[i]); phgMatAddEntries(matC, 1, I + i, M, Ip, C[i]); } } phgVecAssemble(V); n0 = V->map->partition[V->map->rank]; for (k = 0; k < V->map->nlocal; k++) phgMatAddGlobalEntry(S, k + n0, k + n0, s / V->data[k]); phgVecDestroy(&V); phgMatAssemble(S); phgMatSetupDiagonal(S); phgMatAssemble(matA); phgMatAssemble(matM); phgMatAssemble(matC); }
void phgNSBuildPc(NSSolver *ns) { GRID *g = ns->g; SIMPLEX *e; FLOAT *dt = ns->dt; int i, j, q, s, k, l; FLOAT Theta = _nsp->Theta, nu = _nsp->nu, Thet1, nu0 = 0; DOF *tmp_u1 = phgDofNew(g, _nsp->utype, Dim, "tmp u1", func_u); int viscosity_type = ns->viscosity_type; LTYPE ltype = ns->ltype; #if STEADY_STATE assert(fabs(Theta - 1) < 1e-12); Thet1 = 0; Unused(Thet1); Unused(dt); #else Thet1 = 1 - Theta; #endif /* STEADY_STATE */ ForAllElements(g, e) { int M = ns->u[1]->type->nbas; /* num of bases of Velocity */ int N = ns->p[1]->type->nbas; /* num of bases of Pressure */ int order = 2 * DofTypeOrder(ns->p[1], e) + DofTypeOrder(ns->u[1], e) - 1; /* highest order term (u \nabla p, psi) */ FLOAT Ap[N][N], Fp[N][N], Qp[N][N], bufp[N], rhs1 = 1; FLOAT F[M*Dim][M*Dim], B[N][M*Dim], Bt[M*Dim][N]; INT Ip[N]; QUAD *quad; FLOAT vol, det; const FLOAT *w, *p, *vw, *gu, *vTe; quad = phgQuadGetQuad3D(order); vw = phgQuadGetDofValues(e, ns->wind, quad); /* value wind */ gu = phgQuadGetDofValues(e, ns->gradu[1], quad); /* grad u^{n+1} */ if (ns_params->noniter_temp) vTe = phgQuadGetDofValues(e, ns->T[1], quad); /* value temp */ else vTe = phgQuadGetDofValues(e, ns->T[0], quad); /* value temp */ vol = 0; Bzero(Ap); Bzero(Fp); Bzero(Qp); Bzero(F); Bzero(Bt); Bzero(B); Bzero(bufp); p = quad->points; w = quad->weights; for (q = 0; q < quad->npoints; q++) { phgGeomGetCurvedJacobianAtLambda(g, e, p, &det); vol = fabs(det / 6.); for (i = 0; i < N; i++) { const FLOAT *gi = phgQuadGetBasisValues(e, ns->p[1], i, quad) + q; /* phi_i */ const FLOAT *ggi = phgQuadGetBasisCurvedGradient(e, ns->p[1], i, quad, q); /* grad phi_i */ for (j = 0; j < N; j++) { const FLOAT *gj = phgQuadGetBasisValues(e, ns->p[1], j, quad) + q; /* phi_j */ const FLOAT *ggj = phgQuadGetBasisCurvedGradient(e, ns->p[1], j, quad, q); /* grad phi_i */ nu = get_effective_viscosity(gu, *vTe, 0, viscosity_type); if (i == 0 && j == 0) nu0 += nu; #if ICE_BENCH_TEST || \ ESIMINT_TEST || \ HEINO_TEST || \ TEST_CASE == ICE_EXACT || \ TEST_CASE == ICE_GREEN_LAND Unused(dt); /* Note: B Q^-1 Bt ~ Ap(nu=1), * Fp(nu varies) is very different to Ap */ Ap[i][j] += vol*(*w) * INNER_PRODUCT(ggj, ggi); # if USE_QP_ONLY //Qp[i][j] += vol*(*w) * LEN_SCALING * PRES_SCALING /(nu) * (*gj) * (*gi); Qp[i][j] += vol*(*w) * 1. /(EQU_SCALING * nu) * (*gj) * (*gi); /* if (i < NVert && j < NVert) { */ /* Qp[i][j] += vol*(*w) * LEN_SCALING * PRES_SCALING / (nu) * (*gj) * (*gi); */ /* } else if (i == NVert && j == NVert) { */ /* Qp[i][j] += vol*(*w) * LEN_SCALING * PRES_SCALING / (nu) * (*gj) * (*gi); */ /* } */ # else Qp[i][j] += vol*(*w) * (*gj) * (*gi); # endif Fp[i][j] += vol*(*w) * (EQU_SCALING * nu * INNER_PRODUCT(ggj, ggi) ); #elif STEADY_STATE Ap[i][j] += vol*(*w) * INNER_PRODUCT(ggj, ggi); Qp[i][j] += vol*(*w) * (*gj) * (*gi); Fp[i][j] += vol*(*w) * (nu * INNER_PRODUCT(ggj, ggi) * EQU_SCALING ); #elif TIME_DEP_NON Ap[i][j] += vol*(*w) * INNER_PRODUCT(ggj, ggi); Qp[i][j] += vol*(*w) * (*gj) * (*gi); Fp[i][j] += vol*(*w) * ((*gj) * (*gi) / dt[0] + Theta * (nu * INNER_PRODUCT(ggj, ggi) ) ); #else TIME_DEP_LINEAR_ENTRY; /* Unavailable */ #endif /* STEADY_STATE */ } } vw += Dim; gu += DDim; vTe++; w++; p += Dim+1; } /* Map: Element -> system */ for (i = 0; i < N; i++) Ip[i] = phgMapE2L(_pcd->matFp->cmap, 0, e, i); /* * PCD boundary setup I: * Automaticly decide inflow boundary condition using wind direction. * * NOT active. * */ if (FALSE && !_nsp->pin_node) { for (i = 0; i < N; i++) { BOOLEAN flag_inflow = FALSE; for (s = 0; s < NFace; s++) { SHORT bases[NbasFace(ns->p[1])]; FLOAT *coord, vw_[3]; const FLOAT *lam, *normal; if (!(e->bound_type[s] & BDRY_MASK)) //if (!(e->bound_type[s] & INFLOW)) continue; /* boundary face */ phgDofGetBasesOnFace(ns->p[1], e, s, bases); for (j = 0; j < NbasFace(ns->p[1]); j++) if (i == bases[j]) { normal = phgGeomGetFaceOutNormal(g, e, s); coord = phgDofGetElementCoordinates(ns->p[1], e, i); lam = phgGeomXYZ2Lambda(g, e, coord[0], coord[1], coord[2]); phgDofEval(tmp_u1, e, lam, vw_); if (INNER_PRODUCT(vw_, normal) > 1e-8) flag_inflow = TRUE; } } if (flag_inflow) { Bzero(bufp); bufp[i] = 1.0; phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp); //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp); phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1); } else { /* interior node Or Neumann */ phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, Ap[i]); phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]); //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]); } phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]); } } /* * PCD boundary setup II: * Enclose flow: use pinnode boundary. * * Qp is pinned, this is different to open flow. * * */ else if (_nsp->pin_node) { for (i = 0; i < N; i++) { if (phgDofDirichletBC(_pcd->pbc, e, i, NULL, bufp, NULL, DOF_PROJ_NONE)) { #if PIN_AT_ROOT if (g->rank != 0) phgError(1, "Pinned node only on rank 0!\n"); if (e->verts[i] != ns->pinned_node_id) phgError(1, "pinned node [%d] & [%d] doesn't coincide when build pc!\n", e->verts[i], ns->pinned_node_id); #else if (GlobalVertex(g, e->verts[i]) != ns->pinned_node_id) phgError(1, "pinned node [%d] & [%d] doesn't coincide when build pc!\n", e->verts[i], ns->pinned_node_id); #endif /* PIN_AT_ROOT */ phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp); phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1); } else { /* interior node Or Neumann */ phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, Ap[i]); phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]); phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]); } } } /* * PCD boundary setup III: * Open flow: there could be varies kinds of combination on seting up * boundary conditon, but Inflow:Robin & Outflow:scaled Dirich is * prefered. See Ref[2]. * * */ else { for (i = 0; i < N; i++) { /*****************/ /* Inflow */ /*****************/ #warning PCD B.C.: Step 2.1. build mat, all neumann, add dirich entries if (FALSE && phgDofDirichletBC(_pcd->dof_inflow, e, i, NULL, bufp, NULL, DOF_PROJ_NONE)) { phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp); phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1); } else if (FALSE && phgDofDirichletBC(_pcd->dof_outflow, e, i, NULL, bufp, NULL, DOF_PROJ_NONE) && !(phgDofGetElementBoundaryType(ns->p[1], e, i) & INFLOW) ) { ERROR_MSG("Fp, Qp"); nu = get_effective_viscosity(NULL, 0, 0, viscosity_type); phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp); bufp[i] *= EQU_SCALING * nu; phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp); phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1); //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp); } else if (FALSE && phgDofDirichletBC(_pcd->pbc, e, i, NULL, bufp, NULL, DOF_PROJ_NONE)) { phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp); phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp); phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1); } else if (FALSE) { /* interior node Or Neumann */ ERROR_MSG("Fp, Qp"); phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, Ap[i]); phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]); //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]); } /******************/ /* No bdry */ /******************/ //phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]); phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]); } } if (0) { /* Special term <[[p_i]], [[p_j]]> */ int face; nu0 /= quad->npoints; for (face = 0; face < NFace; face++) { FLOAT area = phgGeomGetFaceArea(g, e, face); //FLOAT value = {area, -area}; FLOAT values[2] = {vol * 1. /(EQU_SCALING * nu0), -vol * 1. /(EQU_SCALING * nu0)}; SIMPLEX *e_neigh; phgMatAddEntries(_pcd->matQp, 1, Ip+NVert, 1, Ip+NVert, values); if ((e_neigh = GetNeighbour(e, face)) != NULL) { INT Ip_neigh = phgMapE2L(_pcd->matFp->cmap, 0, e_neigh, NVert); phgMatAddEntries(_pcd->matQp, 1, Ip+NVert, 1, &Ip_neigh, values + 1); } } } } /* end element */