Hi Kaushik,

1) Yes, this is possible, but tricky: `SolutionTransfer` is not capable of 
this feature, and you need to do it manually with the more versatile class 
`CellDataTransfer`. A way to do it has been discussed in #11132 
<https://github.com/dealii/dealii/pull/11132>.

2) I did not know you were trying to interpolate a FENothing element into a 
FEQ element. This should not be possible, as you can not interpolate 
information from simply 'nothing', and some assertion should be triggered 
while trying to do so. The other way round should be possible, i.e., 
interpolation from FEQ to FENothing, since you will simply 'forget' what 
has been on the old cell.

Did you run your program in debug or release mode? If you ran it in debug 
mode without an assertion being triggered, please tell me. The user should 
be warned that they are not allowed to interpolate from a FENothing object 
to a FEQ element (or any element with nodes).

I think what is happening here is that we initialize some container with 
zeros by default, and accidentally overwrite the corresponding dof values. 
I suspect that if you assign both topmost cells to FEQ and the lower ones 
to FENothing, your solution would look okay due to the way we iterate over 
cells (Z-order). Or in short, we first unpack the lower cells, and then 
unpack the upper ones, which may overwrite dof values. This is undefined 
behavior, and we should warn the user about that.

Best,
Marc

On Sunday, December 27, 2020 at 2:28:50 PM UTC-7 k.d...@gmail.com wrote:

