17 #include <deal.II/base/template_constraints.h> 18 #include <deal.II/lac/vector.h> 19 #include <deal.II/lac/block_vector_base.h> 20 #include <deal.II/lac/block_vector.h> 21 #include <deal.II/lac/petsc_vector.h> 22 #include <deal.II/lac/petsc_block_vector.h> 23 #include <deal.II/lac/trilinos_vector.h> 24 #include <deal.II/lac/trilinos_block_vector.h> 26 #include <deal.II/grid/grid_refinement.h> 27 #include <deal.II/grid/tria_accessor.h> 28 #include <deal.II/grid/tria_iterator.h> 29 #include <deal.II/grid/tria.h> 37 DEAL_II_NAMESPACE_OPEN
44 template <
typename number>
49 return *std::max_element(criteria.
begin(), criteria.
end());
53 template <
typename number>
58 return *std::min_element(criteria.
begin(), criteria.
end());
63 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
65 #ifdef DEAL_II_WITH_PETSC 74 for (
unsigned int i=0; i<criteria.
size(); ++i)
75 m = std::max (m, criteria(i));
87 PetscScalar m = criteria(0);
88 for (
unsigned int i=1; i<criteria.
size(); ++i)
89 m = std::min (m, criteria(i));
95 #ifdef DEAL_II_WITH_TRILINOS 100 TrilinosScalar m = 0;
110 TrilinosScalar m = 0;
121 template <
typename VectorType>
122 typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
123 typename VectorType::value_type>::type
124 min_element (
const VectorType &criteria)
126 return internal::min_element (criteria);
130 template <
typename VectorType>
131 typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
132 typename VectorType::value_type>::type
133 max_element (
const VectorType &criteria)
135 return internal::max_element (criteria);
139 template <
typename VectorType>
140 typename constraint_and_return_value<IsBlockVector<VectorType>::value,
141 typename VectorType::value_type>::type
142 min_element (
const VectorType &criteria)
144 typename VectorType::value_type t = internal::min_element(criteria.block(0));
145 for (
unsigned int b=1; b<criteria.n_blocks(); ++b)
146 t = std::min (t, internal::min_element(criteria.block(b)));
152 template <
typename VectorType>
153 typename constraint_and_return_value<IsBlockVector<VectorType>::value,
154 typename VectorType::value_type>::type
155 max_element (
const VectorType &criteria)
157 typename VectorType::value_type t = internal::max_element(criteria.block(0));
158 for (
unsigned int b=1; b<criteria.n_blocks(); ++b)
159 t = std::max (t, internal::max_element(criteria.block(b)));
175 template <
class VectorType>
176 void qsort_index (
const VectorType &a,
177 std::vector<unsigned int> &ind,
182 typename VectorType::value_type v;
196 while ((a(ind[i])>v) && (i<r));
201 while ((a(ind[j])<v) && (j>0));
204 std::swap (ind[i], ind[j]);
206 std::swap (ind[i], ind[r]);
209 qsort_index(a,ind,l,i-1);
210 qsort_index(a,ind,i+1,r);
217 template <
int dim,
class VectorType,
int spacedim>
219 const VectorType &criteria,
220 const double threshold,
221 const unsigned int max_to_mark)
225 Assert (criteria.is_non_negative (), ExcNegativeCriteria());
230 if (criteria.all_zero())
233 const unsigned int n_cells = criteria.size();
237 double new_threshold=threshold;
241 if (new_threshold==0)
243 new_threshold = criteria(0);
244 for (
unsigned int index=1; index<n_cells; ++index)
245 if (criteria(index)>0
246 && (criteria(index)<new_threshold))
247 new_threshold=criteria(index);
250 unsigned int marked=0;
252 cell != tria.
end(); ++cell)
253 if (std::fabs(criteria(cell->active_cell_index())) >= new_threshold)
258 cell->set_refine_flag();
264 template <
int dim,
class VectorType,
int spacedim>
266 const VectorType &criteria,
267 const double threshold)
271 Assert (criteria.is_non_negative (), ExcNegativeCriteria());
274 cell != tria.
end(); ++cell)
275 if (std::fabs(criteria(cell->active_cell_index())) <= threshold)
276 if (!cell->refine_flag_set())
277 cell->set_coarsen_flag();
281 std::pair<double, double>
283 const unsigned int max_n_cells,
284 const double top_fraction,
285 const double bottom_fraction)
287 Assert (top_fraction>=0, ExcInvalidParameterValue());
288 Assert (top_fraction<=1, ExcInvalidParameterValue());
289 Assert (bottom_fraction>=0, ExcInvalidParameterValue());
290 Assert (bottom_fraction<=1, ExcInvalidParameterValue());
291 Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
293 double refine_cells = current_n_cells * top_fraction;
294 double coarsen_cells = current_n_cells * bottom_fraction;
299 std::pair<double, double> adjusted_fractions(top_fraction, bottom_fraction);
303 if (current_n_cells >= max_n_cells)
315 adjusted_fractions.first = 0;
316 coarsen_cells = (current_n_cells - max_n_cells) /
317 cell_decrease_on_coarsen;
318 adjusted_fractions.second = std::min(coarsen_cells/current_n_cells, 1.0);
333 else if (static_cast<unsigned int>
335 + refine_cells * cell_increase_on_refine
336 - coarsen_cells * cell_decrease_on_coarsen)
351 (max_n_cells - current_n_cells)
353 (refine_cells * cell_increase_on_refine
354 - coarsen_cells * cell_decrease_on_coarsen);
356 adjusted_fractions.first = alpha * top_fraction;
357 adjusted_fractions.second = alpha * bottom_fraction;
359 return (adjusted_fractions);
362 template <
int dim,
class VectorType,
int spacedim>
365 const VectorType &criteria,
366 const double top_fraction,
367 const double bottom_fraction,
368 const unsigned int max_n_cells)
372 Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
373 Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
374 Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
375 Assert (criteria.is_non_negative (), ExcNegativeCriteria());
377 const std::pair<double, double> adjusted_fractions =
378 adjust_refine_and_coarsen_number_fraction<dim> (criteria.size(),
383 const int refine_cells =
static_cast<int>(adjusted_fractions.first * criteria.size());
384 const int coarsen_cells =
static_cast<int>(adjusted_fractions.second * criteria.size());
386 if (refine_cells || coarsen_cells)
391 std::nth_element (tmp.
begin(), tmp.
begin()+refine_cells,
393 std::greater<double>());
394 refine (tria, criteria, *(tmp.
begin() + refine_cells));
399 std::nth_element (tmp.
begin(), tmp.
begin()+tmp.
size()-coarsen_cells,
401 std::greater<double>());
403 *(tmp.
begin() + tmp.
size() - coarsen_cells));
410 template <
int dim,
typename VectorType,
int spacedim>
413 const VectorType &criteria,
414 const double top_fraction,
415 const double bottom_fraction,
416 const unsigned int max_n_cells)
420 Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
421 Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
422 Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
423 Assert (criteria.is_non_negative (), ExcNegativeCriteria());
431 const double total_error = tmp.
l1_norm();
435 std::sort (tmp.begin(), tmp.end(), std::greater<double>());
441 (sum<top_fraction*total_error) && (pp!=(tmp.end()-1));
444 double top_threshold = ( pp != tmp.
begin () ?
450 (sum<bottom_fraction*total_error) && (qq!=tmp.
begin());
453 double bottom_threshold = ( qq != (tmp.end()-1) ?
470 const unsigned int refine_cells = pp - tmp.
begin(),
471 coarsen_cells = tmp.end() - qq;
473 if (static_cast<unsigned int>
482 refine_and_coarsen_fixed_number (tria,
484 1.*refine_cells/criteria.size(),
485 1.*coarsen_cells/criteria.size(),
523 if ((top_threshold == max_element(criteria)) &&
525 top_threshold *= 0.999;
527 if (bottom_threshold>=top_threshold)
528 bottom_threshold = 0.999*top_threshold;
531 if (top_threshold < max_element(criteria))
532 refine (tria, criteria, top_threshold, pp - tmp.
begin());
534 if (bottom_threshold > min_element(criteria))
535 coarsen (tria, criteria, bottom_threshold);
540 template <
int dim,
typename VectorType,
int spacedim>
543 const VectorType &criteria,
544 const unsigned int order)
548 Assert (criteria.is_non_negative (), ExcNegativeCriteria());
552 std::vector<unsigned int> tmp(criteria.size());
553 for (
unsigned int i=0; i<criteria.size(); ++i)
556 qsort_index (criteria, tmp, 0, criteria.size()-1);
558 double expected_error_reduction = 0;
559 const double original_error = criteria.l1_norm();
561 const unsigned int N = criteria.size();
564 double min_cost = std::numeric_limits<double>::max();
565 unsigned int min_arg = 0;
567 for (
unsigned int M = 0; M<criteria.size(); ++M)
569 expected_error_reduction += (1-std::pow(2.,-1.*order)) * criteria(tmp[M]);
571 const double cost = std::pow(((std::pow(2.,dim)-1)*(1+M)+N),
573 (original_error-expected_error_reduction);
574 if (cost <= min_cost)
581 refine (tria, criteria, criteria(tmp[min_arg]));
586 #include "grid_refinement.inst" 588 DEAL_II_NAMESPACE_CLOSE
unsigned int n_active_cells() const
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const VectorType &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
static const unsigned int invalid_unsigned_int
real_type l1_norm() const
active_cell_iterator begin_active(const unsigned int level=0) const
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
cell_iterator end() const
#define Assert(cond, exc)
void coarsen(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold)
const Epetra_MultiVector & trilinos_vector() const
std::pair< double, double > adjust_refine_and_coarsen_number_fraction(const unsigned int current_n_cells, const unsigned int max_n_cells, const double top_fraction_of_cells, const double bottom_fraction_of_cells)
void refine(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void refine_and_coarsen_optimize(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const unsigned int order=2)