> > We need to update function v based on solution u in the following loop. > > for ( x in vector_of_nodal_points ) > v(x) = f(x, u(x)) > > where v(x) and u(x) are the nodal point value for functions v and u, > respectively, and f() is some function depending on function u as well as > the location of the nodal point. So in this case, we need to access nodal > point coordinates and nodal point values of the solution in the same loop. >
You can just loop over all degrees of freedom and do for (unsigned int i=0; i<n_dofs; ++i) v(i) = f(map_dof_to_support_points[i], u(i)) or loop over the cells yourself and do Quadrature q(fe.get_unit_support_points()); FEValues fe_values (..., q, update_q_points); for (const auto& cell) ... points = fe_values.get_quadrature_points(); for (unsigned int i=0; i<n_dofs_per_cell; ++i) v(local_dof_indices[i]) = f(points[i], u(local_dof_indices)); (also see https://github.com/dealii/dealii/wiki/Frequently-Asked-Questions#how-to-get-the-mapped-position-of-support-points-of-my-element). > Well, the above function calculates the gradients of a finite element at > the quadrature points of a cell, not at the nodal points of a cell. > Such a need arises in the following situation. > > for ( x in vector_of_nodal_points ) > v(x) = g(x, u(x), grad u(x)) > You can do similarly, Quadrature q(fe.get_unit_support_points()); FEValues fe_values (..., q, update_q_points); for (const auto& cell) ... points = fe_values.get_quadrature_points(); fe_values.get_function_values(values); fe_values.get_function_gradients(gradients); for (unsigned int i=0; i<n_dofs_per_cell; ++i) v(local_dof_indices[i]) = f(points[i], values(i), gradients(i)); Best, Daniel -- The deal.II project is located at http://www.dealii.org/ For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en --- You received this message because you are subscribed to the Google Groups "deal.II User Group" group. To unsubscribe from this group and stop receiving emails from it, send an email to dealii+unsubscr...@googlegroups.com. To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/CAOYDWbL%2Brkbnqd0vbec33eDzSZdnVCzArZghU1fkU9uyhaXW1w%40mail.gmail.com.