Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
scratch_data.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2019 - 2022 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.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_meshworker_scratch_data_h
17#define dealii_meshworker_scratch_data_h
18
19#include <deal.II/base/config.h>
20
22
24
27
31
32#include <boost/any.hpp>
33
34#include <algorithm>
35#include <map>
36#include <vector>
37
39
40namespace MeshWorker
41{
222 template <int dim, int spacedim = dim>
224 {
225 public:
248 const Quadrature<dim> & quadrature,
249 const UpdateFlags & update_flags,
252
272 const Quadrature<dim> & quadrature,
273 const UpdateFlags & update_flags,
274 const UpdateFlags & neighbor_update_flags,
278
293 const Quadrature<dim> & quadrature,
294 const UpdateFlags & update_flags,
297
314 const Quadrature<dim> & quadrature,
315 const UpdateFlags & update_flags,
316 const UpdateFlags & neighbor_update_flags,
320
347
373
392
415
420 // CurrentCellMethods
425
434 reinit(
436
446 const unsigned int face_no);
447
460 const unsigned int face_no,
461 const unsigned int subface_no);
462
479 const unsigned int face_no,
480 const unsigned int sub_face_no,
482 & cell_neighbor,
483 const unsigned int face_no_neighbor,
484 const unsigned int sub_face_no_neighbor);
485
497 get_current_fe_values() const;
498
504
508 const std::vector<Point<spacedim>> &
509 get_quadrature_points() const;
510
514 const std::vector<double> &
515 get_JxW_values() const;
516
520 const std::vector<Tensor<1, spacedim>> &
521 get_normal_vectors() const;
522
527 const std::vector<types::global_dof_index> &
528 get_local_dof_indices() const;
529
534 unsigned int
535 n_dofs_per_cell() const;
536 // CurrentCellMethods
538 // NeighborCellMethods
543
554
565 const unsigned int face_no);
566
580 const unsigned int face_no,
581 const unsigned int subface_no);
582
596
600 const std::vector<double> &
602
606 const std::vector<Tensor<1, spacedim>> &
608
613 const std::vector<types::global_dof_index> &
615
620 unsigned int
622 // NeighborCellMethods
624 // hpCellMethods
629
648 & neighbor_cell,
649 const unsigned int face_no);
650
672 & neighbor_cell,
673 const unsigned int face_no,
674 const unsigned int subface_no);
675
695 & neighbor_cell,
696 const unsigned int face_no);
697
720 & neighbor_cell,
721 const unsigned int face_no,
722 const unsigned int subface_no);
723 // hpCellMethods
725
733
739 const GeneralDataStorage &
741
746 get_mapping() const;
747
752 get_fe() const;
753
757 const Quadrature<dim> &
758 get_cell_quadrature() const;
759
763 const Quadrature<dim - 1> &
764 get_face_quadrature() const;
765
771
776 get_fe_collection() const;
777
783
787 const hp::QCollection<dim - 1> &
789
794 bool
795 has_hp_capabilities() const;
796
801 get_cell_update_flags() const;
802
808
813 get_face_update_flags() const;
814
820 // CurrentCellEvaluation
825
855 template <typename VectorType, typename Number = double>
856 void
857 extract_local_dof_values(const std::string &global_vector_name,
858 const VectorType & input_vector,
859 const Number dummy = Number(0));
860
869 template <typename Number = double>
870 const std::vector<Number> &
871 get_local_dof_values(const std::string &global_vector_name,
872 Number dummy = Number(0)) const;
873
895 template <typename Extractor, typename Number = double>
896 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
897 template solution_value_type<Number>> &
898 get_values(const std::string &global_vector_name,
899 const Extractor & variable,
900 const Number dummy = Number(0));
901
923 template <typename Extractor, typename Number = double>
924 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
925 template solution_gradient_type<Number>> &
926 get_gradients(const std::string &global_vector_name,
927 const Extractor & variable,
928 const Number dummy = Number(0));
929
952 template <typename Extractor, typename Number = double>
953 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
954 template solution_symmetric_gradient_type<Number>> &
955 get_symmetric_gradients(const std::string &global_vector_name,
956 const Extractor & variable,
957 const Number dummy = Number(0));
958
980 template <typename Extractor, typename Number = double>
981 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
982 template solution_divergence_type<Number>> &
983 get_divergences(const std::string &global_vector_name,
984 const Extractor & variable,
985 const Number dummy = Number(0));
986
1008 template <typename Extractor, typename Number = double>
1009 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1010 template solution_curl_type<Number>> &
1011 get_curls(const std::string &global_vector_name,
1012 const Extractor & variable,
1013 const Number dummy = Number(0));
1014
1036 template <typename Extractor, typename Number = double>
1037 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1038 template solution_hessian_type<Number>> &
1039 get_hessians(const std::string &global_vector_name,
1040 const Extractor & variable,
1041 const Number dummy = Number(0));
1042
1064 template <typename Extractor, typename Number = double>
1065 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1066 template solution_laplacian_type<Number>> &
1067 get_laplacians(const std::string &global_vector_name,
1068 const Extractor & variable,
1069 const Number dummy = Number(0));
1070
1092 template <typename Extractor, typename Number = double>
1093 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1094 template solution_third_derivative_type<Number>> &
1095 get_third_derivatives(const std::string &global_vector_name,
1096 const Extractor & variable,
1097 const Number dummy = Number(0));
1098 // CurrentCellEvaluation
1100 // CurrentInterfaceJumpEvaluation
1105
1128 template <typename Extractor, typename Number = double>
1129 const std::vector<
1131 template solution_value_type<Number>> &
1132 get_jumps_in_values(const std::string &global_vector_name,
1133 const Extractor & variable,
1134 const Number dummy = Number(0));
1135
1158 template <typename Extractor, typename Number = double>
1159 const std::vector<
1161 template solution_gradient_type<Number>> &
1162 get_jumps_in_gradients(const std::string &global_vector_name,
1163 const Extractor & variable,
1164 const Number dummy = Number(0));
1165
1188 template <typename Extractor, typename Number = double>
1189 const std::vector<
1191 template solution_hessian_type<Number>> &
1192 get_jumps_in_hessians(const std::string &global_vector_name,
1193 const Extractor & variable,
1194 const Number dummy = Number(0));
1195
1218 template <typename Extractor, typename Number = double>
1219 const std::vector<
1221 template solution_third_derivative_type<Number>> &
1222 get_jumps_in_third_derivatives(const std::string &global_vector_name,
1223 const Extractor & variable,
1224 const Number dummy = Number(0));
1225 // CurrentInterfaceJumpEvaluation
1227 // CurrentInterfaceAverageEvaluation
1232
1255 template <typename Extractor, typename Number = double>
1256 const std::vector<
1258 template solution_value_type<Number>> &
1259 get_averages_of_values(const std::string &global_vector_name,
1260 const Extractor & variable,
1261 const Number dummy = Number(0));
1262
1285 template <typename Extractor, typename Number = double>
1286 const std::vector<
1288 template solution_gradient_type<Number>> &
1289 get_averages_of_gradients(const std::string &global_vector_name,
1290 const Extractor & variable,
1291 const Number dummy = Number(0));
1292
1315 template <typename Extractor, typename Number = double>
1316 const std::vector<
1318 template solution_hessian_type<Number>> &
1319 get_averages_of_hessians(const std::string &global_vector_name,
1320 const Extractor & variable,
1321 const Number dummy = Number(0));
1322 // CurrentInterfaceAverageEvaluation
1324
1325 private:
1330 template <typename Extractor, typename Number = double>
1331 std::string
1332 get_unique_name(const std::string &global_vector_name,
1333 const Extractor & variable,
1334 const std::string &object_type,
1335 const unsigned int size,
1336 const Number & exemplar_number) const;
1337
1341 template <typename Number = double>
1342 std::string
1343 get_unique_dofs_name(const std::string &global_vector_name,
1344 const unsigned int size,
1345 const Number & exemplar_number) const;
1346 // non-hp data
1351
1357
1363
1369
1375
1379 std::unique_ptr<FEValues<dim, spacedim>> fe_values;
1380
1384 std::unique_ptr<FEFaceValues<dim, spacedim>> fe_face_values;
1385
1389 std::unique_ptr<FESubfaceValues<dim, spacedim>> fe_subface_values;
1390
1394 std::unique_ptr<FEValues<dim, spacedim>> neighbor_fe_values;
1395
1399 std::unique_ptr<FEFaceValues<dim, spacedim>> neighbor_fe_face_values;
1400
1404 std::unique_ptr<FESubfaceValues<dim, spacedim>> neighbor_fe_subface_values;
1405
1409 std::unique_ptr<FEInterfaceValues<dim, spacedim>> interface_fe_values;
1410 // non-hp data
1412 // hp data
1417
1423
1429
1435
1441
1447
1451 std::unique_ptr<hp::FEValues<dim, spacedim>> hp_fe_values;
1452
1456 std::unique_ptr<hp::FEFaceValues<dim, spacedim>> hp_fe_face_values;
1457
1461 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>> hp_fe_subface_values;
1462
1466 std::unique_ptr<hp::FEValues<dim, spacedim>> neighbor_hp_fe_values;
1467
1471 std::unique_ptr<hp::FEFaceValues<dim, spacedim>> neighbor_hp_fe_face_values;
1472
1476 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>>
1478
1486 "The current function doesn't make sense when used with a "
1487 "ScratchData object with hp-capabilities.");
1488
1496 "The current function doesn't make sense when used with a "
1497 "ScratchData object without hp-capabilities.");
1498 // hp data
1500
1505
1510
1516
1522
1526 std::vector<types::global_dof_index> local_dof_indices;
1527
1531 std::vector<types::global_dof_index> neighbor_dof_indices;
1532
1537
1542
1548
1554 };
1555
1556#ifndef DOXYGEN
1557 template <int dim, int spacedim>
1558 template <typename Extractor, typename Number>
1559 std::string
1561 const std::string &global_vector_name,
1562 const Extractor & variable,
1563 const std::string &object_type,
1564 const unsigned int size,
1565 const Number & exemplar_number) const
1566 {
1567 return global_vector_name + "_" + variable.get_name() + "_" + object_type +
1568 "_" + Utilities::int_to_string(size) + "_" +
1569 Utilities::type_to_string(exemplar_number);
1570 }
1571
1572
1573
1574 template <int dim, int spacedim>
1575 template <typename Number>
1576 std::string
1578 const std::string &global_vector_name,
1579 const unsigned int size,
1580 const Number & exemplar_number) const
1581 {
1582 return global_vector_name + "_independent_local_dofs_" +
1583 Utilities::int_to_string(size) + "_" +
1584 Utilities::type_to_string(exemplar_number);
1585 }
1586
1587
1588
1589 template <int dim, int spacedim>
1590 template <typename VectorType, typename Number>
1591 void
1593 const std::string &global_vector_name,
1594 const VectorType & input_vector,
1595 const Number dummy)
1596 {
1597 const unsigned int n_dofs = local_dof_indices.size();
1598
1599 const std::string name =
1600 get_unique_dofs_name(global_vector_name, n_dofs, dummy);
1601
1602 auto &independent_local_dofs =
1603 internal_data_storage
1604 .template get_or_add_object_with_name<std::vector<Number>>(name,
1605 n_dofs);
1606
1607 AssertDimension(independent_local_dofs.size(), n_dofs);
1608
1610 for (unsigned int i = 0; i < n_dofs; ++i)
1612 input_vector(local_dof_indices[i]),
1613 i,
1614 n_dofs,
1615 independent_local_dofs[i]);
1616 else
1617 for (unsigned int i = 0; i < n_dofs; ++i)
1618 independent_local_dofs[i] = input_vector(local_dof_indices[i]);
1619 }
1620
1621
1622
1623 template <int dim, int spacedim>
1624 template <typename Number>
1625 const std::vector<Number> &
1627 const std::string &global_vector_name,
1628 Number dummy) const
1629 {
1630 const unsigned int n_dofs = local_dof_indices.size();
1631
1632 const std::string dofs_name =
1633 get_unique_dofs_name(global_vector_name, n_dofs, dummy);
1634
1635 Assert(
1636 internal_data_storage.stores_object_with_name(dofs_name),
1637 ExcMessage(
1638 "You did not call yet extract_local_dof_values with the right types!"));
1639
1640 return internal_data_storage
1641 .template get_object_with_name<std::vector<Number>>(dofs_name);
1642 }
1643
1644
1645
1646 template <int dim, int spacedim>
1647 template <typename Extractor, typename Number>
1648 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1649 template solution_value_type<Number>> &
1650 ScratchData<dim, spacedim>::get_values(const std::string &global_vector_name,
1651 const Extractor & variable,
1652 const Number dummy)
1653 {
1654 const std::vector<Number> &independent_local_dofs =
1655 get_local_dof_values(global_vector_name, dummy);
1656
1657 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1658
1659 const unsigned int n_q_points = fev.n_quadrature_points;
1660
1661 const std::string name = get_unique_name(
1662 global_vector_name, variable, "_values_q", n_q_points, dummy);
1663
1664 // Now build the return type
1665 using RetType =
1666 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1667 template solution_value_type<Number>>;
1668
1669 RetType &ret =
1670 internal_data_storage.template get_or_add_object_with_name<RetType>(
1671 name, n_q_points);
1672
1673 AssertDimension(ret.size(), n_q_points);
1674
1675 fev[variable].get_function_values_from_local_dof_values(
1676 independent_local_dofs, ret);
1677 return ret;
1678 }
1679
1680
1681
1682 template <int dim, int spacedim>
1683 template <typename Extractor, typename Number>
1684 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1685 template solution_gradient_type<Number>> &
1687 const std::string &global_vector_name,
1688 const Extractor & variable,
1689 const Number dummy)
1690 {
1691 const std::vector<Number> &independent_local_dofs =
1692 get_local_dof_values(global_vector_name, dummy);
1693
1694 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1695
1696 const unsigned int n_q_points = fev.n_quadrature_points;
1697
1698 const std::string name = get_unique_name(
1699 global_vector_name, variable, "_gradients_q", n_q_points, dummy);
1700
1701 // Now build the return type
1702 using RetType =
1703 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1704 template solution_gradient_type<Number>>;
1705
1706 RetType &ret =
1707 internal_data_storage.template get_or_add_object_with_name<RetType>(
1708 name, n_q_points);
1709
1710 AssertDimension(ret.size(), n_q_points);
1711
1712 fev[variable].get_function_gradients_from_local_dof_values(
1713 independent_local_dofs, ret);
1714 return ret;
1715 }
1716
1717
1718
1719 template <int dim, int spacedim>
1720 template <typename Extractor, typename Number>
1721 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1722 template solution_hessian_type<Number>> &
1724 const std::string &global_vector_name,
1725 const Extractor & variable,
1726 const Number dummy)
1727 {
1728 const std::vector<Number> &independent_local_dofs =
1729 get_local_dof_values(global_vector_name, dummy);
1730
1731 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1732
1733 const unsigned int n_q_points = fev.n_quadrature_points;
1734
1735 const std::string name = get_unique_name(
1736 global_vector_name, variable, "_hessians_q", n_q_points, dummy);
1737
1738 // Now build the return type
1739 using RetType =
1740 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1741 template solution_hessian_type<Number>>;
1742
1743 RetType &ret =
1744 internal_data_storage.template get_or_add_object_with_name<RetType>(
1745 name, n_q_points);
1746
1747
1748 AssertDimension(ret.size(), n_q_points);
1749
1750 fev[variable].get_function_hessians_from_local_dof_values(
1751 independent_local_dofs, ret);
1752 return ret;
1753 }
1754
1755
1756
1757 template <int dim, int spacedim>
1758 template <typename Extractor, typename Number>
1759 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1760 template solution_laplacian_type<Number>> &
1762 const std::string &global_vector_name,
1763 const Extractor & variable,
1764 const Number dummy)
1765 {
1766 const std::vector<Number> &independent_local_dofs =
1767 get_local_dof_values(global_vector_name, dummy);
1768
1769 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1770
1771 const unsigned int n_q_points = fev.n_quadrature_points;
1772
1773 const std::string name = get_unique_name(
1774 global_vector_name, variable, "_laplacians_q", n_q_points, dummy);
1775
1776 // Now build the return type
1777 using RetType =
1778 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1779 template solution_laplacian_type<Number>>;
1780
1781 RetType &ret =
1782 internal_data_storage.template get_or_add_object_with_name<RetType>(
1783 name, n_q_points);
1784
1785
1786 AssertDimension(ret.size(), n_q_points);
1787
1788 fev[variable].get_function_laplacians_from_local_dof_values(
1789 independent_local_dofs, ret);
1790 return ret;
1791 }
1792
1793
1794
1795 template <int dim, int spacedim>
1796 template <typename Extractor, typename Number>
1797 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1798 template solution_third_derivative_type<Number>> &
1800 const std::string &global_vector_name,
1801 const Extractor & variable,
1802 const Number dummy)
1803 {
1804 const std::vector<Number> &independent_local_dofs =
1805 get_local_dof_values(global_vector_name, dummy);
1806
1807 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1808
1809 const unsigned int n_q_points = fev.n_quadrature_points;
1810
1811 const std::string name = get_unique_name(
1812 global_vector_name, variable, "_third_derivatives_q", n_q_points, dummy);
1813
1814 // Now build the return type
1815 using RetType =
1816 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1817 template solution_third_derivative_type<Number>>;
1818
1819 RetType &ret =
1820 internal_data_storage.template get_or_add_object_with_name<RetType>(
1821 name, n_q_points);
1822
1823
1824 AssertDimension(ret.size(), n_q_points);
1825
1826 fev[variable].get_function_third_derivatives_from_local_dof_values(
1827 independent_local_dofs, ret);
1828 return ret;
1829 }
1830
1831
1832
1833 template <int dim, int spacedim>
1834 template <typename Extractor, typename Number>
1835 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1836 template solution_symmetric_gradient_type<Number>> &
1838 const std::string &global_vector_name,
1839 const Extractor & variable,
1840 const Number dummy)
1841 {
1842 const std::vector<Number> &independent_local_dofs =
1843 get_local_dof_values(global_vector_name, dummy);
1844
1845 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1846
1847 const unsigned int n_q_points = fev.n_quadrature_points;
1848
1849 const std::string name = get_unique_name(
1850 global_vector_name, variable, "_symmetric_gradient_q", n_q_points, dummy);
1851
1852
1853 // Now build the return type
1854 using RetType =
1855 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1856 template solution_symmetric_gradient_type<Number>>;
1857
1858 RetType &ret =
1859 internal_data_storage.template get_or_add_object_with_name<RetType>(
1860 name, n_q_points);
1861
1862
1863 AssertDimension(ret.size(), n_q_points);
1864
1865 fev[variable].get_function_symmetric_gradients_from_local_dof_values(
1866 independent_local_dofs, ret);
1867 return ret;
1868 }
1869
1870
1871 template <int dim, int spacedim>
1872 template <typename Extractor, typename Number>
1873 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1874 template solution_divergence_type<Number>> &
1876 const std::string &global_vector_name,
1877 const Extractor & variable,
1878 const Number dummy)
1879 {
1880 const std::vector<Number> &independent_local_dofs =
1881 get_local_dof_values(global_vector_name, dummy);
1882
1883 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1884
1885 const unsigned int n_q_points = fev.n_quadrature_points;
1886
1887 const std::string name = get_unique_name(
1888 global_vector_name, variable, "_divergence_q", n_q_points, dummy);
1889
1890 // Now build the return type
1891 using RetType =
1892 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1893 template solution_divergence_type<Number>>;
1894
1895 RetType &ret =
1896 internal_data_storage.template get_or_add_object_with_name<RetType>(
1897 name, n_q_points);
1898
1899
1900 AssertDimension(ret.size(), n_q_points);
1901
1902 fev[variable].get_function_divergences_from_local_dof_values(
1903 independent_local_dofs, ret);
1904 return ret;
1905 }
1906
1907
1908
1909 template <int dim, int spacedim>
1910 template <typename Extractor, typename Number>
1911 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1912 template solution_curl_type<Number>> &
1913 ScratchData<dim, spacedim>::get_curls(const std::string &global_vector_name,
1914 const Extractor & variable,
1915 const Number dummy)
1916 {
1917 const std::vector<Number> &independent_local_dofs =
1918 get_local_dof_values(global_vector_name, dummy);
1919
1920 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1921
1922 const unsigned int n_q_points = fev.n_quadrature_points;
1923
1924 const std::string name = get_unique_name(
1925 global_vector_name, variable, "_curl_q", n_q_points, dummy);
1926
1927 // Now build the return type
1928 using RetType =
1929 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1930 template solution_curl_type<Number>>;
1931
1932 RetType &ret =
1933 internal_data_storage.template get_or_add_object_with_name<RetType>(
1934 name, n_q_points);
1935
1936 AssertDimension(ret.size(), n_q_points);
1937
1938 fev[variable].get_function_curls_from_local_dof_values(
1939 independent_local_dofs, ret);
1940 return ret;
1941 }
1942
1943
1944
1945 template <int dim, int spacedim>
1946 template <typename Extractor, typename Number>
1947 const std::vector<typename FEInterfaceViews::View<dim, spacedim, Extractor>::
1948 template solution_value_type<Number>> &
1950 const std::string &global_vector_name,
1951 const Extractor & variable,
1952 const Number dummy)
1953 {
1954 const std::vector<Number> &independent_local_dofs =
1955 get_local_dof_values(global_vector_name, dummy);
1956
1958 get_current_interface_fe_values();
1959
1960 AssertDimension(independent_local_dofs.size(),
1962
1963 const unsigned int n_q_points = feiv.n_quadrature_points;
1964
1965 const std::string name = get_unique_name(
1966 global_vector_name, variable, "_jump_values_q", n_q_points, dummy);
1967
1968 // Now build the return type
1969 using RetType =
1970 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1971 template solution_value_type<Number>>;
1972
1973 RetType &ret =
1974 internal_data_storage.template get_or_add_object_with_name<RetType>(
1975 name, n_q_points);
1976
1977 AssertDimension(ret.size(), n_q_points);
1978
1979 feiv[variable].get_jump_in_function_values_from_local_dof_values(
1980 independent_local_dofs, ret);
1981 return ret;
1982 }
1983
1984
1985
1986 template <int dim, int spacedim>
1987 template <typename Extractor, typename Number>
1988 const std::vector<typename FEInterfaceViews::View<dim, spacedim, Extractor>::
1989 template solution_gradient_type<Number>> &
1991 const std::string &global_vector_name,
1992 const Extractor & variable,
1993 const Number dummy)
1994 {
1995 const std::vector<Number> &independent_local_dofs =
1996 get_local_dof_values(global_vector_name, dummy);
1997
1999 get_current_interface_fe_values();
2000
2001 AssertDimension(independent_local_dofs.size(),
2003
2004 const unsigned int n_q_points = feiv.n_quadrature_points;
2005
2006 const std::string name = get_unique_name(
2007 global_vector_name, variable, "_jump_gradients_q", n_q_points, dummy);
2008
2009 // Now build the return type
2010 using RetType =
2011 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
2012 template solution_gradient_type<Number>>;
2013
2014 RetType &ret =
2015 internal_data_storage.template get_or_add_object_with_name<RetType>(
2016 name, n_q_points);
2017
2018 AssertDimension(ret.size(), n_q_points);
2019
2020 feiv[variable].get_jump_in_function_gradients_from_local_dof_values(
2021 independent_local_dofs, ret);
2022 return ret;
2023 }
2024
2025
2026
2027 template <int dim, int spacedim>
2028 template <typename Extractor, typename Number>
2029 const std::vector<typename FEInterfaceViews::View<dim, spacedim, Extractor>::
2030 template solution_hessian_type<Number>> &
2032 const std::string &global_vector_name,
2033 const Extractor & variable,
2034 const Number dummy)
2035 {
2036 const std::vector<Number> &independent_local_dofs =
2037 get_local_dof_values(global_vector_name, dummy);
2038
2040 get_current_interface_fe_values();
2041
2042 AssertDimension(independent_local_dofs.size(),
2044
2045 const unsigned int n_q_points = feiv.n_quadrature_points;
2046
2047 const std::string name = get_unique_name(
2048 global_vector_name, variable, "_jump_hessians_q", n_q_points, dummy);
2049
2050 // Now build the return type
2051 using RetType =
2052 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
2053 template solution_hessian_type<Number>>;
2054
2055 RetType &ret =
2056 internal_data_storage.template get_or_add_object_with_name<RetType>(
2057 name, n_q_points);
2058
2059 AssertDimension(ret.size(), n_q_points);
2060
2061 feiv[variable].get_jump_in_function_hessians_from_local_dof_values(
2062 independent_local_dofs, ret);
2063 return ret;
2064 }
2065
2066
2067
2068 template <int dim, int spacedim>
2069 template <typename Extractor, typename Number>
2070 const std::vector<typename FEInterfaceViews::View<dim, spacedim, Extractor>::
2071 template solution_third_derivative_type<Number>> &
2073 const std::string &global_vector_name,
2074 const Extractor & variable,
2075 const Number dummy)
2076 {
2077 const std::vector<Number> &independent_local_dofs =
2078 get_local_dof_values(global_vector_name, dummy);
2079
2081 get_current_interface_fe_values();
2082
2083 AssertDimension(independent_local_dofs.size(),
2085
2086 const unsigned int n_q_points = feiv.n_quadrature_points;
2087
2088 const std::string name = get_unique_name(global_vector_name,
2089 variable,
2090 "_jump_third_derivatives_q",
2091 n_q_points,
2092 dummy);
2093
2094 // Now build the return type
2095 using RetType =
2096 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
2097 template solution_third_derivative_type<Number>>;
2098
2099 RetType &ret =
2100 internal_data_storage.template get_or_add_object_with_name<RetType>(
2101 name, n_q_points);
2102
2103 AssertDimension(ret.size(), n_q_points);
2104
2105 feiv[variable].get_jump_in_function_third_derivatives_from_local_dof_values(
2106 independent_local_dofs, ret);
2107 return ret;
2108 }
2109
2110
2111
2112 template <int dim, int spacedim>
2113 template <typename Extractor, typename Number>
2114 const std::vector<typename FEInterfaceViews::View<dim, spacedim, Extractor>::
2115 template solution_value_type<Number>> &
2117 const std::string &global_vector_name,
2118 const Extractor & variable,
2119 const Number dummy)
2120 {
2121 const std::vector<Number> &independent_local_dofs =
2122 get_local_dof_values(global_vector_name, dummy);
2123
2125 get_current_interface_fe_values();
2126
2127 AssertDimension(independent_local_dofs.size(),
2129
2130 const unsigned int n_q_points = feiv.n_quadrature_points;
2131
2132 const std::string name = get_unique_name(
2133 global_vector_name, variable, "_average_values_q", n_q_points, dummy);
2134
2135 // Now build the return type
2136 using RetType =
2137 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
2138 template solution_value_type<Number>>;
2139
2140 RetType &ret =
2141 internal_data_storage.template get_or_add_object_with_name<RetType>(
2142 name, n_q_points);
2143
2144 AssertDimension(ret.size(), n_q_points);
2145
2146 feiv[variable].get_average_of_function_values_from_local_dof_values(
2147 independent_local_dofs, ret);
2148 return ret;
2149 }
2150
2151
2152
2153 template <int dim, int spacedim>
2154 template <typename Extractor, typename Number>
2155 const std::vector<typename FEInterfaceViews::View<dim, spacedim, Extractor>::
2156 template solution_gradient_type<Number>> &
2158 const std::string &global_vector_name,
2159 const Extractor & variable,
2160 const Number dummy)
2161 {
2162 const std::vector<Number> &independent_local_dofs =
2163 get_local_dof_values(global_vector_name, dummy);
2164
2166 get_current_interface_fe_values();
2167
2168 AssertDimension(independent_local_dofs.size(),
2170
2171 const unsigned int n_q_points = feiv.n_quadrature_points;
2172
2173 const std::string name = get_unique_name(
2174 global_vector_name, variable, "_average_gradients_q", n_q_points, dummy);
2175
2176 // Now build the return type
2177 using RetType =
2178 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
2179 template solution_gradient_type<Number>>;
2180
2181 RetType &ret =
2182 internal_data_storage.template get_or_add_object_with_name<RetType>(
2183 name, n_q_points);
2184
2185 AssertDimension(ret.size(), n_q_points);
2186
2187 feiv[variable].get_average_of_function_gradients_from_local_dof_values(
2188 independent_local_dofs, ret);
2189 return ret;
2190 }
2191
2192
2193
2194 template <int dim, int spacedim>
2195 template <typename Extractor, typename Number>
2196 const std::vector<typename FEInterfaceViews::View<dim, spacedim, Extractor>::
2197 template solution_hessian_type<Number>> &
2199 const std::string &global_vector_name,
2200 const Extractor & variable,
2201 const Number dummy)
2202 {
2203 const std::vector<Number> &independent_local_dofs =
2204 get_local_dof_values(global_vector_name, dummy);
2205
2207 get_current_interface_fe_values();
2208
2209 AssertDimension(independent_local_dofs.size(),
2211
2212 const unsigned int n_q_points = feiv.n_quadrature_points;
2213
2214 const std::string name = get_unique_name(
2215 global_vector_name, variable, "_average_hessians_q", n_q_points, dummy);
2216
2217 // Now build the return type
2218 using RetType =
2219 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
2220 template solution_hessian_type<Number>>;
2221
2222 RetType &ret =
2223 internal_data_storage.template get_or_add_object_with_name<RetType>(
2224 name, n_q_points);
2225
2226 AssertDimension(ret.size(), n_q_points);
2227
2228 feiv[variable].get_average_of_function_hessians_from_local_dof_values(
2229 independent_local_dofs, ret);
2230 return ret;
2231 }
2232
2233
2234#endif
2235
2236} // namespace MeshWorker
2237
2239
2240#endif
const unsigned int n_quadrature_points
unsigned n_current_interface_dofs() const
const unsigned int n_quadrature_points
Definition: fe_values.h:2432
Abstract base class for mapping classes.
Definition: mapping.h:311
GeneralDataStorage user_data_storage
hp::QCollection< dim - 1 > face_quadrature_collection
const std::vector< typename FEInterfaceViews::View< dim, spacedim, Extractor >::template solution_gradient_type< Number > > & get_jumps_in_gradients(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
bool has_hp_capabilities() const
const hp::MappingCollection< dim, spacedim > & get_mapping_collection() const
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_value_type< Number > > & get_values(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const FEValuesBase< dim, spacedim > & get_current_fe_values() const
GeneralDataStorage & get_general_data_storage()
const Quadrature< dim - 1 > & get_face_quadrature() const
std::string get_unique_name(const std::string &global_vector_name, const Extractor &variable, const std::string &object_type, const unsigned int size, const Number &exemplar_number) const
const std::vector< typename FEInterfaceViews::View< dim, spacedim, Extractor >::template solution_value_type< Number > > & get_jumps_in_values(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_divergence_type< Number > > & get_divergences(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const std::vector< Tensor< 1, spacedim > > & get_neighbor_normal_vectors()
unsigned int n_dofs_per_cell() const
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > hp_fe_subface_values
unsigned int n_neighbor_dofs_per_cell() const
std::unique_ptr< FEValues< dim, spacedim > > neighbor_fe_values
std::unique_ptr< FEFaceValues< dim, spacedim > > fe_face_values
std::string get_unique_dofs_name(const std::string &global_vector_name, const unsigned int size, const Number &exemplar_number) const
void extract_local_dof_values(const std::string &global_vector_name, const VectorType &input_vector, const Number dummy=Number(0))
std::unique_ptr< FEValues< dim, spacedim > > fe_values
UpdateFlags cell_update_flags
const std::vector< double > & get_JxW_values() const
const Mapping< dim, spacedim > & get_mapping() const
std::unique_ptr< FEInterfaceValues< dim, spacedim > > interface_fe_values
const std::vector< double > & get_neighbor_JxW_values() const
Quadrature< dim > cell_quadrature
const std::vector< typename FEInterfaceViews::View< dim, spacedim, Extractor >::template solution_value_type< Number > > & get_averages_of_values(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const FiniteElement< dim, spacedim > & get_fe() const
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
std::unique_ptr< FEFaceValues< dim, spacedim > > neighbor_fe_face_values
const std::vector< typename FEInterfaceViews::View< dim, spacedim, Extractor >::template solution_gradient_type< Number > > & get_averages_of_gradients(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const hp::QCollection< dim > & get_cell_quadrature_collection() const
std::vector< types::global_dof_index > neighbor_dof_indices
SmartPointer< const FiniteElement< dim, spacedim > > fe
UpdateFlags face_update_flags
const std::vector< types::global_dof_index > & get_neighbor_dof_indices() const
UpdateFlags get_face_update_flags() const
const std::vector< typename FEInterfaceViews::View< dim, spacedim, Extractor >::template solution_hessian_type< Number > > & get_jumps_in_hessians(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_laplacian_type< Number > > & get_laplacians(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
std::unique_ptr< hp::FEValues< dim, spacedim > > hp_fe_values
SmartPointer< const hp::MappingCollection< dim, spacedim > > mapping_collection
UpdateFlags neighbor_cell_update_flags
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_third_derivative_type< Number > > & get_third_derivatives(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
UpdateFlags get_neighbor_cell_update_flags() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_curl_type< Number > > & get_curls(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const FEValues< dim, spacedim > & reinit_neighbor(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > neighbor_hp_fe_face_values
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_hessian_type< Number > > & get_hessians(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
std::vector< types::global_dof_index > local_dof_indices
SmartPointer< const Mapping< dim, spacedim > > mapping
const std::vector< typename FEInterfaceViews::View< dim, spacedim, Extractor >::template solution_third_derivative_type< Number > > & get_jumps_in_third_derivatives(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
SmartPointer< const FEValuesBase< dim, spacedim > > current_fe_values
const std::vector< Number > & get_local_dof_values(const std::string &global_vector_name, Number dummy=Number(0)) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
std::unique_ptr< hp::FEValues< dim, spacedim > > neighbor_hp_fe_values
const FEValuesBase< dim, spacedim > & get_current_neighbor_fe_values() const
const std::vector< types::global_dof_index > & get_local_dof_indices() const
SmartPointer< const FEValuesBase< dim, spacedim > > current_neighbor_fe_values
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_symmetric_gradient_type< Number > > & get_symmetric_gradients(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
Quadrature< dim - 1 > face_quadrature
const std::vector< typename FEInterfaceViews::View< dim, spacedim, Extractor >::template solution_hessian_type< Number > > & get_averages_of_hessians(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const FEInterfaceValues< dim, spacedim > & get_current_interface_fe_values() const
hp::QCollection< dim > cell_quadrature_collection
UpdateFlags get_cell_update_flags() const
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_gradient_type< Number > > & get_gradients(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
UpdateFlags neighbor_face_update_flags
const Quadrature< dim > & get_cell_quadrature() const
const hp::QCollection< dim - 1 > & get_face_quadrature_collection() const
std::unique_ptr< FESubfaceValues< dim, spacedim > > neighbor_fe_subface_values
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > neighbor_hp_fe_subface_values
std::unique_ptr< FESubfaceValues< dim, spacedim > > fe_subface_values
GeneralDataStorage internal_data_storage
UpdateFlags get_neighbor_face_update_flags() const
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > hp_fe_face_values
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
UpdateFlags
@ update_default
No update.
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
static ::ExceptionBase & ExcOnlyAvailableWithoutHP()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcOnlyAvailableWithHP()
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:438
typename ::internal::FEInterfaceViews::ViewType< dim, spacedim, Extractor >::type View
std::string type_to_string(const T &t)
Definition: utilities.h:1144
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:473