<https://lh3.googleusercontent.com/-xs6Y69oHfuY/WZ86zsk0w1I/AAAAAAAABE8/FVYbjIPSRrAZoi9PMw11AEofuNlOsYsBgCLcBGAs/s1600/Picture1.png>
Hi, all.
I have a questions on an algorithm in step-9 that does some iteration over
the cell with its neighbors.
I tried to read explanation carefully, but some part is unclear yet.
In step-9, one class estimates a gradient of solution in each cell by
calculating estimated value from its neighbor cells.
Special care might be needed when our domain consist of different level of
refinement in order to include all neighbor cells, and I think this parts
from step-9 does this... but I am not understanding this..
I made some example, in the simple mesh above,
*case 1: *
when we iterate over on cell-4 in my pic, it has four faces and each face
has its neighbor, (2,3,5,6)
for cell, 2 and 3, it has same refinement level and I assume that 2 and 3
become automatically active.
and I wonder what happens when it comes to 5 and 6. where the neighbor cell
has one coarse level?
The explanation in step 9 tutorial tells us that..
*If we are not in 1d, we collect all neighbor children `behind' the
subfaces of the current face*
for (unsigned int subface_no=0; subface_no<face->n_children(); ++subface_no)
active_neighbors.push_back (
std::get<0>(*cell)->neighbor_child_on_subface(face_no,subface_no));
}
}
and I think I don't understand the blue line properly....and especially the
word subfaces means.
*case 2: *
when we iterative over on cell-5 in my pic, I am assuming that the cell 2,4
and 7 are automatically activated neighbor cell...
and we don't need to anything special..
I would appreciate if someone can help me.
Thanks!!
Jaekwang
template <int dim>
void
GradientEstimation::estimate_cell (const
SynchronousIterators<std_cxx11::tuple<typename
DoFHandler<dim>::active_cell_iterator,
Vector<float>::iterator> > &cell,
EstimateScratchData<dim>
&scratch_data,
const EstimateCopyData &)
{
Tensor<2,dim> Y;
std::vector<typename DoFHandler<dim>::active_cell_iterator>
active_neighbors;
active_neighbors.reserve (GeometryInfo<dim>::faces_per_cell *
GeometryInfo<dim>::max_children_per_face);
typename DoFHandler<dim>::active_cell_iterator
cell_it(std_cxx11::get<0>(cell.iterators));
scratch_data.fe_midpoint_value.reinit (cell_it);
Tensor<1,dim> projected_gradient;
active_neighbors.clear ();
for (unsigned int face_no=0;
face_no<GeometryInfo<dim>::faces_per_cell; ++face_no)
if (! std_cxx11::get<0>(cell.iterators)->at_boundary(face_no))
{
//just definition of abbreviation for the iterator to the
face and the neighbor...
const typename DoFHandler<dim>::face_iterator
face = std_cxx11::get<0>(cell.iterators)->face(face_no);
const typename DoFHandler<dim>::cell_iterator
neighbor = std_cxx11::get<0
>(cell.iterators)->neighbor(face_no);
if (neighbor->active())
active_neighbors.push_back (neighbor); //increase
storage
else
{
// neighbor cell can be one level of coasered
// You should find all neighbor cell
// find a child of neighbor
if (dim == 1) // special care is necessary for 1D case
{
typename DoFHandler<dim>::cell_iterator
neighbor_child = neighbor;
while (neighbor_child->has_children())
neighbor_child = neighbor_child->child
(face_no==0 ? 1 : 0);
Assert (neighbor_child->neighbor(face_no==0 ? 1 : 0)
==std_cxx11::get<0
>(cell.iterators),ExcInternalError());
active_neighbors.push_back (neighbor_child);
}
else // Most cases will fall under this case (2D and
3D Problems..)
for (unsigned int subface_no=0; subface_no <
face->n_children(); ++subface_no)
active_neighbors.push_back (std_cxx11::get<0
>(cell.iterators)->neighbor_child_on_subface(face_no,subface_no));
}
}
--
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.