Dear Professor Bangerth,
thanks again for your reply, I will start answering your questions and give
you more insight into what I try to do at the moment:
That's exactly the question you have to answer first. What do you *want* to
> happen on these cells? If you know what you want to happen, we can talk
> about
> implementing it :-)
>
In my current Cut FEM approach I'm interested in the solution on the
"physical domain", in my case that's the solution on the (moving) flow
domain including all the boundaries (also on the solid moving ball). Since
I don't adjust the grid to the problem I have also a kind of "flow domain"
inside the solid, which is nonphysically. I don't want to see any effect of
these "nonphysical cells" on my "physical cells", so the natural approach
to me was using FE_Nothing - and this is perfect for the static case and
does the job very well (I don't have a fluid domain inside the solid and
additionally I save some DoFs inside the solid, which is a nice side
effect, but not crucial to me). Therefore I would like to also use this
possibility in the dynamic case, so I would like to go for the 1.) approach
in my first post.
Well, you have to define how you want the solution that didn't exist on the
> empty cell to be interpreted if you actually have something on that cell
> now.
> I think that more than implementation, it's a question of what you *want*
> to
> happen.
>
That's indeed the crucial point that I have to answer. What I definitely
need, are more DoFs, when the solid moves outside a FE_Nothing cell. The
question is what the value of my old solution on these cells shall be. I
guess there are several possibilities e.g. a "fancy one" like extrapolating
the solution from the physical domain from the step before to the now
activated cell. But this is not what I want to do in this first step now. I
just want to set the solution on these cell to "0" (knowing that probably
an extrapolation would be better as a first guess, but since it is just a
starting value for my solver, I guess the "0"-approach is good enough).
This behavior is of course not implemented in deal.ii at the moment. What I
did to reflect the above answer is that I improved my naive approach, which
I described in my first post here, a little bit and do the following now:
1. initialize a parallel::distributed::SolutionTransfer object and add
my old (ghosted) solution vectors via prepare_for_coarsening_and_refinement
2. check if we have a cell that get's "activated"(change from
FE_Nothing to my FE_System for my flow problem) or "deactivated" (other
way around). If this is the case, I do:
1. update the fe_indices
2. call execute_coarsening_and_refinement() on the triangulation
3. distribute the new DoFs
4. create a new sparsity pattern and reinitialize my system matrix
with it
5. create a new solution vector with the new (locally owned) DoFs
6. interpolate the solution to this new solution vector using
SolutionTransfer::interpolate()
This works again in the static case, also if I omit the check in step 2 and
do the whole process even if it's not necessary. In the dynamic case step 6
fails:
SolutionTransfer::interpolate() leads to the
SolutionTransfer::unpack_callback() function called by the tria, which
calls set_dof_values_by_interpolation() on the cell.
This calls set_dof_values() and the following assertion fails:
Assert(values.size() == this->get_dof_handler().n_dofs(),
> typename DoFCellAccessor::ExcVectorDoesNotMatch());
>
which is of course understandable. At the moment I'm trying to understand
the internals of the unpack process and look for a nice way to implement
what I described above - so any help or advice is greatly welcomed :)
Greetings,
Mathias
--
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/4723d890-9dc6-4134-a538-932b3cabb9f7%40googlegroups.com.