/usr/include/deal.II/hp/fe_values.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 | // ---------------------------------------------------------------------
//
// Copyright (C) 2003 - 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__hp_fe_values_h
#define dealii__hp_fe_values_h
#include <deal.II/base/config.h>
#include <deal.II/fe/fe.h>
#include <deal.II/hp/fe_collection.h>
#include <deal.II/hp/q_collection.h>
#include <deal.II/hp/mapping_collection.h>
#include <deal.II/fe/fe_values.h>
#include <map>
#include <deal.II/base/std_cxx11/shared_ptr.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim> class FiniteElement;
namespace internal
{
namespace hp
{
/**
* Base class for the <tt>hp::FE*Values</tt> classes, storing the data
* that is common to them. The main task of this class is to provide a
* table where for every combination of finite element, mapping, and
* quadrature object from their corresponding collection objects there is
* a matching ::FEValues, ::FEFaceValues, or ::FESubfaceValues object. To
* make things more efficient, however, these FE*Values objects are only
* created once requested (lazy allocation).
*
* The first template parameter denotes the space dimension we are in, the
* second the dimensionality of the object that we integrate on, i.e. for
* usual @p hp::FEValues it is equal to the first one, while for face
* integration it is one less. The third template parameter indicates the
* type of underlying non-hp FE*Values base type, i.e. it could either be
* dealii::FEValues, dealii::FEFaceValues, or dealii::FESubfaceValues.
*
* @ingroup hp
*
* @author Wolfgang Bangerth, 2003
*/
template <int dim, int q_dim, class FEValuesType>
class FEValuesBase
{
public:
/**
* Constructor. Set the fields of this class to the values indicated by
* the parameters to the constructor.
*/
FEValuesBase (const dealii::hp::MappingCollection<dim,FEValuesType::space_dimension> &mapping_collection,
const dealii::hp::FECollection<dim,FEValuesType::space_dimension> &fe_collection,
const dealii::hp::QCollection<q_dim> &q_collection,
const dealii::UpdateFlags update_flags);
/**
* Constructor. This constructor is equivalent to the other one except
* that it makes the object use a $Q_1$ mapping (i.e., an object of type
* MappingQGeneric(1)) implicitly.
*/
FEValuesBase (const dealii::hp::FECollection<dim,FEValuesType::space_dimension> &fe_collection,
const dealii::hp::QCollection<q_dim> &q_collection,
const UpdateFlags update_flags);
/**
* Get a reference to the collection of finite element objects used
* here.
*/
const dealii::hp::FECollection<dim,FEValuesType::space_dimension> &
get_fe_collection () const;
/**
* Get a reference to the collection of mapping objects used here.
*/
const dealii::hp::MappingCollection<dim,FEValuesType::space_dimension> &
get_mapping_collection () const;
/**
* Get a reference to the collection of quadrature objects used here.
*/
const dealii::hp::QCollection<q_dim> &
get_quadrature_collection () const;
/**
* Get the underlying update flags.
*/
UpdateFlags get_update_flags() const;
/**
* Return a reference to the @p FEValues object selected by the last
* call to select_fe_values(). select_fe_values() in turn is called when
* you called the @p reinit function of the <tt>hp::FE*Values</tt> class
* the last time.
*/
const FEValuesType &get_present_fe_values () const;
protected:
/**
* Select a FEValues object suitable for the given FE, quadrature, and
* mapping indices. If such an object doesn't yet exist, create one.
*
* The function returns a writable reference so that derived classes can
* also reinit() the selected FEValues object.
*/
FEValuesType &
select_fe_values (const unsigned int fe_index,
const unsigned int mapping_index,
const unsigned int q_index);
protected:
/**
* A pointer to the collection of finite elements to be used.
*/
const SmartPointer<const dealii::hp::FECollection<dim,FEValuesType::space_dimension>,
FEValuesBase<dim,q_dim,FEValuesType> > fe_collection;
/**
* A pointer to the collection of mappings to be used.
*/
const SmartPointer<const dealii::hp::MappingCollection<dim, FEValuesType::space_dimension>,
FEValuesBase<dim,q_dim,FEValuesType> > mapping_collection;
/**
* Copy of the quadrature collection object provided to the constructor.
*/
const dealii::hp::QCollection<q_dim> q_collection;
private:
/**
* A table in which we store pointers to fe_values objects for different
* finite element, mapping, and quadrature objects from our collection.
* The first index indicates the index of the finite element within the
* fe_collection, the second the index of the mapping within the mapping
* collection, and the last one the index of the quadrature formula
* within the q_collection.
*
* Initially, all entries have zero pointers, and we will allocate them
* lazily as needed in select_fe_values().
*/
dealii::Table<3,std_cxx11::shared_ptr<FEValuesType> > fe_values_table;
/**
* Set of indices pointing at the fe_values object selected last time
* the select_fe_value() function was called.
*/
TableIndices<3> present_fe_values_index;
/**
* Values of the update flags as given to the constructor.
*/
const UpdateFlags update_flags;
};
}
}
namespace hp
{
/**
* An hp equivalent of the ::FEValues class. See the step-27 tutorial
* program for examples of use.
*
* The idea of this class is as follows: when one assembled matrices in the
* hp finite element method, there may be different finite elements on
* different cells, and consequently one may also want to use different
* quadrature formulas for different cells. On the other hand, the
* ::FEValues efficiently handles pre-evaluating whatever information is
* necessary for a single finite element and quadrature object. This class
* brings these concepts together: it provides a "collection" of ::FEValues
* objects.
*
* Upon construction, one passes not one finite element and quadrature
* object (and possible a mapping), but a whole collection of type
* hp::FECollection and hp::QCollection. Later on, when one sits on a
* concrete cell, one would call the reinit() function for this particular
* cell, just as one does for a regular ::FEValues object. The difference is
* that this time, the reinit() function looks up the active_fe_index of
* that cell, if necessary creates a ::FEValues object that matches the
* finite element and quadrature formulas with that particular index in
* their collections, and then re-initializes it for the current cell. The
* ::FEValues object that then fits the finite element and quadrature
* formula for the current cell can then be accessed using the
* get_present_fe_values() function, and one would work with it just like
* with any ::FEValues object for non-hp DoF handler objects.
*
* The reinit() functions have additional arguments with default values. If
* not specified, the function takes the index into the hp::FECollection,
* hp::QCollection, and hp::MappingCollection objects from the
* active_fe_index of the cell, as explained above. However, one can also
* select different indices for a current cell. For example, by specifying a
* different index into the hp::QCollection class, one does not need to sort
* the quadrature objects in the quadrature collection so that they match
* one-to-one the order of finite element objects in the FE collection (even
* though choosing such an order is certainly convenient).
*
* Note that ::FEValues objects are created on the fly, i.e. only as they
* are needed. This ensures that we do not create objects for every
* combination of finite element, quadrature formula and mapping, but only
* those that will actually be needed.
*
* This class has not yet been implemented for the use in the codimension
* one case (<tt>spacedim != dim </tt>).
*
* @ingroup hp hpcollection
* @author Wolfgang Bangerth, 2003
*/
template <int dim, int spacedim=dim>
class FEValues : public dealii::internal::hp::FEValuesBase<dim,dim,dealii::FEValues<dim,spacedim> >
{
public:
static const unsigned int dimension = dim;
static const unsigned int space_dimension = spacedim;
/**
* Constructor. Initialize this object with the given parameters.
*
* The finite element collection parameter is actually ignored, but is in
* the signature of this function to make it compatible with the signature
* of the respective constructor of the usual FEValues object, with the
* respective parameter in that function also being the return value of
* the DoFHandler::get_fe() function.
*/
FEValues (const dealii::hp::MappingCollection<dim,spacedim> &mapping_collection,
const dealii::hp::FECollection<dim,spacedim> &fe_collection,
const dealii::hp::QCollection<dim> &q_collection,
const UpdateFlags update_flags);
/**
* Constructor. This constructor is equivalent to the other one except
* that it makes the object use a $Q_1$ mapping (i.e., an object of type
* MappingQGeneric(1)) implicitly.
*
* The finite element collection parameter is actually ignored, but is in
* the signature of this function to make it compatible with the signature
* of the respective constructor of the usual FEValues object, with the
* respective parameter in that function also being the return value of
* the DoFHandler::get_fe() function.
*/
FEValues (const hp::FECollection<dim,spacedim> &fe_collection,
const hp::QCollection<dim> &q_collection,
const UpdateFlags update_flags);
/**
* Reinitialize the object for the given cell.
*
* After the call, you can get an FEValues object using the
* get_present_fe_values() function that corresponds to the present cell.
* For this FEValues object, we use the additional arguments described
* below to determine which finite element, mapping, and quadrature
* formula to use. They are order in such a way that the arguments one may
* want to change most frequently come first. The rules for these
* arguments are as follows:
*
* If the @p fe_index argument to this function is left at its default
* value, then we use that finite element within the hp::FECollection
* passed to the constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>. Consequently, the
* hp::FECollection argument given to this object should really be the
* same as that used in the construction of the hp::DofHandler associated
* with the present cell. On the other hand, if a value is given for this
* argument, it overrides the choice of
* <code>cell-@>active_fe_index()</code>.
*
* If the @p q_index argument is left at its default value, then we use
* that quadrature formula within the hp::QCollection passed to the
* constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>, i.e. the same index as that of
* the finite element. In this case, there should be a corresponding
* quadrature formula for each finite element in the hp::FECollection. As
* a special case, if the quadrature collection contains only a single
* element (a frequent case if one wants to use the same quadrature object
* for all finite elements in an hp discretization, even if that may not
* be the most efficient), then this single quadrature is used unless a
* different value for this argument is specified. On the other hand, if a
* value is given for this argument, it overrides the choice of
* <code>cell-@>active_fe_index()</code> or the choice for the single
* quadrature.
*
* If the @p mapping_index argument is left at its default value, then we
* use that mapping object within the hp::MappingCollection passed to the
* constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>, i.e. the same index as that of
* the finite element. As above, if the mapping collection contains only a
* single element (a frequent case if one wants to use a $Q_1$ mapping for
* all finite elements in an hp discretization), then this single mapping
* is used unless a different value for this argument is specified.
*/
template <typename DoFHandlerType, bool lda>
void
reinit (const TriaIterator<DoFCellAccessor<DoFHandlerType,lda> > cell,
const unsigned int q_index = numbers::invalid_unsigned_int,
const unsigned int mapping_index = numbers::invalid_unsigned_int,
const unsigned int fe_index = numbers::invalid_unsigned_int);
/**
* Like the previous function, but for non-hp iterators. The reason this
* (and the other non-hp iterator) function exists is so that one can use
* hp::FEValues not only for hp::DoFhandler objects, but for all sorts of
* DoFHandler objects, and triangulations not associated with DoFHandlers
* in general.
*
* Since <code>cell-@>active_fe_index()</code> doesn't make sense for
* triangulation iterators, this function chooses the zero-th finite
* element, mapping, and quadrature object from the relevant constructions
* passed to the constructor of this object. The only exception is if you
* specify a value different from the default value for any of these last
* three arguments.
*/
void
reinit (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int q_index = numbers::invalid_unsigned_int,
const unsigned int mapping_index = numbers::invalid_unsigned_int,
const unsigned int fe_index = numbers::invalid_unsigned_int);
};
/**
* This is the equivalent of the hp::FEValues class but for face
* integrations, i.e. it is to hp::FEValues what ::FEFaceValues is to
* ::FEValues.
*
* The same comments apply as in the documentation of the hp::FEValues
* class. However, it is important to note that it is here more common that
* one would want to explicitly specify an index to a particular quadrature
* formula in the reinit() functions. This is because the default index
* corresponds to the finite element index on the current function. On the
* other hand, integration on faces will typically have to happen with a
* quadrature formula that is adjusted to the finite elements used on both
* sides of a face. If one sorts the elements of the hp::FECollection with
* ascending polynomial degree, and matches these finite elements with
* corresponding quadrature formulas in the hp::QCollection passed to the
* constructor, then the quadrature index passed to the reinit() function
* should typically be something like <code>std::max
* (cell-@>active_fe_index(), neighbor-@>active_fe_index()</code> to ensure
* that a quadrature formula is chosen that is sufficiently accurate for
* <em>both</em> finite elements.
*
* @ingroup hp hpcollection
* @author Wolfgang Bangerth, 2003
*/
template <int dim, int spacedim=dim>
class FEFaceValues : public dealii::internal::hp::FEValuesBase<dim,dim-1,dealii::FEFaceValues<dim,spacedim> >
{
public:
/**
* Constructor. Initialize this object with the given parameters.
*
* The finite element collection parameter is actually ignored, but is in
* the signature of this function to make it compatible with the signature
* of the respective constructor of the usual FEValues object, with the
* respective parameter in that function also being the return value of
* the <tt>DoFHandler::get_fe()</tt> function.
*/
FEFaceValues (const hp::MappingCollection<dim,spacedim> &mapping_collection,
const hp::FECollection<dim,spacedim> &fe_collection,
const hp::QCollection<dim-1> &q_collection,
const UpdateFlags update_flags);
/**
* Constructor. This constructor is equivalent to the other one except
* that it makes the object use a $Q_1$ mapping (i.e., an object of type
* MappingQGeneric(1)) implicitly.
*
* The finite element collection parameter is actually ignored, but is in
* the signature of this function to make it compatible with the signature
* of the respective constructor of the usual FEValues object, with the
* respective parameter in that function also being the return value of
* the <tt>DoFHandler::get_fe()</tt> function.
*/
FEFaceValues (const hp::FECollection<dim,spacedim> &fe_collection,
const hp::QCollection<dim-1> &q_collection,
const UpdateFlags update_flags);
/**
* Reinitialize the object for the given cell and face.
*
* After the call, you can get an FEFaceValues object using the
* get_present_fe_values() function that corresponds to the present cell.
* For this FEFaceValues object, we use the additional arguments described
* below to determine which finite element, mapping, and quadrature
* formula to use. They are order in such a way that the arguments one may
* want to change most frequently come first. The rules for these
* arguments are as follows:
*
* If the @p fe_index argument to this function is left at its default
* value, then we use that finite element within the hp::FECollection
* passed to the constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>. Consequently, the
* hp::FECollection argument given to this object should really be the
* same as that used in the construction of the hp::DofHandler associated
* with the present cell. On the other hand, if a value is given for this
* argument, it overrides the choice of
* <code>cell-@>active_fe_index()</code>.
*
* If the @p q_index argument is left at its default value, then we use
* that quadrature formula within the hp::QCollection passed to the
* constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>, i.e. the same index as that of
* the finite element. In this case, there should be a corresponding
* quadrature formula for each finite element in the hp::FECollection. As
* a special case, if the quadrature collection contains only a single
* element (a frequent case if one wants to use the same quadrature object
* for all finite elements in an hp discretization, even if that may not
* be the most efficient), then this single quadrature is used unless a
* different value for this argument is specified. On the other hand, if a
* value is given for this argument, it overrides the choice of
* <code>cell-@>active_fe_index()</code> or the choice for the single
* quadrature.
*
* If the @p mapping_index argument is left at its default value, then we
* use that mapping object within the hp::MappingCollection passed to the
* constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>, i.e. the same index as that of
* the finite element. As above, if the mapping collection contains only a
* single element (a frequent case if one wants to use a $Q_1$ mapping for
* all finite elements in an hp discretization), then this single mapping
* is used unless a different value for this argument is specified.
*/
template <typename DoFHandlerType, bool lda>
void
reinit (const TriaIterator<DoFCellAccessor<DoFHandlerType,lda> > cell,
const unsigned int face_no,
const unsigned int q_index = numbers::invalid_unsigned_int,
const unsigned int mapping_index = numbers::invalid_unsigned_int,
const unsigned int fe_index = numbers::invalid_unsigned_int);
/**
* Like the previous function, but for non-hp iterators. The reason this
* (and the other non-hp iterator) function exists is so that one can use
* hp::FEValues not only for hp::DoFhandler objects, but for all sorts of
* DoFHandler objects, and triangulations not associated with DoFHandlers
* in general.
*
* Since <code>cell-@>active_fe_index()</code> doesn't make sense for
* triangulation iterators, this function chooses the zero-th finite
* element, mapping, and quadrature object from the relevant constructions
* passed to the constructor of this object. The only exception is if you
* specify a value different from the default value for any of these last
* three arguments.
*/
void
reinit (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const unsigned int q_index = numbers::invalid_unsigned_int,
const unsigned int mapping_index = numbers::invalid_unsigned_int,
const unsigned int fe_index = numbers::invalid_unsigned_int);
};
/**
* This class implements for subfaces what hp::FEFaceValues does for faces.
* See there for further documentation.
*
* @ingroup hp hpcollection
* @author Wolfgang Bangerth, 2003
*/
template <int dim, int spacedim=dim>
class FESubfaceValues : public dealii::internal::hp::FEValuesBase<dim,dim-1,dealii::FESubfaceValues<dim,spacedim> >
{
public:
/**
* Constructor. Initialize this object with the given parameters.
*
* The finite element collection parameter is actually ignored, but is in
* the signature of this function to make it compatible with the signature
* of the respective constructor of the usual FEValues object, with the
* respective parameter in that function also being the return value of
* the <tt>DoFHandler::get_fe()</tt> function.
*/
FESubfaceValues (const hp::MappingCollection<dim,spacedim> &mapping_collection,
const hp::FECollection<dim,spacedim> &fe_collection,
const hp::QCollection<dim-1> &q_collection,
const UpdateFlags update_flags);
/**
* Constructor. This constructor is equivalent to the other one except
* that it makes the object use a $Q_1$ mapping (i.e., an object of type
* MappingQGeneric(1)) implicitly.
*
* The finite element collection parameter is actually ignored, but is in
* the signature of this function to make it compatible with the signature
* of the respective constructor of the usual FEValues object, with the
* respective parameter in that function also being the return value of
* the <tt>DoFHandler::get_fe()</tt> function.
*/
FESubfaceValues (const hp::FECollection<dim,spacedim> &fe_collection,
const hp::QCollection<dim-1> &q_collection,
const UpdateFlags update_flags);
/**
* Reinitialize the object for the given cell, face, and subface.
*
* After the call, you can get an FESubfaceValues object using the
* get_present_fe_values() function that corresponds to the present cell.
* For this FESubfaceValues object, we use the additional arguments
* described below to determine which finite element, mapping, and
* quadrature formula to use. They are order in such a way that the
* arguments one may want to change most frequently come first. The rules
* for these arguments are as follows:
*
* If the @p q_index argument is left at its default value, then we use
* that quadrature formula within the hp::QCollection passed to the
* constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>, i.e. the same index as that of
* the finite element. In this case, there should be a corresponding
* quadrature formula for each finite element in the hp::FECollection. As
* a special case, if the quadrature collection contains only a single
* element (a frequent case if one wants to use the same quadrature object
* for all finite elements in an hp discretization, even if that may not
* be the most efficient), then this single quadrature is used unless a
* different value for this argument is specified. On the other hand, if a
* value is given for this argument, it overrides the choice of
* <code>cell-@>active_fe_index()</code> or the choice for the single
* quadrature.
*
* If the @p mapping_index argument is left at its default value, then we
* use that mapping object within the hp::MappingCollection passed to the
* constructor of this class with index given by
* <code>cell-@>active_fe_index()</code>, i.e. the same index as that of
* the finite element. As above, if the mapping collection contains only a
* single element (a frequent case if one wants to use a $Q_1$ mapping for
* all finite elements in an hp discretization), then this single mapping
* is used unless a different value for this argument is specified.
*/
template <typename DoFHandlerType, bool lda>
void
reinit (const TriaIterator<DoFCellAccessor<DoFHandlerType,lda> > cell,
const unsigned int face_no,
const unsigned int subface_no,
const unsigned int q_index = numbers::invalid_unsigned_int,
const unsigned int mapping_index = numbers::invalid_unsigned_int,
const unsigned int fe_index = numbers::invalid_unsigned_int);
/**
* Like the previous function, but for non-hp iterators. The reason this
* (and the other non-hp iterator) function exists is so that one can use
* hp::FEValues not only for hp::DoFhandler objects, but for all sorts of
* DoFHandler objects, and triangulations not associated with DoFHandlers
* in general.
*
* Since <code>cell-@>active_fe_index()</code> doesn't make sense for
* triangulation iterators, this function chooses the zero-th finite
* element, mapping, and quadrature object from the relevant constructions
* passed to the constructor of this object. The only exception is if you
* specify a value different from the default value for any of these last
* three arguments.
*/
void
reinit (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const unsigned int subface_no,
const unsigned int q_index = numbers::invalid_unsigned_int,
const unsigned int mapping_index = numbers::invalid_unsigned_int,
const unsigned int fe_index = numbers::invalid_unsigned_int);
};
}
// -------------- inline and template functions --------------
namespace internal
{
namespace hp
{
template <int dim, int q_dim, class FEValuesType>
inline
const FEValuesType &
FEValuesBase<dim,q_dim,FEValuesType>::get_present_fe_values () const
{
return *fe_values_table(present_fe_values_index);
}
template <int dim, int q_dim, class FEValuesType>
inline
const dealii::hp::FECollection<dim,FEValuesType::space_dimension> &
FEValuesBase<dim,q_dim,FEValuesType>::get_fe_collection () const
{
return *fe_collection;
}
template <int dim, int q_dim, class FEValuesType>
inline
const dealii::hp::MappingCollection<dim,FEValuesType::space_dimension> &
FEValuesBase<dim,q_dim,FEValuesType>::get_mapping_collection () const
{
return *mapping_collection;
}
template <int dim, int q_dim, class FEValuesType>
inline
const dealii::hp::QCollection<q_dim> &
FEValuesBase<dim,q_dim,FEValuesType>::get_quadrature_collection () const
{
return q_collection;
}
template <int dim, int q_dim, class FEValuesType>
inline
dealii::UpdateFlags
FEValuesBase<dim,q_dim,FEValuesType>::get_update_flags () const
{
return update_flags;
}
}
}
DEAL_II_NAMESPACE_CLOSE
#endif
|