Time dependent hp refinement with FE_Nothing elements

428 views
Skip to first unread message

Claire

unread,
Jul 31, 2015, 11:20:10 AM7/31/15
to deal.II User Group
Dear all,

I am currently intenting to solve the heat equation on a geometry varying with time (to simulate the additive manufacturing of a part).
My simulations also include a local moving heat source responsible of sharp gradients requiring local mesh refinement.

To simulate the additions of the material, I am using an "element birth" technique, which in deal.II translates to :
- Deactivate the cells affecting them a FE_Nothing element
- Activate them using a user-defined function affecting them FE_Q elements

Since my problem involves time dependent hp-refinement since the FE evolves throughout the simulation, I adapted and mixed methods from both step 26 and 27.

The problem that I currently encounter lies in the "refine_mesh" method.
Here is the error message that I get :

--------------------------------------------------------
An error occurred in line <199> of file </home/claire/dealii/include/deal.II/lac/full_matrix.templates.h> in function
    void dealii::FullMatrix<number>::vmult(dealii::Vector<OtherNumber>&, const dealii::Vector<OtherNumber>&, bool) const [with number2 = double; number = double]
The violated condition was:
    !this->empty()
The name and call sequence of the exception was:
    ExcEmptyMatrix()
Additional Information:
(none)

Stacktrace:
-----------
#0  /home/claire/FEM/dealii_update/lib/libdeal_II.g.so.8.4.pre: void dealii::FullMatrix<double>::vmult<double>(dealii::Vector<double>&, dealii::Vector<double> const&, bool) const
#1  /home/claire/FEM/dealii_update/lib/libdeal_II.g.so.8.4.pre: void dealii::DoFCellAccessor<dealii::hp::DoFHandler<2, 2>, false>::get_interpolated_dof_values<dealii::Vector<double>, double>(dealii::Vector<double> const&, dealii::Vector<double>&, unsigned int) const
#2  /home/claire/FEM/dealii_update/lib/libdeal_II.g.so.8.4.pre: void dealii::DoFCellAccessor<dealii::hp::DoFHandler<2, 2>, false>::get_interpolated_dof_values<dealii::Vector<double>, double>(dealii::Vector<double> const&, dealii::Vector<double>&, unsigned int) const
#3  /home/claire/FEM/dealii_update/lib/libdeal_II.g.so.8.4.pre: dealii::SolutionTransfer<2, dealii::Vector<double>, dealii::hp::DoFHandler<2, 2> >::prepare_for_coarsening_and_refinement(std::vector<dealii::Vector<double>, std::allocator<dealii::Vector<double> > > const&)
#4  /home/claire/FEM/dealii_update/lib/libdeal_II.g.so.8.4.pre: dealii::SolutionTransfer<2, dealii::Vector<double>, dealii::hp::DoFHandler<2, 2> >::prepare_for_coarsening_and_refinement(dealii::Vector<double> const&)
#5  /home/claire/FEM/dealii_update/project_apsc/v7/apsc_v7: apsc_v7::HeatEquation<2>::refine_mesh(unsigned int, unsigned int)
#6  /home/claire/FEM/dealii_update/project_apsc/v7/apsc_v7: apsc_v7::HeatEquation<2>::run()
#7  /home/claire/FEM/dealii_update/project_apsc/v7/apsc_v7: main
--------------------------------------------------------


I am using the latest version of deal.II from the git repository.

This error message appears at the first time step where some cells have been tagged for coarsening.
(The first three steps run fine but only involve refinement)

I think that the problem comes from the fact that I use FE_Nothing elements, because when I use only FE_Q (of different orders) elements the simulations run normally.
However I got from this mailing list that the issue has been dealt with before :
https://groups.google.com/forum/#!topic/dealii/pCj05PmL4Z0
Actually the test case proposed by K. Bzowski is running just fine.
That is why I feel a little bit clueless about the origin of this error message.

Thank you very much by advance for all the help, clues, suggestions that I may receive.

Regards,

Claire

Wolfgang Bangerth

unread,
Jul 31, 2015, 12:34:18 PM7/31/15
to dea...@googlegroups.com

Claire,
I think I have an idea of where the problem comes from, but it's always
better to have a testcase. Do you think you could come up with a small
program that shows the problem by essentially stripping everything you
don't need to show the problem from your current code?

