Re: [deal.II] Assertion error when using FE_Nothing with FESystems with different number of base elements

2024-04-08 Thread Wolfgang Bangerth

On 4/6/24 16:24, Mohamad Ghadban wrote:

*
*
But when I changed the first FESystem into:
*FESystem[FE_Q(1)^2-**FE_Q(1)^1**]*
it worked fine because the number of base elements in both FESystems is now 
the same. My question is then, is this assertion necessary?


Julian's suggestion works, but if you have two variables that are defined 
everywhere and one that is defined in only part of the domain, then what you 
are suggesting (working with only three components) makes sense to me.


You've already found how to make that work without getting the assertion. As 
to the question of whether the assertion is necessary: Perhaps not, but the 
function in question (FESystem::compare_for_domination()) might be 
substantially more complicated to implement. The way it is implemented is to 
assume that the two elements you are comparing have the same number of base 
elements and multiplicities, and in that case determining which element 
dominates is a question of determining which base element dominates. If you 
wanted to remove the assertion (patches as always welcome!) you have to break 
things into smaller components and compare those for domination. Here, if you have

  FESystem[FE_Q(1)^3]
and
  FESystem[FE_Q(1)^2-FE_Nothing^1]
that would not be too difficult. I bet that could be done in 25 extra lines of 
code (again: patches are very much welcome!) but the difficult case is you have

  FESystem[FE_Q(1)^3]
and
  FESystem[FESystem[FE_Q(1),FE_Q(1)]-FE_Nothing^1]
for example.

Best
 Wolfgang

--

Wolfgang Bangerth  email: bange...@colostate.edu
   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 dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/32538aef-9787-40ea-af0e-f230177e0c45%40colostate.edu.


Re: [deal.II] Re: Assertion error when using FE_Nothing with FESystems with different number of base elements

2024-04-08 Thread Mohamad Ghadban
I think it makes sense. What you are saying is that the number of
components would be summed across all domains regardless of whether they
represent the same variables or not, is that correct?
So if, hypothetically speaking, I am solving the same variables u1, u2, and
u3 in five different domains, that means that the total number of
components will be 15, right?
If that's the case, do you foresee any troubles with the workaround I used,
i.e.,

*FESystem[FE_Q(1)^2-**FE_Q(1)^1**]*  for the first part of the domain
*FESystem[FE_Q(1)^2-FE_Nothing^1]* for the second part of the domain

given that it worked fine where the solution I got was as expected? I'm
just unsure if that would cause issues in more complicated cases.

Thanks,
Mohamad

On Mon, Apr 8, 2024 at 3:18 PM julian...@gmail.com 
wrote:

