Hello Prof. Bangerth,
Yes I think I am setting them in a way that is exactly same throughout all
processes.
Since I can check the n_active_cells() before the refinement, and also the
refine_flags and coarsen_flags for each process. They are the same for each
process. But after the refinement, the number of active_cells suddenly
becomes different.
I attached the refine_mesh() code.
Best,
Yiyang
Vector<float> estimated_error_per_cell(this->tria.n_active_cells());
KellyErrorEstimator<dim>::estimate(this->dof_handler,
QGauss<dim - 1>(Q_POINTS),
typename FunctionMap<dim>::type(),
locally_relevant_sln,
estimated_error_per_cell,
ComponentMask(),
0,
MultithreadInfo::n_threads(),
Utilities::MPI::this_mpi_process(this->mpi_communicator));
IndexSet local_active_cells(this->tria.n_active_cells());
local_active_cells.clear();
for(auto cell = this->tria.begin_active(); cell != this->
tria.end(); ++cell){
if(cell->is_locally_owned()){
local_active_cells.add_index(cell->
active_cell_index());
}
}
TrilinosWrappers::MPI::Vector distributed_error_per_cell(
local_active_cells, this->mpi_communicator);
for(auto cell = this->tria.begin_active(); cell != this->
tria.end(); ++cell){
if(cell->is_locally_owned()){
distributed_error_per_cell(cell->
active_cell_index()) = estimated_error_per_cell(cell->active_cell_index());
}
}
distributed_error_per_cell.compress(VectorOperation::insert
);
estimated_error_per_cell = distributed_error_per_cell;
GridRefinement::refine_and_coarsen_fixed_number(this->tria,
estimated_error_per_cell,
param[0],
param[1],
max_cells);
//restrict grid levels if these two parameters are given.
if(max_grid_level > 0){
if (this->tria.n_levels() > max_grid_level) {
for (typename Triangulation<dim>::
active_cell_iterator
cell = this->tria.begin_active(
max_grid_level);
cell != this->tria.end(
max_grid_level); ++cell) {
cell->clear_refine_flag();
}
}
}
if(min_grid_level > 0){
for (typename Triangulation<dim>::
active_cell_iterator
cell = this->tria.begin_active(
min_grid_level);
cell != this->tria.end_active(
min_grid_level); ++cell) {
cell->clear_coarsen_flag();
}
}
int refine_count = 0;
int coarsen_count = 0;
for(auto cell = this->tria.begin_active(); cell != this->
tria.end(); ++cell){
if(cell->refine_flag_set()) refine_count++;
if(cell->coarsen_flag_set()) coarsen_count++;
}
std::cout<<"Process No. "<<Utilities::MPI::this_mpi_process(
this->mpi_communicator) <<". Refine Count = "<<refine_count<<std::endl;
std::cout<<"Process No. "<<Utilities::MPI::this_mpi_process(
this->mpi_communicator) <<". Coarsen Count = "<<coarsen_count<<std::endl;
parallelSolutionTransfer<dim> sol_trans(this->dof_handler);
this->tria.prepare_coarsening_and_refinement();
sol_trans.prepare_for_coarsening_and_refinement(
locally_relevant_sln);
std::cout<<"Process No. "<<Utilities::MPI::this_mpi_process(
this->mpi_communicator) <<". Pre_refine, n_active_cells = "<<this->tria.
n_active_cells()<<std::endl;
this->run_mesh_refinement();
std::cout<<"Process No. "<<Utilities::MPI::this_mpi_process(
this->mpi_communicator) <<". Post_refine, n_active_cells = "<<this->tria.
n_active_cells()<<std::endl;
--
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.