On Sunday, February 12, 2017 at 8:49:34 AM UTC-6, Timo Heister wrote: > > >> > >> > 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? > > I am confused what you you do: > a) loop over cells with get_function_values() > b) loop over locally owned indices and do vec1(idx) = vec2(idx); > c) assign vectors: velocity_lr = solution_lr.block(0) > > If you do a) it should work if the vector is ghosted and you only look > at owned cells. > If you do b) it should work if you only look at locally owned indices > of vec2. You can write anywhere but only read from available entries. > If you do c) I have no idea what will happen. >
I believe that I am currently doing option c). I have two member functions on my velocity class (with the block system) which form the public interface for my velocity to be used in the transport equation. They are get_dof_handler_velocity() and get_velocity() which are defined as follows: 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)); // 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; } and LinearAlgebra::Vector & get_velocity(){return solution_lr.block(0);} These are then used as input for constructing the fe_values and pulling out the local velocity values. It is precisely when using the fe_values.get_function_values(velocity_model->get_velocity(), local_velocity_values) that the error 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). was triggered. Is there a better way to implement the get_velocity() member function? I could fairly easily change the public interface to return std::shared_ptr<LinearAlgebra::Vector> if that makes is simpler as well... > > 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> ? > > The vectors as explained above contain duplicates and ghost indices, > that you should probably remove. They also don't give you a direct > mapping from idx to idx. You could construct a map but you can also > use the vectors as is. > > > Does the DoFHandler with it's attached FiniteElement change anything in > the > > FiniteElement when we renumber? > > No, a DoFHandler does not modify the FiniteElement (look at > distribute_dofs(), it gets a const ref of the FE so it won't modify > it). One of the important concepts of deal.II is we try to be modular: > you can have as many Triangulation, FE, DoFHandler, FEValues, ... > objects at the same time as you want. > > I thought that you had built deal.II in this way, but there were a few comments in the documentation which gave me the idea that there might be side effects with reordering to the finite element being used. I see now that that was a misunderstanding on my part. > -- > 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.