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 shink_ratio{};
31
32 template <class Archive> void serialize(Archive& ar)
33 {
39 }
40 };
41
52 template <ModelType Model> class ParticlesContainer
53 {
54 public:
55 using UsedModel = Model;
56
60
62 std::size_t n_particle,
63 std::size_t _n_samples);
64 ParticlesContainer(); //=default;
72
77
78 // NOLINTBEGIN(cppcoreguidelines-non-private-member-variables-in-classes)
79 Model::SelfParticle model;
85 // NOLINTEND(cppcoreguidelines-non-private-member-variables-in-classes)
86
87 KOKKOS_INLINE_FUNCTION auto sample(const std::size_t idx,
88 const std::size_t i_sample) const
89 {
90 return random(idx, i_sample);
91 }
92
93 KOKKOS_INLINE_FUNCTION auto samples()
94 {
95 return random;
96 }
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
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
135 [[nodiscard]] KOKKOS_INLINE_FUNCTION bool
136 handle_division(const MC::pool_type& random_pool, std::size_t idx1) const;
137
141 [[nodiscard]] KOKKOS_INLINE_FUNCTION double
142 get_weight(std::size_t idx) const;
143
144 // HOST
145
151
156
160 template <class Archive> void save(Archive& ar) const;
161
165 template <class Archive> void load(Archive& ar);
166
170 [[nodiscard]] double get_allocation_factor() const noexcept;
171
172 void change_nsample(const std::size_t new_n_sample)
173 {
174 if (new_n_sample != n_samples)
175 {
176 n_samples = new_n_sample;
177 Kokkos::resize(random, n_allocated_elements, n_samples);
178 }
179 }
180
185 [[nodiscard]] std::size_t capacity() const noexcept;
186
190 [[nodiscard]] std::size_t get_inactive() const noexcept;
191
197 void update_and_remove_inactive(std::size_t out, std::size_t dead);
198
205 [[nodiscard]] KOKKOS_INLINE_FUNCTION std::size_t n_particles() const;
206
207 void change_runtime(RuntimeParameters&& parameters) noexcept;
208
209 void _sort(size_t n_c);
210
228 void remove_inactive_particles(std::size_t to_remove);
229
230// FIXME: used only in unit test
231#ifndef NDEBUG
232 [[maybe_unused]] [[nodiscard]] auto get_buffer_index() const;
233#endif
234
235 private:
236 Model::SelfParticle buffer_model;
238 Kokkos::View<uint64_t, Kokkos::SharedSpace> buffer_index;
241 std::size_t inactive_counter;
242
244 void _resize(std::size_t new_size, bool force = false);
245 std::size_t n_samples;
247
248 int begin;
249 int end;
250 };
251
252} // namespace MC
253
254/*
255 * IMPLEMENTATION
256 */
257namespace MC
258{
259 namespace
260 {
261
262 using TeamPolicy = Kokkos::TeamPolicy<ComputeSpace>;
263 using TeamMember = TeamPolicy::member_type;
264
303 template <ModelType M> struct CompactParticlesFunctor
304 {
305
306 CompactParticlesFunctor(MC::ParticleStatus _status,
307 M::SelfParticle _model,
308 MC::ParticlePositions _position,
309 MC::ParticleAges _ages,
310 std::size_t _to_remove,
311 std::size_t _last_used_index)
312 : status(std::move(_status)), model(std::move(_model)),
313 position(std::move(_position)), ages(std::move(_ages)),
314 offset("offset"), to_remove(_to_remove),
315 last_used_index(_last_used_index)
316 {
317
318 Kokkos::deep_copy(offset, 0);
319 }
320
321 KOKKOS_INLINE_FUNCTION void
322 operator()(const int i, std::size_t& update, const bool final) const
323 {
324
331
332 // Check if the particle is non-idle.
333 const bool is_inactive = (status(i) != MC::Status::Idle);
334
335 // Capture the current prefix sum of inactive particles.
336 const std::size_t scan_index = update;
337
338 // Increment the shared counter if the particle is inactive.
339 update += is_inactive ? 1 : 0;
340
341 // Final pass: Compact the container by replacing inactive particles.
342 if (final && is_inactive && scan_index < to_remove)
343 {
344
345 const auto inactive_slot = i; // Index of the "gap" to fill
346
347 // Atomically find a replacement particle from the end of the
348 // container. The replacement must be alive (idle) and not the same as
349 // the particle being removed.
350 // Index of the particle that will "fill" the inactive slot
351 auto replacement_index =
352 last_used_index - Kokkos::atomic_fetch_add(&offset(), 1);
353 while (status(replacement_index) != MC::Status::Idle ||
354 replacement_index == static_cast<std::size_t>(inactive_slot))
355 {
356 replacement_index =
357 last_used_index - Kokkos::atomic_fetch_add(&offset(), 1);
358 }
359
360 // Mark the removed particle's position as idle.
361 status(inactive_slot) = MC::Status::Idle;
362
363 // Merge position EZ
364 // No need atomic: final executed exactly once
365 Kokkos::atomic_exchange(&position(inactive_slot),
366 position(replacement_index));
367
368 // TODO Use hierachical parallism here, thread range is likely to work
369 for (std::size_t i_properties = 0; i_properties < M::n_var;
370 ++i_properties)
371 {
372 model(inactive_slot, i_properties) =
373 model(replacement_index, i_properties);
374 }
375 ages(inactive_slot, 0) = ages(replacement_index, 0);
376 ages(inactive_slot, 1) = ages(replacement_index, 1);
377 }
378 }
379
380 MC::ParticleStatus status;
381 M::SelfParticle model;
382 MC::ParticlePositions position;
383 MC::ParticleAges ages;
384 Kokkos::View<std::size_t, ComputeSpace> offset;
385 std::size_t to_remove;
386 std::size_t last_used_index;
387 };
388
405 template <ModelType M> struct InsertFunctor
406 {
407 InsertFunctor(std::size_t _original_size,
408 M::SelfParticle _model,
409 MC::ParticlePositions _position,
410 MC::ParticleAges _ages,
411 M::SelfParticle _buffer_model,
412 MC::ParticlePositions _buffer_position)
413 : original_size(_original_size), model(std::move(_model)),
414 ages(std::move(_ages)), position(std::move(_position)),
415 buffer_model(std::move(_buffer_model)),
416 buffer_position(std::move(_buffer_position))
417 {
418 }
419
420 // KOKKOS_INLINE_FUNCTION void operator()(const int i) const
421 // {
422 // position(original_size + i) = buffer_position(i);
423 // for (std::size_t j = 0; j < M::n_var; ++j)
424 // {
425 // model(original_size + i, j) = buffer_model(i, j);
426 // }
427 // }
428
429 // TODO try this functor
430 // TODO find a way to organise data to not copy non needed data (like
431 // contribs). Split model in two arrays?
432
433 KOKKOS_INLINE_FUNCTION
434 void operator()(const TeamMember& team) const
435 {
436 auto range = M::n_var;
437 const int i = team.league_rank();
438
439 Kokkos::parallel_for(
440 Kokkos::TeamVectorRange(team, range),
441 [&](const int& j)
442 { model(original_size + i, j) = buffer_model(i, j); });
443 position(original_size + i) = buffer_position(i);
444
445 // Actually needs buffer to store mother's hydraulic time
446 // But set new hydraulic time to 0 to not create new buffer a save
447 // memory usage
448 ages(original_size + i, 0) = 0;
449 ages(original_size + i, 1) = 0;
450 }
451
452 std::size_t original_size;
453 M::SelfParticle model;
454 MC::ParticleAges ages;
455 MC::ParticlePositions position;
456 M::SelfParticle buffer_model;
457 MC::ParticlePositions buffer_position;
458 };
459
460 }; // namespace
461
462 template <ModelType Model>
463 [[nodiscard]] KOKKOS_INLINE_FUNCTION std::size_t
468
469 template <ModelType Model>
470 [[nodiscard]] std::size_t
472 {
473 return inactive_counter;
474 }
475
476 template <ModelType Model> void ParticlesContainer<Model>::force_remove_dead()
477 {
479 }
480#ifndef NDEBUG
481
482 template <ModelType Model>
483 [[maybe_unused]] [[nodiscard]] auto
488#endif
489
490 template <ModelType Model>
491 [[nodiscard]] double
493 {
494 return rt_params.allocation_factor;
495 }
496
497 template <ModelType Model>
498 [[nodiscard]] std::size_t ParticlesContainer<Model>::capacity() const noexcept
499 {
501 }
502
503 template <ModelType Model>
504 template <class Archive>
506 {
507 // Basically store everything that is usefull
509
510 // Allocation is done in the deserialize
511 deserialize_view(ar, weights);
512 deserialize_view(ar, position);
513 deserialize_view(ar, status);
514 deserialize_view(ar, model);
515 deserialize_view(ar, ages);
516#ifndef NDEBUG
517 Kokkos::printf("ParticlesContainer::load: Check if load_tuning_constant "
518 "works with different value");
519#endif
520
521 __allocate_buffer__(); // Dont forget to allocate buffer
522 }
523
524 template <ModelType Model>
525 template <class Archive>
526 void ParticlesContainer<Model>::save(Archive& ar) const
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
541 const std::size_t dead)
542 {
543 // Actually out particle and dead ones are just inactive
544 inactive_counter += out;
545 inactive_counter += dead;
546
547 const auto _threshold = std::max(
548 rt_params.minimum_dead_particle_removal,
549 static_cast<uint64_t>(static_cast<double>(n_used_elements) *
550 rt_params.dead_particle_ratio_threshold));
551
552 // TODO: May change threshold as container is now cleaned before exporting,
553 if (inactive_counter > _threshold)
554 {
556 }
557 }
558
559 template <ModelType Model>
560 KOKKOS_INLINE_FUNCTION bool
562 std::size_t idx1) const
563 {
564 if (Kokkos::atomic_load(&buffer_index()) < buffer_model.extent(0))
565 {
566 const auto idx2 = Kokkos::atomic_fetch_add(&buffer_index(), 1);
567 Model::division(random_pool, idx1, idx2, model, buffer_model);
568 buffer_position(idx2) = position(idx1);
569 ages(idx1, 1) = 0;
570 return true;
571 }
572 return false;
573 }
574
575 template <ModelType Model> void ParticlesContainer<Model>::merge_buffer()
576 {
577 PROFILE_SECTION("ParticlesContainer::merge_buffer")
578 const auto original_size = n_used_elements;
579 const auto n_add_item = buffer_index();
580 if (n_add_item == 0)
581 {
582 return;
583 }
584 _resize(original_size + n_add_item);
585 Kokkos::parallel_for("insert_merge",
586 TeamPolicy(n_add_item, Kokkos::AUTO, Model::n_var),
587 InsertFunctor<Model>(original_size,
588 model,
589 position,
590 ages,
593
594 buffer_index() = 0;
595 n_used_elements += n_add_item;
597 }
598
599 template <ModelType Model>
600 void ParticlesContainer<Model>::_resize(std::size_t new_size, bool force)
601 {
602 PROFILE_SECTION("ParticlesContainer::_resize")
603
604 // Ensure new_size is greater than zero
605 if (new_size > 0)
606 {
607 // Determine if resizing is necessary based on the condition
608 if (new_size > n_allocated_elements || force)
609 {
610 // Calculate the new allocated size
611 const auto new_allocated_size = static_cast<std::size_t>(std::ceil(
612 static_cast<double>(new_size) * rt_params.allocation_factor));
613
614 // Update the allocated size
615 n_allocated_elements = new_allocated_size;
616
617 // Perform the resizing on all relevant data containers
618 Kokkos::resize(position, n_allocated_elements);
619 Kokkos::resize(model,
621 Model::n_var); // use 2nd dim resize if dynamic
622 Kokkos::resize(status, n_allocated_elements);
623 Kokkos::resize(ages, n_allocated_elements);
624
625 if (n_samples != 0)
626 {
627 Kokkos::resize(random, n_allocated_elements, n_samples);
628 }
629
630 // Handle resizing for weights based on model type
631 if constexpr (ConstWeightModelType<Model>)
632 {
633 Kokkos::resize(weights, 1); // Fixed size for ConstWeightModelType
634 }
635 else
636 {
637 Kokkos::resize(weights, n_allocated_elements);
638 }
639 }
640 }
641 }
642
643 template <ModelType Model>
645 {
646 PROFILE_SECTION("ParticlesContainer::__allocate_buffer__")
647 auto buffer_size = buffer_position.extent(0);
648 const auto buffer_ratio = rt_params.buffer_ratio;
649 if (static_cast<double>(buffer_size) /
650 static_cast<double>(n_allocated_elements) <
651 buffer_ratio)
652 {
653 buffer_size = static_cast<std::size_t>(
654 std::ceil(static_cast<double>(n_allocated_elements) * buffer_ratio));
655
656 // Realloc because not needed to keep buffer as it has been copied
657 Kokkos::realloc(buffer_position, buffer_size);
658 Kokkos::realloc(buffer_model,
659 buffer_size,
660 Model::n_var); // use 2nd dim resize if dynamic
661 buffer_index() = 0;
662 }
663 }
664
665 template <ModelType M>
667 std::size_t n_particle,
668 std::size_t _n_samples)
669 : model(Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_model"),
670 0,
671 0),
672 position(Kokkos::view_alloc(Kokkos::WithoutInitializing,
673 "particle_position"),
674 0),
675 status(
676 Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_status"),
677 0),
678 weights(
679 Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_weigth"),
680 0),
681 ages(Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_age"),
682 0),
683 random(
684 Kokkos::view_alloc(Kokkos::WithoutInitializing, "particle_random"),
685 0,
686 0),
687 buffer_model("buffer_particle_model", 0),
688 buffer_position("buffer_particle_position", 0),
689 buffer_index("buffer_index"), n_allocated_elements(0),
690 n_used_elements(n_particle), inactive_counter(0), n_samples(_n_samples),
691 rt_params(rt_param), begin(0), end(0)
692 {
693
694 // load_tuning_constant();
695
696 if (n_particle != 0)
697 {
698 _resize(n_particle);
700 }
701
702 auto bounds = M::get_bounds();
703 begin = bounds.begin;
704 end = bounds.end;
705 }
706
707 template <ModelType M>
712
713 template <ModelType M>
715 {
716
717 PROFILE_SECTION("ParticlesContainer::remove_inactive_particles")
718 if (to_remove == 0)
719 {
720 return;
721 }
722
723 if (to_remove == n_used_elements)
724 {
725 _resize(0, true);
726 n_used_elements = 0;
728 }
729 else if (to_remove > n_used_elements)
730 {
731 throw std::runtime_error(
732 "remove_inactive_particles: Error in kernel cannot remove more "
733 "element than existing");
734 }
735 else
736 {
737
738 const auto new_used_item = n_used_elements - to_remove;
739
740 const auto last_used_index = n_used_elements - 1;
741 Kokkos::parallel_scan(
742 "find_and_fill_gap",
743 Kokkos::RangePolicy<ComputeSpace>(0, n_used_elements),
744 CompactParticlesFunctor<M>(
745 status, model, position, ages, to_remove, last_used_index));
746
747 Kokkos::fence();
748 KOKKOS_ASSERT(this->position.extent(0) == n_allocated_elements);
749 KOKKOS_ASSERT(this->model.extent(0) == n_allocated_elements);
750 KOKKOS_ASSERT(this->status.extent(0) == n_allocated_elements);
751 if (n_samples != 0)
752 {
753 KOKKOS_ASSERT(this->random.extent(0) == n_allocated_elements);
754 KOKKOS_ASSERT(this->random.extent(1) == n_samples);
755 }
756
757 n_used_elements = new_used_item;
758 if (static_cast<double>(n_used_elements) /
759 static_cast<double>(n_allocated_elements) <=
760 rt_params.shink_ratio)
761 {
762 _resize(n_used_elements * rt_params.allocation_factor, false);
763 }
765 };
766 }
767
768 template <ModelType M>
769 [[nodiscard]] KOKKOS_INLINE_FUNCTION double
771 {
772 if constexpr (ConstWeightModelType<M>)
773 {
774 return weights(0);
775 }
776 else
777 {
778 return weights(idx);
779 }
780 }
781
782 template <ModelType M>
783 void
785 {
786 this->rt_params = parameters;
787 }
788
789 template <ModelType M> void ParticlesContainer<M>::_sort(std::size_t n_c)
790 {
791 (void)n_c;
792 }
793
794} // namespace MC
795
796#endif
RuntimeParameters rt_params
Definition particles_container.hpp:246
void update_and_remove_inactive(std::size_t out, std::size_t dead)
Definition particles_container.hpp:540
void change_nsample(const std::size_t new_n_sample)
Definition particles_container.hpp:172
double get_allocation_factor() const noexcept
Get the allocation factor.
Definition particles_container.hpp:492
~ParticlesContainer()=default
Default destructor.
void force_remove_dead()
Clean all the non-idle particle even if number smaller than threshold.
Definition particles_container.hpp:476
KOKKOS_INLINE_FUNCTION auto sample(const std::size_t idx, const std::size_t i_sample) const
Definition particles_container.hpp:87
Model::SelfParticle model
Definition particles_container.hpp:79
KOKKOS_INLINE_FUNCTION std::size_t n_particles() const
Definition particles_container.hpp:464
ParticleWeigths weights
Definition particles_container.hpp:82
int end
Definition particles_container.hpp:249
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:561
ParticleSamples random
Definition particles_container.hpp:84
Kokkos::View< uint64_t, Kokkos::SharedSpace > buffer_index
Definition particles_container.hpp:238
void _resize(std::size_t new_size, bool force=false)
Definition particles_container.hpp:600
void remove_inactive_particles(std::size_t to_remove)
Definition particles_container.hpp:714
void __allocate_buffer__()
Definition particles_container.hpp:644
ParticlesContainer()
Definition particles_container.hpp:708
ParticlesContainer & operator=(ParticlesContainer &&)=default
std::size_t get_inactive() const noexcept
Definition particles_container.hpp:471
std::size_t capacity() const noexcept
Returns the total number of elements the container can hold without reallocating.
Definition particles_container.hpp:498
uint64_t n_used_elements
Definition particles_container.hpp:240
auto get_buffer_index() const
Definition particles_container.hpp:484
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:575
Model::SelfParticle buffer_model
Definition particles_container.hpp:236
KOKKOS_INLINE_FUNCTION double get_weight(std::size_t idx) const
Return the particle weight.
Definition particles_container.hpp:770
std::size_t n_samples
Definition particles_container.hpp:245
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:666
ParticlePositions buffer_position
Definition particles_container.hpp:237
MC::ParticlePositions position
Definition particles_container.hpp:80
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:241
ParticleAges ages
Definition particles_container.hpp:83
Model UsedModel
Definition particles_container.hpp:55
ParticlesContainer & operator=(const ParticlesContainer &)=default
void _sort(size_t n_c)
Definition particles_container.hpp:789
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:505
ParticlesContainer(ParticlesContainer &&)=default
std::size_t n_allocated_elements
Definition particles_container.hpp:239
void change_runtime(RuntimeParameters &&parameters) noexcept
Definition particles_container.hpp:784
int begin
Definition particles_container.hpp:248
KOKKOS_INLINE_FUNCTION auto samples()
Definition particles_container.hpp:93
MC::ParticleStatus status
Definition particles_container.hpp:81
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:14
Kokkos::View< double *, ComputeSpace > ParticleWeigths
Definition alias.hpp:72
Kokkos::View< Status *, ComputeSpace > ParticleStatus
Definition alias.hpp:71
Kokkos::View< uint64_t *, ComputeSpace > ParticlePositions
Definition alias.hpp:70
@ Idle
Definition alias.hpp:58
ParticleAgesBase< ComputeSpace > ParticleAges
Definition alias.hpp:73
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:38
Kokkos:: View< Kokkos::Experimental::half_t **, Kokkos::LayoutRight, ComputeSpace > ParticleSamples
Definition alias.hpp:75
Definition particles_container.hpp:25
double shink_ratio
Definition particles_container.hpp:29
void serialize(Archive &ar)
Definition particles_container.hpp:32
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 allocation_factor
Definition particles_container.hpp:28