@@ -54,14 +54,16 @@ namespace aspect
54
54
{
55
55
CitationInfo::add (" particles" );
56
56
if (particle_load_balancing & ParticleLoadBalancing::repartition)
57
- #if DEAL_II_VERSION_GTE(9,4,0)
58
57
this ->get_triangulation ().signals .weight .connect (
58
+ #if DEAL_II_VERSION_GTE(9,6,0)
59
+ [&] (const typename parallel::distributed::Triangulation<dim>::cell_iterator &cell,
60
+ const CellStatus status)
61
+ -> unsigned int
59
62
#else
60
- this ->get_triangulation ().signals .cell_weight .connect (
61
- #endif
62
63
[&] (const typename parallel::distributed::Triangulation<dim>::cell_iterator &cell,
63
64
const typename parallel::distributed::Triangulation<dim>::CellStatus status)
64
65
-> unsigned int
66
+ #endif
65
67
{
66
68
return this ->cell_weight (cell, status);
67
69
});
@@ -390,7 +392,12 @@ namespace aspect
390
392
template <int dim>
391
393
unsigned int
392
394
World<dim>::cell_weight(const typename parallel::distributed::Triangulation<dim>::cell_iterator &cell,
393
- const typename parallel::distributed::Triangulation<dim>::CellStatus status)
395
+ #if DEAL_II_VERSION_GTE(9,6,0)
396
+ const CellStatus status
397
+ #else
398
+ const typename parallel::distributed::Triangulation<dim>::CellStatus status
399
+ #endif
400
+ )
394
401
{
395
402
if (cell->is_active () && !cell->is_locally_owned ())
396
403
return 0 ;
@@ -400,6 +407,20 @@ namespace aspect
400
407
unsigned int n_particles_in_cell = 0 ;
401
408
switch (status)
402
409
{
410
+ #if DEAL_II_VERSION_GTE(9,6,0)
411
+ case CellStatus::cell_will_persist:
412
+ case CellStatus::cell_will_be_refined:
413
+ n_particles_in_cell = particle_handler->n_particles_in_cell (cell);
414
+ break ;
415
+
416
+ case CellStatus::cell_invalid:
417
+ break ;
418
+
419
+ case CellStatus::children_will_be_coarsened:
420
+ for (const auto &child : cell->child_iterators ())
421
+ n_particles_in_cell += particle_handler->n_particles_in_cell (child);
422
+ break ;
423
+ #else
403
424
case parallel::distributed::Triangulation<dim>::CELL_PERSIST:
404
425
case parallel::distributed::Triangulation<dim>::CELL_REFINE:
405
426
n_particles_in_cell = particle_handler->n_particles_in_cell (cell);
@@ -412,7 +433,7 @@ namespace aspect
412
433
for (const auto &child : cell->child_iterators ())
413
434
n_particles_in_cell += particle_handler->n_particles_in_cell (child);
414
435
break ;
415
-
436
+ # endif
416
437
default :
417
438
Assert (false , ExcInternalError ());
418
439
break ;
0 commit comments