> Hi Mohamad,
>
> Yes, you understood correctly.
> In your example you should have 5 components.
> Basically with FENothing you have placeholders for the components that are
> not used on the current domain.
> In your example you have the solution components:
> - in domain 1: u_1^1, u_2^1, u_3^1, _, _
> - in domain 2: _, _, _, u_1^2, u_2^2.
> In your code the _ will be replaced by FENothing components and I used
> ^1/^2 to show on which domain the component is defined.
> I hope this helps.
>  Kind regards,
>
> Julian
>
>
> On Monday, April 8, 2024 at 11:07:52 PM UTC+2 mgha...@ualberta.ca wrote:
>
>> Hi Julian
>>
>> Thanks for the tips. My code is actually based on Step 46, but the
>> difference is that in my case I am not solving for a completely new set of
>> variables in each part of the domain.
>> Based on what you are suggesting, wouldn't that mean that the total
>> number of variables/components that I would be solving for is 5 (by summing
>> the multiplicities) instead of 3? That's at least what I inferred from Step
>> 46.
>> Also, is there a more intuitive approach to follow? Because your
>> suggestion means that variables u1 and u2 belong to the first base element
>> in this FESystem
>>
>> FESystem[FE_Q(1)^3-FE_Nothing^2]
>>
>> but they belong to a completely different base element in the other
>> FESystem
>>
>> FESystem[FE_Nothing^3-FE_Q(1)^2]
>>
>> I think the way the logic is conveyed in Step 46 is more intuitive
>> because Stokes variables which exist in the first base element in one
>> FESystem still exist in the first base base element in the other FESystem
>> but are represented by FE_Nothing instead of FE_Q elements.
>>
>> Thanks,
>>  Mohamad
>>
>> On Sun, Apr 7, 2024 at 11:52 PM julian...@gmail.com 
>> wrote:
>>
>>> Hi Mohamad,
>>>
>>> You need to create the two FESystems as follows:
>>> FESystem[FE_Q(1)^3-FE_Nothing^2] and FESystem[FE_Nothing^3-FE_Q(1)^2].
>>> This way the first three variables u1, u2 and u3 will be only defined in
>>> the first domain and in the second domain we define u1 and u2 as the last
>>> two variables.
>>> For more information on FE_Nothing look at step 46 of the deal.II
>>> tutorials.
>>> Additionally, I also have some code online that uses FE_Nothing:
>>> https://github.com/mathmerizing/TemporalMultirateFEM/blob/main/src%2Fmonolithic_heatwave%2Fmain.cc
>>> Kind regards,
>>>
>>> Julian
>>> On Sunday, April 7, 2024 at 12:24:53 AM UTC+2 mgha...@ualberta.ca wrote:
>>>
 Hello,

 I am solving a problem in which I use the FE_Nothing class to represent
 variables that do not exist in parts of the physical domain. My domain is
 split into two parts: In the first part, I solve for all three variables
 (I'll call them u1, u2, and u3) and in the second part of the domain I
 solve for u1 and u2 only.
 For my problem, I use two FESystem objects in the parameter file as
 follows:
 *FESystem[FE_Q(1)^3]* for the first part of the domain
 *FESystem[FE_Q(1)^2-FE_Nothing^1]* for the second part of the domain.

 However, this results in the following assertion error:













 *An error occurred in line <2478> of file
 <../deal.II/src/source/fe/fe_system.cc> in function
 dealii::FiniteElementDomination::Domination dealii::FESystem>>> spacedim>::compare_for_domination(const dealii::FiniteElement>>> spacedim>&, unsigned int) const [with int dim = 2; int spacedim = 2]The
 violated condition was: this->n_base_elements() ==
 fe_sys_other->n_base_elements()Additional information: You are trying
 to use functionality in deal.II that is currently notimplemented. In
 many cases, this indicates that there simply didn'tappear much of a
 need for it, or that the author of the original codedid not have the
 time to implement a particular case. If you hit thisexception, it is
 therefore worth the time to look into the code tofind out whether you
 may be able to implement the missingfunctionality. If you do, please
 consider providing a patch to thedeal.II development sources (see the
 deal.II website on how to

Re: [deal.II] Re: Assertion error when using FE_Nothing with FESystems with different number of base elements

2024-04-08 Thread julian...@gmail.com
Hi Mohamad,

Yes, you understood correctly.
In your example you should have 5 components.
Basically with FENothing you have placeholders for the components that are 
not used on the current domain.
In your example you have the solution components:
- in domain 1: u_1^1, u_2^1, u_3^1, _, _
- in domain 2: _, _, _, u_1^2, u_2^2.
In your code the _ will be replaced by FENothing components and I used 
^1/^2 to show on which domain the component is defined.
I hope this helps.
 Kind regards,

Julian


On Monday, April 8, 2024 at 11:07:52 PM UTC+2 mgha...@ualberta.ca wrote:

> Hi Julian
>
> Thanks for the tips. My code is actually based on Step 46, but the 
> difference is that in my case I am not solving for a completely new set of 
> variables in each part of the domain.
> Based on what you are suggesting, wouldn't that mean that the total number 
> of variables/components that I would be solving for is 5 (by summing the 
> multiplicities) instead of 3? That's at least what I inferred from Step 46.
> Also, is there a more intuitive approach to follow? Because your 
> suggestion means that variables u1 and u2 belong to the first base element 
> in this FESystem
>
> FESystem[FE_Q(1)^3-FE_Nothing^2]
>
> but they belong to a completely different base element in the other 
> FESystem
>  
> FESystem[FE_Nothing^3-FE_Q(1)^2]
>
> I think the way the logic is conveyed in Step 46 is more intuitive because 
> Stokes variables which exist in the first base element in one FESystem 
> still exist in the first base base element in the other FESystem but are 
> represented by FE_Nothing instead of FE_Q elements.
>
> Thanks,
>  Mohamad
>
> On Sun, Apr 7, 2024 at 11:52 PM julian...@gmail.com  
> wrote:
>
>> Hi Mohamad,
>>
>> You need to create the two FESystems as follows:
>> FESystem[FE_Q(1)^3-FE_Nothing^2] and FESystem[FE_Nothing^3-FE_Q(1)^2].
>> This way the first three variables u1, u2 and u3 will be only defined in 
>> the first domain and in the second domain we define u1 and u2 as the last 
>> two variables.
>> For more information on FE_Nothing look at step 46 of the deal.II 
>> tutorials.
>> Additionally, I also have some code online that uses FE_Nothing: 
>> https://github.com/mathmerizing/TemporalMultirateFEM/blob/main/src%2Fmonolithic_heatwave%2Fmain.cc
>> Kind regards,
>>
>> Julian
>> On Sunday, April 7, 2024 at 12:24:53 AM UTC+2 mgha...@ualberta.ca wrote:
>>
>>> Hello,
>>>
>>> I am solving a problem in which I use the FE_Nothing class to represent 
>>> variables that do not exist in parts of the physical domain. My domain is 
>>> split into two parts: In the first part, I solve for all three variables 
>>> (I'll call them u1, u2, and u3) and in the second part of the domain I 
>>> solve for u1 and u2 only. 
>>> For my problem, I use two FESystem objects in the parameter file as 
>>> follows:
>>> *FESystem[FE_Q(1)^3]* for the first part of the domain
>>> *FESystem[FE_Q(1)^2-FE_Nothing^1]* for the second part of the domain.
>>>
>>> However, this results in the following assertion error:
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>> *An error occurred in line <2478> of file 
>>> <../deal.II/src/source/fe/fe_system.cc> in function
>>> dealii::FiniteElementDomination::Domination dealii::FESystem>> spacedim>::compare_for_domination(const dealii::FiniteElement>> spacedim>&, unsigned int) const [with int dim = 2; int spacedim = 2]The 
>>> violated condition was: this->n_base_elements() == 
>>> fe_sys_other->n_base_elements()Additional information: You are trying 
>>> to use functionality in deal.II that is currently notimplemented. In 
>>> many cases, this indicates that there simply didn'tappear much of a 
>>> need for it, or that the author of the original codedid not have the 
>>> time to implement a particular case. If you hit thisexception, it is 
>>> therefore worth the time to look into the code tofind out whether you 
>>> may be able to implement the missingfunctionality. If you do, please 
>>> consider providing a patch to thedeal.II development sources (see the 
>>> deal.II website on how tocontribute).*
>>>
>>> But when I changed the first FESystem into:
>>> *FESystem[FE_Q(1)^2-**FE_Q(1)^1**]* 
>>> it worked fine because the number of base elements in both FESystems is 
>>> now the same. My question is then, is this assertion necessary?
>>>
>>> Thank you,
>>> Mohamad
>>>
>> -- 
>> 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/csAJHjYuK88/unsubscribe.
>> To unsubscribe from this group and all its topics, send an email to 
>> dealii+un...@googlegroups.com.
>> To view this discussion on the web visit 
>> 

Re: [deal.II] Re: Assertion error when using FE_Nothing with FESystems with different number of base elements

2024-04-08 Thread Mohamad Ghadban
Hi Julian

Thanks for the tips. My code is actually based on Step 46, but the
difference is that in my case I am not solving for a completely new set of
variables in each part of the domain.
Based on what you are suggesting, wouldn't that mean that the total number
of variables/components that I would be solving for is 5 (by summing the
multiplicities) instead of 3? That's at least what I inferred from Step 46.
Also, is there a more intuitive approach to follow? Because your suggestion
means that variables u1 and u2 belong to the first base element in this
FESystem

FESystem[FE_Q(1)^3-FE_Nothing^2]

but they belong to a completely different base element in the other FESystem

FESystem[FE_Nothing^3-FE_Q(1)^2]

I think the way the logic is conveyed in Step 46 is more intuitive because
Stokes variables which exist in the first base element in one FESystem
still exist in the first base base element in the other FESystem but are
represented by FE_Nothing instead of FE_Q elements.

Thanks,
 Mohamad

On Sun, Apr 7, 2024 at 11:52 PM julian...@gmail.com 
wrote:

> Hi Mohamad,
>
> You need to create the two FESystems as follows:
> FESystem[FE_Q(1)^3-FE_Nothing^2] and FESystem[FE_Nothing^3-FE_Q(1)^2].
> This way the first three variables u1, u2 and u3 will be only defined in
> the first domain and in the second domain we define u1 and u2 as the last
> two variables.
> For more information on FE_Nothing look at step 46 of the deal.II
> tutorials.
> Additionally, I also have some code online that uses FE_Nothing:
> https://github.com/mathmerizing/TemporalMultirateFEM/blob/main/src%2Fmonolithic_heatwave%2Fmain.cc
> Kind regards,
>
> Julian
> On Sunday, April 7, 2024 at 12:24:53 AM UTC+2 mgha...@ualberta.ca wrote:
>
>> Hello,
>>
>> I am solving a problem in which I use the FE_Nothing class to represent
>> variables that do not exist in parts of the physical domain. My domain is
>> split into two parts: In the first part, I solve for all three variables
>> (I'll call them u1, u2, and u3) and in the second part of the domain I
>> solve for u1 and u2 only.
>> For my problem, I use two FESystem objects in the parameter file as
>> follows:
>> *FESystem[FE_Q(1)^3]* for the first part of the domain
>> *FESystem[FE_Q(1)^2-FE_Nothing^1]* for the second part of the domain.
>>
>> However, this results in the following assertion error:
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>> *An error occurred in line <2478> of file
>> <../deal.II/src/source/fe/fe_system.cc> in function
>> dealii::FiniteElementDomination::Domination dealii::FESystem> spacedim>::compare_for_domination(const dealii::FiniteElement> spacedim>&, unsigned int) const [with int dim = 2; int spacedim = 2]The
>> violated condition was: this->n_base_elements() ==
>> fe_sys_other->n_base_elements()Additional information: You are trying
>> to use functionality in deal.II that is currently notimplemented. In
>> many cases, this indicates that there simply didn'tappear much of a
>> need for it, or that the author of the original codedid not have the
>> time to implement a particular case. If you hit thisexception, it is
>> therefore worth the time to look into the code tofind out whether you
>> may be able to implement the missingfunctionality. If you do, please
>> consider providing a patch to thedeal.II development sources (see the
>> deal.II website on how tocontribute).*
>>
>> But when I changed the first FESystem into:
>> *FESystem[FE_Q(1)^2-**FE_Q(1)^1**]*
>> it worked fine because the number of base elements in both FESystems is
>> now the same. My question is then, is this assertion necessary?
>>
>> Thank you,
>> Mohamad
>>
> --
> 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/csAJHjYuK88/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/0f0b1e5b-5e10-40b9-a5ed-bf6eba6f4d4en%40googlegroups.com
> 
> .
>

-- 
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/CAP2_OwEcTkCjE54g7QP54PijhcGci4fZfuCFhLH3Ah1Ngw2-2w%40mail.gmail.com.


[deal.II] deal.II Newsletter #280

2024-04-08 Thread 'Rene Gassmoeller' via deal.II User Group
Hello everyone!

This is deal.II newsletter #280.
It automatically reports recently merged features and discussions about the 
deal.II finite element library.


## Below you find a list of recently proposed or merged features:

#16871: fix maybe-unused warning in grid_tools_dofhandler (proposed by tjhei) 
https://github.com/dealii/dealii/pull/16871

#16870: Update documentation for GridTools::distort_random(). (proposed by 
bangerth) https://github.com/dealii/dealii/pull/16870

#16869: Add DoFTools::extract_level_constant_modes() (proposed by peterrum) 
https://github.com/dealii/dealii/pull/16869

#16868: Remove redundant defaulted copy assignment operator for 
VectorizedArrayIterator (proposed by masterleinad) 
https://github.com/dealii/dealii/pull/16868

#16866: Small fixes in step-89 (proposed by peterrum; merged) 
https://github.com/dealii/dealii/pull/16866

#16865: Vectorize class Tensor without requiring alignment or padding. 
(proposed by bangerth) https://github.com/dealii/dealii/pull/16865

#16864: Simplify some functions by just deferring to member functions. 
(proposed by bangerth; merged) https://github.com/dealii/dealii/pull/16864

#16863: Unclutter step-40's run() function. (proposed by bangerth; merged) 
https://github.com/dealii/dealii/pull/16863

#16862: In step-40, move timer section into the function being timed. (proposed 
by bangerth; merged) https://github.com/dealii/dealii/pull/16862

#16859: Remove debug output (proposed by peterrum; merged) 
https://github.com/dealii/dealii/pull/16859

#16858: Test: cgal/cgal_mesh_criteria: disable floating point exceptions 
(proposed by tamiko; merged) https://github.com/dealii/dealii/pull/16858

#16857: Boost 1.80 (proposed by bangerth) 
https://github.com/dealii/dealii/pull/16857

#16856: Fix MG transfer for Tetrahedrons  (proposed by dominiktassilostill) 
https://github.com/dealii/dealii/pull/16856

#16855: Avoid std::pow in BarycentricPolynomials (proposed by kronbichler; 
merged) https://github.com/dealii/dealii/pull/16855

#16854:  NonMatching::MappingInfo: Avoid memory allocations in reinit() 
(proposed by kronbichler) https://github.com/dealii/dealii/pull/16854

#16853: Add Quadrature::initialize() function taking ArrayView (proposed by 
kronbichler; merged) https://github.com/dealii/dealii/pull/16853

