Re: [deal.II] Re: Parallel distributed hp solution transfer with FE_nothing

2021-01-09 Thread Wolfgang Bangerth



Thanks to all of you, I can write a small test that seems to be doing what I 
want following Marc's suggestion. But I have a few questions (the test code is 
attached)


  *   cell->set_dof_values(local_data, m_completely_distributed_solution); 

Does set_dof_values set entries in m_completely_distributed_solution that are 
not owned by the current rank?


No. You determine up front upon construction of the vector elements are 
locally owned, which are ghosted, and which are stored elsewhere. Writing into 
entries does not change this decision.




  * m_completely_distributed_solution.compress(VectorOperation::insert);

After set_dof_values, if I call compress(VectorOperation::insert) on that 
PETSc vector, now, does that mean that I now have entries in that distributed 
vector that are not owned by the current rank?


Same here.



  * m_locally_relevant_solution = m_completely_distributed_solution;

Will this assignment of a distributed vector to a locally relevant vector 
correctly set the ghost values?


Yes. That's the point :-)


If this test code is correct, it will be great to add this as a test to 
deal.ii. My current work will be dependent on this part working correctly, and 
it will be nice to have a test that can alert us if something changes in the 
future.


We're always happy to add more tests to the test suite!

Best
 W.

--

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/800a858e-1294-f4d7-cf7a-9deed2de02c9%40colostate.edu.


Re: [deal.II] FEInterfaceValues

2021-01-09 Thread Timo Heister
> If this is correct, do you see a way to circumvent the TRILINOS problem I am 
> facing?

If you need entries i,j in the sparsity pattern and none of the
make_*_sparsity_pattern() function apply, you will need to add them
manually. This should be as simple as following the same logic you are
using in the assembly, where you identify DoF indices (you called them
local_interface_dof_indices) and you add them using

dsp.add(i,j)

where i and j are all indices in local_interface_dof_indices.


-- 
Timo Heister
http://www.math.clemson.edu/~heister/

-- 
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/CAMRj59HonM1P7Vm0oDc8PqJ-t7tkoYafcBVdBL-iuMGSr7-quQ%40mail.gmail.com.


Re: [deal.II] Re: Parallel distributed hp solution transfer with FE_nothing

2021-01-09 Thread Kaushik Das
Hello Prof. Bangerth and others,

Thanks to all of you, I can write a small test that seems to be doing what
I want following Marc's suggestion. But I have a few questions (the test
code is attached)

   -  cell->set_dof_values(local_data, m_completely_distributed_solution);

Does set_dof_values set entries in m_completely_distributed_solution that
are not owned by the current rank?

   - m_completely_distributed_solution.compress(VectorOperation::insert);

After  set_dof_values, if I call  compress(VectorOperation::insert) on that
PETSc vector, now, does that mean that I now have entries in that
distributed vector that are not owned by the current rank?

   - m_locally_relevant_solution = m_completely_distributed_solution;

Will this assignment of a distributed vector to a locally relevant vector
correctly set the ghost values?

If this test code is correct, it will be great to add this as a test to
deal.ii. My current work will be dependent on this part working correctly,
and it will be nice to have a test that can alert us if something changes
in the future.

Thank you very much,
Kaushik


On Mon, Jan 4, 2021 at 4:54 PM Wolfgang Bangerth 
wrote:

>
> Kaushik
> Marc and others have already answered the technical details, so just one
> overall comment:
>
> > Let me explain what I am trying to do and why.
> > I want to solve a transient heat transfer problem of the additive
> > manufacturing (AM) process. In AM processes, metal powder is deposited
> in
> > layers, and then a laser source scans each layer and melts and bonds the
> > powder to the layer underneath it. To simulate this layer by layer
> process, I
> > want to start with a grid that covers all the layers, but initially,
> only the
> > bottom-most layer is active and all other layers are inactive, and then
> as
> > time progresses, I want to activate one layer at a time. I read on this
> > mailing list that cell "birth" or "activation" can be done by changing a
> cell
> > from FE_Nothing to FE_Q using p-refinement. I was trying to keep all
> cells of
> > the grid initially to FE_Nothing except the bottom-most layer. And then
> > convert one layer at a time to FE_Q. My questions are:
> > 1. Does this make sense?
>
> Yes, this is how I would do things as well. It is quite possible that
> nobody
> has ever tried this, and that some of the steps don't work without further
> modification -- but whatever doesn't work should be treated as either a
> missing feature, or a bug. The general approach is sound.
>
> When I try to build a code with a work flow that is not quite standard, I
> (like you) frequently run into things that don't quite work. My usual
> approach
> is then to write a small and self-contained test case that illustrates the
> issue without the overhead of the actual application. Most of the time,
> one
> can show the issue with <100 lines of code. This then goes into a github
> issue
> so that one can write a fix for this particular problem without having to
> understand the overall code architecture, the application, etc. I would
> encourage you to follow this kind of work-flow!
>
> Best
>   W.
>
> --
> 
> 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 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+unsubscr...@googlegroups.com.
> To view this discussion on the web visit
> https://groups.google.com/d/msgid/dealii/5ebefd8a-5dde-603e-76d1-6967c7783504%40colostate.edu
> .
>

-- 
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/CAC-fs6sJ-GLAumx6YqiXMDbG3kQeg3Xzy_mQzZ4t%3DLhkM_xdFQ%40mail.gmail.com.
#include 
#include 
#include 
#include 
#include 
#include 
#include 

