Dear Daniel
My name is Pasha. Thank you for your response. I want to extract sum of
values from system RHS ( which is a
dealii:TrilinosWrappers::MPI::BlockVector) correspond to all degrees of
freedom which are at the boundary. I am working with a
parallel::distributed::Triangulation. I wrote the following code to extract
value for component 1 at boundary 3:
IndexSet in_set;
std::vector<bool> component_select (n_components, false);
component_select[1] = true;
std::set< types::boundary_id > boundary_ids;
boundary_ids.insert(3);
DoFTools::extract_boundary_dofs(dof_handler, component_select,
in_set, boundary_ids);
typename DoFHandler<dim>::active_cell_iterator
cell = dof_handler.begin_active(),
endc = dof_handler.end();
double local_load = 0.0 ;
for (; cell!=endc; ++cell)
if (cell->is_locally_owned())
{
cell->get_dof_indices(local_dof_indices);
dealii::IndexSet::ElementIterator
index = in_set.begin(),
endind = in_set.end();
for (; index!=endind; ++index)
{
dealii::types::global_dof_index gdi = *index;
unsigned int pos = find(local_dof_indices.begin(),
local_dof_indices.end(),
gdi) - local_dof_indices.begin();
if(pos >= local_dof_indices.size()) {
//not found
continue;
}
// now DoF is in this cell
local_load -= saved_residual(gdi);
}
}
double load_value = Utilities::MPI::sum(local_load, mpi_com);
is the above method correct? I am confusing because degrees of freedom on
the boundary that belong to the locally relevant set.
>From deal II online documentation for DoFTools::extract_boundary_dofs():
If the DofHandler object is indeed defined on a
parallel::distributed::Triangulation, then the selected_dofs index set will
contain only those degrees of freedom on the boundary that belong to the
locally relevant set.
*Locally relevant degrees of freedom*
This concept identifies a subset of all degrees of freedom when using
distributed meshes, see the Parallel computing with multiple processors
using distributed memory
<https://www.dealii.org/8.4.0/doxygen/deal.II/group__distributed.html> module.
Locally relevant degrees of freedom are those that live on locally owned or
ghost cells. Consequently, they may be owned by different processors.
Thank you again,
Pasha
On Thursday, June 30, 2016 at 9:26:54 PM UTC+4:30, Daniel Arndt wrote:
> Dear ?,
>
> dealii::IndexSet::ElementIterator::operator*() returns the stored
> global_dof_index [1].
> Therefore, all you have to do is to dereference index:
>
> dealii::types::global_dof_index gdi = *index
>
> Best,
> Daniel
>
> [1]
> https://www.dealii.org/8.4.1/doxygen/deal.II/classIndexSet_1_1ElementIterator.html#a6218ab9e8e76699fc998834bf8e006ee
>
> Am Donnerstag, 30. Juni 2016 15:12:00 UTC+2 schrieb [email protected]:
>>
>> Hi
>>
>> I am trying to extract all degrees of freedom which are at the boundary
>> and belong to specified components using DoFTools::extract_boundary_dofs()
>> function and IndexSet for a parallel distributed triangulation. Now, how I
>> can iterate over the IndexSet and how I can access to the Dof index in it?
>>
>> Something like this:
>>
>> dealii::IndexSet::ElementIterator index = in_set.begin(), endind
>> = in_set.end();
>>
>> for (; index!=endind; ++index)
>> {
>> how to access to dof index here?
>> }
>>
>> Thank you
>>
>
--
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 [email protected].
For more options, visit https://groups.google.com/d/optout.