On Friday, February 10, 2017 at 12:15:10 PM UTC-6, Timo Heister wrote:
>
> I wanted to give you some info on your original question in case your 
> want to still use PETSc. 
>
> > It appears that for petsc, the assumption that the locally owned dofs 
> Index 
> > Sets are contiguous is really throwing a wrench in our plans for the non 
> > block system approach.  I have seen the other discussions where everyone 
> > says essentially it is not currently possible to get petsc and petsc 
> > wrappers to a point where it could use non contiguous locally owned 
> index 
> > sets. 
>
> Yes, you need the locally owned range of indices to be contiguous. 
> This means, that you can split your matrices/vectors into blocks 
> however you want (for example by component), but they need to be 
> contiguous within a block. 
>
>
The blockwise reordering works just fine for block matrices when using the 
iterative solvers since they boil down to matrix vector operations but in 
the case of the direct solvers we don't have a block interface so it is 
necessary to use a PETScWrappers::MPI::Vector but then the contiguous-ness 
comes back into play if we do blockwise renumbering...  
 

>
> > fe_values.get_function_values() on the velocity vector. 
> >     ExcIndexNotPresent (n) 
>
> get_function_values() needs to see the ghost values of your vector 
> (all DoFs on a cell). 
>

yes,  that is what is being passed to it.  I use the post flags _lo and _lr 
for my vector names to distinguish between the completely distributed 
(locally owned) and ghosted (locally relevant) vectors.  I am not positive, 
but I believe the issue is that we are using  velocity_lr = 
solution_lr.block{0} where solution_lr has not been blockwise renumbered 
and so the locally owned IndexSet corresponding to the first block is not a 
single interval and in fact contains global indices larger than the number 
of dofs on the block. The locally relevant index set has the same problem. 
 Could it be that there is a difference in dof index between the global 
indexset of solution_lr and the dof index on the block and the use of the 
incorrect index is what triggers the ExcIndexNotPresent() error?
 

>
> > So, in this new case, what would be a good way to accomplish the desired 
> > result? 
>
> how about: 
> ... 
> DoFRenumbering::hierarchical (dof_handler); 
>
> std::vector<unsigned int> velocity_dof_indices1; 
> std::vector<unsigned int> velocity_dof_indices2; 
>
> for each active cell: 
>   get_dof_indices(...) 
>   for (unsigned int i;i<n_local_dofs;}}i) 
>      { 
>        unsigned int component = fe.system_to_component_index(i).first; 
>        if (component<dim) // velocity 
>              velocity_dof_indices1.push_back(local_dof_indices[i]); 
> } 
>
> DoFRenumbering::component_wise(dof_handler); 
>
> for each active cell: 
>   get_dof_indices(...) 
>   for (unsigned int i;i<n_local_dofs;}}i) 
>      { 
>        unsigned int component = fe.system_to_component_index(i).first; 
>        if (component<dim) // velocity 
>              velocity_dof_indices2.push_back(local_dof_indices[i]); 
> } 
>
> Now you have a mapping between the two different numberings and you 
> can use that to convert solution vectors. Note that you probably only 
> want to consider locally owned DoFs and use a more efficient data 
> structure that actually maps between the pairs of indices directly. 
>
>
Yes I thought there would be an approach like this. Thanks for laying it 
out.   As I am only needing to convert one direction, the vectors as 
constructed above would likely be sufficient but are you thinking of a data 
structure like std::map<global_dof_index, global_dof_index>  ?

Does the DoFHandler with it's attached FiniteElement change anything in the 
FiniteElement when we renumber?  In other words,  do I need to renumber the 
dofs before I extract the dofhandler for the first block as per the 
following code?

    std::shared_ptr<DoFHandler<dim> > get_dof_handler_velocity()

      {

        std::shared_ptr<DoFHandler<dim> > dof_velocity_ptr 
(new DoFHandler<dim>(triangulation));

        
dof_velocity_ptr->distribute_dofs(fe_system_ptr->base_element(0)); // give 
the velocity block fe


        // order same as dof_handler_ptr in setup_system

        DoFRenumbering::hierarchical (*(dof_velocity_ptr));

        DoFRenumbering::component_wise (*(dof_velocity_ptr),

                                        std::vector<unsigned int> (dim,0)); 
// this might not actually do anything

                                                                            
// since it is essentially blockwise

        return dof_velocity_ptr;

      }


or do I just need to construct the new reordered solution vector?

Otherwise, I suppose I need the old dofhandler which will continue to be 
used and the new dofhandler reordered and then i can use the above code 
with the new reordered dofhandler.  Would I need a new Finite Element 
object as well or could I use the same one for both the full dofhandlers?

Does that make sense what I am asking?

Spencer


 

> -- 
> 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.
For more options, visit https://groups.google.com/d/optout.

Reply via email to