That would make finding and fixing the problem much simpler!

Best
W.


--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@math.tamu.edu
www: http://www.math.tamu.edu/~bangerth/

Claire

unread,
Aug 3, 2015, 6:35:57 AM8/3/15
to deal.II User Group
Dear Pr. Bangerth,

First of all, thank you very much for your answer.

Here is a piece of code generating the same error message that I get.

It is actually the same test case as the one written by K. Bzowski to which I added a second step.
This step is actually just the coarsening of the cells previously refined.

I hope that it is this little program is actually what you need,

kind regards,

Claire
test_case.cc

Claire

unread,
Sep 22, 2015, 10:09:59 AM9/22/15
to deal.II User Group
Dear all,

I would just like to know if there have been any improvement on the issue that I highlighted.
I would also like to know if the test case that I sent is suitable for the debugging.
It case it is not, feel free to tell me what is missing and I will be glad to improve it.

Thank you very much by advance,

Regards,

Claire


Wolfgang Bangerth

unread,
Sep 23, 2015, 8:25:58 AM9/23/15
to dea...@googlegroups.com

> I would just like to know if there have been any improvement on the issue that
> I highlighted.
> I would also like to know if the test case that I sent is suitable for the
> debugging.
> It case it is not, feel free to tell me what is missing and I will be glad to
> improve it.

Like I mentioned in my previous answer, I think I know where the problem may
lie, but it would be very useful to have a small program that doesn't do
anything useful other than demonstrate the problem you are seeing. It's very
difficult to debug problems if you don't have a program that shows the problem.

Claire

unread,
Sep 23, 2015, 9:42:30 AM9/23/15
to deal.II User Group
Dear Pr. Bangerth,

Is the test case attached to the message that I posted on August 3rd suitable?
In case there was a problem and you could not get it, I enclose it to the present post.
It actually returns the same error message that I get when doing the full simulations.
Please, let me know if it is the kind of program you expect.


Thank you very much by advance,

Kind regards,

Claire
test_case.cc

Wolfgang Bangerth

unread,
Sep 23, 2015, 2:29:57 PM9/23/15
to dea...@googlegroups.com
On 09/23/2015 08:42 AM, Claire wrote:
>
> Is the test case attached to the message that I posted on August 3rd
> suitable?
> In case there was a problem and you could not get it, I enclose it to
> the present post.
> It actually returns the same error message that I get when doing the
> full simulations.
> Please, let me know if it is the kind of program you expect.

Claire -- I must have missed the attachment when you sent it the first
time around. It is actually an excellent testcase, short and concise and
shows exactly what the problem is.

I've opened a bug report here:
https://github.com/dealii/dealii/issues/1658
The fix is here:
https://github.com/dealii/dealii/pull/1660

I would like to credit you in part with this fix, but I only know your
first name. Would you offer your full name so I can list it in the
changelog?

Best
Wolfgang

claire br

unread,
Sep 24, 2015, 2:37:59 AM9/24/15
to dea...@googlegroups.com
Dear Pr. Bangerth,

Thank you very much for having taken care of this issue.
I am glad if I was helpful in some way to the deal.ii library.

My full name is Claire Bruna-Rosso.

Regards,

Claire


                                www: http://www.math.tamu.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 a topic in the Google Groups "deal.II User Group" group.
To unsubscribe from this topic, visit https://groups.google.com/d/topic/dealii/ETuJ6V4f5oM/unsubscribe.
To unsubscribe from this group and all its topics, send an email to dealii+un...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Wolfgang Bangerth

unread,
Sep 24, 2015, 7:29:07 AM9/24/15
to dea...@googlegroups.com
On 09/24/2015 01:37 AM, claire br wrote:
>
> Thank you very much for having taken care of this issue.
> I am glad if I was helpful in some way to the deal.ii library.

Definitely. Small testcases make it so much simpler to find and fix bugs.
Yours was excellent in this regard. You deserve the credit to be mentioned in
the file that lists changes since the last release.


> My full name is Claire Bruna-Rosso.

Thanks -- I've sent a patch to update your name:
https://github.com/dealii/dealii/pull/1665/files

Best
W.
Reply all
Reply to author
Forward
0 new messages