#16852: MappingQ: Remove outdated documentation (proposed by kronbichler; 
merged) https://github.com/dealii/dealii/pull/16852

#16851: Import BOOST 1.83. (proposed by bangerth) 
https://github.com/dealii/dealii/pull/16851

#16850: Import BOOST 1.85 beta 1 (proposed by bangerth) 
https://github.com/dealii/dealii/pull/16850

#16849: Test fe/fe_values_view_{07|10|13|16}: make Assert more robust (proposed 
by tamiko; merged) https://github.com/dealii/dealii/pull/16849

#16848: Use vector intrinsics for tensor operations. (proposed by bangerth) 
https://github.com/dealii/dealii/pull/16848

#16847: MappingQ: Favor get_vertices() over compute_mapping_support_points() if 
possible (proposed by kronbichler; merged) 
https://github.com/dealii/dealii/pull/16847

#16845: Pass some SmartPointer variables by reference (proposed by kronbichler; 
merged) https://github.com/dealii/dealii/pull/16845

#16844: Use ArrayView for interface to tensor-product point evaluator (proposed 
by kronbichler; merged) https://github.com/dealii/dealii/pull/16844

#16843: VectorTools::interpolate: Make common case of vector problem faster 
(proposed by kronbichler; merged) https://github.com/dealii/dealii/pull/16843

