HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
LEPNC_StefanPME.hpp
Go to the documentation of this file.
1// Implementation of the NC Poly scheme in 2D for the Stefan and PME models
2//
3// { u - div(K \grad(zeta(u))) = f, inside Omega
4// { K \grad(zeta(u)) . nTF = g, on GammaN
5// { u = g, on GammaD
6//
7// At the moment, only pure Neumann or pure Dirichlet
8//
9// Author: Jerome Droniou (jerome.droniou@monash.edu)
10//
11
12#ifndef _LEPNC_STEFANPME_HPP
13#define _LEPNC_STEFANPME_HPP
14
15#include <functional>
16#include <utility>
17#include <iostream>
18
19#include <boost/timer/timer.hpp>
20
21// Matrices and linear solvers
22#include <Eigen/Sparse>
23#include <Eigen/Dense>
24//#include "Eigen/MA41.cpp"
25
26#include "mesh.hpp"
27#include "lepnccore.hpp"
28#include "quad2d.hpp"
32
37namespace HArDCore2D {
38
43// ----------------------------------------------------------------------------
44// Class definition
45// ----------------------------------------------------------------------------
46/* The LEPNC_StefanPME class provides an implementation of the LEPNC method for the stationnary Stefan and PME problems.
47
48* If using this code in a scientific publication, please cite the reference for the LEPNC:
49*
50* Non-conforming finite elements on polytopal mesh,
51* J. Droniou, R. Eymard, T. Gallouët and R. Herbin
52* url:
53*
54*/
56
58
60
61// Types
62public:
63 using solution_function_type = std::function<double(double,double)>;
64 using source_function_type = std::function<double(double,double,Cell*)>;
65 using grad_function_type = std::function<VectorRd(double,double,Cell*)>;
66 using tensor_function_type = std::function<Eigen::Matrix2d(double,double,Cell*)>;
67
70 LEPNCCore &nc,
72 size_t deg_kappa,
75 solution_function_type exact_solution,
76 grad_function_type grad_exact_solution,
78 double weight,
79 std::string solver_type,
80 std::ostream & output = std::cout
81 );
82
84 Eigen::VectorXd solve();
85
87 Eigen::VectorXd apply_nonlinearity(const Eigen::VectorXd& Y, const std::string type) const;
88
90 double L2_MassLumped(const Eigen::VectorXd Xh) const;
91
93 double EnergyNorm(const Eigen::VectorXd Xh) const;
94
95 double get_assembly_time() const;
96 double get_solving_time() const;
97 double get_solving_error() const;
98 double get_itime(size_t idx) const;
99
100private:
102 Eigen::MatrixXd diffusion_operator(const size_t iT) const;
103
105 Eigen::VectorXd load_operator(const size_t iT) const;
106 Eigen::VectorXd load_operator_ml(const size_t iT) const;
107
110 Eigen::VectorXd residual_scheme(const Eigen::VectorXd& Xh) const;
111
112 const LEPNCCore& nc;
113 const tensor_function_type kappa;
114 size_t _deg_kappa;
115 const source_function_type source;
116 const BoundaryConditions m_BC;
117 const solution_function_type exact_solution;
118 const grad_function_type grad_exact_solution;
120 const double _weight;
121 const std::string solver_type;
122 std::ostream & m_output;
123
124 // To store local diffusion, mass and source terms
125 std::vector<Eigen::MatrixXd> DiffT;
126 std::vector<Eigen::MatrixXd> MassT;
127 std::vector<Eigen::VectorXd> SourceT;
128
129 // Computation statistics
130 size_t _assembly_time;
131 size_t _solving_time;
132 double _solving_error;
133 mutable std::vector<size_t> _itime = std::vector<size_t>(10, 0);
134
135};
136
137LEPNC_StefanPME::LEPNC_StefanPME(LEPNCCore &nc, tensor_function_type kappa, size_t deg_kappa, source_function_type source, BoundaryConditions BC, solution_function_type exact_solution, grad_function_type grad_exact_solution, TestCaseNonLinearity::nonlinearity_function_type zeta, double weight, std::string solver_type, std::ostream & output)
138 : nc(nc),
139 kappa(kappa),
140 _deg_kappa(deg_kappa),
141 source(source),
142 m_BC(BC),
143 exact_solution(exact_solution),
144 grad_exact_solution(grad_exact_solution),
145 zeta(zeta),
146 _weight(weight),
147 solver_type(solver_type),
148 m_output(output) {
149
150 //---- Compute diffusion, mass matrix and source (do not change during Newton) ------//
151 const Mesh* mesh = nc.get_mesh();
152 DiffT.resize(mesh->n_cells());
153 MassT.resize(mesh->n_cells());
154 SourceT.resize(mesh->n_cells());
155 for (size_t iT = 0; iT < mesh->n_cells(); iT++) {
156 size_t n_edgesT = mesh->cell(iT)->n_edges();
157 size_t n_local_dofs = DimPoly<Cell>(1) + n_edgesT;
158 double measT = mesh->cell(iT)->measure();
159 DiffT[iT] = diffusion_operator(iT);
160 MassT[iT] = Eigen::MatrixXd::Zero(n_local_dofs, n_local_dofs);
161 MassT[iT].topLeftCorner(DimPoly<Cell>(1), DimPoly<Cell>(1)) = (1-_weight) * (measT / DimPoly<Cell>(1)) * Eigen::MatrixXd::Identity(DimPoly<Cell>(1), DimPoly<Cell>(1));
162 MassT[iT].bottomRightCorner(n_edgesT, n_edgesT) = _weight * (measT / n_edgesT) * Eigen::MatrixXd::Identity(n_edgesT, n_edgesT);
163 SourceT[iT] = load_operator_ml(iT);
164 }
165}
166
167// ASSEMBLE AND SOLVE THE SCHEME
168Eigen::VectorXd LEPNC_StefanPME::solve() {
169
170 boost::timer::cpu_timer timer;
171 boost::timer::cpu_timer timerint;
172 const auto mesh = nc.get_mesh();
173 size_t n_edges_dofs = mesh->n_edges();
174 size_t n_cell_dofs = mesh->n_cells() * DimPoly<Cell>(1);
175
176 //----- PARAMETERS FOR NEWTON ------------//
177 constexpr double tol = 1e-6;
178 constexpr size_t maxiter = 400;
179 // relaxation parameter
180 double relax = 1;
181 size_t iter = 0;
182 Eigen::VectorXd Xhprev = Eigen::VectorXd::Zero(n_cell_dofs + n_edges_dofs);
183 Eigen::VectorXd Xh = Eigen::VectorXd::Zero(n_cell_dofs + n_edges_dofs);
184
185 // Vector of Dirichlet BC, either on u (if weight>0) or zeta(u) (if weight=0)
186 Eigen::VectorXd DirBC = Eigen::VectorXd::Zero(n_cell_dofs + n_edges_dofs);
187 for (size_t ibF = 0; ibF < mesh->n_b_edges(); ibF++){
188 Edge* edge = mesh->b_edge(ibF);
189 if (m_BC.type(*edge)=="dir"){
190 size_t iF = edge->global_index();
192 for (QuadratureNode quadrule : quadF){
193 if (_weight == 0){
194 DirBC(n_cell_dofs + iF) += quadrule.w * zeta(exact_solution(quadrule.x, quadrule.y), "fct");
195 }else{
196 DirBC(n_cell_dofs + iF) += quadrule.w * exact_solution(quadrule.x, quadrule.y);
197 }
198 }
199 // Take average
200 DirBC(n_cell_dofs + iF) /= mesh->b_edge(ibF)->measure();
201 }
202 }
203 Xhprev = DirBC;
204
205
206 // ---- NEWTON ITERATONS ------ //
207 Eigen::VectorXd RES = residual_scheme(Xhprev); // Scheme is F(X)=C, then RES = C - F(u_prev)
208 std::vector<double> residual(maxiter, 1.0);
209 residual[iter] = RES.norm();
210
211 while ( (iter < maxiter) && (residual[iter] > tol) ){
212 iter++;
213
214 // System matrix and initialisation of RHS (only on the edge dofs)
215 Eigen::SparseMatrix<double> GlobMat(n_edges_dofs, n_edges_dofs);
216 std::vector<Eigen::Triplet<double>> triplets_GlobMat;
217 Eigen::VectorXd GlobRHS = RES.tail(n_edges_dofs);
218
219 // Static condensation: matrix and source to recover cell unknowns
220 Eigen::SparseMatrix<double> ScMat(n_cell_dofs, n_edges_dofs);
221 std::vector<Eigen::Triplet<double>> triplets_ScMat;
222 Eigen::VectorXd ScRHS = Eigen::VectorXd::Zero(n_cell_dofs);
223
224 // Assemble local contributions
225 for (size_t iT = 0; iT < mesh->n_cells(); iT++) {
226 Cell* cell = mesh->cell(iT);
227 size_t n_edgesT = cell->n_edges();
228
229 // Local static condensation of element unknowns
230 Eigen::VectorXd XTprev = nc.nc_restr(Xhprev, iT);
231 Eigen::VectorXd ZetaprimeXTprev = apply_nonlinearity(XTprev, "der");
232 Eigen::MatrixXd MatZetaprimeXTprev = Eigen::MatrixXd::Zero(DimPoly<Cell>(1)+n_edgesT, DimPoly<Cell>(1)+n_edgesT);
233 for (size_t i = 0; i < DimPoly<Cell>(1) + n_edgesT; i++){
234 MatZetaprimeXTprev(i, i) = ZetaprimeXTprev(i);
235 }
236 Eigen::MatrixXd MatT = DiffT[iT] * MatZetaprimeXTprev + MassT[iT];
237 Eigen::MatrixXd ATT = MatT.topLeftCorner(DimPoly<Cell>(1), DimPoly<Cell>(1));
238 Eigen::MatrixXd ATF = MatT.topRightCorner(DimPoly<Cell>(1), n_edgesT);
239 Eigen::MatrixXd AFT = MatT.bottomLeftCorner(n_edgesT, DimPoly<Cell>(1));
240 Eigen::MatrixXd AFF = MatT.bottomRightCorner(n_edgesT, n_edgesT);
241
242 Eigen::PartialPivLU<Eigen::MatrixXd> invATT;
243 invATT.compute(ATT);
244
245 Eigen::MatrixXd invATT_ATF = invATT.solve(ATF);
246 Eigen::VectorXd RES_cell = RES.segment(iT * DimPoly<Cell>(1), DimPoly<Cell>(1));
247 Eigen::VectorXd invATT_RES_cell = invATT.solve(RES_cell);
248
249 // Local matrix and right-hand side on the face unknowns
250 Eigen::MatrixXd MatF = Eigen::MatrixXd::Zero(n_edgesT, n_edgesT);
251 Eigen::VectorXd bF = Eigen::VectorXd::Zero(n_edgesT);
252 MatF = AFF - AFT * invATT_ATF;
253 bF = - AFT * invATT_RES_cell; // only additional RHS coming from static condensation, beyond RES
254
255 // Assemble static condensation operator
256 ScRHS.segment(iT * DimPoly<Cell>(1), DimPoly<Cell>(1)) = invATT_RES_cell;
257 for (size_t i = 0; i < DimPoly<Cell>(1); i++){
258 size_t iGlobal = iT * DimPoly<Cell>(1) + i;
259 for (size_t j = 0; j < n_edgesT; j++){
260 size_t jGlobal = cell->edge(j)->global_index();
261 triplets_ScMat.emplace_back(iGlobal, jGlobal, invATT_ATF(i, j));
262 }
263 }
264
265 // GLOBAL MATRIX and RHS on the edge unknowns, using the matrices MatF obtained after static condensation
266 for (size_t i = 0; i < n_edgesT; i++){
267 size_t iGlobal = cell->edge(i)->global_index();
268 for (size_t j = 0; j < n_edgesT; j++){
269 size_t jGlobal = cell->edge(j)->global_index();
270 triplets_GlobMat.emplace_back(iGlobal, jGlobal, MatF(i, j));
271 }
272 GlobRHS(iGlobal) += bF(i);
273 }
274 }
275
276 // Assemble the global linear system (without BC), and matrix to recover statically-condensed cell dofs
277 GlobMat.setFromTriplets(std::begin(triplets_GlobMat), std::end(triplets_GlobMat));
278 ScMat.setFromTriplets(std::begin(triplets_ScMat), std::end(triplets_ScMat));
279
280 // Dirichlet boundary conditions are trivial. Newton requires to solve
281 // F'(X^k) (X^{k+1} - X^k) = B - F(X^k)
282 // Both X^k and X^{k+1} have fixed Dirichlet values on boundary edges, so the system we solve
283 // is Dirichlet homogeneous, and we just select the non-dirichlet edges (first ones in the list)
284 size_t n_edge_unknowns = mesh->n_edges() - m_BC.n_dir_edges();
285 Eigen::SparseMatrix<double> A = GlobMat.topLeftCorner(n_edge_unknowns, n_edge_unknowns);
286//std::cout << "Number non-zero terms in matrix: " << A.nonZeros() << "\n";
287 Eigen::VectorXd B = GlobRHS.head(n_edge_unknowns);
288
289 // Solve condensed system and recover cell unknowns. dX = X^{k+1}-X^k//
290 Eigen::BiCGSTAB<Eigen::SparseMatrix<double> > solver;
291 solver.compute(A);
292 Eigen::VectorXd dX_edge_unknowns = solver.solve(B);
293
294 Eigen::VectorXd dX = Eigen::VectorXd::Zero(n_cell_dofs + n_edges_dofs);
295 dX.segment(n_cell_dofs, n_edge_unknowns) = dX_edge_unknowns;
296 dX.head(n_cell_dofs) = ScRHS - ScMat * dX.tail(n_edges_dofs);
297
298 // Recover the fixed boundary values and cell unknowns (from static condensation and Newton).
299 // We start by putting the Dirichlet BC at the end.
300 Xh = DirBC;
301 Xh.head(n_cell_dofs + n_edge_unknowns) = relax*dX.head(n_cell_dofs + n_edge_unknowns) + Xhprev.head(n_cell_dofs + n_edge_unknowns);
302
303 // Compute new residual. We might not accept new Xh but start over reducing relax. Otherwise, accept Xh and increase relax
304 // to a max of 1
305 Eigen::VectorXd RES_temp = residual_scheme(Xh);
306 double residual_temp = RES_temp.norm();
307 if (iter > 1 && residual_temp > 10*residual[iter-1]){
308 // We don't update the residual and Xh, but we will do another iteration with relax reduced
309 relax = relax/2.0;
310 iter--;
311 } else {
312 // Increase relax
313 relax = std::min(1.0, 1.2*relax);
314 // Store Xh in Xhprev, compute new residual
315 Xhprev = Xh;
316 RES = RES_temp;
317 residual[iter] = residual_temp;
318 }
319
320 m_output << " ...Iteration: " << iter << "; residual = " << residual[iter] << "; relax = "<< relax << "\n";
321
322 } //-------- END NEWTON ITERATIONS -----------//
323
324
325 return Xh;
326}
327
328//****************************************
329// apply nonlinearity zeta to a vector
330//****************************************
331
332Eigen::VectorXd LEPNC_StefanPME::apply_nonlinearity(const Eigen::VectorXd& Y, const std::string type) const {
333 // Type="fct", "nlin" or "der" as per the nonlinear function zeta
334 // The function applied to each coefficient of Y depends on its nature. If it's a coefficient appearing in
335 // the mass-lumping (which depends on _weight), the unknown represents u and we apply the full nonlinear function zeta.
336 // Otherwise, the unknown is zeta(u) itself and the function we apply is the identity s->s (with nlin/der: s-> 1)
337 //
338
339 Eigen::VectorXd ZetaY;
340 // Initialisation depends on "type". We initialise as if the nonlinearity was the identity, we'll change the
341 // appropriate coefficients afterwards
342 if (type == "fct"){
343 ZetaY = Y;
344 }else if (type == "der" || type == "nlin"){
345 ZetaY = Eigen::VectorXd::Ones(Y.size());
346 }else{
347 m_output << "type unknown in apply_nonlinearity: " << type << "\n";
348 exit(EXIT_FAILURE);
349 }
350
351 const Mesh* mesh = nc.get_mesh();
352 // We check if Y is a local vector or a global one, this determines the number of cell and edge unknowns
353 size_t n_cell_dofs = 0;
354 if (size_t(Y.size()) == DimPoly<Cell>(1)*mesh->n_cells() + mesh->n_edges()){
355 n_cell_dofs = DimPoly<Cell>(1)*mesh->n_cells();
356 }else{
357 n_cell_dofs = DimPoly<Cell>(1);
358 }
359 size_t n_edges_dofs = Y.size() - n_cell_dofs;
360
361 // The type of nonlinearity we apply on each coefficient depends on _weight
362 if (_weight == 0){
363 for (size_t i = 0; i < n_cell_dofs; i++){
364 ZetaY(i) = zeta(Y(i), type);
365 }
366 }else if (_weight == 1){
367 for (size_t i = n_cell_dofs; i < n_cell_dofs + n_edges_dofs; i++){
368 ZetaY(i) = zeta(Y(i), type);
369 }
370 }else{
371 for (size_t i = 0; i < n_cell_dofs + n_edges_dofs; i++){
372 ZetaY(i) = zeta(Y(i), type);
373 }
374 }
375
376
377 return ZetaY;
378}
379
380
381//********************************
382// local diffusion matrix
383//********************************
384
385Eigen::MatrixXd LEPNC_StefanPME::diffusion_operator(const size_t iT) const {
386
387 const auto mesh = nc.get_mesh();
388 Cell* cell = mesh->cell(iT);
389 const size_t n_edgesT = cell->n_edges();
390
391 // Total number of degrees of freedom local to this cell
392 size_t local_dofs = n_edgesT + DimPoly<Cell>(1);
393
394 // Diffusion in the cell.
395 std::function<Eigen::Matrix2d(double,double)> kappaT = [&](double x, double y){
396 return kappa(x, y, cell);
397 };
398
399 // Compute cell quadrature nodes and values of cell basis functions and gradients at these nodes
400 size_t doe = 2 + _deg_kappa;
401 QuadratureRule quadT = generate_quadrature_rule(*cell, doe, true);
402 std::vector<Eigen::ArrayXXd> Dnc_phiT_quadT = nc.grad_nc_basis_quad(iT, quadT);
403 // Diffusion tensor at the quadrature nodes
404 std::vector<Eigen::Matrix2d> kappaT_quadT(quadT.size());
405 std::transform(quadT.begin(), quadT.end(), kappaT_quadT.begin(),
406 [&kappaT](QuadratureNode qr) -> Eigen::MatrixXd { return kappaT(qr.x, qr.y); });
407
408 return nc.gram_matrix(Dnc_phiT_quadT, Dnc_phiT_quadT, local_dofs, local_dofs, quadT, true, kappaT_quadT);
409}
410
411
412//********************************
413// local load term
414//********************************
415
416Eigen::VectorXd LEPNC_StefanPME::load_operator(const size_t iT) const {
417 // Load for the cell DOFs (first indices) and face DOFs (last indices)
418 const auto mesh = nc.get_mesh();
419 Cell* cell = mesh->cell(iT);
420 size_t n_edgesT = cell->n_edges();
421 size_t local_dofs = DimPoly<Cell>(1) + n_edgesT;
422 Eigen::VectorXd b = Eigen::VectorXd::Zero(local_dofs);
423
424 // Quadrature nodes and values of cell basis functions at these nodes
425 size_t doe = 5;
426 QuadratureRule quadT = generate_quadrature_rule(*cell, doe, true);
427 size_t nbq = quadT.size();
428 std::vector<Eigen::ArrayXd> nc_phiT_quadT = nc.nc_basis_quad(iT, quadT);
429
430 // Value of source times quadrature weights at the quadrature nodes
431 Eigen::ArrayXd weight_source_quad = Eigen::ArrayXd::Zero(nbq);
432 for (size_t iqn = 0; iqn < nbq; iqn++){
433 weight_source_quad(iqn) = quadT[iqn].w * source(quadT[iqn].x, quadT[iqn].y, cell);
434 }
435
436 for (size_t i=0; i < local_dofs; i++){
437 b(i) = (weight_source_quad * nc_phiT_quadT[i]).sum();
438 }
439
440 // Boundary values, if we have a boundary cell and Neumann edges
441 if (cell->is_boundary()){
442 // Boundary values only on boundary Neumann edges
443 for (size_t ilF = 0; ilF < cell->n_edges(); ilF++) {
444 Edge* edge = cell->edge(ilF);
445 if (m_BC.type(*edge)=="neu"){
446 // BC on boundary faces
447 if (edge->is_boundary()){
448 const auto& nTF = cell->edge_normal(ilF);
449 const auto& phi_i = nc.nc_basis(iT, DimPoly<Cell>(1) + ilF);
451 std::function<double(VectorRd)> Kgrad_n = [&](VectorRd p){
452 return zeta(exact_solution(p.x(), p.y()), "der") * nTF.dot(kappa(p.x(),p.y(),cell) * grad_exact_solution(p.x(),p.y(),cell)) * phi_i(p.x(),p.y());
453 };
454 for (QuadratureNode qr : quadF){
455 b(DimPoly<Cell>(1) + ilF) += qr.w * Kgrad_n(qr.vector());
456 }
457 }
458 }
459 }
460 }
461
462 return b;
463}
464
465Eigen::VectorXd LEPNC_StefanPME::load_operator_ml(const size_t iT) const {
466 // Load for the cell DOFs (first indices) and face DOFs (last indices, will be 0)
467 const auto mesh = nc.get_mesh();
468 Cell* cell = mesh->cell(iT);
469 size_t n_edgesT = cell->n_edges();
470 size_t local_dofs = DimPoly<Cell>(1) + n_edgesT;
471
472 // Beware, MassT[iT] must have been created
473 if (MassT.size() < iT || size_t(MassT[iT].rows()) != local_dofs){
474 m_output << "Called load_operator_ml without creating MassT[iT]\n";
475 exit(EXIT_FAILURE);
476 }
477 Eigen::VectorXd fvec = Eigen::VectorXd::Zero(DimPoly<Cell>(1) + n_edgesT);
478 for (size_t i = 0; i < local_dofs; i++){
479 VectorRd node;
480 if (i < DimPoly<Cell>(1)){
481 node = nc.ml_node(iT, i);
482 }else{
483 node = cell->edge(i-DimPoly<Cell>(1))->center_mass();
484 }
485 fvec(i) = source(node.x(), node.y(), cell);
486 }
487 // Contribution of source to loading term
488 Eigen::VectorXd load = MassT[iT] * fvec;
489
490 // Adding Neumann boundary conditions
491 // The boundary conditions coming from Neumann edges are computed based on the piecewise constant
492 // discrete trace: if g is the value of the outer flux,
493 // \int g v -> \int g Tv
494 // where (Tv)_e = average of v on e, for e Neumann edge.
495 // Only edge-based basis functions therefore have a contribution, which will simply be \int_e g.
496 if (cell->is_boundary()){
497 for (size_t ilF = 0; ilF < cell->n_edges(); ilF++) {
498 Edge* edge = cell->edge(ilF);
499 if (m_BC.type(*edge)=="neu"){
500 if (edge->is_boundary()){
501 const auto& nTF = cell->edge_normal(ilF);
503 std::function<double(VectorRd)> Kgrad_n = [&](VectorRd p){
504 return zeta(exact_solution(p.x(), p.y()), "der") * nTF.dot(kappa(p.x(),p.y(),cell) * grad_exact_solution(p.x(),p.y(),cell));
505 };
506 for (QuadratureNode qr : quadF){
507 load(DimPoly<Cell>(1) + ilF) += qr.w * Kgrad_n(qr.vector());
508 }
509 }
510 }
511 }
512 }
513
514 return load;
515}
516
517//****************************
518// Residual non-linear scheme
519//****************************
520
521Eigen::VectorXd LEPNC_StefanPME::residual_scheme(const Eigen::VectorXd& Xh) const{
522
523 const Mesh* mesh = nc.get_mesh();
524 Eigen::VectorXd RES = Eigen::VectorXd::Zero(mesh->n_cells() * DimPoly<Cell>(1) + mesh->n_edges());
525
526 // Compute the residual: C - (DIFF * ZetaXh + MASS * Xh), where DIFF is the diffusion matrix, MASS is the mass matrix and C is the source term
527 size_t n_cell_dofs = mesh->n_cells() * DimPoly<Cell>(1);
528 for (size_t iT = 0; iT < mesh->n_cells(); iT++){
529 Cell* cell = mesh->cell(iT);
530 size_t n_edgesT = cell->n_edges();
531 // Local unknown in the cell, and zeta of these unknowns
532 Eigen::VectorXd XT = nc.nc_restr(Xh, iT);
533 Eigen::VectorXd ZetaXT = apply_nonlinearity(XT, "fct");
534
535 Eigen::VectorXd REST = SourceT[iT] - (DiffT[iT] * ZetaXT + MassT[iT] * XT);
536 RES.segment(iT*DimPoly<Cell>(1), DimPoly<Cell>(1)) += REST.head(DimPoly<Cell>(1));
537 for (size_t ilE = 0; ilE < n_edgesT; ilE++){
538 size_t iE = cell->edge(ilE)->global_index();
539 RES(n_cell_dofs + iE) += REST(DimPoly<Cell>(1)+ilE);
540 }
541 }
542
543 // No residual on Dirichlet edges (do not correspond to equations of the system)
544 RES.tail(m_BC.n_dir_edges()) = Eigen::VectorXd::Zero(m_BC.n_dir_edges());
545
546 return RES;
547}
548
549//*******************************************************
550// Norms: L2 mass lumped, Energy (pure diffusion)
551//*******************************************************
552
553double LEPNC_StefanPME::L2_MassLumped(const Eigen::VectorXd Xh) const {
554 const Mesh* mesh = nc.get_mesh();
555 double value = 0.0;
556
557 for (size_t iT = 0; iT < mesh-> n_cells(); iT++) {
558 Eigen::VectorXd XTF = nc.nc_restr(Xh, iT);
559 value += XTF.transpose() * MassT[iT] * XTF;
560 }
561
562 return sqrt(value);
563}
564
565double LEPNC_StefanPME::EnergyNorm(const Eigen::VectorXd Xh) const {
566 const Mesh* mesh = nc.get_mesh();
567 double value = 0.0;
568
569 for (size_t iT = 0; iT < mesh-> n_cells(); iT++) {
570 Eigen::VectorXd XTF = nc.nc_restr(Xh, iT);
571 value += XTF.transpose() * DiffT[iT] * XTF;
572 }
573
574 return sqrt(value);
575}
576
577
579 return double(_assembly_time) * pow(10, -9);
580}
581
583 return double(_solving_time) * pow(10, -9);
584}
585
586double LEPNC_StefanPME::get_itime(size_t idx) const {
587 return double(_itime[idx]) * pow(10, -9);
588}
589
591 return _solving_error;
592}
593
595} // end of namespace HArDCore2D
596
597#endif //_LEPNC_STEFANPME_HPP
The BoundaryConditions class provides definition of boundary conditions.
Definition BoundaryConditions.hpp:45
const std::string type(const Edge &edge) const
Test the boundary condition of an edge.
Definition BoundaryConditions.cpp:41
const size_t n_dir_edges() const
Returns the number of Dirichlet edges.
Definition BoundaryConditions.hpp:70
Definition lepnccore.hpp:33
LEPNC scheme for diffusion equation .
Definition LEPNC_StefanPME.hpp:59
Definition Mesh2D.hpp:26
std::function< double(double, std::string)> nonlinearity_function_type
type for nonlinear function
Definition TestCaseNonLinearity.hpp:48
Compute max and min eigenvalues of all matrices for i
Definition compute_eigs.m:5
Eigen::Vector2d VectorRd
Definition basis.hpp:55
const size_t DimPoly< Cell >(const int m)
Compute the size of the basis of 2-variate polynomials up to degree m.
Definition hybridcore.hpp:63
const Mesh * get_mesh() const
Returns a pointer to the mesh.
Definition hybridcore.hpp:201
double L2_MassLumped(const Eigen::VectorXd Xh) const
Mass-lumped L2 norm of a function given by a vector.
Definition LEPNC_StefanPME.hpp:553
Eigen::MatrixXd gram_matrix(const std::vector< Eigen::ArrayXd > &f_quad, const std::vector< Eigen::ArrayXd > &g_quad, const size_t &nrows, const size_t &ncols, const QuadratureRule &quad, const bool &sym, std::vector< double > L2weight={}) const
Definition lepnccore.hpp:517
const std::vector< Eigen::ArrayXd > nc_basis_quad(const size_t iT, const QuadratureRule quad) const
Computes non-conforming basis functions at the given quadrature nodes.
Definition lepnccore.hpp:432
Eigen::VectorXd nc_restr(const Eigen::VectorXd &Xh, size_t iT) const
Extract from a global vector Xh of unknowns the non-conforming unknowns corresponding to cell iT.
Definition lepnccore.hpp:678
double get_solving_error() const
residual after solving the scheme
Definition LEPNC_StefanPME.hpp:590
LEPNC_StefanPME(LEPNCCore &nc, tensor_function_type kappa, size_t deg_kappa, source_function_type source, BoundaryConditions BC, solution_function_type exact_solution, grad_function_type grad_exact_solution, TestCaseNonLinearity::nonlinearity_function_type zeta, double weight, std::string solver_type, std::ostream &output=std::cout)
Constructor of the class.
Definition LEPNC_StefanPME.hpp:137
const std::vector< Eigen::ArrayXXd > grad_nc_basis_quad(const size_t iT, const QuadratureRule quad) const
Compute at the quadrature nodes, where are the cell basis functions.
Definition lepnccore.hpp:471
std::function< double(double, double)> solution_function_type
type for solution
Definition LEPNC_StefanPME.hpp:63
std::function< Eigen::Matrix2d(double, double, Cell *)> tensor_function_type
type for diffusion tensor
Definition LEPNC_StefanPME.hpp:66
double EnergyNorm(const Eigen::VectorXd Xh) const
Discrete energy norm (associated to the diffusion operator)
Definition LEPNC_StefanPME.hpp:565
Eigen::VectorXd apply_nonlinearity(const Eigen::VectorXd &Y, const std::string type) const
Compute non-linearity on vector (depends if weight=0, weight=1 or weight\in (0,1) )
Definition LEPNC_StefanPME.hpp:332
double get_assembly_time() const
cpu time to assemble the scheme
Definition LEPNC_StefanPME.hpp:578
const VectorRd ml_node(size_t iT, size_t i) const
Return the i'th nodal point in cell iT (for mass lumping)
Definition lepnccore.hpp:422
std::function< double(double, double, Cell *)> source_function_type
type for source
Definition LEPNC_StefanPME.hpp:64
double get_solving_time() const
cpu time to solve the scheme
Definition LEPNC_StefanPME.hpp:582
const cell_basis_type & nc_basis(size_t iT, size_t i) const
Return a reference to the i'th non-conforming basis function of the cell iT.
Definition lepnccore.hpp:410
std::function< VectorRd(double, double, Cell *)> grad_function_type
type for gradient
Definition LEPNC_StefanPME.hpp:65
Eigen::VectorXd solve()
Assemble and solve the scheme.
Definition LEPNC_StefanPME.hpp:168
double get_itime(size_t idx) const
various intermediate assembly times
Definition LEPNC_StefanPME.hpp:586
double measure() const
Return the Lebesgue measure of the Polytope.
Definition Polytope2D.hpp:76
Polytope< DIMENSION > Cell
Definition Polytope2D.hpp:151
Cell * cell(std::size_t index) const
get a constant pointer to a cell using its global index
Definition Mesh2D.hpp:178
size_t n_edges() const
Return the number of edges of the Polytope.
Definition Polytope2D.hpp:88
Polytope< 1 > Edge
A Face is a Polytope with object_dim = DIMENSION - 1.
Definition Polytope2D.hpp:147
std::size_t n_cells() const
number of cells in the mesh.
Definition Mesh2D.hpp:63
std::size_t n_edges() const
number of edges in the mesh.
Definition Mesh2D.hpp:61
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:55
QuadratureRule generate_quadrature_rule(const Cell &T, const int doe, const bool force_split)
Generate quadrature rule on mesh element.
Definition quadraturerule.cpp:10
A
Definition mmread.m:102
for j
Definition mmread.m:174
Definition ddr-klplate.hpp:27
Description of one node and one weight from a quadrature rule.
Definition quadraturerule.hpp:41