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 DiscreteProblem::assemble_matrix_and_rhs(bool rhsonly) { int i, j, k, l, m, n; bool bnd[4], nat[neq]; EdgePos ep[4]; if (!ndofs) return; warned_order = false; if (!rhsonly) { alloc_matrix_values(); if (!quiet) { verbose("Assembling stiffness matrix..."); begin_time(); } } else { memset(RHS, 0, sizeof(scalar) * ndofs); if (!quiet) { verbose("Assembling RHS..."); begin_time(); } } // create slave pss's for test functions, init quadrature points PrecalcShapeset* spss[neq]; PrecalcShapeset *fu, *fv; for (i = 0; i < 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); } // initialize buffer buffer = NULL; mat_size = 0; get_matrix_buffer(9); // initialize assembly lists, refmap AsmList al[neq], *am, *an; RefMap* refmap = new RefMap[neq]; for (i = 0; i < neq; i++) refmap[i].set_quad_2d(&g_quad_2d_std); for (i = 0; i < num_extern; i++) extern_fns[i]->set_quad_2d(&g_quad_2d_std); // init multi-mesh traversal int nm = neq + num_extern; Mesh* meshes[nm]; Transformable* fn[nm]; for (i = 0; i < neq; i++) meshes[i] = spaces[i]->get_mesh(); memcpy(fn, pss, neq * sizeof(Transformable*)); for (i = 0; i < num_extern; i++) { meshes[neq+i] = extern_fns[i]->get_mesh(); fn[neq+i] = extern_fns[i]; } // todo: kdyz maji nektere slozky stejnou sit, at sdili i refmapy // - ale to bysme potrebovali slave RefMap // loop through all elements Element** e; Traverse trav; trav.begin(nm, meshes, fn); while ((e = trav.get_next_state(bnd, ep)) != NULL) { // set maximum integration order for use in integrals, see limit_order() update_limit_table(e[0]->get_mode()); // obtain assembly lists for the element at all spaces, set appropriate mode for each pss for (i = 0; i < neq; i++) { spaces[i]->get_element_assembly_list(e[i], al + i); // todo: neziskavat znova, pokud se element nezmenil if (is_equi) for (j = 0; j < al[i].cnt; j++) if (al[i].dof[j] >= 0) al[i].coef[j] /= equi[al[i].dof[j]]; spss[i]->set_active_element(e[i]); spss[i]->set_master_transform(); refmap[i].set_active_element(e[i]); refmap[i].force_transform(pss[i]->get_transform(), pss[i]->get_ctm()); } // go through all equation-blocks of the element stiffness matrix, assemble volume integrals for (m = 0, am = al; m < neq; m++, am++) { fv = spss[m]; if (!rhsonly) { for (n = 0, an = al; n < neq; n++, an++) { fu = pss[n]; BiForm* bf = biform[m] + n; if (!bf->sym && !bf->unsym) continue; if (bf->unsym == BF_SYM || bf->unsym == BF_ANTISYM) continue; bool tra = (biform[n][m].unsym == BF_SYM || biform[n][m].unsym == BF_ANTISYM); // assemble the (m,n)-block of the stiffness matrix scalar sy, un, **mat = get_matrix_buffer(std::max(am->cnt, an->cnt)); for (i = 0; i < am->cnt; i++) { if (!tra && (k = am->dof[i]) < 0) continue; fv->set_active_shape(am->idx[i]); // unsymmetric block if (!bf->sym) { for (j = 0; j < an->cnt; j++) { fu->set_active_shape(an->idx[j]); un = bf->unsym(fu, fv, refmap+n, refmap+m) * an->coef[j] * am->coef[i]; if (an->dof[j] < 0) Dir[k] -= un; else mat[i][j] = un; } } // symmetric block else { for (j = 0; j < an->cnt; j++) { scalar coef = an->coef[j] * am->coef[i]; if (an->dof[j] < 0) { fu->set_active_shape(an->idx[j]); un = bf->unsym ? bf->unsym(fu, fv, refmap+n, refmap+m) * coef : 0.0; sy = bf->sym(fu, fv, refmap+n, refmap+m) * coef; Dir[k] -= (un + sy); } if (j >= i) { fu->set_active_shape(an->idx[j]); un = bf->unsym ? bf->unsym(fu, fv, refmap+n, refmap+m) * coef : 0.0; mat[j][i] = sy = bf->sym (fu, fv, refmap+n, refmap+m) * coef; mat[i][j] = (un + sy); } else if (bf->unsym) { fu->set_active_shape(an->idx[j]); mat[i][j] += bf->unsym(fu, fv, refmap+n, refmap+m) * coef; } } } } // insert the local stiffness matrix into the global one insert_matrix(mat, am->dof, an->dof, am->cnt, an->cnt); // insert also the off-diagonal (anti-)symmetric block, if required if (tra) { if (biform[n][m].unsym == BF_ANTISYM) chsgn(mat, am->cnt, an->cnt); transpose(mat, am->cnt, an->cnt); insert_matrix(mat, an->dof, am->dof, an->cnt, am->cnt); // we also need to take care of the RHS... for (j = 0; j < am->cnt; j++) if (am->dof[j] < 0) for (i = 0; i < an->cnt; i++) if (an->dof[i] >= 0) Dir[an->dof[i]] -= mat[i][j]; } } } // assemble rhs (linear form) if (!liform[m].lf) continue; for (i = 0; i < am->cnt; i++) { if (am->dof[i] < 0) continue; fv->set_active_shape(am->idx[i]); RHS[am->dof[i]] += liform[m].lf(fv, refmap+m) * am->coef[i]; } } // assemble surface integrals now: loop through boundary edges of the element if (rhsonly) continue; // fixme for (int edge = 0; edge < e[0]->nvert; edge++) { if (!bnd[edge]) continue; // obtain the list of shape functions which are nonzero on this edge for (i = 0; i < neq; i++) if ((nat[i] = (spaces[i]->bc_type_callback(ep[edge].marker) == BC_NATURAL))) spaces[i]->get_edge_assembly_list(e[i], edge, al + i); // loop through the equation-blocks for (m = 0, am = al; m < neq; m++, am++) { if (!nat[m]) continue; fv = spss[m]; ep[edge].base = trav.get_base(); ep[edge].space_v = spaces[m]; for (n = 0, an = al; n < neq; n++, an++) { if (!nat[n]) continue; BiForm* bf = biform[m] + n; if (!bf->surf) continue; fu = pss[n]; ep[edge].space_u = spaces[n]; // assemble the surface part of the bilinear form scalar bi, **mat = get_matrix_buffer(std::max(am->cnt, an->cnt)); for (i = 0; i < am->cnt; i++) { if ((k = am->dof[i]) < 0) continue; fv->set_active_shape(am->idx[i]); for (j = 0; j < an->cnt; j++) { fu->set_active_shape(an->idx[j]); bi = bf->surf(fu, fv, refmap+n, refmap+m, ep+edge) * an->coef[j] * am->coef[i]; if (an->dof[j] >= 0) mat[i][j] = bi; else Dir[k] -= bi; } } insert_matrix(mat, am->dof, an->dof, am->cnt, an->cnt); } // assemble the surface part of the linear form if (!liform[m].surf) continue; for (i = 0; i < am->cnt; i++) { if (am->dof[i] < 0) continue; fv->set_active_shape(am->idx[i]); RHS[am->dof[i]] += liform[m].surf(fv, refmap+m, ep+edge) * am->coef[i]; } } } } trav.finish(); for (i = 0; i < ndofs; i++) RHS[i] += Dir[i]; if (!quiet) verbose(" (time: %g sec)", end_time()); for (i = 0; i < neq; i++) delete spss[i]; delete [] buffer; delete [] refmap; }