#16841: Augment documentation of Utilities::type_to_string(). (proposed by 
bangerth; merged) https://github.com/dealii/dealii/pull/16841

#16840: Vectorization type trait (proposed by tamiko; merged) 
https://github.com/dealii/dealii/pull/16840

#16839: Make tensor*tensor implementation easier to read. (proposed by 
bangerth; merged) https://github.com/dealii/dealii/pull/16839

#16838: Mark DEAL_II_ALWAYS_INLINE functions as 'inline'. (proposed by 
bangerth) https://github.com/dealii/dealii/pull/16838

#16837: Provide an implementation for VectorizedArray::dot_product(). (proposed 
by bangerth; merged) https://github.com/dealii/dealii/pull/16837

#16836: Minor updates to VectorizedArray. (proposed by bangerth; merged) 
https://github.com/dealii/dealii/pull/16836

#16835: Fix ensure_single_trailing_newline() on macs. (proposed by drwells; 
merged) https://github.com/dealii/dealii/pull/16835

#16834: Export compile definitions for dependent targets (proposed by 
gassmoeller) https://github.com/dealii/dealii/pull/16834

#16833: Switch SolverFGMRES to new delayed classical Gram-Schmidt (proposed by 
kronbichler; merged) https://github.com/dealii/dealii/pull/16833

