Reference documentation for deal.II version 8.4.2
fe_abf.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/utilities.h>
18 #include <deal.II/base/quadrature.h>
19 #include <deal.II/base/quadrature_lib.h>
20 #include <deal.II/base/qprojector.h>
21 #include <deal.II/base/table.h>
22 #include <deal.II/grid/tria.h>
23 #include <deal.II/grid/tria_iterator.h>
24 #include <deal.II/dofs/dof_accessor.h>
25 #include <deal.II/fe/fe.h>
26 #include <deal.II/fe/mapping.h>
27 #include <deal.II/fe/fe_abf.h>
28 #include <deal.II/fe/fe_values.h>
29 #include <deal.II/fe/fe_tools.h>
30 
31 #include <sstream>
32 #include <iostream>
33 
34 //TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
35 //adjust_line_dof_index_for_line_orientation_table fields, and write tests
36 //similar to bits/face_orientation_and_fe_q_*
37 
38 
39 DEAL_II_NAMESPACE_OPEN
40 
41 
42 template <int dim>
43 FE_ABF<dim>::FE_ABF (const unsigned int deg)
44  :
45  FE_PolyTensor<PolynomialsABF<dim>, dim> (
46  deg,
47  FiniteElementData<dim>(get_dpo_vector(deg),
48  dim, deg+1, FiniteElementData<dim>::Hdiv),
49  std::vector<bool>(PolynomialsABF<dim>::compute_n_pols(deg), true),
50  std::vector<ComponentMask>(PolynomialsABF<dim>::compute_n_pols(deg),
51  std::vector<bool>(dim,true))),
52  rt_order(deg)
53 {
54  Assert (dim >= 2, ExcImpossibleInDim(dim));
55  const unsigned int n_dofs = this->dofs_per_cell;
56 
58  // First, initialize the
59  // generalized support points and
60  // quadrature weights, since they
61  // are required for interpolation.
63  // Now compute the inverse node
64  //matrix, generating the correct
65  //basis functions from the raw
66  //ones.
67  FullMatrix<double> M(n_dofs, n_dofs);
69 
70  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
71  this->inverse_node_matrix.invert(M);
72  // From now on, the shape functions
73  // will be the correct ones, not
74  // the raw shape functions anymore.
75 
76  // Reinit the vectors of
77  // restriction and prolongation
78  // matrices to the right sizes.
79  // Restriction only for isotropic
80  // refinement
82  // Fill prolongation matrices with embedding operators
84 
86 
87  // TODO[TL]: for anisotropic refinement we will probably need a table of submatrices with an array for each refine case
88  std::vector<FullMatrix<double> >
89  face_embeddings(1<<(dim-1), FullMatrix<double>(this->dofs_per_face,
90  this->dofs_per_face));
91  // TODO: Something goes wrong there. The error of the least squares fit
92  // is to large ...
93  // FETools::compute_face_embedding_matrices(*this, &face_embeddings[0], 0, 0);
94  this->interface_constraints.reinit((1<<(dim-1)) * this->dofs_per_face,
95  this->dofs_per_face);
96  unsigned int target_row=0;
97  for (unsigned int d=0; d<face_embeddings.size(); ++d)
98  for (unsigned int i=0; i<face_embeddings[d].m(); ++i)
99  {
100  for (unsigned int j=0; j<face_embeddings[d].n(); ++j)
101  this->interface_constraints(target_row,j) = face_embeddings[d](i,j);
102  ++target_row;
103  }
104 }
105 
106 
107 
108 template <int dim>
109 std::string
111 {
112  // note that the
113  // FETools::get_fe_from_name
114  // function depends on the
115  // particular format of the string
116  // this function returns, so they
117  // have to be kept in synch
118 
119  std::ostringstream namebuf;
120 
121  namebuf << "FE_ABF<" << dim << ">(" << rt_order << ")";
122 
123  return namebuf.str();
124 }
125 
126 
127 
128 template <int dim>
131 {
132  return new FE_ABF<dim>(rt_order);
133 }
134 
135 
136 //---------------------------------------------------------------------------
137 // Auxiliary and internal functions
138 //---------------------------------------------------------------------------
139 
140 
141 
142 // Version for 2d and higher. See above for 1d version
143 template <int dim>
144 void
146 {
147  QGauss<dim> cell_quadrature(deg+2);
148  const unsigned int n_interior_points = cell_quadrature.size();
149 
150  unsigned int n_face_points = (dim>1) ? 1 : 0;
151  // compute (deg+1)^(dim-1)
152  for (unsigned int d=1; d<dim; ++d)
153  n_face_points *= deg+1;
154 
156  + n_interior_points);
157  this->generalized_face_support_points.resize (n_face_points);
158 
159 
160  // These might be required when the faces contribution is computed
161  // Therefore they will be initialised at this point.
162  std::vector<AnisotropicPolynomials<dim>* > polynomials_abf(dim);
163 
164  // Generate x_1^{i} x_2^{r+1} ...
165  for (unsigned int dd=0; dd<dim; ++dd)
166  {
167  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
168  for (unsigned int d=0; d<dim; ++d)
169  poly[d].push_back (Polynomials::Monomial<double> (deg+1));
171 
172  polynomials_abf[dd] = new AnisotropicPolynomials<dim>(poly);
173  }
174 
175  // Number of the point being entered
176  unsigned int current = 0;
177 
178  if (dim>1)
179  {
180  QGauss<dim-1> face_points (deg+1);
181  TensorProductPolynomials<dim-1> legendre =
183 
184  boundary_weights.reinit(n_face_points, legendre.n());
185 
186 // Assert (face_points.size() == this->dofs_per_face,
187 // ExcInternalError());
188 
189  for (unsigned int k=0; k<n_face_points; ++k)
190  {
191  this->generalized_face_support_points[k] = face_points.point(k);
192  // Compute its quadrature
193  // contribution for each
194  // moment.
195  for (unsigned int i=0; i<legendre.n(); ++i)
196  {
197  boundary_weights(k, i)
198  = face_points.weight(k)
199  * legendre.compute_value(i, face_points.point(k));
200  }
201  }
202 
204  for (; current<GeometryInfo<dim>::faces_per_cell*n_face_points;
205  ++current)
206  {
207  // Enter the support point
208  // into the vector
209  this->generalized_support_points[current] = faces.point(current);
210  }
211 
212 
213  // Now initialise edge interior weights for the ABF elements.
214  // These are completely independent from the usual edge moments. They
215  // stem from applying the Gauss theorem to the nodal values, which
216  // was necessary to cast the ABF elements into the deal.II framework
217  // for vector valued elements.
218  boundary_weights_abf.reinit(faces.size(), polynomials_abf[0]->n() * dim);
219  for (unsigned int k=0; k < faces.size(); ++k)
220  {
221  for (unsigned int i=0; i<polynomials_abf[0]->n() * dim; ++i)
222  {
223  boundary_weights_abf(k,i) = polynomials_abf[i%dim]->
224  compute_value(i / dim, faces.point(k)) * faces.weight(k);
225  }
226  }
227  }
228 
229  // Create Legendre basis for the
230  // space D_xi Q_k
231  if (deg>0)
232  {
233  std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
234 
235  for (unsigned int dd=0; dd<dim; ++dd)
236  {
237  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
238  for (unsigned int d=0; d<dim; ++d)
241 
242  polynomials[dd] = new AnisotropicPolynomials<dim>(poly);
243  }
244 
245  interior_weights.reinit(TableIndices<3>(n_interior_points, polynomials[0]->n(), dim));
246 
247  for (unsigned int k=0; k<cell_quadrature.size(); ++k)
248  {
249  for (unsigned int i=0; i<polynomials[0]->n(); ++i)
250  for (unsigned int d=0; d<dim; ++d)
251  interior_weights(k,i,d) = cell_quadrature.weight(k)
252  * polynomials[d]->compute_value(i,cell_quadrature.point(k));
253  }
254 
255  for (unsigned int d=0; d<dim; ++d)
256  delete polynomials[d];
257  }
258 
259 
260  // Decouple the creation of the generalized support points
261  // from computation of interior weights.
262  for (unsigned int k=0; k<cell_quadrature.size(); ++k)
263  this->generalized_support_points[current++] = cell_quadrature.point(k);
264 
265  // Additional functionality for the ABF elements
266  // TODO: Here the canonical extension of the principle
267  // behind the ABF elements is implemented. It is unclear,
268  // if this really leads to the ABF spaces in 3D!
270  polynomials_abf[0]->n() * dim, dim));
271  Tensor<1, dim> poly_grad;
272 
273  for (unsigned int k=0; k<cell_quadrature.size(); ++k)
274  {
275  for (unsigned int i=0; i<polynomials_abf[0]->n() * dim; ++i)
276  {
277  poly_grad = polynomials_abf[i%dim]->compute_grad(i / dim,cell_quadrature.point(k))
278  * cell_quadrature.weight(k);
279  // The minus sign comes from the use of the Gauss theorem to replace the divergence.
280  for (unsigned int d=0; d<dim; ++d)
281  interior_weights_abf(k,i,d) = -poly_grad[d];
282  }
283  }
284 
285  for (unsigned int d=0; d<dim; ++d)
286  delete polynomials_abf[d];
287 
288  Assert (current == this->generalized_support_points.size(),
289  ExcInternalError());
290 }
291 
292 
293 
294 // This function is the same Raviart-Thomas interpolation performed by
295 // interpolate. Still, we cannot use interpolate, since it was written
296 // for smooth functions. The functions interpolated here are not
297 // smooth, maybe even not continuous. Therefore, we must double the
298 // number of quadrature points in each direction in order to integrate
299 // only smooth functions.
300 
301 // Then again, the interpolated function is chosen such that the
302 // moments coincide with the function to be interpolated.
303 
304 template <int dim>
305 void
307 {
308  if (dim==1)
309  {
311  for (unsigned int i=0; i<GeometryInfo<dim>::max_children_per_cell; ++i)
312  this->restriction[iso][i].reinit(0,0);
313  return;
314  }
316  QGauss<dim-1> q_base (rt_order+1);
317  const unsigned int n_face_points = q_base.size();
318  // First, compute interpolation on
319  // subfaces
320  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
321  {
322  // The shape functions of the
323  // child cell are evaluated
324  // in the quadrature points
325  // of a full face.
326  Quadrature<dim> q_face
327  = QProjector<dim>::project_to_face(q_base, face);
328  // Store shape values, since the
329  // evaluation suffers if not
330  // ordered by point
332  for (unsigned int k=0; k<q_face.size(); ++k)
333  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
334  cached_values(i,k)
335  = this->shape_value_component(i, q_face.point(k),
337 
338  for (unsigned int sub=0; sub<GeometryInfo<dim>::max_children_per_face; ++sub)
339  {
340  // The weight functions for
341  // the coarse face are
342  // evaluated on the subface
343  // only.
344  Quadrature<dim> q_sub
345  = QProjector<dim>::project_to_subface(q_base, face, sub);
346  const unsigned int child
349 
350  // On a certain face, we must
351  // compute the moments of ALL
352  // fine level functions with
353  // the coarse level weight
354  // functions belonging to
355  // that face. Due to the
356  // orthogonalization process
357  // when building the shape
358  // functions, these weights
359  // are equal to the
360  // corresponding shape
361  // functions.
362  for (unsigned int k=0; k<n_face_points; ++k)
363  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
364  for (unsigned int i_face = 0; i_face < this->dofs_per_face; ++i_face)
365  {
366  // The quadrature
367  // weights on the
368  // subcell are NOT
369  // transformed, so we
370  // have to do it here.
371  this->restriction[iso][child](face*this->dofs_per_face+i_face,
372  i_child)
373  += Utilities::fixed_power<dim-1>(.5) * q_sub.weight(k)
374  * cached_values(i_child, k)
375  * this->shape_value_component(face*this->dofs_per_face+i_face,
376  q_sub.point(k),
378  }
379  }
380  }
381 
382  if (rt_order==0) return;
383 
384  // Create Legendre basis for the
385  // space D_xi Q_k. Here, we cannot
386  // use the shape functions
387  std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
388  for (unsigned int dd=0; dd<dim; ++dd)
389  {
390  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
391  for (unsigned int d=0; d<dim; ++d)
394 
395  polynomials[dd] = new AnisotropicPolynomials<dim>(poly);
396  }
397 
398  QGauss<dim> q_cell(rt_order+1);
399  const unsigned int start_cell_dofs
401 
402  // Store shape values, since the
403  // evaluation suffers if not
404  // ordered by point
405  Table<3,double> cached_values(this->dofs_per_cell, q_cell.size(), dim);
406  for (unsigned int k=0; k<q_cell.size(); ++k)
407  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
408  for (unsigned int d=0; d<dim; ++d)
409  cached_values(i,k,d) = this->shape_value_component(i, q_cell.point(k), d);
410 
411  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell; ++child)
412  {
413  Quadrature<dim> q_sub = QProjector<dim>::project_to_child(q_cell, child);
414 
415  for (unsigned int k=0; k<q_sub.size(); ++k)
416  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
417  for (unsigned int d=0; d<dim; ++d)
418  for (unsigned int i_weight=0; i_weight<polynomials[d]->n(); ++i_weight)
419  {
420  this->restriction[iso][child](start_cell_dofs+i_weight*dim+d,
421  i_child)
422  += q_sub.weight(k)
423  * cached_values(i_child, k, d)
424  * polynomials[d]->compute_value(i_weight, q_sub.point(k));
425  }
426  }
427 
428  for (unsigned int d=0; d<dim; ++d)
429  delete polynomials[d];
430 }
431 
432 
433 
434 template <int dim>
435 std::vector<unsigned int>
437 {
438  if (dim == 1)
439  {
440  Assert (false, ExcImpossibleInDim(1));
441  return std::vector<unsigned int>();
442  }
443 
444  // the element is face-based (not
445  // to be confused with George
446  // W. Bush's Faith Based
447  // Initiative...), and we have
448  // (rt_order+1)^(dim-1) DoFs per face
449  unsigned int dofs_per_face = 1;
450  for (int d=0; d<dim-1; ++d)
451  dofs_per_face *= rt_order+1;
452 
453  // and then there are interior dofs
454  const unsigned int
455  interior_dofs = dim*(rt_order+1)*dofs_per_face;
456 
457  std::vector<unsigned int> dpo(dim+1);
458  dpo[dim-1] = dofs_per_face;
459  dpo[dim] = interior_dofs;
460 
461  return dpo;
462 }
463 
464 //---------------------------------------------------------------------------
465 // Data field initialization
466 //---------------------------------------------------------------------------
467 
468 template <int dim>
469 bool
470 FE_ABF<dim>::has_support_on_face (const unsigned int shape_index,
471  const unsigned int face_index) const
472 {
473  Assert (shape_index < this->dofs_per_cell,
474  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
476  ExcIndexRange (face_index, 0, GeometryInfo<dim>::faces_per_cell));
477 
478  // Return computed values if we
479  // know them easily. Otherwise, it
480  // is always safe to return true.
481  switch (rt_order)
482  {
483  case 0:
484  {
485  switch (dim)
486  {
487  case 2:
488  {
489  // only on the one
490  // non-adjacent face
491  // are the values
492  // actually zero. list
493  // these in a table
494  return (face_index != GeometryInfo<dim>::opposite_face[shape_index]);
495  }
496 
497  default:
498  return true;
499  };
500  };
501 
502  default: // other rt_order
503  return true;
504  };
505 
506  return true;
507 }
508 
509 
510 
511 template <int dim>
512 void
514  std::vector<double> &,
515  const std::vector<double> &) const
516 {
517  Assert(false, ExcNotImplemented());
518 }
519 
520 
521 
522 template <int dim>
523 void
525  std::vector<double> &local_dofs,
526  const std::vector<Vector<double> > &values,
527  unsigned int offset) const
528 {
529  Assert (values.size() == this->generalized_support_points.size(),
530  ExcDimensionMismatch(values.size(), this->generalized_support_points.size()));
531  Assert (local_dofs.size() == this->dofs_per_cell,
532  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
533  Assert (values[0].size() >= offset+this->n_components(),
534  ExcDimensionMismatch(values[0].size(),offset+this->n_components()));
535 
536  std::fill(local_dofs.begin(), local_dofs.end(), 0.);
537 
538  const unsigned int n_face_points = boundary_weights.size(0);
539  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
540  for (unsigned int k=0; k<n_face_points; ++k)
541  for (unsigned int i=0; i<boundary_weights.size(1); ++i)
542  {
543  local_dofs[i+face*this->dofs_per_face] += boundary_weights(k,i)
544  * values[face*n_face_points+k](GeometryInfo<dim>::unit_normal_direction[face]+offset);
545  }
546 
547  const unsigned int start_cell_dofs = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
548  const unsigned int start_cell_points = GeometryInfo<dim>::faces_per_cell*n_face_points;
549 
550  for (unsigned int k=0; k<interior_weights.size(0); ++k)
551  for (unsigned int i=0; i<interior_weights.size(1); ++i)
552  for (unsigned int d=0; d<dim; ++d)
553  local_dofs[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * values[k+start_cell_points](d+offset);
554 
555  const unsigned int start_abf_dofs = start_cell_dofs + interior_weights.size(1) * dim;
556 
557  // Cell integral of ABF terms
558  for (unsigned int k=0; k<interior_weights_abf.size(0); ++k)
559  for (unsigned int i=0; i<interior_weights_abf.size(1); ++i)
560  for (unsigned int d=0; d<dim; ++d)
561  local_dofs[start_abf_dofs+i] += interior_weights_abf(k,i,d) * values[k+start_cell_points](d+offset);
562 
563  // Face integral of ABF terms
564  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
565  {
566  double n_orient = (double) GeometryInfo<dim>::unit_normal_orientation[face];
567  for (unsigned int fp=0; fp < n_face_points; ++fp)
568  {
569  // TODO: Check what the face_orientation, face_flip and face_rotation have to be in 3D
570  unsigned int k = QProjector<dim>::DataSetDescriptor::face (face, false, false, false, n_face_points);
571  for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
572  local_dofs[start_abf_dofs+i] += n_orient * boundary_weights_abf(k + fp, i)
573  * values[k + fp](GeometryInfo<dim>::unit_normal_direction[face]+offset);
574  }
575  }
576 
577  // TODO: Check if this "correction" can be removed.
578  for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
579  if (std::fabs (local_dofs[start_abf_dofs+i]) < 1.0e-16)
580  local_dofs[start_abf_dofs+i] = 0.0;
581 }
582 
583 template <int dim>
584 void
586  std::vector<double> &local_dofs,
587  const VectorSlice<const std::vector<std::vector<double> > > &values) const
588 {
589  Assert (values.size() == this->n_components(),
590  ExcDimensionMismatch(values.size(), this->n_components()));
591  Assert (values[0].size() == this->generalized_support_points.size(),
592  ExcDimensionMismatch(values[0].size(), this->generalized_support_points.size()));
593  Assert (local_dofs.size() == this->dofs_per_cell,
594  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
595 
596  std::fill(local_dofs.begin(), local_dofs.end(), 0.);
597 
598  const unsigned int n_face_points = boundary_weights.size(0);
599  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
600  for (unsigned int k=0; k<n_face_points; ++k)
601  for (unsigned int i=0; i<boundary_weights.size(1); ++i)
602  {
603  local_dofs[i+face*this->dofs_per_face] += boundary_weights(k,i)
604  * values[GeometryInfo<dim>::unit_normal_direction[face]][face*n_face_points+k];
605  }
606 
607  const unsigned int start_cell_dofs = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
608  const unsigned int start_cell_points = GeometryInfo<dim>::faces_per_cell*n_face_points;
609 
610  for (unsigned int k=0; k<interior_weights.size(0); ++k)
611  for (unsigned int i=0; i<interior_weights.size(1); ++i)
612  for (unsigned int d=0; d<dim; ++d)
613  local_dofs[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * values[d][k+start_cell_points];
614 
615  const unsigned int start_abf_dofs = start_cell_dofs + interior_weights.size(1) * dim;
616 
617  // Cell integral of ABF terms
618  for (unsigned int k=0; k<interior_weights_abf.size(0); ++k)
619  for (unsigned int i=0; i<interior_weights_abf.size(1); ++i)
620  for (unsigned int d=0; d<dim; ++d)
621  local_dofs[start_abf_dofs+i] += interior_weights_abf(k,i,d) * values[d][k+start_cell_points];
622 
623  // Face integral of ABF terms
624  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
625  {
626  double n_orient = (double) GeometryInfo<dim>::unit_normal_orientation[face];
627  for (unsigned int fp=0; fp < n_face_points; ++fp)
628  {
629  // TODO: Check what the face_orientation, face_flip and face_rotation have to be in 3D
630  unsigned int k = QProjector<dim>::DataSetDescriptor::face (face, false, false, false, n_face_points);
631  for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
632  local_dofs[start_abf_dofs+i] += n_orient * boundary_weights_abf(k + fp, i)
633  * values[GeometryInfo<dim>::unit_normal_direction[face]][k + fp];
634  }
635  }
636 
637  // TODO: Check if this "correction" can be removed.
638  for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
639  if (std::fabs (local_dofs[start_abf_dofs+i]) < 1.0e-16)
640  local_dofs[start_abf_dofs+i] = 0.0;
641 }
642 
643 
644 template <int dim>
645 std::size_t
647 {
648  Assert (false, ExcNotImplemented ());
649  return 0;
650 }
651 
652 
653 
654 /*-------------- Explicit Instantiations -------------------------------*/
655 #include "fe_abf.inst"
656 
657 DEAL_II_NAMESPACE_CLOSE
void initialize_restriction()
Definition: fe_abf.cc:306
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1003
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2106
std::vector< Point< dim > > generalized_support_points
Definition: fe.h:2157
Table< 3, double > interior_weights
Definition: fe_abf.h:222
FullMatrix< double > interface_constraints
Definition: fe.h:2132
Table< 2, double > boundary_weights
Definition: fe_abf.h:215
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:1009
virtual std::string get_name() const
Definition: fe_abf.cc:110
const Point< dim > & point(const unsigned int i) const
virtual std::size_t memory_consumption() const
Definition: fe_abf.cc:646
friend class FE_ABF
Definition: fe_abf.h:245
void invert(const FullMatrix< number2 > &M)
STL namespace.
void compute_embedding_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
Definition: fe_tools.cc:891
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_abf.cc:470
Definition: fe_abf.h:103
Table< 2, double > boundary_weights_abf
Definition: fe_abf.h:232
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
Definition: fe_abf.cc:513
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2120
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim-1 > &face_refinement_case=RefinementCase< dim-1 >::isotropic_refinement)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
Table< 3, double > interior_weights_abf
Definition: fe_abf.h:239
#define Assert(cond, exc)
Definition: exceptions.h:294
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim-1 > &ref_case=RefinementCase< dim-1 >::isotropic_refinement)
static std::vector< Polynomial< number > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:680
virtual FiniteElement< dim > * clone() const
Definition: fe_abf.cc:130
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:283
std::vector< Point< dim-1 > > generalized_face_support_points
Definition: fe.h:2163
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
std::vector< Tensor< 1, dim > > cached_values
const unsigned int rt_order
Definition: fe_abf.h:145
unsigned int n_components() const
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
const unsigned int dofs_per_face
Definition: fe_base.h:276
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_abf.cc:436
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:278
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
unsigned int size(const unsigned int i) const
void initialize_support_points(const unsigned int rt_degree)
Definition: fe_abf.cc:145
void compute_node_matrix(FullMatrix< double > &M, const FiniteElement< dim, spacedim > &fe)
Definition: fe_tools.cc:633
double weight(const unsigned int i) const
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: quadrature.cc:1081