libMesh
Public Member Functions | Private Attributes | List of all members
libMesh::BoundaryProjectSolution Class Reference

This class implements projecting an arbitrary boundary function to the current mesh. More...

Public Member Functions

 BoundaryProjectSolution (const std::set< boundary_id_type > &b_in, const std::vector< unsigned int > &variables_in, const System &system_in, FunctionBase< Number > *f_in, FunctionBase< Gradient > *g_in, const Parameters &parameters_in, NumericVector< Number > &new_v_in)
 
 BoundaryProjectSolution (const BoundaryProjectSolution &in)
 
void operator() (const ConstElemRange &range) const
 

Private Attributes

const std::set< boundary_id_type > & b
 
const std::vector< unsigned int > & variables
 
const Systemsystem
 
std::unique_ptr< FunctionBase< Number > > f
 
std::unique_ptr< FunctionBase< Gradient > > g
 
const Parametersparameters
 
NumericVector< Number > & new_vector
 

Detailed Description

This class implements projecting an arbitrary boundary function to the current mesh.

This may be executed in parallel on multiple threads.

Definition at line 192 of file system_projection.C.

Constructor & Destructor Documentation

◆ BoundaryProjectSolution() [1/2]

libMesh::BoundaryProjectSolution::BoundaryProjectSolution ( const std::set< boundary_id_type > &  b_in,
const std::vector< unsigned int > &  variables_in,
const System system_in,
FunctionBase< Number > *  f_in,
FunctionBase< Gradient > *  g_in,
const Parameters parameters_in,
NumericVector< Number > &  new_v_in 
)
inline

Definition at line 204 of file system_projection.C.

References libMesh::libmesh_assert().

210  :
211  b(b_in),
212  variables(variables_in),
213  system(system_in),
214  f(f_in ? f_in->clone() : std::unique_ptr<FunctionBase<Number>>()),
215  g(g_in ? g_in->clone() : std::unique_ptr<FunctionBase<Gradient>>()),
216  parameters(parameters_in),
217  new_vector(new_v_in)
218  {
219  libmesh_assert(f.get());
220  f->init();
221  if (g.get())
222  g->init();
223  }
std::unique_ptr< FunctionBase< Number > > f
const std::vector< unsigned int > & variables
NumericVector< Number > & new_vector
libmesh_assert(ctx)
virtual std::unique_ptr< FunctionBase< Output > > clone() const =0
const std::set< boundary_id_type > & b
std::unique_ptr< FunctionBase< Gradient > > g

◆ BoundaryProjectSolution() [2/2]

libMesh::BoundaryProjectSolution::BoundaryProjectSolution ( const BoundaryProjectSolution in)
inline

Definition at line 225 of file system_projection.C.

References libMesh::libmesh_assert().

225  :
226  b(in.b),
227  variables(in.variables),
228  system(in.system),
229  f(in.f.get() ? in.f->clone() : std::unique_ptr<FunctionBase<Number>>()),
230  g(in.g.get() ? in.g->clone() : std::unique_ptr<FunctionBase<Gradient>>()),
231  parameters(in.parameters),
232  new_vector(in.new_vector)
233  {
234  libmesh_assert(f.get());
235  f->init();
236  if (g.get())
237  g->init();
238  }
std::unique_ptr< FunctionBase< Number > > f
const std::vector< unsigned int > & variables
NumericVector< Number > & new_vector
libmesh_assert(ctx)
const std::set< boundary_id_type > & b
std::unique_ptr< FunctionBase< Gradient > > g

Member Function Documentation

◆ operator()()

void libMesh::BoundaryProjectSolution::operator() ( const ConstElemRange range) const

This method projects an arbitrary boundary solution to the current mesh. The input function f gives the arbitrary solution, while the new_vector (which should already be correctly sized) gives the solution (to be computed) on the current mesh.

Definition at line 1473 of file system_projection.C.