#16832: Fix naming of output files. (proposed by bangerth; merged) 
https://github.com/dealii/dealii/pull/16832

#16831: Fix two places where TriaAccessor::n_vertices is called repeatedly 
(proposed by kronbichler; merged) https://github.com/dealii/dealii/pull/16831

#16830: MappingCartesian: Some 

Re: [deal.II] gridin.read_exodusii to import exodusii mesh in dealii 9.5.1 fails

2024-04-08 Thread Sina Tajfirooz
Thanks!

On Friday, March 29, 2024 at 2:04:08 PM UTC+1 d.arnd...@gmail.com wrote:

> Sina,
>
> You'll have to provide the path directly as a function argument, i.e.,
>
> gridin.read_exodusii(case_path + "mesh.exo");
>
> Best,
> Daniel
>
>
> On Fri, Mar 29, 2024 at 7:32 AM Sina Tajfirooz  
> wrote:
> >
> > Dear developers/users,
> >
> > How can I import an exodusii mesh in dealii version 9.5.1, which is 
> configured with Trilinos?
> >
> > doing
> >
> > "
> > gridin.attach_triangulation(triangulation);
> > std::ifstream f(case_path + "mesh.exo");
> > gridin.read_exodusii(f);
> >
> > "
> >
> > raises the error:
> > "
> > error: no matching function for call to 
> ‘dealii::GridIn<3>::read_exodusii(std::ifstream&)’
> > gridin.read_exodusii(f);
> > "
> >
> > How can I resolve this?
> >
> > Best,
> > Sina
> >
> > --
> > 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+un...@googlegroups.com.
> > To view this discussion on the web visit 
> https://groups.google.com/d/msgid/dealii/2ea642fb-b3fb-480f-88e7-2b2c96a264e2n%40googlegroups.com
> .
>

-- 
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/c2a5f6ea-7df1-43f5-9b36-338fe9075304n%40googlegroups.com.