> Hello Marc,
> Thank you very much. 
> I have modified my test code as you suggested and it is working fine now. 
> That the code is attached. I have a few more questions that I added to the 
> attached PNG file below along with the results from the test code. 
> 1. Is it possible to specify an initial value on dofs that are activated 
> then a FE_Nothing cell becomes a FE_Q cell?
> 2. What happens to the dofs that were shared between a FE_Nothing cell and 
> a FE_Q cell before the p-refinement, and are shared between two FE_Q 
> elements after p-refinement? Are those always set to zeros after the 
> refinement? 
>
> [image: image.png]
> Thank you,
> Kaushik 
>
> On Wed, Dec 23, 2020 at 5:35 PM Marc Fehling <mafe...@gmail.com> wrote:
>
>> Hi Kaushik,
>>
>> Be careful on what you are doing here: You prepare your solution to be 
>> transferred on refinement, but at the point where you wish to interpolate 
>> your solution the mesh differs from the one your SolutionTransfer object 
>> expects to encounter. That is because you changed the assigned finite 
>> element between executing the refinement and interpolating the old solution 
>> to the new mesh.
>>
>> You are basically doing two steps here at once: You perform h-refinement 
>> as your first step, and then alter your function space by assigning 
>> different finite elements (p-adaptation). I would suggest you split your 
>> intentions into two separate steps.
>>
>> First, only perform h-refinement on that one cell and interpolate the 
>> entire solution on the new grid. Next, tell your DoFHandler object that you 
>> intend to change the finite element on some of your cells by assigning a 
>> corresponding future finite element index to them (see here 
>> <https://www.dealii.org/developer/doxygen/deal.II/classDoFCellAccessor.html#ae4d4d8562cb47b70b797369b8872b04d>),
>>  
>> prepare your solution for refinement using a SolutionTransfer object once 
>> more, and finally execute refinement again. This should accomplish what you 
>> whish to intend.
>>
>> In addition, there have been issues using p::d::SolutionTransfer objects 
>> with FENothing elements which have been fixed in #10592 
>> <https://github.com/dealii/dealii/pull/10592>. Please incorporate these 
>> upstream fixes into your deal.II library by building it on the current 
>> master branch.
>>
>> Hope this helps!
>>
>> Best,
>> Marc
>>
>> On Wednesday, December 23, 2020 at 2:33:09 PM UTC-7 k.d...@gmail.com 
>> wrote:
>>
>>> Hi Marc:
>>> Thank you again for your help.
>>> I have another problem. 
>>> A small test code is attached. 
>>>
>>> I have one cell of FEQ element. I refine that into four cells and then 
>>> assign FE_q to two of them and FE_nothing to the other two child cells. 
>>> Then when I try to transfer the solution, the code aborts. 
>>>
>>> Is this a limitation? 
>>>
>>> Thank you very much,
>>> Kaushik 
>>>
>>>
>>> On Wed, Dec 9, 2020 at 7:16 PM Kaushik Das <k.d...@gmail.com> wrote:
>>>
>>>> Thank you, Mark. I just built dealii from the source 
>>>> (deal.II-9.3.0-pre). And my little test is passing now. 
>>>> Thanks for the help.
>>>> -Kaushik 
>>>>
>>>> On Wed, Dec 9, 2020 at 4:36 PM Kaushik Das <k.d...@gmail.com> wrote:
>>>>
>>>>> Thank you Mark. 
>>>>> I am using the dealii lib that I got from apt-get from  
>>>>> deal.ii-9.2.0-backports.
>>>>> I used PETSc and the abort was on even 1 cpus. I tried 2, 3, 6 cpus 
>>>>> and all aborted similarly. 
>>>>>
>>>>> I will get the latest master branch and build that. 
>>>>>
>>>>>
>>>>> Thanks,
>>>>> Kaushik
>>>>>
>>>>> On Wed, Dec 9, 2020 at 4:23 PM Marc Fehling <mafe...@gmail.com> wrote:
>>>>>
>>>>>> From your stacktrace I can see you are using PETSc and deal.II 9.2.0 
>>>>>> which already incorporates the specified patch. Would you try to build 
>>>>>> the 
>>>>>> actual master branch anyways?
>>>>>> On Wednesday, December 9, 2020 at 2:11:59 PM UTC-7 Marc Fehling wrote:
>>>>>>
>>>>>>> Hi Kaushik,
>>>>>>>
>>>>>>> I am unable to reproduce your problem with the code you provided on 
>>>>>>> the latest build of deal.II and Trilinos.
>>>>>>>
>>>>>>>    - On how many processes did you run your program?
>>>>>>>    - Did you use PETSc or Trilinos?
>>>>>>>    - Could you try to build deal.II on the latest master branch? 
>>>>>>>    There is a chance that your issue has been solved upstream. Chances 
>>>>>>> are 
>>>>>>>    high that fix #8860 <https://github.com/dealii/dealii/pull/8860> 
>>>>>>>    and the changes made to `get_interpolated_dof_values()` will solve 
>>>>>>> your 
>>>>>>>    problem.
>>>>>>>    
>>>>>>> Marc
>>>>>>> On Wednesday, December 9, 2020 at 7:14:57 AM UTC-7 k.d...@gmail.com 
>>>>>>> wrote:
>>>>>>>
>>>>>>>> Hi Marc and Bruno,
>>>>>>>> I was able to reproduce this abort on an even simpler test. Please 
>>>>>>>> see the attached file. 
>>>>>>>>
>>>>>>>> Initial grid:
>>>>>>>>  /*
>>>>>>>> * -----------
>>>>>>>> * |  0 |  0 |
>>>>>>>> * -----------
>>>>>>>> * |  1 |  1 | 0 - FEQ, 1 - FE_Nothing
>>>>>>>> * -----------
>>>>>>>> */
>>>>>>>>
>>>>>>>> /* Set refine flags:
>>>>>>>> * -----------
>>>>>>>> * |  R |  R |    FEQ
>>>>>>>> * -----------
>>>>>>>> * |      |      |    FE_Nothing
>>>>>>>> * -----------
>>>>>>>> */
>>>>>>>>
>>>>>>>> Then refine and solution trans. During the  
>>>>>>>> execute_coarsening_and_refinement, it aborts. 
>>>>>>>>
>>>>>>>> Here is a stack trace:
>>>>>>>> --------------------------------------------------------
>>>>>>>>
>>>>>>>> An error occurred in line <1167> of file 
>>>>>>>> </build/deal.ii-vFp8uU/deal.ii-9.2.0/include/deal.II/lac/vector.h> in 
>>>>>>>> function
>>>>>>>>     Number& 
>>>>>>>> dealii::Vector<Number>::operator()(dealii::Vector<Number>::size_type) 
>>>>>>>> [with 
>>>>>>>> Number = double; dealii::Vector<Number>::size_type = unsigned int]
>>>>>>>> The violated condition was: 
>>>>>>>>     static_cast<typename ::dealii::internal::argument_type<void( 
>>>>>>>> typename std::common_type<decltype(i), decltype(size())>::type)>:: 
>>>>>>>> type>(i) 
>>>>>>>> < static_cast<typename ::dealii::internal::argument_type<void( 
>>>>>>>> typename 
>>>>>>>> std::common_type<decltype(i), decltype(size())>::type)>:: type>(size())
>>>>>>>> Additional information: 
>>>>>>>>     Index 0 is not in the half-open range [0,0). In the current 
>>>>>>>> case, this half-open range is in fact empty, suggesting that you are 
>>>>>>>> accessing an element of an empty collection such as a vector that has 
>>>>>>>> not 
>>>>>>>> been set to the correct size.
>>>>>>>>
>>>>>>>> Stacktrace:
>>>>>>>> -----------
>>>>>>>> #0  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> dealii::Vector<double>::operator()(unsigned int)
>>>>>>>> #1  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> #2  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> dealii::parallel::distributed::SolutionTransfer<2, 
>>>>>>>> dealii::PETScWrappers::MPI::Vector, dealii::hp::DoFHandler<2, 2> 
>>>>>>>> >::pack_callback(dealii::TriaIterator<dealii::CellAccessor<2, 2> > 
>>>>>>>> >const&, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus)
>>>>>>>> #3  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> dealii::parallel::distributed::SolutionTransfer<2, 
>>>>>>>> dealii::PETScWrappers::MPI::Vector, dealii::hp::DoFHandler<2, 2> 
>>>>>>>> >::register_data_attach()::{lambda(dealii::TriaIterator<dealii::CellAccessor<2,
>>>>>>>> > 
>>>>>>>> 2> > const&, dealii::Triangulation<2, 
>>>>>>>> 2>::CellStatus)#1}::operator()(dealii::TriaIterator<dealii::CellAccessor<2,
>>>>>>>>  
>>>>>>>> 2> > const&, dealii::Triangulation<2, 2>::CellStatus) const
>>>>>>>> #4  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> std::_Function_handler<std::vector<char, std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> > const&, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus), 
>>>>>>>> dealii::parallel::distributed::SolutionTransfer<2, 
>>>>>>>> dealii::PETScWrappers::MPI::Vector, dealii::hp::DoFHandler<2, 2> 
>>>>>>>> >::register_data_attach()::{lambda(dealii::TriaIterator<dealii::CellAccessor<2,
>>>>>>>> > 
>>>>>>>> 2> > const&, dealii::Triangulation<2, 
>>>>>>>> 2>::CellStatus)#1}>::_M_invoke(std::_Any_data const&, 
>>>>>>>> dealii::TriaIterator<dealii::CellAccessor<2, 2> > const&, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus&&)
>>>>>>>> #5  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> std::function<std::vector<char, std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> > const&, 
>>>>>>>> dealii::Triangulation<2, 
>>>>>>>> 2>::CellStatus)>::operator()(dealii::TriaIterator<dealii::CellAccessor<2,
>>>>>>>>  
>>>>>>>> 2> > const&, dealii::Triangulation<2, 2>::CellStatus) const
>>>>>>>> #6  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> std::_Function_handler<std::vector<char, std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> >, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus), 
>>>>>>>> std::function<std::vector<char, 
>>>>>>>> std::allocator<char> > (dealii::TriaIterator<dealii::CellAccessor<2, 
>>>>>>>> 2> > 
>>>>>>>> const&, dealii::Triangulation<2, 2>::CellStatus)> 
>>>>>>>> >::_M_invoke(std::_Any_data const&, 
>>>>>>>> dealii::TriaIterator<dealii::CellAccessor<2, 2> >&&, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus&&)
>>>>>>>> #7  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> std::function<std::vector<char, std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> >, 
>>>>>>>> dealii::Triangulation<2, 
>>>>>>>> 2>::CellStatus)>::operator()(dealii::TriaIterator<dealii::CellAccessor<2,
>>>>>>>>  
>>>>>>>> 2> >, dealii::Triangulation<2, 2>::CellStatus) const
>>>>>>>> #8  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> dealii::parallel::distributed::Triangulation<2, 
>>>>>>>> 2>::DataTransfer::pack_data(std::vector<std::tuple<p4est_quadrant*, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus, 
>>>>>>>> dealii::TriaIterator<dealii::CellAccessor<2, 2> > >, 
>>>>>>>> std::allocator<std::tuple<p4est_quadrant*, dealii::Triangulation<2, 
>>>>>>>> 2>::CellStatus, dealii::TriaIterator<dealii::CellAccessor<2, 2> > > > 
>>>>>>>> > 
>>>>>>>> const&, std::vector<std::function<std::vector<char, 
>>>>>>>> std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> >, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus)>, 
>>>>>>>> std::allocator<std::function<std::vector<char, std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> >, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus)> > > const&, 
>>>>>>>> std::vector<std::function<std::vector<char, std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> >, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus)>, 
>>>>>>>> std::allocator<std::function<std::vector<char, std::allocator<char> > 
>>>>>>>> (dealii::TriaIterator<dealii::CellAccessor<2, 2> >, 
>>>>>>>> dealii::Triangulation<2, 2>::CellStatus)> > > const&)
>>>>>>>> #9  /lib/x86_64-linux-gnu/libdeal.ii.g.so.9.2.0: 
>>>>>>>> dealii::parallel::distributed::Triangulation<2, 
>>>>>>>> 2>::execute_coarsening_and_refinement()
>>>>>>>> #10  /home/kdas/FE_NothingTest/build/test_distributed: main
>>>>>>>>
>>>>>>>> Thank you very much,
>>>>>>>> Kaushik 
>>>>>>>>
>>>>>>>> On Tue, Dec 8, 2020 at 6:25 PM Marc Fehling <mafe...@gmail.com> 
>>>>>>>> wrote:
>>>>>>>>
>>>>>>>>> Hi Kaushik,
>>>>>>>>>
>>>>>>>>> the `p::d::SolutionTransfer` class should be able to deal with 
>>>>>>>>> `FENothing` elements in your example. The tricky cases are when 
>>>>>>>>> you're 
>>>>>>>>> coarsening a `FENothing` element with others as Bruno already pointed 
>>>>>>>>> out 
>>>>>>>>> (h-coarsening), or if you change a `FENothing` element to a different 
>>>>>>>>> element in the process (p-adaptation). But with your recent 
>>>>>>>>> modification, 
>>>>>>>>> you avoid these cases.
>>>>>>>>>
>>>>>>>>> I suspect that something else causes the error in your program. 
>>>>>>>>> Could you run a debugger on this and give us the backtrace for the 
>>>>>>>>> exception? This will give us more clues to figure out what goes wrong!
>>>>>>>>>
>>>>>>>>> Best,
>>>>>>>>> Marc
>>>>>>>>> On Tuesday, December 8, 2020 at 7:13:29 AM UTC-7 k.d...@gmail.com 
>>>>>>>>> wrote:
>>>>>>>>>
>>>>>>>>>> Hi Bruno:
>>>>>>>>>> Thanks for pointing that out. 
>>>>>>>>>> I tried to not refine FE_nothing cells by modifying the refine 
>>>>>>>>>> loop:
>>>>>>>>>> (The modified test is attached).
>>>>>>>>>>
>>>>>>>>>> for (auto &cell : dgq_dof_handler.active_cell_iterators()) 
>>>>>>>>>>       if (cell->is_locally_owned() && cell->active_fe_index () 
>>>>>>>>>> != 0)
>>>>>>>>>>         {
>>>>>>>>>>           if (counter > ((dim == 2) ? 4 : 8))
>>>>>>>>>>             cell->set_coarsen_flag();
>>>>>>>>>>           else
>>>>>>>>>>             cell->set_refine_flag();
>>>>>>>>>>         }
>>>>>>>>>>
>>>>>>>>>> But this still aborts. 
>>>>>>>>>> Kaushik 
>>>>>>>>>>
>>>>>>>>>> On Tue, Dec 8, 2020 at 8:36 AM Bruno Turcksin <
>>>>>>>>>> bruno.t...@gmail.com> wrote:
>>>>>>>>>>
>>>>>>>>>>> Hi,
>>>>>>>>>>>
>>>>>>>>>>> Are you sure that your test makes sense? You randomly assign FE 
>>>>>>>>>>> indices to cells then you refine and coarsen cells. But what does 
>>>>>>>>>>> it mean 
>>>>>>>>>>> to coarsen 4 cells together when one of them is FE_Nothing? What 
>>>>>>>>>>> would you 
>>>>>>>>>>> expect to happen?
>>>>>>>>>>>
>>>>>>>>>>> Best,
>>>>>>>>>>>
>>>>>>>>>>> Bruno
>>>>>>>>>>>
>>>>>>>>>>> On Monday, December 7, 2020 at 5:54:10 PM UTC-5 k.d...@gmail.com 
>>>>>>>>>>> wrote:
>>>>>>>>>>>
>>>>>>>>>>>> Hi all:
>>>>>>>>>>>>
>>>>>>>>>>>> I modified the test  tests/mpi/solution_transfer_05.cc to add a 
>>>>>>>>>>>> FE_Nohting element to the FECollection. I also modified the other 
>>>>>>>>>>>> elements 
>>>>>>>>>>>> to FE_Q. 
>>>>>>>>>>>>
>>>>>>>>>>>> When I run the test, it's aborting in solution transfer. 
>>>>>>>>>>>> Is there any limitations in using FE_Nothing with parallel 
>>>>>>>>>>>> solution transfer? 
>>>>>>>>>>>> The modified test is attached.
>>>>>>>>>>>> Thank you very much.
>>>>>>>>>>>> Kaushik 
>>>>>>>>>>>>
>>>>>>>>>>>> ---- 
>>>>>>>>>>>> An error occurred in line <1167> of file 
>>>>>>>>>>>> </build/deal.ii-vFp8uU/deal.ii-9.2.0/include/deal.II/lac/vector.h> 
>>>>>>>>>>>> in 
>>>>>>>>>>>> function
>>>>>>>>>>>>     Number& 
>>>>>>>>>>>> dealii::Vector<Number>::operator()(dealii::Vector<Number>::size_type)
>>>>>>>>>>>>  [with 
>>>>>>>>>>>> Number = double; dealii::Vector<Number>::size_type = unsigned int]
>>>>>>>>>>>> The violated condition was:
>>>>>>>>>>>>     static_cast<typename 
>>>>>>>>>>>> ::dealii::internal::argument_type<void( typename 
>>>>>>>>>>>> std::common_type<decltype(i), decltype(size())>::type)>:: type>(i) 
>>>>>>>>>>>> < 
>>>>>>>>>>>> static_cast<typename ::dealii::internal::argument_type<void( 
>>>>>>>>>>>> typename 
>>>>>>>>>>>> std::common_type<decltype(i), decltype(size())>::type)>:: 
>>>>>>>>>>>> type>(size())
>>>>>>>>>>>> Additional information:
>>>>>>>>>>>>     Index 0 is not in the half-open range [0,0). In the current 
>>>>>>>>>>>> case, this half-open range is in fact empty, suggesting that you 
>>>>>>>>>>>> are 
>>>>>>>>>>>> accessing an element of an empty collection such as a vector that 
>>>>>>>>>>>> has not 
>>>>>>>>>>>> been set to the correct size.
>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>> -- 
>>>>>>>>>>> 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/ssEva6C2PU8/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 
>>>>>>>>>>> https://groups.google.com/d/msgid/dealii/609b0457-ba90-4aa2-a7d1-5b798d5349ebn%40googlegroups.com
>>>>>>>>>>>  
>>>>>>>>>>> <https://groups.google.com/d/msgid/dealii/609b0457-ba90-4aa2-a7d1-5b798d5349ebn%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/ssEva6C2PU8/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 
>>>>>>>>> https://groups.google.com/d/msgid/dealii/4cd4582b-3168-4b8b-8195-e316b049cfadn%40googlegroups.com
>>>>>>>>>  
>>>>>>>>> <https://groups.google.com/d/msgid/dealii/4cd4582b-3168-4b8b-8195-e316b049cfadn%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/ssEva6C2PU8/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 
>>>>>> https://groups.google.com/d/msgid/dealii/ada069fe-e4cb-455a-bb4b-06a776d4490bn%40googlegroups.com
>>>>>>  
>>>>>> <https://groups.google.com/d/msgid/dealii/ada069fe-e4cb-455a-bb4b-06a776d4490bn%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/ssEva6C2PU8/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 
>> https://groups.google.com/d/msgid/dealii/7f41a733-77ea-4a6e-b8d9-412de5d9e963n%40googlegroups.com
>>  
>> <https://groups.google.com/d/msgid/dealii/7f41a733-77ea-4a6e-b8d9-412de5d9e963n%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 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/b47d4e5d-4d94-40db-b7f8-d8b28f0e5bf2n%40googlegroups.com.

Reply via email to