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