#include 
#include 
#include 

#include 

#include 
#include 
#include 
#include 
#include 
#include 

#include 
#include 
#include 

#include 
#include 
#include 
#include 

//
//  Including cell_data_transfer.h did not work. Had to include 
cell_data_transfer.templates.h why?
//
//#include 
#include 

namespace LA
{
#if defined(DEAL_II_WITH_PETSC) && !defined(DEAL_II_PETSC_WITH_COMPLEX) && \

Re: [deal.II] FEInterfaceValues

2021-01-09 Thread Alberto Salvadori
Thank you, Wolfgang. The explanation was clear as usual, I understood the
principle but overlooked the DoFTools class. This is now set.
Unfortunately, it seems that it does not solve my issue and I guess I
figured out why.

In brief, I am building an
interface matrix and an interface rhs in the very same way cell matrix and
cell rhs have been built in step 18 - as depicted already (I won't repeat
it here).

Issues arise when I attempt at distributing the matrices into system
matrix:  I coded this

  local_interface_dof_indices = fiv.get_interface_dof_indices() ;


this->hanging_node_constraints.distribute_local_to_global(interface_matrix,
interface_rhs,  local_interface_dof_indices, this->system_matrix,
this->system_rhs);
using the very same hanging_node_constraints I used for the usual continuum
elements.
This last statement causes an error FOR TRILINOS (AND NOT FOR PETSc, which
works OK).

According to your remark "make_flux_sparsity_pattern() also adds entries
for DoFs i,j that are located
on cells that are face neighbors of each other" the error may arise because
in defining the interfaces I am:
1 - duplicating the nodes in order to separate the initial triangulation
into two
2 - building a new data structure that provide connectivity for the dofs of
the separated elements, exploiting the class FEInterfaceValues
I wonder if by doing so, the notion of "cells that are face neighbors of
each other" may not apply anymore and hence, as you suspected, I am
computing terms for the interface matrix
that are not part of the sparsity pattern of the Trilinos matrix.

If this is correct, do you see a way to circumvent the TRILINOS problem I
am facing?

Thank you,
Alberto


*Alberto Salvadori* Dipartimento di Ingegneria Meccanica e Industriale
(DIMI)
 Universita` di Brescia, via Branze 43, 25123 Brescia
 Italy
 tel 030 3715426

e-mail:
 alberto.salvad...@unibs.it
web-page:
 http://m4lab.unibs.it/faculty.html



On Fri, Jan 8, 2021 at 9:10 PM Wolfgang Bangerth 
wrote:

> On 1/8/21 10:19 AM, Alberto Salvadori wrote:
> >
> > thank you for the hint. I looked at setp 47 sparsity pattern. If I
> understood
> > it right,
> > the main difference stands in the locally_relevant_dofs, which are no
> longer
> > extracted and used.
> > Rather, the sparsity pattern is based on the whole dof_handler, as per
> the
> > instruction
> >
> > DoFTools::make_flux_sparsity_pattern( this->dof_handler, dsp,
> > this->hanging_node_constraints, true );
>
> No, you mistook what the difference is. The difference is that
> make_sparsity_pattern() only adds entries for the matrix that correspond
> to
> DoFs i and j that are located on the same cell.
>
> make_flux_sparsity_pattern() also adds entries for DoFs i,j that are
> located
> on cells that are face neighbors of each other.
>
> In both cases, you can use the same locally_relevant_set when you set up
> the
> sparsity pattern object.
>
> Best
>   W.
>
>
> --
> 
> 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/a7c6b38b-8537-26e1-5ea6-72b06cc5b8f3%40colostate.edu
> .
>

-- 


Informativa sulla Privacy: http://www.unibs.it/node/8155 


-- 
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/CABcATpfwhHhyDhCENUshUraRm3md1Uvnh5hhzqAFGXhXt7TZkw%40mail.gmail.com.


Re: [deal.II] Stabilized FEM implementation with bubble function

2021-01-09 Thread Lixing Zhu
Hi Wolfgang,

Many thanks! The issue is resolved and the code is working.

Regards,
Lixing
On Saturday, January 9, 2021 at 4:07:27 AM UTC+8 Wolfgang Bangerth wrote:

> On 1/8/21 12:15 PM, Lixing Zhu wrote:
> > 
> > This way I can freely use bubble functions and only construct local K 
> and F of 
> > FE_1. However, one issue is that, apparently, triangulation.cell does 
> not 
> > support cell->get_dof_indices(local_dof_indices). Is there any way to 
> call the 
> > DOFhandler by the cell in triangulation and then access the 
> local_dof_indices?
>
> You want a statement like this (from step-43):
>
> auto cell = darcy_dof_handler.begin_active();
> const auto endc = darcy_dof_handler.end();
> auto saturation_cell = saturation_dof_handler.begin_active();
>
> for (; cell != endc; ++cell, ++saturation_cell)
>
> Here, we're walking iterators over two DoFHandler objects in parallel. In 
> your 
> case, one will be the DoFHandler iterator and the other one the 
> Triangulation 
> iterator. Or, you could do
> fe_values_b.reinit
> (static_cast::active_cell_iterator>
> (cell));
> to cast the DoFHandler iterator to a Triangulation iterator.
>
> Best
> W.
>
> -- 
> 
> Wolfgang Bangerth email: bang...@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/d85fad4e-baa7-4984-bde2-b1f1f5093bd0n%40googlegroups.com.