On 02/26/2017 07:35 AM, 'Seyed Ali Mohseni' via deal.II User Group wrote:
Thank you everyone. I understand your suggestions. From what I learned now in
deal.II is, that each cell is owned by a specific processor which means I
cannot access the information within a locally owned cell. That's why I have
to loop over all cells and work with distributed data without using locally
owned entities.

Seyed -- no, that's exactly the opposite of what you can do. "Locally owned" means precisely that this is the kind of data the current processor has access to. You cannot access data that is *not locally owned* (with the exception of ghost elements).


My approach using local_dof_indices is somehow impossible. It only works for 1
core according to the below code:

|
std::vector<double> displacement_norms, loaded_nodes, loaded_elements;
double max_displacement;
int max_loaded_node, max_loaded_element;

typename DoFHandler<dim>::active_cell_iterator cell =
dof_handler.begin_active(), endc = dof_handler.end();
std::vector<unsigned int> local_dof_indices(fe.dofs_per_cell);
unsigned int cell_number = 0;

for (; cell != endc; ++cell, ++cell_number)
{
     pcout << "CELL ID: " << cell_number << std::endl;     // this works
surprisingly.

     cell->get_dof_indices(local_dof_indices);    // Output or accessing this
won't work I assume.

It totally should. But it only will if cell->locally_owned()==true. You do not have this information for cells that are not locally owned.


     for (unsigned int i = 0, j = 1; i < fe.dofs_per_cell; i += 2, j += 2)
     {
  double displacement_x = locally_relevant_solution(local_dof_indices[i]);
  double displacement_y = locally_relevant_solution(local_dof_indices[j]);

          double displacement_norm = sqrt(displacement_x * displacement_x +
displacement_y * displacement_y);

          loaded_nodes.push_back(cell->vertex_index(v));
          loaded_elements.push_back(cell_number);
          displacement_norms.push_back(displacement_norm);
     }
}

max_displacement = *max_element(displacement_norms.begin(),
displacement_norms.end());

double max_index = distance(displacement_norms.begin(),
max_element(displacement_norms.begin(), displacement_norms.end()));

// Node and element number of maximal displacements
max_loaded_node = loaded_nodes[max_index];
max_loaded_element = loaded_elements[max_index];
|

At this point, with the
  if (cell->is_locally_owned())
statement added to the code above, every processor has now computed the maximal displacement among all nodes it knows about. You then need to find the maximum among all processors, for example using Utilities::MPI::max().

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

Reply via email to