BioCMAMC-ST
particles_container.hpp
1#ifndef __PARTICLES_CONTAINER_HPP__
2#define __PARTICLES_CONTAINER_HPP__
3
4#include "Kokkos_Macros.hpp"
5#include "common/alg.hpp"
6#include <Kokkos_Core.hpp>
7#include <biocma_cst_config.hpp>
8#include <cmath>
9#include <common/common.hpp>
10#include <common/env_var.hpp>
11#include <common/has_serialize.hpp>
12#include <cstdint>
13#include <mc/alias.hpp>
14#include <mc/prng/prng.hpp>
15#include <mc/traits.hpp>
16
17#include <Kokkos_Sort.hpp>
18#include <sorting/Kokkos_BinSortPublicAPI.hpp>
19#include <sorting/impl/Kokkos_CopyOpsForBinSortImpl.hpp>
20#include <sorting/impl/Kokkos_SortByKeyImpl.hpp>
21#include <utility>
22namespace MC
23{
24
26 {
28 double buffer_ratio{};
29 double allocation_factor = {};
30 double shink_ratio{};
32
33 template <class Archive>
34 void
43 };
44
56 template <ModelType Model> class ParticlesContainer
57 {
58 public:
59 using UsedModel = Model;
60
64
66 std::size_t n_particle,
67 std::size_t _n_samples);
68 ParticlesContainer(); //=default;
76
81
82 // NOLINTBEGIN(cppcoreguidelines-non-private-member-variables-in-classes)
83 Model::SelfParticle model;
89 // NOLINTEND(cppcoreguidelines-non-private-member-variables-in-classes)
90
92 KOKKOS_INLINE_FUNCTION auto sample(std::size_t idx,
93 std::size_t i_sample) const;
94
96 KOKKOS_INLINE_FUNCTION auto samples();
97
102 template <typename CviewType>
103 KOKKOS_INLINE_FUNCTION void
104 get_contributions(const std::size_t idx,
105 const CviewType& contributions) const
106 {
107 if (begin < end)
108 {
109 static_assert(ConstWeightModelType<Model>,
110 "ModelType: Constapply_weight()");
111
112 const double weight = get_weight(idx);
113 const auto pos = position(idx);
114 auto access = contributions.access();
115
116 for (int i = begin; i < end; ++i)
117 {
118 const int rel = i - begin;
119 access(rel, pos) += weight * model(idx, i);
120 }
121 }
122 }
123
137 [[nodiscard]] KOKKOS_INLINE_FUNCTION bool
138 handle_division(const MC::pool_type& random_pool, std::size_t idx1) const;
139
143 [[nodiscard]] KOKKOS_INLINE_FUNCTION double
144 get_weight(std::size_t idx) const;
145
146 // HOST
147
153
158
162 template <class Archive> void save(Archive& ar) const;
163
167 template <class Archive> void load(Archive& ar);
168
172 [[nodiscard]] double get_allocation_factor() const noexcept;
173
174 void change_nsample(std::size_t new_n_sample);
179 [[nodiscard]] std::size_t capacity() const noexcept;
180
184 [[nodiscard]] std::size_t get_inactive() const noexcept;
185
191 void update_and_remove_inactive(std::size_t out, std::size_t dead);
192
199 [[nodiscard]] KOKKOS_INLINE_FUNCTION std::size_t n_particles() const;
200
201 void change_runtime(RuntimeParameters&& parameters) noexcept;
202
203 void _sort(size_t n_c);
204
222 void remove_inactive_particles(std::size_t to_remove);
223
224// FIXME: used only in unit test
225#ifndef NDEBUG
226 [[maybe_unused]] [[nodiscard]] auto get_buffer_index() const;
227#endif
228
229 private:
230 Model::SelfParticle buffer_model;
232 Kokkos::View<uint64_t, Kokkos::SharedSpace> buffer_index;
235 std::size_t inactive_counter;
236
238 void _resize(std::size_t new_size, bool force = false);
239 std::size_t n_samples;
241
242 int begin;
243 int end;
244 };
245
246} // namespace MC
247
248/*
249 * IMPLEMENTATION
250 */
251namespace MC
252{
253 namespace
254 {
255
256 using TeamPolicy = Kokkos::TeamPolicy<ComputeSpace>;
257 using TeamMember = TeamPolicy::member_type;
258
297 template <ModelType M> struct CompactParticlesFunctor
298 {
299
300 CompactParticlesFunctor(MC::ParticleStatus _status,
301 M::SelfParticle _model,
302 MC::ParticlePositions _position,
303 MC::ParticleAges _ages,
304 std::size_t _to_remove,
305 std::size_t _last_used_index)
306 : status(std::move(_status)), model(std::move(_model)),
307 position(std::move(_position)), ages(std::move(_ages)),
308 offset("offset"), to_remove(_to_remove),
309 last_used_index(_last_used_index)
310 {
311
312 Kokkos::deep_copy(offset, 0);
313 }
314
315 KOKKOS_INLINE_FUNCTION void
316 operator()(const int i, std::size_t& update, const bool final) const
317 {
318
325
326 // Check if the particle is non-idle.
327 const bool is_inactive = (status(i) != MC::Status::Idle);
328
329 // Capture the current prefix sum of inactive particles.
330 const std::size_t scan_index = update;
331
332 // Increment the shared counter if the particle is inactive.
333 update += is_inactive ? 1 : 0;
334
335 // Final pass: Compact the container by replacing inactive particles.
336 if (final && is_inactive && scan_index < to_remove)
337 {
338
339 const auto inactive_slot = i; // Index of the "gap" to fill
340
341 // Atomically find a replacement particle from the end of the
342 // container. The replacement must be alive (idle) and not the same as
343 // the particle being removed.
344 // Index of the particle that will "fill" the inactive slot
345 auto replacement_index
346 = last_used_index - Kokkos::atomic_fetch_add(&offset(), 1);
347 while (status(replacement_index) != MC::Status::Idle
348 || replacement_index
349 == static_cast<std::size_t>(inactive_slot))
350 {
351 replacement_index
352 = last_used_index - Kokkos::atomic_fetch_add(&offset(), 1);
353 }
354
355 // Mark the removed particle's position as idle.
356 status(inactive_slot) = MC::Status::Idle;
357
358 // Merge position EZ
359 // No need atomic: final executed exactly once
360 Kokkos::atomic_exchange(&position(inactive_slot),
361 position(replacement_index));
362
363 // TODO Use hierachical parallism here, thread range is likely to work
364 for (std::size_t i_properties = 0; i_properties < M::n_var;
365 ++i_properties)
366 {
367 model(inactive_slot, i_properties)
368 = model(replacement_index, i_properties);
369 }
370 ages(inactive_slot, 0) = ages(replacement_index, 0);
371 ages(inactive_slot, 1) = ages(replacement_index, 1);
372 }
373 }
374
375 MC::ParticleStatus status;
376 M::SelfParticle model;
377 MC::ParticlePositions position;
378 MC::ParticleAges ages;
379 Kokkos::View<std::size_t, ComputeSpace> offset;
380 std::size_t to_remove;
381 std::size_t last_used_index;
382 };
383
400 template <ModelType M> struct InsertFunctor
401 {
402 InsertFunctor(std::size_t _original_size,
403 M::SelfParticle _model,
404 MC::ParticlePositions _position,
405 MC::ParticleAges _ages,
406 M::SelfParticle _buffer_model,
407 MC::ParticlePositions _buffer_position)
408 : original_size(_original_size), model(std::move(_model)),
409 ages(std::move(_ages)), position(std::move(_position)),
410 buffer_model(std::move(_buffer_model)),
411 buffer_position(std::move(_buffer_position))
412 {
413 }
414
415 // KOKKOS_INLINE_FUNCTION void operator()(const int i) const
416 // {
417 // position(original_size + i) = buffer_position(i);
418 // for (std::size_t j = 0; j < M::n_var; ++j)
419 // {
420 // model(original_size + i, j) = buffer_model(i, j);
421 // }
422 // }
423
424 // TODO try this functor
425 // TODO find a way to organise data to not copy non needed data (like
426 // contribs). Split model in two arrays?
427
428 KOKKOS_INLINE_FUNCTION
429 void
430 operator()(const TeamMember& team) const
431 {
432 auto range = M::n_var;
433 const int i = team.league_rank();
434
435 Kokkos::parallel_for(
436 Kokkos::TeamVectorRange(team, range),
437 [&](const int& j)
438 { model(original_size + i, j) = buffer_model(i, j); });
439 position(original_size + i) = buffer_position(i);
440
441 // Actually needs buffer to store mother's hydraulic time
442 // But set new hydraulic time to 0 to not create new buffer a save
443 // memory usage
444 ages(original_size + i, 0) = 0;
445 ages(original_size + i, 1) = 0;
446 }
447
448 std::size_t original_size;
449 M::SelfParticle model;
450 MC::ParticleAges ages;
451 MC::ParticlePositions position;
452 M::SelfParticle buffer_model;
453 MC::ParticlePositions buffer_position;
454 };
455
456 }; // namespace
457
458 template <ModelType Model>
459 [[nodiscard]] KOKKOS_INLINE_FUNCTION std::size_t
464
465 template <ModelType Model>
466 [[nodiscard]] std::size_t
468 {
469 return inactive_counter;
470 }
471
472 template <ModelType Model>
473 void
478#ifndef NDEBUG
479
480 template <ModelType Model>
481 [[maybe_unused]] [[nodiscard]] auto
486#endif
487
488 template <ModelType Model>
489 [[nodiscard]] double
491 {
492 return rt_params.allocation_factor;
493 }
494
495 template <ModelType Model>
496 [[nodiscard]] std::size_t
498 {
500 }
501
502 template <ModelType Model>
503 KOKKOS_INLINE_FUNCTION auto
508 template <ModelType Model>
509 KOKKOS_INLINE_FUNCTION auto
511 const std::size_t i_sample) const
512 {
513 return random(idx, i_sample);
514 }
515
516 template <ModelType Model>
517 template <class Archive>
518 void
520 {
521 // Basically store everything that is usefull
523
524 // Allocation is done in the deserialize
525 deserialize_view(ar, weights);
526 deserialize_view(ar, position);
527 deserialize_view(ar, status);
528 deserialize_view(ar, model);
529 deserialize_view(ar, ages);
530#ifndef NDEBUG
531 Kokkos::printf("ParticlesContainer::load: Check if load_tuning_constant "
532 "works with different value");
533#endif
534
535 __allocate_buffer__(); // Dont forget to allocate buffer
536 }
537
538 template <ModelType Model>
539 template <class Archive>
540 void
542 {
543
544 // Basically store everything that is usefull
546 serialize_view(ar, weights);
547 serialize_view(ar, position);
548 serialize_view(ar, status);
549 serialize_view(ar, model);
550 serialize_view(ar, ages);
551 }
552
553 template <ModelType Model>
554 void
555 ParticlesContainer<Model>::change_nsample(const std::size_t new_n_sample)
556 {
557 if (new_n_sample != n_samples)
558 {
559 n_samples = new_n_sample;
560 Kokkos::resize(random, n_allocated_elements, n_samples);
561 }
562 }
563
564 template <ModelType Model>
565 void
567 const std::size_t dead)
568 {
569 // Actually out particle and dead ones are just inactive
570 inactive_counter += out;
571 inactive_counter += dead;
572
573 const auto _threshold = std::max(
574 rt_params.minimum_dead_particle_removal,
575 static_cast<uint64_t>(static_cast<double>(n_used_elements)
576 * rt_params.dead_particle_ratio_threshold));
577
578 // TODO: May change threshold as container is now cleaned before exporting,
579 if (inactive_counter > _threshold)
580 {
582 }
583 }
584
585 template <ModelType Model>
586 KOKKOS_INLINE_FUNCTION bool
588 std::size_t idx1) const
589 {
590 if (Kokkos::atomic_load(&buffer_index()) < buffer_model.extent(0))
591 {
592 const auto idx2 = Kokkos::atomic_fetch_add(&buffer_index(), 1);
593 Model::division(random_pool, idx1, idx2, model, buffer_model);
594 buffer_position(idx2) = position(idx1);
595 ages(idx1, 1) = 0;
596 return true;
597 }
598 return false;
599 }
600
601 template <ModelType Model>
602 void
604 {
605 PROFILE_SECTION("ParticlesContainer::merge_buffer")
606 const auto original_size = n_used_elements;
607 const auto n_add_item = buffer_index();
608 if (n_add_item == 0)
609 {
610 return;
611 }
612 _resize(original_size + n_add_item);
613 Kokkos::parallel_for("insert_merge",
614 TeamPolicy(n_add_item, Kokkos::AUTO, Model::n_var),
615 InsertFunctor<Model>(original_size,
616 model,
617 position,
618 ages,
621
622 buffer_index() = 0;
623 n_used_elements += n_add_item;
625 }
626
627 template <ModelType Model>
628 void
629 ParticlesContainer<Model>::_resize(std::size_t new_size, bool force)
630 {
631 PROFILE_SECTION("ParticlesContainer::_resize")
632
633 // Ensure new_size is greater than zero
634 if (new_size > 0)
635 {
636 // Determine if resizing is necessary based on the condition
637 if (new_size > n_allocated_elements || force)
638 {
639 // Calculate the new allocated size
640 const auto new_allocated_size = static_cast<std::size_t>(std::ceil(
641 static_cast<double>(new_size) * rt_params.allocation_factor));
642
643 // Update the allocated size
644 n_allocated_elements = new_allocated_size;
645
646 // Perform the resizing on all relevant data containers
647 Kokkos::resize(position, n_allocated_elements);
648 Kokkos::resize(model,
650 Model::n_var); // use 2nd dim resize if dynamic
651 Kokkos::resize(status, n_allocated_elements);
652 Kokkos::resize(ages, n_allocated_elements);
653
654 if (n_samples != 0)
655 {
656 Kokkos::resize(random, n_allocated_elements, n_samples);
657 }
658
659 // Handle resizing for weights based on model type
660 if constexpr (ConstWeightModelType<Model>)
661 {
662 Kokkos::resize(weights,
663 1); // Fixed size for ConstWeightModelType
664 }
665 else
666 {
667 Kokkos::resize(weights, n_allocated_elements);
668 }
669 }
670 }
671 }
672
673 template <ModelType Model>
674 void
676 {
677 PROFILE_SECTION("ParticlesContainer::__allocate_buffer__")
678 auto buffer_size = buffer_position.extent(0);
679 const auto buffer_ratio = rt_params.buffer_ratio;
680 if (static_cast<double>(buffer_size)
681 / static_cast<double>(n_allocated_elements)
682 < buffer_ratio)
683 {
684 buffer_size = static_cast<std::size_t>(
685 std::ceil(static_cast<double>(n_allocated_elements) * buffer_ratio));
686
687 // Realloc because not needed to keep buffer as it has been copied
688 Kokkos::realloc(buffer_position, buffer_size);
689 Kokkos::realloc(buffer_model,
690 buffer_size,
691 Model::n_var); // use 2nd dim resize if dynamic
692 buffer_index() = 0;
693 }
694 }
695
696 template <ModelType M>
698 std::size_t n_particle,
699 std::size_t _n_samples)
700 : model(Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_model"),
701 0,
702 0),
703 position(Kokkos::view_alloc(Kokkos::WithoutInitializing,
704 "particle_position"),
705 0),
706 status(
707 Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_status"),
708 0),
709 weights(
710 Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_weigth"),
711 0),
712 ages(Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_age"),
713 0),
714 random(
715 Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_random"),
716 0,
717 0),
718 buffer_model("buffer_particle_model", 0),
719 buffer_position("buffer_particle_position", 0),
720 buffer_index("buffer_index"), n_allocated_elements(0),
721 n_used_elements(n_particle), inactive_counter(0), n_samples(_n_samples),
722 rt_params(rt_param), begin(0), end(0)
723 {
724
725 // load_tuning_constant();
726
727 if (n_particle != 0)
728 {
729 _resize(n_particle);
731 }
732
733 const auto bounds = M::get_bounds();
734 begin = bounds.begin;
735 end = bounds.end;
736 }
737
738 template <ModelType M>
743
744 template <ModelType M>
745 void
747 {
748
749 PROFILE_SECTION("ParticlesContainer::remove_inactive_particles")
750 if (to_remove == 0)
751 {
752 return;
753 }
754
755 if (to_remove == n_used_elements)
756 {
757 _resize(0, true);
758 n_used_elements = 0;
760 }
761 else if (to_remove > n_used_elements)
762 {
763 throw std::runtime_error(
764 "remove_inactive_particles: Error in kernel cannot remove more "
765 "element than existing");
766 }
767 else
768 {
769
770 const auto new_used_item = n_used_elements - to_remove;
771
772 const auto last_used_index = n_used_elements - 1;
773 Kokkos::parallel_scan(
774 "find_and_fill_gap",
775 Kokkos::RangePolicy<ComputeSpace>(0, n_used_elements),
776 CompactParticlesFunctor<M>(
777 status, model, position, ages, to_remove, last_used_index));
778
779 Kokkos::fence();
780 KOKKOS_ASSERT(this->position.extent(0) == n_allocated_elements);
781 KOKKOS_ASSERT(this->model.extent(0) == n_allocated_elements);
782 KOKKOS_ASSERT(this->status.extent(0) == n_allocated_elements);
783 if (n_samples != 0)
784 {
785 KOKKOS_ASSERT(this->random.extent(0) == n_allocated_elements);
786 KOKKOS_ASSERT(this->random.extent(1) == n_samples);
787 }
788
789 n_used_elements = new_used_item;
790 if (static_cast<double>(n_used_elements)
791 / static_cast<double>(n_allocated_elements)
792 <= rt_params.shink_ratio)
793 {
794 // Kokkos::printf("SHRINK \r\n");
795 // force to true if we want to shrink
796 _resize(n_used_elements * rt_params.allocation_factor, true);
797 }
799 };
800 }
801
802 template <ModelType M>
803 [[nodiscard]] KOKKOS_INLINE_FUNCTION double
805 {
806 if constexpr (ConstWeightModelType<M>)
807 {
808 return weights(0);
809 }
810 else
811 {
812 return weights(idx);
813 }
814 }
815
816 template <ModelType M>
817 void
819 {
820 this->rt_params = parameters;
821 }
822
823 template <ModelType M>
824 void
826 {
828 n_c,
831 position, // ref
832 position,
833 model,
834 status);
835 }
836
837} // namespace MC
838
839#endif
RuntimeParameters rt_params
Definition particles_container.hpp:240
void update_and_remove_inactive(std::size_t out, std::size_t dead)
Definition particles_container.hpp:566
double get_allocation_factor() const noexcept
Get the allocation factor.
Definition particles_container.hpp:490
~ParticlesContainer()=default
Default destructor.
void force_remove_dead()
Clean all the non-idle particle even if number smaller than threshold.
Definition particles_container.hpp:474
Model::SelfParticle model
Definition particles_container.hpp:83
KOKKOS_INLINE_FUNCTION std::size_t n_particles() const
Definition particles_container.hpp:460
ParticleWeigths weights
Definition particles_container.hpp:86
int end
Definition particles_container.hpp:243
KOKKOS_INLINE_FUNCTION bool handle_division(const MC::pool_type &random_pool, std::size_t idx1) const
Attempts to spawn a new particle by performing division on the specified particle.
Definition particles_container.hpp:587
ParticleSamples random
Definition particles_container.hpp:88
Kokkos::View< uint64_t, Kokkos::SharedSpace > buffer_index
Definition particles_container.hpp:232
void _resize(std::size_t new_size, bool force=false)
Definition particles_container.hpp:629
void change_nsample(std::size_t new_n_sample)
Definition particles_container.hpp:555
void remove_inactive_particles(std::size_t to_remove)
Definition particles_container.hpp:746
void __allocate_buffer__()
Definition particles_container.hpp:675
ParticlesContainer()
Definition particles_container.hpp:739
ParticlesContainer & operator=(ParticlesContainer &&)=default
std::size_t get_inactive() const noexcept
Definition particles_container.hpp:467
std::size_t capacity() const noexcept
Definition particles_container.hpp:497
uint64_t n_used_elements
Definition particles_container.hpp:234
auto get_buffer_index() const
Definition particles_container.hpp:482
KOKKOS_INLINE_FUNCTION void get_contributions(const std::size_t idx, const CviewType &contributions) const
Get the contribution if particle at index idx.
Definition particles_container.hpp:104
void merge_buffer()
Insert particle buffer into the main container.
Definition particles_container.hpp:603
Model::SelfParticle buffer_model
Definition particles_container.hpp:230
KOKKOS_INLINE_FUNCTION double get_weight(std::size_t idx) const
Return the particle weight.
Definition particles_container.hpp:804
std::size_t n_samples
Definition particles_container.hpp:239
KOKKOS_INLINE_FUNCTION auto sample(std::size_t idx, std::size_t i_sample) const
get the sample at idx/i_sample position
Definition particles_container.hpp:510
ParticlesContainer(RuntimeParameters rt_param, std::size_t n_particle, std::size_t _n_samples)
Alias for the model used by the container.
Definition particles_container.hpp:697
ParticlePositions buffer_position
Definition particles_container.hpp:231
MC::ParticlePositions position
Definition particles_container.hpp:84
void save(Archive &ar) const
Save data into ar for serialization.
Definition particles_container.hpp:541
std::size_t inactive_counter
Definition particles_container.hpp:235
ParticleAges ages
Definition particles_container.hpp:87
Model UsedModel
Definition particles_container.hpp:59
ParticlesContainer & operator=(const ParticlesContainer &)=default
void _sort(size_t n_c)
Definition particles_container.hpp:825
ParticlesContainer(const ParticlesContainer &)=default
Default copy and move constructors and assignment operators.
void load(Archive &ar)
Load data from ar for deserialization.
Definition particles_container.hpp:519
ParticlesContainer(ParticlesContainer &&)=default
std::size_t n_allocated_elements
Definition particles_container.hpp:233
void change_runtime(RuntimeParameters &&parameters) noexcept
Definition particles_container.hpp:818
int begin
Definition particles_container.hpp:242
KOKKOS_INLINE_FUNCTION auto samples()
get the sample view
Definition particles_container.hpp:504
MC::ParticleStatus status
Definition particles_container.hpp:85
Concept to check if a model type has uniform_weight
Definition traits.hpp:204
Namespace that contains classes and structures related to Monte Carlo (MC) simulations.
Definition alias.hpp:15
Kokkos::View< double *, ComputeSpace > ParticleWeigths
Definition alias.hpp:75
Kokkos::View< Status *, ComputeSpace > ParticleStatus
Definition alias.hpp:74
Kokkos::View< uint64_t *, ComputeSpace > ParticlePositions
Definition alias.hpp:73
@ Idle
Definition alias.hpp:59
ParticleAgesBase< ComputeSpace > ParticleAges
Definition alias.hpp:76
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:39
Kokkos:: View< Kokkos::Experimental::half_t **, Kokkos::LayoutRight, ComputeSpace > ParticleSamples
Definition alias.hpp:78
void sort_soa(std::size_t n0, std::size_t nn, std::size_t n_e, std::size_t n_max, RefView ref_sorting, ViewType... views)
Definition alg.hpp:65
Definition particles_container.hpp:26
double shink_ratio
Definition particles_container.hpp:30
void serialize(Archive &ar)
Definition particles_container.hpp:35
double buffer_ratio
Definition particles_container.hpp:28
double dead_particle_ratio_threshold
Definition particles_container.hpp:31
uint64_t minimum_dead_particle_removal
Definition particles_container.hpp:27
double allocation_factor
Definition particles_container.hpp:29