On 10/16/22 16:45, Corbin Foucart wrote:
* I'm looking to add a finite element solution of one polynomial order to
another (say, adding a p=2 field to a p=3 field, without loss of generality)
* I can do this by hand the usual way
o declare a nodal quadrature rule (at p=3)
o iterate over the active cells in DoFHandler, using
FEValues::get_function_values to add the results appropriately, either
manually or with MeshWorker
* However, I suspect that deal.ii must have a more idiomatic way of doing
so, since this must be a common operation
Is there a library function to do operations like this? Any guidance would be
appreciated.
There may be a function in VectorTools that interpolates from one finite
element space (=DoFHandler) to another, I don't recall.
But if there isn't: You can do this on a single cell. All you need is the
interpolation matrix for the two finite element spaces on a single cell. Take
a look at FiniteElement::get_interpolation_matrix(). Then all you have to do is do
cell_1->get_dof_values()
for the original FE field, multiply it by that matrix to get the corresponding
values with regard to the other FE field, and then call
cell_2->set_dof_values()
Best
W.
--
------------------------------------------------------------------------
Wolfgang Bangerth email: [email protected]
www: http://www.math.colostate.edu/~bangerth/
--
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 [email protected].
To view this discussion on the web visit
https://groups.google.com/d/msgid/dealii/a88c4289-78a8-a4b1-2a25-55d0b476aed1%40colostate.edu.