203 const unsigned int dim =
mesh.mesh_dimension();
213 for (
const auto & elem : range)
217 for (
unsigned int var=0; var<
n_vars; var++)
219 const FEType & fe_type = dof_map.variable_type(var);
220 const Order element_order = fe_type.order + elem->p_level();
224 fe->attach_quadrature_rule(qrule.get());
226 const std::vector<Real> & JxW = fe->get_JxW();
227 const std::vector<Point> & q_point = fe->get_xyz();
228 const std::vector<std::vector<Real>> * phi = &(fe->get_phi());
230 std::vector<dof_id_type> dof_indices;
232 unsigned int matsize = element_order + 1;
235 matsize *= (element_order + 2);
240 matsize *= (element_order + 3);
244 DenseMatrix<Number> Kp(matsize, matsize);
245 DenseVector<Number> F;
246 DenseVector<Number> Pu_h;
249 Pu_h.resize(matsize);
254 dof_map.dof_indices(elem, dof_indices, var);
255 libmesh_assert_equal_to(dof_indices.size(), phi->size());
257 const unsigned int n_dofs = cast_int<unsigned int>(dof_indices.size());
258 const unsigned int n_qp = qrule->n_points();
260 for (
unsigned int qp = 0; qp < n_qp; qp++)
262 std::vector<Real> psi =
legepoly(
dim, element_order, q_point[qp], matsize);
263 const unsigned int psi_size = cast_int<unsigned int>(psi.size());
265 for (
unsigned int i = 0; i < matsize; i++)
266 for (
unsigned int j = 0; j < matsize; j++)
267 Kp(i, j) += JxW[qp] * psi[i] * psi[j];
270 for (
unsigned int i = 0; i < n_dofs; i++)
273 for (
unsigned int i = 0; i < psi_size; i++)
274 F(i) += JxW[qp] * u_h * psi[i];
277 Kp.lu_solve(F, Pu_h);
280 std::vector<unsigned int> total_degree_per_index;
282 for (
unsigned int poly_deg = 0; poly_deg <= static_cast<unsigned int>(element_order); ++poly_deg)
287 for (
int i = poly_deg; i >= 0; --i)
288 for (
int j = poly_deg - i; j >= 0; --j)
290 int k = poly_deg - i - j;
291 total_degree_per_index.push_back(i + j + k);
296 for (
int i = poly_deg; i >= 0; --i)
298 int j = poly_deg - i;
299 total_degree_per_index.push_back(i + j);
304 total_degree_per_index.push_back(poly_deg);
308 libmesh_error_msg(
"Invalid dimension dim " <<
dim);
313 std::map<unsigned int, std::vector<Number>> coeff_by_degree;
315 for (
unsigned int i = 0; i < Pu_h.size(); ++i)
317 unsigned int degree = total_degree_per_index[i];
321 coeff_by_degree[degree].push_back(Pu_h(i));
325 std::vector<unsigned int> degrees;
326 std::vector<double> norms;
328 for (
const auto & pair : coeff_by_degree)
330 unsigned int deg = pair.first;
331 const std::vector<Number> & coeffs = pair.second;
334 for (
const auto & c : coeffs)
339 degrees.push_back(deg);
340 norms.push_back(
norm);
345 std::vector<double> log_norms, log_degrees;
347 for (
unsigned int i = 0; i < degrees.size(); ++i)
350 if (norms[i] > 1e-12)
352 log_degrees.push_back(std::log(static_cast<double>(degrees[i])));
353 log_norms.push_back(std::log(norms[i]));
359 double Sx = 0., Sy = 0., Sxx = 0., Sxy = 0.;
360 const size_t N = log_degrees.size();
362 for (
size_t i = 0; i < N; ++i)
364 Sx += log_degrees[i];
366 Sxx += log_degrees[i] * log_degrees[i];
367 Sxy += log_degrees[i] * log_norms[i];
370 const double regularity = -
compute_slope(N, Sx, Sy, Sxx, Sxy);
Order
defines an enum for polynomial orders.
ErrorVector & smoothness_per_cell
static std::vector< Real > legepoly(const unsigned int dim, const Order order, const Point p, const unsigned int matsize)
int _extra_order
Extra order to use for quadrature rule.
const MeshBase & get_mesh() const
Number current_solution(const dof_id_type global_dof_number) const
static Real compute_slope(int N, Real Sx, Real Sy, Real Sxx, Real Sxy)
Computes slop in a linear regression.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
const SmoothnessEstimator & smoothness_estimator
auto norm(const libMesh::TypeVector< T > &vector) -> decltype(std::norm(T()))
unsigned int n_vars() const
const DofMap & get_dof_map() const
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.