DOF * phgDofCopy_(DOF *src, DOF **dest_ptr, DOF_TYPE *newtype, const char *name, const char *srcfile, int srcline) /* returns a new DOF whose type is 'newtype', and whose values are evaluated * using 'src' (a copy of src). */ { GRID *g = src->g; SIMPLEX *e; FLOAT w = 1.0, *basvalues = NULL, *a, *d, *buffer = NULL; int i, k, dim, dim_new, count = 0, nvalues, nd; INT n; DOF *wgts = NULL; DOF *dest = (dest_ptr == NULL ? NULL : *dest_ptr); BYTE *flags0, *flags; char *auto_name; DOF_TYPE *oldtype; BOOLEAN basflag = FALSE; MagicCheck(DOF, dest) if (dest != NULL && newtype != NULL && dest->type != newtype) { phgDofFree(&dest); dest = NULL; } dim = DofDim(src); /* the name of dest */ if ((auto_name = (void *)name) == NULL) { #if 0 auto_name = phgAlloc(strlen(src->name) + 8 + 1); sprintf(auto_name, "copy of %s", src->name); #else auto_name = strdup(src->name); #endif } if (dest == NULL) { if (newtype == NULL) newtype = src->type; dim_new = (newtype == NULL ? DofTypeDim(src) : newtype->dim); assert(dim % dim_new == 0); dest = phgDofNew_(g, newtype, newtype == NULL ? src->hp : NULL, dim / dim_new, auto_name, DofNoAction, srcfile, srcline); if (dest_ptr != NULL) *dest_ptr = dest; } else { assert(dim == DofDim(dest)); phgFree(dest->name); dest->name = strdup(auto_name); dest->srcfile = srcfile; dest->srcline = srcline; phgFree(dest->cache_func); dest->cache_func = NULL; newtype = dest->type; if (!SpecialDofType(newtype)) memset(dest->data, 0, DofGetDataCount(dest) * sizeof(*dest->data)); } if (auto_name != name) phgFree(auto_name); phgDofClearCache(NULL, dest, NULL, NULL, FALSE); dest->DB_mask = src->DB_mask; if (src->DB_masks != NULL) { if (dest->DB_masks == NULL) dest->DB_masks = phgAlloc(dest->dim * sizeof(*dest->DB_masks)); memcpy(dest->DB_masks, src->DB_masks, dest->dim * sizeof(*dest->DB_masks)); } dest->invariant = src->invariant; if (SpecialDofType(newtype)) { assert(newtype == src->type); dest->userfunc = src->userfunc; dest->userfunc_lambda = src->userfunc_lambda; if (newtype == DOF_CONSTANT) memcpy(dest->data, src->data, dim * sizeof(*dest->data)); return dest; } if ((newtype != NULL && newtype == src->type) || (newtype == NULL && src->hp == dest->hp)) { /* simply duplicate the data */ size_t size = DofGetDataCount(dest); if (size > 0) memcpy(dest->data, src->data, sizeof(FLOAT) * size); dest->userfunc = src->userfunc; dest->userfunc_lambda = src->userfunc_lambda; return dest; } if (src->type == DOF_ANALYTIC) { if (src->userfunc_lambda != NULL) phgDofSetDataByLambdaFunction(dest, src->userfunc_lambda); else phgDofSetDataByFunction(dest, src->userfunc); return dest; } if (newtype != NULL && newtype->BasFuncs == NULL) phgError(1, "phgDofCopy: basis funcs for DOF type \"%s\" undefined.\n", newtype->name); dest->userfunc = src->userfunc; dest->userfunc_lambda = src->userfunc_lambda; oldtype = src->type; if (oldtype == NULL) oldtype = src->hp->info->types[src->hp->info->min_order]; if (!SpecialDofType(oldtype) && newtype != NULL && newtype->points != NULL && !DofIsHP(src)) { basflag = TRUE; count = oldtype->nbas * oldtype->dim; basvalues = phgAlloc(newtype->nbas * count * sizeof(*basvalues)); if (oldtype->invariant == TRUE) get_bas_funcs(src, dest, src->g->roots, basvalues); } if (newtype == NULL) newtype = dest->hp->info->types[dest->hp->max_order]; flags0 = phgCalloc((newtype->np_vert > 0 ? g->nvert : 0) + (newtype->np_edge > 0 ? g->nedge : 0) + (newtype->np_face > 0 ? g->nface : 0), sizeof(*flags0)); if (!SpecialDofType(oldtype) && oldtype->continuity < 0) { static DOF_TYPE DOF_WGTS = {DofCache, "Weights", NULL, NULL, NULL, NULL, NULL, NULL, NULL, phgDofInitFuncPoint, NULL, NULL, NULL, FE_None, FALSE, FALSE, -1, 0, 0, -1, 1, 0, 0, 0, 0 }; DOF_WGTS.np_vert = (newtype->np_vert > 0) ? 1 : 0; DOF_WGTS.np_edge = (newtype->np_edge > 0) ? 1 : 0; DOF_WGTS.np_face = (newtype->np_face > 0) ? 1 : 0; DOF_WGTS.nbas = DOF_WGTS.np_vert * NVert + DOF_WGTS.np_edge * NEdge + DOF_WGTS.np_face * NFace; if (DOF_WGTS.nbas > 0) { /* Other cases will be implemented later when really needed */ wgts = phgDofNew(g, &DOF_WGTS, 1, "weights", DofNoAction); phgDofSetDataByValue(wgts, 0.0); phgDofSetDataByValue(dest, 0.0); /* allocate buffer for storing weighted vertex/edge/face data */ if ((n = DofGetDataCount(dest) - DofGetElementDataCount(dest)) > 0) buffer = phgCalloc(n, sizeof(*buffer)); } } nvalues = dest->dim; cache_dof = src; bas_count = count; ForAllElements(g, e) { if (wgts != NULL) { #if 0 /* use element volume as weight */ w = phgGeomGetVolume(g, e); #else /* use 1 as weight */ w = 1.0; #endif } if (basflag && oldtype->invariant == FALSE) get_bas_funcs(src, dest, e, basvalues); bas = basvalues; flags = flags0; if (DofIsHP(dest)) newtype = dest->hp->info->types[dest->hp->elem_order[e->index]]; if (newtype->np_vert > 0) { nd = nvalues * newtype->np_vert; for (k = 0; k < NVert; k++) { if (flags[n = e->verts[k]] && wgts == NULL) { /* note: count==0 and bas==NULL for variable order DOF */ bas += count * newtype->np_vert; continue; } flags[n] = TRUE; a = DofVertexData(dest, n); newtype->InitFunc(dest, e, VERTEX, k, NULL, func, NULL, a, NULL); if (wgts != NULL) { d = buffer + (a - DofData(dest)); for (i = 0; i < nd; i++) *(d++) += *(a++) * w; *DofVertexData(wgts, n) += w; } } flags += g->nvert; } if (newtype->np_edge > 0) { nd = nvalues * newtype->np_edge; for (k = 0; k < NEdge; k++) { if (flags[n = e->edges[k]] && wgts == NULL) { /* note: count==0 and bas==NULL for variable order DOF */ bas += count * newtype->np_edge; continue; } flags[n] = TRUE; a = DofEdgeData(dest, n); newtype->InitFunc(dest, e, EDGE, k, NULL, func, NULL, a, NULL); if (wgts != NULL) { d = buffer + (a - DofData(dest)); if (DofIsHP(dest)) nd = dest->dim * (dest->hp->edge_index[n + 1] - dest->hp->edge_index[n]); for (i = 0; i < nd; i++) *(d++) += *(a++) * w; *DofEdgeData(wgts, n) += w; } } flags += g->nedge; } if (newtype->np_face > 0) { nd = nvalues * newtype->np_face; for (k = 0; k < NFace; k++) { if (flags[n = e->faces[k]] && wgts == NULL) { /* note: count==0 and bas==NULL for variable order DOF */ bas += count * newtype->np_face; continue; } flags[n] = TRUE; a = DofFaceData(dest, n); newtype->InitFunc(dest, e, FACE, k, NULL, func, NULL, a, NULL); if (wgts != NULL) { d = buffer + (a - DofData(dest)); if (DofIsHP(dest)) nd = dest->dim * (dest->hp->face_index[n + 1] - dest->hp->face_index[n]); for (i = 0; i < nd; i++) *(d++) += *(a++) * w; *DofFaceData(wgts, n) += w; } } } if (newtype->np_elem > 0) { a = DofElementData(dest, e->index); newtype->InitFunc(dest, e, ELEMENT, 0, NULL, func, NULL, a, NULL); } } phgFree(basvalues); phgFree(flags0); if (wgts == NULL) return dest; if ((n = DofGetDataCount(dest) - DofGetElementDataCount(dest)) > 0) { memcpy(DofData(dest), buffer, n * sizeof(*buffer)); phgFree(buffer); } if (DofIsHP(dest)) newtype = dest->hp->info->types[dest->hp->max_order]; a = DofData(dest); d = DofData(wgts); #if USE_MPI /* FIXME: directly interacting with the vector assembly code in solver.c is more efficient */ if (g->nprocs > 1) { SOLVER *solver = phgSolverCreate(SOLVER_BUILTIN, dest, NULL); INT K = 0, I; int j; if (newtype->np_vert > 0) { nd = dest->count_vert; for (n = 0; n < g->nvert; n++, d++) { if (g->types_vert[n] == UNREFERENCED) { K += nd; a += nd; continue; } for (j = 0; j < nd; j++, K++, a++) { I = phgSolverMapD2L(solver, 0, K); assert(I >= 0 && I < solver->mat->rmap->localsize); phgSolverAddRHSEntry(solver, I, *a); phgSolverAddMatrixEntry(solver, I, I, *d); } } } if (newtype->np_edge > 0) { nd = nvalues * newtype->np_edge; for (n = 0; n < g->nedge; n++, d++) { if (DofIsHP(dest)) nd = dest->dim * (dest->hp->edge_index[n + 1] - dest->hp->edge_index[n]); if (g->types_edge[n] == UNREFERENCED) { K += nd; a += nd; continue; } for (j = 0; j < nd; j++, K++, a++) { I = phgSolverMapD2L(solver, 0, K); assert(I >= 0 && I < solver->mat->rmap->localsize); phgSolverAddRHSEntry(solver, I, *a); phgSolverAddMatrixEntry(solver, I, I, *d); } } } if (newtype->np_face > 0) { nd = nvalues * newtype->np_face; for (n = 0; n < g->nface; n++, d++) { if (DofIsHP(dest)) nd = dest->dim * (dest->hp->face_index[n + 1] - dest->hp->face_index[n]); if (g->types_face[n] == UNREFERENCED) { K += nd; a += nd; continue; } for (j = 0; j < nd; j++, K++, a++) { I = phgSolverMapD2L(solver, 0, K); assert(I >= 0 && I < solver->mat->rmap->localsize); phgSolverAddRHSEntry(solver, I, *a); phgSolverAddMatrixEntry(solver, I, I, *d); } } } if (newtype->np_elem > 0) { nd = nvalues * newtype->np_elem; for (n = 0; n < g->nelem; n++) { if (DofIsHP(dest)) nd = dest->dim * (dest->hp->elem_index[n + 1] - dest->hp->elem_index[n]); if (g->types_elem[n] == UNREFERENCED) { K += nd; a += nd; continue; } for (j = 0; j < nd; j++, K++, a++) { I = phgSolverMapD2L(solver, 0, K); assert(I >= 0 && I < solver->mat->rmap->localsize); phgSolverAddRHSEntry(solver, I, *a); phgSolverAddMatrixEntry(solver, I, I, 1.0); } } } phgSolverSolve(solver, TRUE, dest, NULL); phgSolverDestroy(&solver); phgDofFree(&wgts); return dest; } #endif if (newtype->np_vert > 0) { k = nvalues * newtype->np_vert; for (n = 0; n < g->nvert; n++) { if ((w = *(d++)) == 0.0) { a += k; } else if (k == 1) { *(a++) /= w; } else { w = 1.0 / w; for (i = 0; i < k; i++) *(a++) *= w; } } } if (newtype->np_edge > 0) { k = nvalues * newtype->np_edge; for (n = 0; n < g->nedge; n++) { if (DofIsHP(dest)) k = dest->dim * (dest->hp->edge_index[n + 1] - dest->hp->edge_index[n]); if ((w = *(d++)) == 0.0) { a += k; } else if (k == 1) { *(a++) /= w; } else { w = 1.0 / w; for (i = 0; i < k; i++) *(a++) *= w; } } } if (newtype->np_face > 0) { k = nvalues * newtype->np_face; for (n = 0; n < g->nface; n++) { if (DofIsHP(dest)) k = dest->dim * (dest->hp->face_index[n + 1] - dest->hp->face_index[n]); if ((w = *(d++)) == 0.0) { a += k; } else if (k == 1) { *(a++) /= w; } else { w = 1.0 / w; for (i = 0; i < k; i++) *(a++) *= w; } } } phgDofFree(&wgts); return dest; }
static void getMoveDirection(MovingMesh *mmesh) { GRID *g = _m->g; SIMPLEX *e; DOF *logical_move = _m->logical_move; DOF *logical_node = _m->logical_node; DOF *move = _m->move; DOF *logical_node_new, *grad_logical_node; INT i, j, k, l0, l1; SOLVER *solver; DOF *mass_lumping; FLOAT *v_move, *v_mass; FLOAT max_vol = -1e10, min_vol = 1e10; /* Get new monitor */ _m->get_monitor(mmesh); DOF_SCALE(_m->monitor); logical_node_new = phgDofCopy(logical_node, NULL, NULL, "new logical node coordinates"); grad_logical_node = phgDofGradient(logical_node, NULL, NULL, "grad logical node coordinates"); /* if (logical_node_new->DB_mask == UNDEFINED) */ /* phgError(1, "g->types_vert and DB_mask is VAILD in moving mesh method!!!"); */ /* Create move solver */ phgOptionsPush(); #if 0 phgOptionsSetOptions("-default_solver mm_ml " "-mm_ml_solver gmres " "-mm_ml_pc mm_ml " "-mm_ml_sweep0 6 " "-solver_maxit 100 " "-solver_rtol 1e-12 "); #else phgOptionsSetOptions("-default_solver hypre " "-hypre_solver gmres " "-hypre_pc boomeramg " "-solver_maxit 100 " "-solver_rtol 1e-12 "); #endif phgOptionsSetOptions(_mp->move_opts); setMoveConstrain(mmesh, logical_node_new); solver = phgSolverCreate(SOLVER_DEFAULT, logical_node_new, NULL); solver->mat->mv_data = phgAlloc(sizeof(*solver->mat->mv_data)); solver->mat->mv_data[0] = (void *) mmesh; phgOptionsPop(); /* build mat */ ForAllElements(g, e) { int order = 1, q; int N = DofGetNBas(logical_node, e); /* number of bases in the element */ FLOAT A[N][Dim][N][Dim], rhs[N][Dim]; INT I[N][Dim]; QUAD *quad; const FLOAT *w, *lambda; FLOAT vol, mon; vol = phgGeomGetVolume(g, e); max_vol = MAX(max_vol, vol); min_vol = MIN(min_vol, vol); for (i = 0; i < N; i++) for (k = 0; k < Dim; k++) I[i][k] = phgMapE2L(solver->rhs->map, 0, e, i * Dim + k); bzero(A, sizeof(A)); bzero(rhs, sizeof(rhs)); quad = phgQuadGetQuad3D(order); mon = *DofElementData(_m->monitor, e->index); //printf("mon:%e\n", mon); lambda = quad->points; w = quad->weights; assert(quad->npoints == 1); for (q = 0; q < quad->npoints; q++) { for (i = 0; i < N; i++) { for (j = 0; j < N; j++) { const FLOAT *ggi = phgQuadGetBasisGradient(e, _m->move, i, quad) + q*Dim; /* grad phi_i */ const FLOAT *ggj = phgQuadGetBasisGradient(e, _m->move, j, quad) + q*Dim; /* grad phi_j */ FLOAT a = vol*(*w) * mon * INNER_PRODUCT(ggj, ggi); for (k = 0; k < Dim; k++) { A[i][k][j][k] += a; } } } w++; lambda += Dim + 1; } /* end quad point */ for (i = 0; i < N; i++) { INT vert0 = e->verts[i]; if (_mp->fix_bdry || MM_FIXED(g->types_vert[vert0])) { /* no move */ bzero(A[i], sizeof(A[i])); for (k = 0; k < Dim; k++) { A[i][k][i][k] = 1.; rhs[i][k] = logical_node->data[vert0 * Dim + k]; } } else if (g->types_vert[vert0] & MM_CONSTR) { VERT_CONSTRAIN *vert_constr = _m->vert_constr + _m->vert_constr_index[vert0]; FLOAT *trans = vert_constr->Trans; FLOAT *bb = vert_constr->bb; assert(_m->vert_constr_index[vert0] >= 0); assert(vert0 == vert_constr->index); trans_left_ (&A[i][0][0][0], Dim*N, Dim*N, trans); trans_right_(&A[0][0][i][0], Dim*N, Dim*N, trans); if (vert_constr->n_constrain == 2) { bzero(A[i][0], sizeof(A[i][0])); bzero(A[i][1], sizeof(A[i][1])); A[i][0][i][0] = 1.; A[i][1][i][1] = 1.; rhs[i][0] = bb[0]; rhs[i][1] = bb[1]; } else if (vert_constr->n_constrain == 1) { bzero(A[i][0], sizeof(A[i][0])); A[i][0][i][0] = 1.; rhs[i][0] = bb[0]; } else { abort(); } } } for (i = 0; i < N; i++) phgSolverAddMatrixEntries(solver, Dim, &I[i][0], N * Dim, I[0], &A[i][0][0][0]); phgSolverAddRHSEntries(solver, N * Dim, &I[0][0], &rhs[0][0]); }
/*************************************************************************** * 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; }
int main(int argc, char *argv[]) { GRID *g; DOF *u, *v, *w, *u0, *v0, *w0; SOLVER *solver; INT i, j; char *fn = "cube.dat"; phgVerbosity = 0; phgInit(&argc, &argv); /*phgPause(0);*/ g = phgNewGrid(-1); if (!phgImport(g, fn, FALSE)) phgError(1, "can't read file \"%s\".\n", fn); for (i = 0; i < 4; i++) phgRefineAllElements(g, 1); phgBalanceGrid(g, 1.0, 0, NULL, 0.); for (i = 0; i < 8; i++) { phgRefineRandomElements(g, "25%"); phgBalanceGrid(g, 1.2, -1, NULL, 0.); } phgCheckConformity(g); u = phgDofNew(g, DOF_P1, 1, "u", DofNoAction); v = phgDofNew(g, DOF_ND1, 1, "v", DofNoAction); w = phgDofNew(g, DOF_P1, 1, "w", DofNoAction); phgDofSetDataByValue(u, 0.0); phgDofSetDataByValue(v, 0.0); phgDofSetDataByValue(w, 0.0); u0 = phgDofNew(g, u->type, u->dim, "u0", DofNoAction); v0 = phgDofNew(g, v->type, v->dim, "v0", DofNoAction); w0 = phgDofNew(g, w->type, w->dim, "w0", DofNoAction); solver = phgSolverCreate(SOLVER_DEFAULT, u, v, w, NULL); for (i = 0; i < DofGetDataCount(u); i++) { j = phgSolverMapD2L(solver, 0, i); u0->data[i] = (FLOAT)(phgSolverMapL2G(solver, j) % 3); phgSolverAddMatrixEntry(solver, j, j, 1.0); phgSolverAddRHSEntry(solver, j, u0->data[i]); } for (i = 0; i < DofGetDataCount(v); i++) { j = phgSolverMapD2L(solver, 1, i); v0->data[i] = (FLOAT)(phgSolverMapL2G(solver, j) % 5); phgSolverAddMatrixEntry(solver, j, j, 1.0); phgSolverAddRHSEntry(solver, j, v0->data[i]); } for (i = 0; i < DofGetDataCount(w); i++) { j = phgSolverMapD2L(solver, 2, i); w0->data[i] = (FLOAT)(phgSolverMapL2G(solver, j) % 7); phgSolverAddMatrixEntry(solver, j, j, 1.0); phgSolverAddRHSEntry(solver, j, w0->data[i]); } phgSolverSolve(solver, TRUE, u, v, w, NULL); phgSolverDestroy(&solver); phgDofAXPY(-1.0, u0, &u); phgDofAXPY(-1.0, v0, &v); phgDofAXPY(-1.0, w0, &w); phgPrintf("Error: u = %lg, v = %lg, w = %lg\n", (double)phgDofNormInftyVec(u), (double)phgDofNormInftyVec(v), (double)phgDofNormInftyVec(w)); phgDofFree(&u); phgDofFree(&v); phgDofFree(&w); phgDofFree(&u0); phgDofFree(&v0); phgDofFree(&w0); phgFreeGrid(&g); phgFinalize(); return 0; }