References std::abs(), libMesh::Variable::active_on_subdomain(), libMesh::BoundaryInfo::boundary_ids(), libMesh::C_ONE, libMesh::C_ZERO, libMesh::DenseMatrix< T >::cholesky_solve(), libMesh::FEType::default_quadrature_rule(), dim, libMesh::DISCONTINUOUS, libMesh::DofMap::dof_indices(), libMesh::BoundaryInfo::edge_boundary_ids(), libMesh::FEType::family, libMesh::HERMITE, libMesh::index_range(), libMesh::libmesh_assert(), libMesh::make_range(), n_nodes, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::SCALAR, libMesh::Threads::spin_mtx, libMesh::TOLERANCE, libMesh::Variable::type(), libMesh::DofMap::variable(), libMesh::DenseMatrix< T >::zero(), and libMesh::DenseVector< T >::zero().

1474 {
1475  // We need data to project
1476  libmesh_assert(f.get());
1477 
1485  // The dimensionality of the current mesh
1486  const unsigned int dim = system.get_mesh().mesh_dimension();
1487 
1488  // The DofMap for this system
1489  const DofMap & dof_map = system.get_dof_map();
1490 
1491  // Boundary info for the current mesh
1492  const BoundaryInfo & boundary_info =
1494 
1495  // The element matrix and RHS for projections.
1496  // Note that Ke is always real-valued, whereas
1497  // Fe may be complex valued if complex number
1498  // support is enabled
1499  DenseMatrix<Real> Ke;
1500  DenseVector<Number> Fe;
1501  // The new element coefficients
1502  DenseVector<Number> Ue;
1503 
1504 
1505  // Loop over all the variables we've been requested to project
1506  for (auto v : make_range(variables.size()))
1507  {
1508  const unsigned int var = variables[v];
1509 
1510  const Variable & variable = dof_map.variable(var);
1511 
1512  const FEType & fe_type = variable.type();
1513 
1514  if (fe_type.family == SCALAR)
1515  continue;
1516 
1517  const unsigned int var_component =
1519 
1520  // Get FE objects of the appropriate type
1521  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
1522 
1523  // Prepare variables for projection
1524  std::unique_ptr<QBase> qedgerule (fe_type.default_quadrature_rule(1));
1525  std::unique_ptr<QBase> qsiderule (fe_type.default_quadrature_rule(dim-1));
1526 
1527  // The values of the shape functions at the quadrature
1528  // points
1529  const std::vector<std::vector<Real>> & phi = fe->get_phi();
1530 
1531  // The gradients of the shape functions at the quadrature
1532  // points on the child element.
1533  const std::vector<std::vector<RealGradient>> * dphi = nullptr;
1534 
1535  const FEContinuity cont = fe->get_continuity();
1536 
1537  if (cont == C_ONE)
1538  {
1539  // We'll need gradient data for a C1 projection
1540  libmesh_assert(g.get());
1541 
1542  const std::vector<std::vector<RealGradient>> &
1543  ref_dphi = fe->get_dphi();
1544  dphi = &ref_dphi;
1545  }
1546 
1547  // The Jacobian * quadrature weight at the quadrature points
1548  const std::vector<Real> & JxW =
1549  fe->get_JxW();
1550 
1551  // The XYZ locations of the quadrature points
1552  const std::vector<Point> & xyz_values =
1553  fe->get_xyz();
1554 
1555  // The global DOF indices
1556  std::vector<dof_id_type> dof_indices;
1557  // Side/edge DOF indices
1558  std::vector<unsigned int> side_dofs;
1559 
1560  // Container to catch IDs passed back from BoundaryInfo.
1561  std::vector<boundary_id_type> bc_ids;
1562 
1563  // Iterate over all the elements in the range
1564  for (const auto & elem : range)
1565  {
1566  // Per-subdomain variables don't need to be projected on
1567  // elements where they're not active
1568  if (!variable.active_on_subdomain(elem->subdomain_id()))
1569  continue;
1570 
1571  const unsigned short n_nodes = elem->n_nodes();
1572  const unsigned short n_edges = elem->n_edges();
1573  const unsigned short n_sides = elem->n_sides();
1574 
1575  // Find out which nodes, edges and sides are on a requested
1576  // boundary:
1577  std::vector<bool> is_boundary_node(n_nodes, false),
1578  is_boundary_edge(n_edges, false),
1579  is_boundary_side(n_sides, false);
1580 
1581  // We also maintain a separate list of nodeset-based boundary nodes
1582  std::vector<bool> is_boundary_nodeset(n_nodes, false);
1583 
1584  for (unsigned char s=0; s != n_sides; ++s)
1585  {
1586  // First see if this side has been requested
1587  boundary_info.boundary_ids (elem, s, bc_ids);
1588  bool do_this_side = false;
1589  for (const auto & bc_id : bc_ids)
1590  if (b.count(bc_id))
1591  {
1592  do_this_side = true;
1593  break;
1594  }
1595  if (!do_this_side)
1596  continue;
1597 
1598  is_boundary_side[s] = true;
1599 
1600  // Then see what nodes and what edges are on it
1601  for (unsigned int n=0; n != n_nodes; ++n)
1602  if (elem->is_node_on_side(n,s))
1603  is_boundary_node[n] = true;
1604  for (unsigned int e=0; e != n_edges; ++e)
1605  if (elem->is_edge_on_side(e,s))
1606  is_boundary_edge[e] = true;
1607  }
1608 
1609  // We can also project on nodes, so we should also independently
1610  // check whether the nodes have been requested
1611  for (unsigned int n=0; n != n_nodes; ++n)
1612  {
1613  boundary_info.boundary_ids (elem->node_ptr(n), bc_ids);
1614 
1615  for (const auto & bc_id : bc_ids)
1616  if (b.count(bc_id))
1617  {
1618  is_boundary_node[n] = true;
1619  is_boundary_nodeset[n] = true;
1620  }
1621  }
1622 
1623  // We can also project on edges, so we should also independently
1624  // check whether the edges have been requested
1625  for (unsigned short e=0; e != n_edges; ++e)
1626  {
1627  boundary_info.edge_boundary_ids (elem, e, bc_ids);
1628 
1629  for (const auto & bc_id : bc_ids)
1630  if (b.count(bc_id))
1631  is_boundary_edge[e] = true;
1632  }
1633 
1634  // Update the DOF indices for this element based on
1635  // the current mesh
1636  dof_map.dof_indices (elem, dof_indices, var);
1637 
1638  // The number of DOFs on the element
1639  const unsigned int n_dofs =
1640  cast_int<unsigned int>(dof_indices.size());
1641 
1642  // Fixed vs. free DoFs on edge/face projections
1643  std::vector<char> dof_is_fixed(n_dofs, false); // bools
1644  std::vector<int> free_dof(n_dofs, 0);
1645 
1646  // Zero the interpolated values
1647  Ue.resize (n_dofs); Ue.zero();
1648 
1649  // In general, we need a series of
1650  // projections to ensure a unique and continuous
1651  // solution. We start by interpolating boundary nodes, then
1652  // hold those fixed and project boundary edges, then hold
1653  // those fixed and project boundary faces,
1654 
1655  // Interpolate node values first
1656  unsigned int current_dof = 0;
1657  for (unsigned short n = 0; n != n_nodes; ++n)
1658  {
1659  // FIXME: this should go through the DofMap,
1660  // not duplicate dof_indices code badly!
1661 
1662  // This call takes into account elem->p_level() internally.
1663  const unsigned int nc =
1664  FEInterface::n_dofs_at_node (fe_type, elem, n);
1665 
1666  if ((!elem->is_vertex(n) || !is_boundary_node[n]) &&
1667  !is_boundary_nodeset[n])
1668  {
1669  current_dof += nc;
1670  continue;
1671  }
1672  if (cont == DISCONTINUOUS)
1673  {
1674  libmesh_assert_equal_to (nc, 0);
1675  }
1676  // Assume that C_ZERO elements have a single nodal
1677  // value shape function
1678  else if (cont == C_ZERO)
1679  {
1680  libmesh_assert_equal_to (nc, 1);
1681  Ue(current_dof) = f->component(var_component,
1682  elem->point(n),
1683  system.time);
1684  dof_is_fixed[current_dof] = true;
1685  current_dof++;
1686  }
1687  // The hermite element vertex shape functions are weird
1688  else if (fe_type.family == HERMITE)
1689  {
1690  Ue(current_dof) = f->component(var_component,
1691  elem->point(n),
1692  system.time);
1693  dof_is_fixed[current_dof] = true;
1694  current_dof++;
1695  Gradient grad = g->component(var_component,
1696  elem->point(n),
1697  system.time);
1698  // x derivative
1699  Ue(current_dof) = grad(0);
1700  dof_is_fixed[current_dof] = true;
1701  current_dof++;
1702 #if LIBMESH_DIM > 1
1703  if (dim > 1)
1704  {
1705  // We'll finite difference mixed derivatives
1706  Point nxminus = elem->point(n),
1707  nxplus = elem->point(n);
1708  nxminus(0) -= TOLERANCE;
1709  nxplus(0) += TOLERANCE;
1710  Gradient gxminus = g->component(var_component,
1711  nxminus,
1712  system.time);
1713  Gradient gxplus = g->component(var_component,
1714  nxplus,
1715  system.time);
1716  // y derivative
1717  Ue(current_dof) = grad(1);
1718  dof_is_fixed[current_dof] = true;
1719  current_dof++;
1720  // xy derivative
1721  Ue(current_dof) = (gxplus(1) - gxminus(1))
1722  / 2. / TOLERANCE;
1723  dof_is_fixed[current_dof] = true;
1724  current_dof++;
1725 
1726 #if LIBMESH_DIM > 2
1727  if (dim > 2)
1728  {
1729  // z derivative
1730  Ue(current_dof) = grad(2);
1731  dof_is_fixed[current_dof] = true;
1732  current_dof++;
1733  // xz derivative
1734  Ue(current_dof) = (gxplus(2) - gxminus(2))
1735  / 2. / TOLERANCE;
1736  dof_is_fixed[current_dof] = true;
1737  current_dof++;
1738  // We need new points for yz
1739  Point nyminus = elem->point(n),
1740  nyplus = elem->point(n);
1741  nyminus(1) -= TOLERANCE;
1742  nyplus(1) += TOLERANCE;
1743  Gradient gyminus = g->component(var_component,
1744  nyminus,
1745  system.time);
1746  Gradient gyplus = g->component(var_component,
1747  nyplus,
1748  system.time);
1749  // xz derivative
1750  Ue(current_dof) = (gyplus(2) - gyminus(2))
1751  / 2. / TOLERANCE;
1752  dof_is_fixed[current_dof] = true;
1753  current_dof++;
1754  // Getting a 2nd order xyz is more tedious
1755  Point nxmym = elem->point(n),
1756  nxmyp = elem->point(n),
1757  nxpym = elem->point(n),
1758  nxpyp = elem->point(n);
1759  nxmym(0) -= TOLERANCE;
1760  nxmym(1) -= TOLERANCE;
1761  nxmyp(0) -= TOLERANCE;
1762  nxmyp(1) += TOLERANCE;
1763  nxpym(0) += TOLERANCE;
1764  nxpym(1) -= TOLERANCE;
1765  nxpyp(0) += TOLERANCE;
1766  nxpyp(1) += TOLERANCE;
1767  Gradient gxmym = g->component(var_component,
1768  nxmym,
1769  system.time);
1770  Gradient gxmyp = g->component(var_component,
1771  nxmyp,
1772  system.time);
1773  Gradient gxpym = g->component(var_component,
1774  nxpym,
1775  system.time);
1776  Gradient gxpyp = g->component(var_component,
1777  nxpyp,
1778  system.time);
1779  Number gxzplus = (gxpyp(2) - gxmyp(2))
1780  / 2. / TOLERANCE;
1781  Number gxzminus = (gxpym(2) - gxmym(2))
1782  / 2. / TOLERANCE;
1783  // xyz derivative
1784  Ue(current_dof) = (gxzplus - gxzminus)
1785  / 2. / TOLERANCE;
1786  dof_is_fixed[current_dof] = true;
1787  current_dof++;
1788  }
1789 #endif // LIBMESH_DIM > 2
1790  }
1791 #endif // LIBMESH_DIM > 1
1792  }
1793  // Assume that other C_ONE elements have a single nodal
1794  // value shape function and nodal gradient component
1795  // shape functions
1796  else if (cont == C_ONE)
1797  {
1798  libmesh_assert_equal_to (nc, 1 + dim);
1799  Ue(current_dof) = f->component(var_component,
1800  elem->point(n),
1801  system.time);
1802  dof_is_fixed[current_dof] = true;
1803  current_dof++;
1804  Gradient grad = g->component(var_component,
1805  elem->point(n),
1806  system.time);
1807  for (unsigned int i=0; i!= dim; ++i)
1808  {
1809  Ue(current_dof) = grad(i);
1810  dof_is_fixed[current_dof] = true;
1811  current_dof++;
1812  }
1813  }
1814  else
1815  libmesh_error_msg("Unknown continuity " << cont);
1816  }
1817 
1818  // In 3D, project any edge values next
1819  if (dim > 2 && cont != DISCONTINUOUS)
1820  for (unsigned short e = 0; e != n_edges; ++e)
1821  {
1822  if (!is_boundary_edge[e])
1823  continue;
1824 
1825  FEInterface::dofs_on_edge(elem, dim, fe_type, e,
1826  side_dofs);
1827 
1828  const unsigned int n_side_dofs =
1829  cast_int<unsigned int>(side_dofs.size());
1830 
1831  // Some edge dofs are on nodes and already
1832  // fixed, others are free to calculate
1833  unsigned int free_dofs = 0;
1834  for (auto i : make_range(n_side_dofs))
1835  if (!dof_is_fixed[side_dofs[i]])
1836  free_dof[free_dofs++] = i;
1837 
1838  // There may be nothing to project
1839  if (!free_dofs)
1840  continue;
1841 
1842  Ke.resize (free_dofs, free_dofs); Ke.zero();
1843  Fe.resize (free_dofs); Fe.zero();
1844  // The new edge coefficients
1845  DenseVector<Number> Uedge(free_dofs);
1846 
1847  // Initialize FE data on the edge
1848  fe->attach_quadrature_rule (qedgerule.get());
1849  fe->edge_reinit (elem, e);
1850  const unsigned int n_qp = qedgerule->n_points();
1851 
1852  // Loop over the quadrature points
1853  for (unsigned int qp=0; qp<n_qp; qp++)
1854  {
1855  // solution at the quadrature point
1856  Number fineval = f->component(var_component,
1857  xyz_values[qp],
1858  system.time);
1859  // solution grad at the quadrature point
1860  Gradient finegrad;
1861  if (cont == C_ONE)
1862  finegrad = g->component(var_component,
1863  xyz_values[qp],
1864  system.time);
1865 
1866  // Form edge projection matrix
1867  for (unsigned int sidei=0, freei=0;
1868  sidei != n_side_dofs; ++sidei)
1869  {
1870  unsigned int i = side_dofs[sidei];
1871  // fixed DoFs aren't test functions
1872  if (dof_is_fixed[i])
1873  continue;
1874  for (unsigned int sidej=0, freej=0;
1875  sidej != n_side_dofs; ++sidej)
1876  {
1877  unsigned int j = side_dofs[sidej];
1878  if (dof_is_fixed[j])
1879  Fe(freei) -= phi[i][qp] * phi[j][qp] *
1880  JxW[qp] * Ue(j);
1881  else
1882  Ke(freei,freej) += phi[i][qp] *
1883  phi[j][qp] * JxW[qp];
1884  if (cont == C_ONE)
1885  {
1886  if (dof_is_fixed[j])
1887  Fe(freei) -= ((*dphi)[i][qp] *
1888  (*dphi)[j][qp]) *
1889  JxW[qp] * Ue(j);
1890  else
1891  Ke(freei,freej) += ((*dphi)[i][qp] *
1892  (*dphi)[j][qp])
1893  * JxW[qp];
1894  }
1895  if (!dof_is_fixed[j])
1896  freej++;
1897  }
1898  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
1899  if (cont == C_ONE)
1900  Fe(freei) += (finegrad * (*dphi)[i][qp]) *
1901  JxW[qp];
1902  freei++;
1903  }
1904  }
1905 
1906  Ke.cholesky_solve(Fe, Uedge);
1907 
1908  // Transfer new edge solutions to element
1909  for (unsigned int i=0; i != free_dofs; ++i)
1910  {
1911  Number & ui = Ue(side_dofs[free_dof[i]]);
1913  std::abs(ui - Uedge(i)) < TOLERANCE);
1914  ui = Uedge(i);
1915  dof_is_fixed[side_dofs[free_dof[i]]] = true;
1916  }
1917  }
1918 
1919  // Project any side values (edges in 2D, faces in 3D)
1920  if (dim > 1 && cont != DISCONTINUOUS)
1921  for (unsigned short s = 0; s != n_sides; ++s)
1922  {
1923  if (!is_boundary_side[s])
1924  continue;
1925 
1926  FEInterface::dofs_on_side(elem, dim, fe_type, s,
1927  side_dofs);
1928 
1929  // Some side dofs are on nodes/edges and already
1930  // fixed, others are free to calculate
1931  unsigned int free_dofs = 0;
1932  for (auto i : index_range(side_dofs))
1933  if (!dof_is_fixed[side_dofs[i]])
1934  free_dof[free_dofs++] = i;
1935 
1936  // There may be nothing to project
1937  if (!free_dofs)
1938  continue;
1939 
1940  Ke.resize (free_dofs, free_dofs); Ke.zero();
1941  Fe.resize (free_dofs); Fe.zero();
1942  // The new side coefficients
1943  DenseVector<Number> Uside(free_dofs);
1944 
1945  // Initialize FE data on the side
1946  fe->attach_quadrature_rule (qsiderule.get());
1947  fe->reinit (elem, s);
1948  const unsigned int n_qp = qsiderule->n_points();
1949 
1950  const unsigned int n_side_dofs =
1951  cast_int<unsigned int>(side_dofs.size());
1952 
1953  // Loop over the quadrature points
1954  for (unsigned int qp=0; qp<n_qp; qp++)
1955  {
1956  // solution at the quadrature point
1957  Number fineval = f->component(var_component,
1958  xyz_values[qp],
1959  system.time);
1960  // solution grad at the quadrature point
1961  Gradient finegrad;
1962  if (cont == C_ONE)
1963  finegrad = g->component(var_component,
1964  xyz_values[qp],
1965  system.time);
1966 
1967  // Form side projection matrix
1968  for (unsigned int sidei=0, freei=0;
1969  sidei != n_side_dofs; ++sidei)
1970  {
1971  unsigned int i = side_dofs[sidei];
1972  // fixed DoFs aren't test functions
1973  if (dof_is_fixed[i])
1974  continue;
1975  for (unsigned int sidej=0, freej=0;
1976  sidej != n_side_dofs; ++sidej)
1977  {
1978  unsigned int j = side_dofs[sidej];
1979  if (dof_is_fixed[j])
1980  Fe(freei) -= phi[i][qp] * phi[j][qp] *
1981  JxW[qp] * Ue(j);
1982  else
1983  Ke(freei,freej) += phi[i][qp] *
1984  phi[j][qp] * JxW[qp];
1985  if (cont == C_ONE)
1986  {
1987  if (dof_is_fixed[j])
1988  Fe(freei) -= ((*dphi)[i][qp] *
1989  (*dphi)[j][qp]) *
1990  JxW[qp] * Ue(j);
1991  else
1992  Ke(freei,freej) += ((*dphi)[i][qp] *
1993  (*dphi)[j][qp])
1994  * JxW[qp];
1995  }
1996  if (!dof_is_fixed[j])
1997  freej++;
1998  }
1999  Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
2000  if (cont == C_ONE)
2001  Fe(freei) += (finegrad * (*dphi)[i][qp]) *
2002  JxW[qp];
2003  freei++;
2004  }
2005  }
2006 
2007  Ke.cholesky_solve(Fe, Uside);
2008 
2009  // Transfer new side solutions to element
2010  for (unsigned int i=0; i != free_dofs; ++i)
2011  {
2012  Number & ui = Ue(side_dofs[free_dof[i]]);
2014  std::abs(ui - Uside(i)) < TOLERANCE);
2015  ui = Uside(i);
2016  dof_is_fixed[side_dofs[free_dof[i]]] = true;
2017  }
2018  }
2019 
2020  const dof_id_type
2021  first = new_vector.first_local_index(),
2022  last = new_vector.last_local_index();
2023 
2024  // Lock the new_vector since it is shared among threads.
2025  {
2026  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
2027 
2028  for (unsigned int i = 0; i < n_dofs; i++)
2029  if (dof_is_fixed[i] &&
2030  (dof_indices[i] >= first) &&
2031  (dof_indices[i] < last))
2032  {
2033  new_vector.set(dof_indices[i], Ue(i));
2034  }
2035  }
2036  } // end elem loop
2037  } // end variables loop
2038 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1595
unsigned int variable_scalar_number(std::string_view var, unsigned int component) const
Definition: system.h:2408
static constexpr Real TOLERANCE
unsigned int dim
std::unique_ptr< FunctionBase< Number > > f
const std::vector< unsigned int > & variables
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:159
const MeshBase & get_mesh() const
Definition: system.h:2277
static void dofs_on_side(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int s, std::vector< unsigned int > &di, bool add_p_level=true)
Fills the vector di with the local degree of freedom indices associated with side s of element elem A...
Definition: fe_interface.C:839
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:57
const dof_id_type n_nodes
Definition: tecplot_io.C:67
NumericVector< Number > & new_vector
NumberVectorValue Gradient
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libmesh_assert(ctx)
template class LIBMESH_EXPORT DenseMatrix< Real >
Definition: dense_matrix.C:35
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
Definition: fe_interface.C:679
virtual numeric_index_type first_local_index() const =0
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
unsigned int mesh_dimension() const
Definition: mesh_base.C:324
static void dofs_on_edge(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int e, std::vector< unsigned int > &di, bool add_p_level=true)
Fills the vector di with the local degree of freedom indices associated with edge e of element elem A...
Definition: fe_interface.C:853
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
const std::set< boundary_id_type > & b
std::unique_ptr< FunctionBase< Gradient > > g
const DofMap & get_dof_map() const
Definition: system.h:2293
virtual numeric_index_type last_local_index() const =0
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:111
uint8_t dof_id_type
Definition: id_types.h:67
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30

Member Data Documentation

◆ b

const std::set<boundary_id_type>& libMesh::BoundaryProjectSolution::b
private

Definition at line 195 of file system_projection.C.

◆ f

std::unique_ptr<FunctionBase<Number> > libMesh::BoundaryProjectSolution::f
private

Definition at line 198 of file system_projection.C.

◆ g

std::unique_ptr<FunctionBase<Gradient> > libMesh::BoundaryProjectSolution::g
private

Definition at line 199 of file system_projection.C.

◆ new_vector

NumericVector<Number>& libMesh::BoundaryProjectSolution::new_vector
private

Definition at line 201 of file system_projection.C.

◆ parameters

const Parameters& libMesh::BoundaryProjectSolution::parameters
private

Definition at line 200 of file system_projection.C.

◆ system

const System& libMesh::BoundaryProjectSolution::system
private

Definition at line 197 of file system_projection.C.

◆ variables

const std::vector<unsigned int>& libMesh::BoundaryProjectSolution::variables
private

Definition at line 196 of file system_projection.C.


The documentation for this class was generated from the following file: