void FeProblem::create(SparseMatrix* mat) { assert(mat != NULL); if (is_up_to_date()) { verbose("Reusing matrix sparse structure."); mat->zero(); return; } // spaces have changed: create the matrix from scratch mat->free(); int ndof = get_num_dofs(); mat->prealloc(ndof); AUTOLA_CL(AsmList, al, wf->neq); AUTOLA_OR(Mesh*, meshes, wf->neq); bool **blocks = wf->get_blocks(); // init multi-mesh traversal for (int i = 0; i < wf->neq; i++) meshes[i] = spaces[i]->get_mesh(); Traverse trav; trav.begin(wf->neq, meshes); // loop through all elements Element **e; while ((e = trav.get_next_state(NULL, NULL)) != NULL) { // obtain assembly lists for the element at all spaces for (int i = 0; i < wf->neq; i++) // TODO: do not get the assembly list again if the element was not changed if (e[i] != NULL) spaces[i]->get_element_assembly_list(e[i], al + i); // go through all equation-blocks of the local stiffness matrix for (int m = 0; m < wf->neq; m++) for (int n = 0; n < wf->neq; n++) if (blocks[m][n] && e[m] != NULL && e[n] != NULL) { AsmList *am = al + m; AsmList *an = al + n; // pretend assembling of the element stiffness matrix // register nonzero elements for (int i = 0; i < am->cnt; i++) if (am->dof[i] >= 0) for (int j = 0; j < an->cnt; j++) if (an->dof[j] >= 0) mat->pre_add_ij(am->dof[i], an->dof[j]); } } trav.finish(); delete [] blocks; mat->alloc(); // save space seq numbers and weakform seq number, so we can detect their changes for (int i = 0; i < wf->neq; i++) sp_seq[i] = spaces[i]->get_seq(); wf_seq = wf->get_seq(); struct_changed = true; have_matrix = true; }
void FeProblem::assemble(_Vector *rhs, _Matrix *jac, Tuple<Solution*> u_ext) { // Sanity checks. if (!have_spaces) error("You have to call set_spaces() before calling assemble()."); for (int i=0; i<this->wf->neq; i++) { if (this->spaces[i] == NULL) error("A space is NULL in assemble()."); } int k, m, n, marker; AUTOLA_CL(AsmList, al, wf->neq); AsmList *am, *an; bool bnd[4]; AUTOLA_OR(bool, nat, wf->neq); AUTOLA_OR(bool, isempty, wf->neq); EdgePos ep[4]; reset_warn_order(); // create slave pss's for test functions, init quadrature points AUTOLA_OR(PrecalcShapeset*, spss, wf->neq); PrecalcShapeset *fu, *fv; AUTOLA_CL(RefMap, refmap, wf->neq); for (int i = 0; i < wf->neq; i++) { spss[i] = new PrecalcShapeset(pss[i]); pss [i]->set_quad_2d(&g_quad_2d_std); spss[i]->set_quad_2d(&g_quad_2d_std); refmap[i].set_quad_2d(&g_quad_2d_std); } // initialize buffer buffer = NULL; mat_size = 0; get_matrix_buffer(9); // obtain a list of assembling stages std::vector<WeakForm::Stage> stages; wf->get_stages(spaces, NULL, stages, jac == NULL); // Loop through all assembling stages -- the purpose of this is increased performance // in multi-mesh calculations, where, e.g., only the right hand side uses two meshes. // In such a case, the matrix forms are assembled over one mesh, and only the rhs // traverses through the union mesh. On the other hand, if you don't use multi-mesh // at all, there will always be only one stage in which all forms are assembled as usual. Traverse trav; for (unsigned ss = 0; ss < stages.size(); ss++) { WeakForm::Stage* s = &stages[ss]; for (unsigned i = 0; i < s->idx.size(); i++) s->fns[i] = pss[s->idx[i]]; for (unsigned i = 0; i < s->ext.size(); i++) s->ext[i]->set_quad_2d(&g_quad_2d_std); trav.begin(s->meshes.size(), &(s->meshes.front()), &(s->fns.front())); // assemble one stage Element** e; while ((e = trav.get_next_state(bnd, ep)) != NULL) { // find a non-NULL e[i] Element* e0; for (unsigned int i = 0; i < s->idx.size(); i++) if ((e0 = e[i]) != NULL) break; if (e0 == NULL) continue; // set maximum integration order for use in integrals, see limit_order() update_limit_table(e0->get_mode()); // obtain assembly lists for the element at all spaces, set appropriate mode for each pss memset(isempty, 0, sizeof(bool) * wf->neq); for (unsigned int i = 0; i < s->idx.size(); i++) { int j = s->idx[i]; if (e[i] == NULL) { isempty[j] = true; continue; } spaces[j]->get_element_assembly_list(e[i], al+j); spss[j]->set_active_element(e[i]); spss[j]->set_master_transform(); refmap[j].set_active_element(e[i]); refmap[j].force_transform(pss[j]->get_transform(), pss[j]->get_ctm()); u_ext[j]->set_active_element(e[i]); u_ext[j]->force_transform(pss[j]->get_transform(), pss[j]->get_ctm()); } marker = e0->marker; init_cache(); //// assemble volume matrix forms ////////////////////////////////////// if (jac != NULL) { for (unsigned ww = 0; ww < s->mfvol.size(); ww++) { WeakForm::MatrixFormVol* mfv = s->mfvol[ww]; if (isempty[mfv->i] || isempty[mfv->j]) continue; if (mfv->area != H2D_ANY && !wf->is_in_area(marker, mfv->area)) continue; m = mfv->i; fv = spss[m]; am = &al[m]; n = mfv->j; fu = pss[n]; an = &al[n]; bool tra = (m != n) && (mfv->sym != 0); bool sym = (m == n) && (mfv->sym == 1); // assemble the local stiffness matrix for the form mfv scalar bi, **mat = get_matrix_buffer(std::max(am->cnt, an->cnt)); for (int i = 0; i < am->cnt; i++) { if (!tra && (k = am->dof[i]) < 0) continue; fv->set_active_shape(am->idx[i]); if (!sym) // unsymmetric block { for (int j = 0; j < an->cnt; j++) { fu->set_active_shape(an->idx[j]); bi = eval_form(mfv, u_ext, fu, fv, refmap+n, refmap+m) * an->coef[j] * am->coef[i]; if (an->dof[j] >= 0) mat[i][j] = bi; } } else // symmetric block { for (int j = 0; j < an->cnt; j++) { if (j < i && an->dof[j] >= 0) continue; fu->set_active_shape(an->idx[j]); bi = eval_form(mfv, u_ext, fu, fv, refmap+n, refmap+m) * an->coef[j] * am->coef[i]; if (an->dof[j] >= 0) mat[i][j] = mat[j][i] = bi; } } } // insert the local stiffness matrix into the global one jac->add(am->cnt, an->cnt, mat, am->dof, an->dof); // insert also the off-diagonal (anti-)symmetric block, if required if (tra) { if (mfv->sym < 0) chsgn(mat, am->cnt, an->cnt); transpose(mat, am->cnt, an->cnt); jac->add(am->cnt, an->cnt, mat, am->dof, an->dof); } } } //// assemble volume linear forms //////////////////////////////////////// if (rhs != NULL) { for (unsigned int ww = 0; ww < s->vfvol.size(); ww++) { WeakForm::VectorFormVol* vfv = s->vfvol[ww]; if (isempty[vfv->i]) continue; if (vfv->area != H2D_ANY && !wf->is_in_area(marker, vfv->area)) continue; m = vfv->i; fv = spss[m]; am = &al[m]; for (int i = 0; i < am->cnt; i++) { if (am->dof[i] < 0) continue; fv->set_active_shape(am->idx[i]); rhs->add(am->dof[i], eval_form(vfv, u_ext, fv, refmap + m) * am->coef[i]); } } } // assemble surface integrals now: loop through boundary edges of the element for (unsigned int edge = 0; edge < e0->nvert; edge++) { if (!bnd[edge]) continue; marker = ep[edge].marker; // obtain the list of shape functions which are nonzero on this edge for (unsigned int i = 0; i < s->idx.size(); i++) { if (e[i] == NULL) continue; int j = s->idx[i]; if ((nat[j] = (spaces[j]->bc_type_callback(marker) == BC_NATURAL))) spaces[j]->get_edge_assembly_list(e[i], edge, al + j); } // assemble surface matrix forms /////////////////////////////////// if (jac != NULL) { for (unsigned int ww = 0; ww < s->mfsurf.size(); ww++) { WeakForm::MatrixFormSurf* mfs = s->mfsurf[ww]; if (isempty[mfs->i] || isempty[mfs->j]) continue; if (mfs->area != H2D_ANY && !wf->is_in_area(marker, mfs->area)) continue; m = mfs->i; fv = spss[m]; am = &al[m]; n = mfs->j; fu = pss[n]; an = &al[n]; if (!nat[m] || !nat[n]) continue; ep[edge].base = trav.get_base(); ep[edge].space_v = spaces[m]; ep[edge].space_u = spaces[n]; scalar bi, **mat = get_matrix_buffer(std::max(am->cnt, an->cnt)); for (int i = 0; i < am->cnt; i++) { if ((k = am->dof[i]) < 0) continue; fv->set_active_shape(am->idx[i]); for (int j = 0; j < an->cnt; j++) { fu->set_active_shape(an->idx[j]); bi = eval_form(mfs, u_ext, fu, fv, refmap+n, refmap+m, ep+edge) * an->coef[j] * am->coef[i]; if (an->dof[j] >= 0) mat[i][j] = bi; } } jac->add(am->cnt, an->cnt, mat, am->dof, an->dof); } } // assemble surface linear forms ///////////////////////////////////// if (rhs != NULL) { for (unsigned ww = 0; ww < s->vfsurf.size(); ww++) { WeakForm::VectorFormSurf* vfs = s->vfsurf[ww]; if (isempty[vfs->i]) continue; if (vfs->area != H2D_ANY && !wf->is_in_area(marker, vfs->area)) continue; m = vfs->i; fv = spss[m]; am = &al[m]; if (!nat[m]) continue; ep[edge].base = trav.get_base(); ep[edge].space_v = spaces[m]; for (int i = 0; i < am->cnt; i++) { if (am->dof[i] < 0) continue; fv->set_active_shape(am->idx[i]); rhs->add(am->dof[i], eval_form(vfs, u_ext, fv, refmap+m, ep+edge) * am->coef[i]); } } } } delete_cache(); } trav.finish(); } for (int i = 0; i < wf->neq; i++) delete spss[i]; delete [] buffer; buffer = NULL; mat_size = 0; }
void L2OrthoHP::get_optimal_refinement(Element* e, int order, Solution* rsln, int& split, int p[4], bool h_only, bool iso_only, double conv_exp, int max_order) { int i, j, k, n = 0; const int maxcand = 300; order = std::max(get_h_order(order), get_v_order(order)); bool tri = e->is_triangle(); // calculate maximal order of elements // linear elements = 9 // curvilinear elements = depends on iro_cache (how curved they are) if (max_order == -1) max_order = (20 - e->iro_cache)/2 - 2; // default else max_order = std::min( max_order, (20 - e->iro_cache)/2 - 2); // user specified AUTOLA_CL(Cand, cand, maxcand); #define make_p_cand(q) { \ assert(n < maxcand); \ cand[n].split = -1; \ cand[n].p[1] = cand[n].p[2] = cand[n].p[3] = 0; \ cand[n++].p[0] = (q); } #define make_hp_cand(q0, q1, q2, q3) { \ assert(n < maxcand); \ cand[n].split = 0; \ cand[n].p[0] = (q0); \ cand[n].p[1] = (q1); \ cand[n].p[2] = (q2); \ cand[n++].p[3] = (q3); } #define make_ani_cand(q0, q1, iso) { \ assert(n < maxcand); \ cand[n].split = iso; \ cand[n].p[2] = cand[n].p[3] = 0; \ cand[n].p[0] = (q0); \ cand[n++].p[1] = (q1); }\ if (h_only) { make_p_cand(order); make_hp_cand(order, order, order, order); make_ani_cand(order, order, 1); make_ani_cand(order, order, 2); } else { // prepare p-candidates int p0, p1 = std::min(max_order, order+1); for (p0 = order; p0 <= p1; p0++) make_p_cand(p0); // prepare hp-candidates p0 = (order+1) / 2; p1 = std::min(p0 + 3, order); int q0, q1, q2, q3; for (q0 = p0; q0 <= p1; q0++) for (q1 = p0; q1 <= p1; q1++) for (q2 = p0; q2 <= p1; q2++) for (q3 = p0; q3 <= p1; q3++) make_hp_cand(q0, q1, q2, q3); // prepare anisotropic candidates // only for quadrilaterals // too distorted (curved) elements cannot have aniso refinement (produces even worse elements) if ((!tri) && (e->iro_cache < 8) && !iso_only) { p0 = 2 * (order+1) / 3; int p_max = std::min(max_order, order+1); p1 = std::min(p0 + 3, p_max); for (q0 = p0; q0 <= p1; q0++) for (q1 = p0; q1 <= p1; q1++) { if ((q0 < order+1) || (q1 < order+1)) { make_ani_cand(q0, q1, 1); make_ani_cand(q0, q1, 2); } } } } // calculate (partial) projection errors double herr[8][11], perr[11]; calc_projection_errors(e, order, rsln, herr, perr); // evaluate candidates (sum partial projection errors, calculate dofs) double avg = 0.0; double dev = 0.0; for (i = k = 0; i < n; i++) { Cand* c = cand + i; if (c->split == 0) { c->error = 0.0; c->dofs = tri ? 6 : 9; for (j = 0; j < 4; j++) { int o = c->p[j]; c->error += herr[j][o] * 0.25; // spravny vypocet chyby if (tri) { c->dofs += (o-2)*(o-1)/2; if (j < 3) c->dofs += std::min(o, c->p[3])-1 + 2*(o-1); } else { c->dofs += sqr(o)-1; c->dofs += 2 * std::min(o, c->p[j>0 ? j-1 : 3]) - 1; } } } else if (c->split == 1 || c->split == 2) // aniso splits { c->dofs = 6 /* vertex */ + 3*(c->p[0] - 1 + c->p[1] - 1); // edge fns c->dofs += std::min(c->p[0], c->p[1]) - 1; // common edge c->dofs += sqr(c->p[0] - 1) + sqr(c->p[1] - 1); // bubbles for (c->error = 0.0, j = 0; j < 2; j++) c->error += herr[(c->split == 1) ? j+4 : j+6][c->p[j]] * 0.5; // spravny vypocet chyby } else { int o = c->p[0]; c->error = perr[o]; c->dofs = tri ? (o+1)*(o+2)/2 : sqr(o+1); } c->error = sqrt(c->error); //verbose("Cand #%d: Orders %d %d %d %d, Error %g, Dofs %d", i, c->p[0],c->p[1],c->p[2],c->p[3],c->error, c->dofs); if (!i || c->error <= cand[0].error) { avg += log(c->error); dev += sqr(log(c->error)); k++; } } avg /= k; // mean dev /= k; // second moment dev = sqrt(dev - sqr(avg)); // deviation is square root of variance // select an above-average candidate with the steepest error decrease int imax = 0; double score, maxscore = 0.0; for (i = 1; i < n; i++) { if ((log(cand[i].error) < avg + dev) && (cand[i].dofs > cand[0].dofs)) { score = (log(cand[0].error) - log(cand[i].error)) / //(pow(cand[i].dofs, conv_exp) - pow(cand[0].dofs, conv_exp)); pow(cand[i].dofs - cand[0].dofs, conv_exp); if (score > maxscore) { maxscore = score; imax = i; } } } // return result split = cand[imax].split; memcpy(p, cand[imax].p, 4*sizeof(int)); //verbose("Selected Candidate #%d: Orders %d %d %d %d\n", imax, p[0],p[1],p[2],p[3]); }