/usr/include/deal.II/fe/fe_dgp_nonparametric.h is in libdeal.ii-dev 8.4.2-2+b1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 | // ---------------------------------------------------------------------
//
// Copyright (C) 2002 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------
#ifndef dealii__fe_dgp_nonparametric_h
#define dealii__fe_dgp_nonparametric_h
#include <deal.II/base/config.h>
#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/fe/fe.h>
DEAL_II_NAMESPACE_OPEN
template <int dim> class PolynomialSpace;
template <int dim, int spacedim> class MappingQ;
/*!@addtogroup fe */
/*@{*/
/**
* Discontinuous finite elements evaluated at the mapped quadrature points.
*
* Warning: this class does not work properly, yet. Don't use it!
*
* This finite element implements complete polynomial spaces, that is,
* $d$-dimensional polynomials of order $k$.
*
* The polynomials are not mapped. Therefore, they are constant, linear,
* quadratic, etc. on any grid cell.
*
* Since the polynomials are evaluated at the quadrature points of the actual
* grid cell, no grid transfer and interpolation matrices are available.
*
* The purpose of this class is experimental, therefore the implementation
* will remain incomplete.
*
* Besides, this class is not implemented for the codimension one case
* (<tt>spacedim != dim</tt>).
*
*
* <h3>Visualization of shape functions</h3> In 2d, the shape functions of
* this element look as follows.
*
* <h4>$P_0$ element</h4>
*
* <table> <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P1/P1_DGPNonparametric_shape0000.png
* </td>
*
* <td align="center"> </td> </tr> <tr> <td align="center"> $P_0$ element,
* shape function 0 </td>
*
* <td align="center"></tr> </table>
*
* <h4>$P_1$ element</h4>
*
* <table> <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P1/P1_DGPNonparametric_shape0000.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P1/P1_DGPNonparametric_shape0001.png
* </td> </tr> <tr> <td align="center"> $P_1$ element, shape function 0 </td>
*
* <td align="center"> $P_1$ element, shape function 1 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P1/P1_DGPNonparametric_shape0002.png
* </td>
*
* <td align="center"> </td> </tr> <tr> <td align="center"> $P_1$ element,
* shape function 2 </td>
*
* <td align="center"></td> </tr> </table>
*
*
* <h4>$P_2$ element</h4>
*
* <table> <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P2/P2_DGPNonparametric_shape0000.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P2/P2_DGPNonparametric_shape0001.png
* </td> </tr> <tr> <td align="center"> $P_2$ element, shape function 0 </td>
*
* <td align="center"> $P_2$ element, shape function 1 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P2/P2_DGPNonparametric_shape0002.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P2/P2_DGPNonparametric_shape0003.png
* </td> </tr> <tr> <td align="center"> $P_2$ element, shape function 2 </td>
*
* <td align="center"> $P_2$ element, shape function 3 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P2/P2_DGPNonparametric_shape0004.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P2/P2_DGPNonparametric_shape0005.png
* </td> </tr> <tr> <td align="center"> $P_2$ element, shape function 4 </td>
*
* <td align="center"> $P_2$ element, shape function 5 </td> </tr> </table>
*
*
* <h4>$P_3$ element</h4>
*
* <table> <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0000.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0001.png
* </td> </tr> <tr> <td align="center"> $P_3$ element, shape function 0 </td>
*
* <td align="center"> $P_3$ element, shape function 1 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0002.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0003.png
* </td> </tr> <tr> <td align="center"> $P_3$ element, shape function 2 </td>
*
* <td align="center"> $P_3$ element, shape function 3 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0004.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0005.png
* </td> </tr> <tr> <td align="center"> $P_3$ element, shape function 4 </td>
*
* <td align="center"> $P_3$ element, shape function 5 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0006.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0007.png
* </td> </tr> <tr> <td align="center"> $P_3$ element, shape function 6 </td>
*
* <td align="center"> $P_3$ element, shape function 7 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0008.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P3/P3_DGPNonparametric_shape0009.png
* </td> </tr> <tr> <td align="center"> $P_3$ element, shape function 8 </td>
*
* <td align="center"> $P_3$ element, shape function 9 </td> </tr> </table>
*
*
* <h4>$P_4$ element</h4> <table> <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0000.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0001.png
* </td> </tr> <tr> <td align="center"> $P_4$ element, shape function 0 </td>
*
* <td align="center"> $P_4$ element, shape function 1 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0002.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0003.png
* </td> </tr> <tr> <td align="center"> $P_4$ element, shape function 2 </td>
*
* <td align="center"> $P_4$ element, shape function 3 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0004.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0005.png
* </td> </tr> <tr> <td align="center"> $P_4$ element, shape function 4 </td>
*
* <td align="center"> $P_4$ element, shape function 5 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0006.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0007.png
* </td> </tr> <tr> <td align="center"> $P_4$ element, shape function 6 </td>
*
* <td align="center"> $P_4$ element, shape function 7 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0008.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0009.png
* </td> </tr> <tr> <td align="center"> $P_4$ element, shape function 8 </td>
*
* <td align="center"> $P_4$ element, shape function 9 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0010.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0011.png
* </td> </tr> <tr> <td align="center"> $P_4$ element, shape function 10 </td>
*
* <td align="center"> $P_4$ element, shape function 11 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0012.png
* </td>
*
* <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0013.png
* </td> </tr> <tr> <td align="center"> $P_4$ element, shape function 12 </td>
*
* <td align="center"> $P_4$ element, shape function 13 </td> </tr>
*
* <tr> <td align="center">
* @image html http://www.dealii.org/images/shape-functions/DGPNonparametric/P4/P4_DGPNonparametric_shape0014.png
* </td>
*
* <td align="center"> </td> </tr> <tr> <td align="center"> $P_4$ element,
* shape function 14 </td>
*
* <td align="center"></td> </tr> </table>
*
*
* <h3> Implementation details </h3>
*
* This element does not have an InternalData class, unlike all other
* elements, because the InternalData classes are used to store things that
* can be computed once and reused multiple times (such as the values of shape
* functions at quadrature points on the reference cell). However, because the
* element is not mapped, this element has nothing that could be computed on
* the reference cell -- everything needs to be computed on the real cell --
* and consequently there is nothing we'd like to store in such an object. We
* can thus simply use the members already provided by
* FiniteElement::InternalDataBase without adding anything in a derived class
* in this class.
*
* @author Guido Kanschat, 2002
*/
template <int dim, int spacedim=dim>
class FE_DGPNonparametric : public FiniteElement<dim,spacedim>
{
public:
/**
* Constructor for tensor product polynomials of degree @p k.
*/
FE_DGPNonparametric (const unsigned int k);
/**
* Return a string that uniquely identifies a finite element. This class
* returns <tt>FE_DGPNonparametric<dim>(degree)</tt>, with @p dim and @p
* degree replaced by appropriate values.
*/
virtual std::string get_name () const;
// for documentation, see the FiniteElement base class
virtual
UpdateFlags
requires_update_flags (const UpdateFlags update_flags) const;
/**
* This function is intended to return the value of a shape function at a
* point on the reference cell. However, since the current element does not
* implement shape functions by mapping from a reference cell, no shape
* functions exist on the reference cell.
*
* Consequently, as discussed in the corresponding function in the base
* class, FiniteElement::shape_value(), this function throws an exception of
* type FiniteElement::ExcUnitShapeValuesDoNotExist.
*/
virtual double shape_value (const unsigned int i,
const Point<dim> &p) const;
/**
* This function is intended to return the value of a shape function at a
* point on the reference cell. However, since the current element does not
* implement shape functions by mapping from a reference cell, no shape
* functions exist on the reference cell.
*
* Consequently, as discussed in the corresponding function in the base
* class, FiniteElement::shape_value_component(), this function throws an
* exception of type FiniteElement::ExcUnitShapeValuesDoNotExist.
*/
virtual double shape_value_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
/**
* This function is intended to return the gradient of a shape function at a
* point on the reference cell. However, since the current element does not
* implement shape functions by mapping from a reference cell, no shape
* functions exist on the reference cell.
*
* Consequently, as discussed in the corresponding function in the base
* class, FiniteElement::shape_grad(), this function throws an exception of
* type FiniteElement::ExcUnitShapeValuesDoNotExist.
*/
virtual Tensor<1,dim> shape_grad (const unsigned int i,
const Point<dim> &p) const;
/**
* This function is intended to return the gradient of a shape function at a
* point on the reference cell. However, since the current element does not
* implement shape functions by mapping from a reference cell, no shape
* functions exist on the reference cell.
*
* Consequently, as discussed in the corresponding function in the base
* class, FiniteElement::shape_grad_component(), this function throws an
* exception of type FiniteElement::ExcUnitShapeValuesDoNotExist.
*/
virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
/**
* This function is intended to return the Hessian of a shape function at a
* point on the reference cell. However, since the current element does not
* implement shape functions by mapping from a reference cell, no shape
* functions exist on the reference cell.
*
* Consequently, as discussed in the corresponding function in the base
* class, FiniteElement::shape_grad_grad(), this function throws an
* exception of type FiniteElement::ExcUnitShapeValuesDoNotExist.
*/
virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
const Point<dim> &p) const;
/**
* This function is intended to return the Hessian of a shape function at a
* point on the reference cell. However, since the current element does not
* implement shape functions by mapping from a reference cell, no shape
* functions exist on the reference cell.
*
* Consequently, as discussed in the corresponding function in the base
* class, FiniteElement::shape_grad_grad_component(), this function throws
* an exception of type FiniteElement::ExcUnitShapeValuesDoNotExist.
*/
virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
/**
* Return the polynomial degree of this finite element, i.e. the value
* passed to the constructor.
*/
unsigned int get_degree () const;
/**
* Return the matrix interpolating from a face of of one element to the face
* of the neighboring element. The size of the matrix is then
* <tt>source.dofs_per_face</tt> times <tt>this->dofs_per_face</tt>.
*
* Derived elements will have to implement this function. They may only
* provide interpolation matrices for certain source finite elements, for
* example those from the same family. If they don't implement interpolation
* from a given element, then they must throw an exception of type
* FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented.
*/
virtual void
get_face_interpolation_matrix (const FiniteElement<dim,spacedim> &source,
FullMatrix<double> &matrix) const;
/**
* Return the matrix interpolating from a face of of one element to the face
* of the neighboring element. The size of the matrix is then
* <tt>source.dofs_per_face</tt> times <tt>this->dofs_per_face</tt>.
*
* Derived elements will have to implement this function. They may only
* provide interpolation matrices for certain source finite elements, for
* example those from the same family. If they don't implement interpolation
* from a given element, then they must throw an exception of type
* FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented.
*/
virtual void
get_subface_interpolation_matrix (const FiniteElement<dim,spacedim> &source,
const unsigned int subface,
FullMatrix<double> &matrix) const;
/**
* @name Functions to support hp
* @{
*/
/**
* If, on a vertex, several finite elements are active, the hp code first
* assigns the degrees of freedom of each of these FEs different global
* indices. It then calls this function to find out which of them should get
* identical values, and consequently can receive the same global DoF index.
* This function therefore returns a list of identities between DoFs of the
* present finite element object with the DoFs of @p fe_other, which is a
* reference to a finite element object representing one of the other finite
* elements active on this particular vertex. The function computes which of
* the degrees of freedom of the two finite element objects are equivalent,
* both numbered between zero and the corresponding value of dofs_per_vertex
* of the two finite elements. The first index of each pair denotes one of
* the vertex dofs of the present element, whereas the second is the
* corresponding index of the other finite element.
*
* This being a discontinuous element, the set of such constraints is of
* course empty.
*/
virtual
std::vector<std::pair<unsigned int, unsigned int> >
hp_vertex_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Same as hp_vertex_dof_indices(), except that the function treats degrees
* of freedom on lines.
*
* This being a discontinuous element, the set of such constraints is of
* course empty.
*/
virtual
std::vector<std::pair<unsigned int, unsigned int> >
hp_line_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Same as hp_vertex_dof_indices(), except that the function treats degrees
* of freedom on quads.
*
* This being a discontinuous element, the set of such constraints is of
* course empty.
*/
virtual
std::vector<std::pair<unsigned int, unsigned int> >
hp_quad_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Return whether this element implements its hanging node constraints in
* the new way, which has to be used to make elements "hp compatible".
*
* For the FE_DGPNonparametric class the result is always true (independent
* of the degree of the element), as it has no hanging nodes (being a
* discontinuous element).
*/
virtual bool hp_constraints_are_implemented () const;
/**
* Return whether this element dominates the one given as argument when they
* meet at a common face, whether it is the other way around, whether
* neither dominates, or if either could dominate.
*
* For a definition of domination, see FiniteElementBase::Domination and in
* particular the
* @ref hp_paper "hp paper".
*/
virtual
FiniteElementDomination::Domination
compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* @}
*/
/**
* This function returns @p true, if the shape function @p shape_index has
* non-zero function values somewhere on the face @p face_index.
*/
virtual bool has_support_on_face (const unsigned int shape_index,
const unsigned int face_index) const;
/**
* Determine an estimate for the memory consumption (in bytes) of this
* object.
*
* This function is made virtual, since finite element objects are usually
* accessed through pointers to their base class, rather than the class
* itself.
*/
virtual std::size_t memory_consumption () const;
private:
/**
* Declare a nested class which will hold static definitions of various
* matrices such as constraint and embedding matrices. The definition of the
* various static fields are in the files <tt>fe_dgp_[123]d.cc</tt> in the
* source directory.
*/
struct Matrices
{
/**
* Pointers to the embedding matrices, one for each polynomial degree
* starting from constant elements
*/
static const double *const embedding[][GeometryInfo<dim>::max_children_per_cell];
/**
* Number of elements (first index) the above field has. Equals the
* highest polynomial degree plus one for which the embedding matrices
* have been computed.
*/
static const unsigned int n_embedding_matrices;
/**
* As @p embedding but for projection matrices.
*/
static const double *const projection_matrices[][GeometryInfo<dim>::max_children_per_cell];
/**
* As @p n_embedding_matrices but for projection matrices.
*/
static const unsigned int n_projection_matrices;
};
protected:
/**
* @p clone function instead of a copy constructor.
*
* This function is needed by the constructors of @p FESystem.
*/
virtual FiniteElement<dim,spacedim> *clone() const;
/**
* Prepare internal data structures and fill in values independent of the
* cell.
*/
virtual
typename FiniteElement<dim,spacedim>::InternalDataBase *
get_data (const UpdateFlags update_flags,
const Mapping<dim,spacedim> &mapping,
const Quadrature<dim> &quadrature,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;
virtual
void
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const CellSimilarity::Similarity cell_similarity,
const Quadrature<dim> &quadrature,
const Mapping<dim,spacedim> &mapping,
const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;
virtual
void
fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const Quadrature<dim-1> &quadrature,
const Mapping<dim,spacedim> &mapping,
const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;
virtual
void
fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const unsigned int sub_no,
const Quadrature<dim-1> &quadrature,
const Mapping<dim,spacedim> &mapping,
const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;
private:
/**
* Only for internal use. Its full name is @p get_dofs_per_object_vector
* function and it creates the @p dofs_per_object vector that is needed
* within the constructor to be passed to the constructor of @p
* FiniteElementData.
*/
static
std::vector<unsigned int>
get_dpo_vector (const unsigned int degree);
/**
* Degree of the polynomials.
*/
const unsigned int degree;
/**
* Pointer to an object representing the polynomial space used here.
*/
const PolynomialSpace<dim> polynomial_space;
/**
* Allow access from other dimensions.
*/
template <int, int> friend class FE_DGPNonparametric;
};
/*@}*/
#ifndef DOXYGEN
// declaration of explicit specializations of member variables, if the
// compiler allows us to do that (the standard says we must)
#ifndef DEAL_II_MEMBER_VAR_SPECIALIZATION_BUG
template <>
const double *const FE_DGPNonparametric<1,1>::Matrices::embedding[][GeometryInfo<1>::max_children_per_cell];
template <>
const unsigned int FE_DGPNonparametric<1,1>::Matrices::n_embedding_matrices;
template <>
const double *const FE_DGPNonparametric<1,1>::Matrices::projection_matrices[][GeometryInfo<1>::max_children_per_cell];
template <>
const unsigned int FE_DGPNonparametric<1,1>::Matrices::n_projection_matrices;
template <>
const double *const FE_DGPNonparametric<2,2>::Matrices::embedding[][GeometryInfo<2>::max_children_per_cell];
template <>
const unsigned int FE_DGPNonparametric<2,2>::Matrices::n_embedding_matrices;
template <>
const double *const FE_DGPNonparametric<2,2>::Matrices::projection_matrices[][GeometryInfo<2>::max_children_per_cell];
template <>
const unsigned int FE_DGPNonparametric<2,2>::Matrices::n_projection_matrices;
template <>
const double *const FE_DGPNonparametric<3,3>::Matrices::embedding[][GeometryInfo<3>::max_children_per_cell];
template <>
const unsigned int FE_DGPNonparametric<3,3>::Matrices::n_embedding_matrices;
template <>
const double *const FE_DGPNonparametric<3,3>::Matrices::projection_matrices[][GeometryInfo<3>::max_children_per_cell];
template <>
const unsigned int FE_DGPNonparametric<3,3>::Matrices::n_projection_matrices;
#endif
#endif // DOXYGEN
DEAL_II_NAMESPACE_CLOSE
#endif
|