Wasim,

Why do you have two different triangulations? Can you have single
Triangulation and use FE_Nothing ? For example Step-70 solves two different
equations on different parts of the domain using a single Triangulation.

Best,

Bruno

Le lun. 21 nov. 2022 à 15:27, Wasim Niyaz Munshi ce21d400 <
ce21d...@smail.iitm.ac.in> a écrit :

> Thank you, Prof. Bruno.
>
> I used the above line of code and realized where I was going wrong. When
> solving the damage equation, I also call a function (let's call it H). This
> function is needed to find the element stiffness matrix of the damage
> equation.
> This function returns some value for a  given quadrature point, cell in
> the damage mesh and the elasticity solution vector. But, this function H
> computes that quantity by using the fe_elasticity object.
> *The error, according to me, comes in this line in H function :
> fe_values_elastic.reinit(cell);*
> I need to *reinit *the fe_values_elasticity for each cell, but the
> problem is that the *cell *I am passing as an argument lives on the
> traingulation_damage while fe_values_elasticity lives on
> triangulation_elasticity.
>
>
>
> On Mon, Nov 21, 2022 at 9:43 PM Bruno Turcksin <bruno.turck...@gmail.com>
> wrote:
>
>> Wasim,
>>
>> Before the line that errors out can you do
>>
>> std::cout<<cell->get_fe().get_name()<<std::endl;
>>
>> and check that it matches what you have set in FEValues.
>>
>> Best,
>>
>> Bruno
>>
>> Le lun. 21 nov. 2022 à 11:08, Wasim Niyaz Munshi ce21d400 <
>> ce21d...@smail.iitm.ac.in> a écrit :
>>
>>> I am working with 2 DoFHandlers. I checked my code. I am using the
>>> correct DoFHandlers for the 2 equations.
>>>
>>> On Mon, Nov 21, 2022 at 6:58 PM Bruno Turcksin <bruno.turck...@gmail.com>
>>> wrote:
>>>
>>>> Wasim,
>>>>
>>>> It's hard to say without seeing any code. Are you working with one or
>>>> two DoFHandler? If you are working with two DoFHandler, are you sure that
>>>> you are using the correct one?
>>>>
>>>> Best,
>>>>
>>>> Bruno
>>>>
>>>> On Monday, November 21, 2022 at 5:54:22 AM UTC-5
>>>> ce21...@smail.iitm.ac.in wrote:
>>>>
>>>>> Hello everyone.
>>>>>
>>>>> I am trying to solve elasticity and *laplace (damage)* equations in
>>>>> the same program. The rhs of laplace equation depends on the elasticity
>>>>> solution vector. Both equations have already been solved in step3 and 8. I
>>>>> am following the same approach. I am creating my objects (like
>>>>> triangulation_elasticity, fe_elasticity,dof_handler_elasticity) for
>>>>> elasticity as in step8 and solving the elasticity equation. Everything
>>>>> works fine till this point.
>>>>>
>>>>>
>>>>> Then I tried to solve the laplace equation using the elasticity
>>>>> solution. I followed step 3 approach and created the corresponding
>>>>> objects(triangulation_damage, fe_damage,dof_handler_damage).
>>>>>
>>>>> However, I am getting the following error:
>>>>>
>>>>> An error occurred in line <4310> of file
>>>>> </home/wasim/dealii-candi/tmp/unpack/deal.II-v9.4.0/source/fe/fe_values.cc>
>>>>> in function
>>>>>     void dealii::FEValues<dim, spacedim>::reinit(const
>>>>> dealii::TriaIterator<dealii::DoFCellAccessor<dim, spacedim, lda> >&) [with
>>>>> bool level_dof_access = false; int dim = 2; int spacedim = 2]
>>>>> The violated condition was:
>>>>>     static_cast<const FiniteElementData<dim> &>(*this->fe) ==
>>>>> static_cast<const FiniteElementData<dim> &>(cell->get_fe())
>>>>> Additional information:
>>>>>
>>>>>
>>>>> * The FiniteElement you provided to FEValues and the FiniteElement
>>>>> that    belongs to the DoFHandler that provided the cell iterator do not
>>>>> match*.
>>>>>
>>>>> Stacktrace:
>>>>> -----------
>>>>> #0  /home/wasim/dealii-candi/deal.II-v9.4.0/lib/libdeal_II.g.so.9.4.0:
>>>>> void dealii::FEValues<2,
>>>>> 2>::reinit<false>(dealii::TriaIterator<dealii::DoFCellAccessor<2, 2, 
>>>>> false>
>>>>> > const&)
>>>>> #1  ./step-200: float
>>>>> step200::PhaseField::H_plus<dealii::TriaActiveIterator<dealii::DoFCellAccessor<2,
>>>>> 2, false> > >(dealii::Vector<double>,
>>>>> dealii::TriaActiveIterator<dealii::DoFCellAccessor<2, 2, false> >, 
>>>>> unsigned
>>>>> int)
>>>>> #2  ./step-200: step200::PhaseField::assemble_system_damage()
>>>>> #3  ./step-200: step200::PhaseField::damage_mesh()
>>>>> #4  ./step-200: main
>>>>> ---------------------
>>>>>
>>>>> It says that in assemble_damage, the finite element provided to
>>>>> fevalues and that of dof_handler are different.
>>>>> I don't know what this exactly means as I am only creating a single
>>>>> object,* FE_Q<2>          fe_damage* for solving my damage equation.
>>>>>
>>>>> Thanks and regards
>>>>> Wasim
>>>>>
>>>> --
>>>> 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/cV63qxefdTM/unsubscribe.
>>>> To unsubscribe from this group and all its topics, send an email to
>>>> dealii+unsubscr...@googlegroups.com.
>>>> To view this discussion on the web visit
>>>> https://groups.google.com/d/msgid/dealii/5f30d3ae-264a-4ff0-ba85-27e6959ecc92n%40googlegroups.com
>>>> <https://groups.google.com/d/msgid/dealii/5f30d3ae-264a-4ff0-ba85-27e6959ecc92n%40googlegroups.com?utm_medium=email&utm_source=footer>
>>>> .
>>>>
>>> --
>>> 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/cV63qxefdTM/unsubscribe.
>>> To unsubscribe from this group and all its topics, send an email to
>>> dealii+unsubscr...@googlegroups.com.
>>> To view this discussion on the web visit
>>> https://groups.google.com/d/msgid/dealii/CAM8ps5C1GHxcm18SM8Zyr9hMQvExX-j9tj50NYvg7hN%3D9JZ4tw%40mail.gmail.com
>>> <https://groups.google.com/d/msgid/dealii/CAM8ps5C1GHxcm18SM8Zyr9hMQvExX-j9tj50NYvg7hN%3D9JZ4tw%40mail.gmail.com?utm_medium=email&utm_source=footer>
>>> .
>>>
>> --
>> 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/cV63qxefdTM/unsubscribe.
>> To unsubscribe from this group and all its topics, send an email to
>> dealii+unsubscr...@googlegroups.com.
>> To view this discussion on the web visit
>> https://groups.google.com/d/msgid/dealii/CAGVt9eOZCVvceE4c_EHsqqtbWPAAB74n8gifS1XzNKE0qjPyqQ%40mail.gmail.com
>> <https://groups.google.com/d/msgid/dealii/CAGVt9eOZCVvceE4c_EHsqqtbWPAAB74n8gifS1XzNKE0qjPyqQ%40mail.gmail.com?utm_medium=email&utm_source=footer>
>> .
>>
> --
> 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/cV63qxefdTM/unsubscribe.
> To unsubscribe from this group and all its topics, send an email to
> dealii+unsubscr...@googlegroups.com.
> To view this discussion on the web visit
> https://groups.google.com/d/msgid/dealii/CAM8ps5B%3D%2BF3XSCpKcg%3DMEio9JHaNmxsBo-gjKkxZUe8U0r%3D0CQ%40mail.gmail.com
> <https://groups.google.com/d/msgid/dealii/CAM8ps5B%3D%2BF3XSCpKcg%3DMEio9JHaNmxsBo-gjKkxZUe8U0r%3D0CQ%40mail.gmail.com?utm_medium=email&utm_source=footer>
> .
>

-- 
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/CAGVt9ePPz2QShcmv8Nu3YajrO7e-o6nGn3OVSN%3D7hjCt0rAvsg%40mail.gmail.com.

Reply via email to