12#ifndef _LEPNC_STEFANPME_HPP
13#define _LEPNC_STEFANPME_HPP
19#include <boost/timer/timer.hpp>
22#include <Eigen/Sparse>
79 std::string solver_type,
80 std::ostream & output = std::cout
84 Eigen::VectorXd
solve();
87 Eigen::VectorXd
apply_nonlinearity(
const Eigen::VectorXd& Y,
const std::string type)
const;
93 double EnergyNorm(
const Eigen::VectorXd Xh)
const;
102 Eigen::MatrixXd diffusion_operator(
const size_t iT)
const;
105 Eigen::VectorXd load_operator(
const size_t iT)
const;
106 Eigen::VectorXd load_operator_ml(
const size_t iT)
const;
110 Eigen::VectorXd residual_scheme(
const Eigen::VectorXd& Xh)
const;
120 const double _weight;
121 const std::string solver_type;
122 std::ostream & m_output;
125 std::vector<Eigen::MatrixXd> DiffT;
126 std::vector<Eigen::MatrixXd> MassT;
127 std::vector<Eigen::VectorXd> SourceT;
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);
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)
140 _deg_kappa(deg_kappa),
143 exact_solution(exact_solution),
144 grad_exact_solution(grad_exact_solution),
147 solver_type(solver_type),
154 SourceT.resize(mesh->
n_cells());
155 for (
size_t iT = 0; iT < mesh->
n_cells(); iT++) {
159 DiffT[iT] = diffusion_operator(iT);
160 MassT[iT] = Eigen::MatrixXd::Zero(n_local_dofs, n_local_dofs);
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);
170 boost::timer::cpu_timer timer;
171 boost::timer::cpu_timer timerint;
173 size_t n_edges_dofs = mesh->
n_edges();
177 constexpr double tol = 1e-6;
178 constexpr size_t maxiter = 400;
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);
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();
194 DirBC(n_cell_dofs + iF) += quadrule.w * zeta(exact_solution(quadrule.x, quadrule.y),
"fct");
196 DirBC(n_cell_dofs + iF) += quadrule.w * exact_solution(quadrule.x, quadrule.y);
200 DirBC(n_cell_dofs + iF) /= mesh->b_edge(ibF)->measure();
207 Eigen::VectorXd RES = residual_scheme(Xhprev);
208 std::vector<double> residual(maxiter, 1.0);
209 residual[iter] = RES.norm();
211 while ( (iter < maxiter) && (residual[iter] > tol) ){
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);
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);
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();
230 Eigen::VectorXd XTprev = nc.
nc_restr(Xhprev, iT);
233 for (
size_t i = 0; i < DimPoly<Cell>(1) + n_edgesT;
i++){
234 MatZetaprimeXTprev(
i,
i) = ZetaprimeXTprev(
i);
236 Eigen::MatrixXd MatT = DiffT[iT] * MatZetaprimeXTprev + MassT[iT];
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);
242 Eigen::PartialPivLU<Eigen::MatrixXd> invATT;
245 Eigen::MatrixXd invATT_ATF = invATT.solve(ATF);
247 Eigen::VectorXd invATT_RES_cell = invATT.solve(RES_cell);
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;
257 for (
size_t i = 0; i < 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));
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));
272 GlobRHS(iGlobal) += bF(
i);
277 GlobMat.setFromTriplets(std::begin(triplets_GlobMat), std::end(triplets_GlobMat));
278 ScMat.setFromTriplets(std::begin(triplets_ScMat), std::end(triplets_ScMat));
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);
287 Eigen::VectorXd B = GlobRHS.head(n_edge_unknowns);
290 Eigen::BiCGSTAB<Eigen::SparseMatrix<double> > solver;
292 Eigen::VectorXd dX_edge_unknowns = solver.solve(B);
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);
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);
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]){
313 relax = std::min(1.0, 1.2*relax);
317 residual[iter] = residual_temp;
320 m_output <<
" ...Iteration: " << iter <<
"; residual = " << residual[iter] <<
"; relax = "<< relax <<
"\n";
339 Eigen::VectorXd ZetaY;
344 }
else if (type ==
"der" || type ==
"nlin"){
345 ZetaY = Eigen::VectorXd::Ones(Y.size());
347 m_output <<
"type unknown in apply_nonlinearity: " << type <<
"\n";
353 size_t n_cell_dofs = 0;
359 size_t n_edges_dofs = Y.size() - n_cell_dofs;
363 for (
size_t i = 0;
i < n_cell_dofs;
i++){
364 ZetaY(
i) = zeta(Y(
i), type);
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);
371 for (
size_t i = 0;
i < n_cell_dofs + n_edges_dofs;
i++){
372 ZetaY(
i) = zeta(Y(
i), type);
385Eigen::MatrixXd LEPNC_StefanPME::diffusion_operator(
const size_t iT)
const {
388 Cell* cell = mesh->cell(iT);
389 const size_t n_edgesT = cell->n_edges();
395 std::function<Eigen::Matrix2d(
double,
double)> kappaT = [&](
double x,
double y){
396 return kappa(x, y, cell);
400 size_t doe = 2 + _deg_kappa;
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); });
408 return nc.
gram_matrix(Dnc_phiT_quadT, Dnc_phiT_quadT, local_dofs, local_dofs, quadT,
true, kappaT_quadT);
416Eigen::VectorXd LEPNC_StefanPME::load_operator(
const size_t iT)
const {
419 Cell* cell = mesh->cell(iT);
420 size_t n_edgesT = cell->n_edges();
422 Eigen::VectorXd b = Eigen::VectorXd::Zero(local_dofs);
427 size_t nbq = quadT.size();
428 std::vector<Eigen::ArrayXd> nc_phiT_quadT = nc.
nc_basis_quad(iT, quadT);
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);
436 for (
size_t i=0;
i < local_dofs;
i++){
437 b(
i) = (weight_source_quad * nc_phiT_quadT[
i]).sum();
441 if (cell->is_boundary()){
443 for (
size_t ilF = 0; ilF < cell->n_edges(); ilF++) {
444 Edge* edge = cell->edge(ilF);
445 if (m_BC.
type(*edge)==
"neu"){
447 if (edge->is_boundary()){
448 const auto& nTF = cell->edge_normal(ilF);
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());
465Eigen::VectorXd LEPNC_StefanPME::load_operator_ml(
const size_t iT)
const {
468 Cell* cell = mesh->cell(iT);
469 size_t n_edgesT = cell->n_edges();
473 if (MassT.size() < iT ||
size_t(MassT[iT].rows()) != local_dofs){
474 m_output <<
"Called load_operator_ml without creating MassT[iT]\n";
477 Eigen::VectorXd fvec = Eigen::VectorXd::Zero(
DimPoly<Cell>(1) + n_edgesT);
478 for (
size_t i = 0;
i < local_dofs;
i++){
485 fvec(
i) = source(node.x(), node.y(), cell);
488 Eigen::VectorXd load = MassT[iT] * fvec;
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);
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));
521Eigen::VectorXd LEPNC_StefanPME::residual_scheme(
const Eigen::VectorXd& Xh)
const{
528 for (
size_t iT = 0; iT < mesh->
n_cells(); iT++){
530 size_t n_edgesT = cell->
n_edges();
532 Eigen::VectorXd XT = nc.
nc_restr(Xh, iT);
535 Eigen::VectorXd REST = SourceT[iT] - (DiffT[iT] * ZetaXT + MassT[iT] * XT);
537 for (
size_t ilE = 0; ilE < n_edgesT; ilE++){
538 size_t iE = cell->edge(ilE)->global_index();
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;
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;
579 return double(_assembly_time) * pow(10, -9);
583 return double(_solving_time) * pow(10, -9);
587 return double(_itime[idx]) * pow(10, -9);
591 return _solving_error;
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
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
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