Shamrock 2025.10.0
Astrophysical Code
Loading...
Searching...
No Matches
Model.cpp
Go to the documentation of this file.
1// -------------------------------------------------------//
2//
3// SHAMROCK code for hydrodynamics
4// Copyright (c) 2021-2026 Timothée David--Cléris <tim.shamrock@proton.me>
5// SPDX-License-Identifier: CeCILL Free Software License Agreement v2.1
6// Shamrock is licensed under the CeCILL 2.1 License, see LICENSE for more information
7//
8// -------------------------------------------------------//
9
20#include "shambase/memory.hpp"
22#include "shambase/string.hpp"
23#include "shamcomm/logs.hpp"
37#include <functional>
38#include <random>
39#include <stdexcept>
40#include <utility>
41#include <vector>
42
43template<class Tvec, template<class> class SPHKernel>
45 auto tmp = solver.evolve_once_time_expl(t_curr, dt_input);
46 solver.print_timestep_logs();
47 return tmp;
48}
49
50template<class Tvec, template<class> class SPHKernel>
52 return solver.evolve_once();
53}
54
55template<class Tvec, template<class> class SPHKernel>
57
58 if (solver.solver_config.scheduler_conf.split_load_value == 0) {
60 "Scheduler load value should be greater than 0");
61 }
62
63 solver.init_required_fields();
64 ctx.init_sched(
65 solver.solver_config.scheduler_conf.split_load_value,
66 solver.solver_config.scheduler_conf.merge_load_value);
67
68 using namespace shamrock::patch;
69
70 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
71
72 sched.add_root_patch();
73
74 shamlog_debug_ln("Sys", "build local scheduler tables");
77 sched.update_local_load_value([&](shamrock::patch::Patch p) {
78 return sched.patch_data.owned_data.get(p.id_patch).get_obj_cnt();
79 });
80 solver.init_ghost_layout();
81
82 solver.init_solver_graph();
83}
84
85template<class Tvec, template<class> class SPHKernel>
87 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
88 return shamalgs::collective::allreduce_sum(sched.get_rank_count());
89}
90
91template<class Tvec, template<class> class SPHKernel>
93 return totmass / get_total_part_count();
94}
95
96template<class Tvec, template<class> class SPHKernel>
98 StackEntry stack_loc{};
99
100 using namespace shamrock::patch;
101
104
105 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
106
107 sched.for_each_patchdata_nonempty([&](const Patch, PatchDataLayer &pdat) {
108 auto acc = pdat.get_field<Tvec>(0).get_buf().copy_to_stdvec();
109
110 u32 cnt = pdat.get_obj_cnt();
111
112 for (u32 i = 0; i < cnt; i++) {
113 Tvec tmp = acc[i];
114 Tvec dr = tmp - pos;
115 Tscal dist2 = sycl::dot(dr, dr);
116 if (dist2 < best_dist2) {
117 best_dr = dr;
118 best_dist2 = dist2;
119 }
120 }
121 });
122
123 std::vector<Tvec> list_dr{};
124 shamalgs::collective::vector_allgatherv(std::vector<Tvec>{best_dr}, list_dr, MPI_COMM_WORLD);
125
126 // reset distances because if two rank find the same distance the return value won't be the same
127 // this bug took me a whole day to fix, aaaaaaaaaaaaah !!!!!
128 // maybe this should be moved somewhere else to prevent similar issues
129 // TODO (in a year maybe XD )
132
133 for (Tvec tmp : list_dr) {
134 Tvec dr = tmp - pos;
135 Tscal dist2 = sycl::dot(dr, dr);
136 if (dist2 < best_dist2) {
137 best_dr = dr;
138 best_dist2 = dist2;
139 }
140 }
141
142 return pos + best_dr;
143}
144
145template<class Tvec, template<class> class SPHKernel>
146auto shammodels::sph::Model<Tvec, SPHKernel>::get_ideal_fcc_box(Tscal dr, std::pair<Tvec, Tvec> box)
147 -> std::pair<Tvec, Tvec> {
148 StackEntry stack_loc{};
149 auto [a, b] = generic::setup::generators::get_ideal_fcc_box<Tscal>(dr, box);
150 return {a, b};
151}
152
153template<class Tvec, template<class> class SPHKernel>
154void shammodels::sph::Model<Tvec, SPHKernel>::remap_positions(std::function<Tvec(Tvec)> map) {
155 StackEntry stack_loc{};
156
157 using namespace shamrock::patch;
158
159 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
160 sched.for_each_patchdata_nonempty([&](const Patch, PatchDataLayer &pdat) {
161 auto &xyz = pdat.get_field<Tvec>(0).get_buf();
162 auto acc = xyz.copy_to_stdvec();
163
164 u32 cnt = pdat.get_obj_cnt();
165
166 for (u32 i = 0; i < cnt; i++) {
167 acc[i] = map(acc[i]);
168 }
169
170 xyz.copy_from_stdvec(acc);
171 });
172
173 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
174 .update_load_balancing();
175 sched.scheduler_step(false, false);
176
177 {
178 StackEntry stack_loc{};
180 sched.patch_tree, sched.get_sim_box().get_patch_transform<Tvec>());
182 sptree.attach_buf();
183 reatrib.reatribute_patch_objects(sptree, "xyz");
184 sched.check_patchdata_locality_correctness();
185 }
186
187 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
188 .update_load_balancing();
189 sched.scheduler_step(true, true);
190
191 {
192 StackEntry stack_loc{};
194 sched.patch_tree, sched.get_sim_box().get_patch_transform<Tvec>());
195
197 sptree.attach_buf();
198 reatrib.reatribute_patch_objects(sptree, "xyz");
199 sched.check_patchdata_locality_correctness();
200 }
201}
202
203template<class Tvec>
204inline void post_insert_data(PatchScheduler &sched) {
205 StackEntry stack_loc{};
206
207 // logger::raw_ln(sched.dump_status());
208 sched.scheduler_step(false, false);
209
210 /*
211 if(shamcomm::world_rank() == 7){
212 logger::raw_ln(sched.dump_status());
213 }
214 */
215
216 auto [m, M] = sched.get_box_tranform<Tvec>();
217
218 {
219 StackEntry stack_loc{};
221 sched.patch_tree, sched.get_sim_box().get_patch_transform<Tvec>());
223 sptree.attach_buf();
224 reatrib.reatribute_patch_objects(sptree, "xyz");
225 sched.check_patchdata_locality_correctness();
226 }
227
228 sched.scheduler_step(true, true);
229
230 {
231 StackEntry stack_loc{};
233 sched.patch_tree, sched.get_sim_box().get_patch_transform<Tvec>());
234
236 sptree.attach_buf();
237 reatrib.reatribute_patch_objects(sptree, "xyz");
238 sched.check_patchdata_locality_correctness();
239 }
240
241 std::string log = "";
242
243 using namespace shamrock::patch;
244
245 u32 smallest_count = u32_max;
246 u32 largest_count = 0;
247
248 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
249 u32 tmp = pdat.get_obj_cnt();
250 smallest_count = sham::min(tmp, smallest_count);
251 largest_count = sham::max(tmp, largest_count);
252 });
253
254 smallest_count = shamalgs::collective::allreduce_min(smallest_count);
255 largest_count = shamalgs::collective::allreduce_max(largest_count);
256
257 if (shamcomm::world_rank() == 0) {
258 logger::info_ln(
259 "Model", "current particle counts : min = ", smallest_count, "max = ", largest_count);
260 }
261
262 // sched.for_each_local_patchdata([&](const Patch p, PatchData &pdat) {
263 // log += shambase::format(
264 // "\n patch id={}, N={} particles", p.id_patch, pdat.get_obj_cnt());
265 // });
266 //
267 // std::string log_gathered = "";
268 // shamalgs::collective::gather_str(log, log_gathered);
269 //
270 // if (shamcomm::world_rank() == 0)
271 // logger::info_ln("Model", "current particle counts : ", log_gathered);
272}
273
274template<class Tvec, template<class> class SPHKernel>
276 std::vector<Tvec> &part_pos_insert,
277 std::vector<Tscal> &part_hpart_insert,
278 std::vector<Tscal> &part_u_insert) {
279 StackEntry stack_loc{};
280
281 using namespace shamrock::patch;
282
283 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
284
285 std::string log = "";
286
287 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
288 PatchCoordTransform<Tvec> ptransf = sched.get_sim_box().get_patch_transform<Tvec>();
289
290 shammath::CoordRange<Tvec> patch_coord = ptransf.to_obj_coord(p);
291
292 std::vector<Tvec> vec_acc;
293 std::vector<Tscal> hpart_acc;
294 std::vector<Tscal> u_acc;
295 for (u32 i = 0; i < part_pos_insert.size(); i++) {
296 Tvec r = part_pos_insert[i];
297 Tscal u = part_u_insert[i];
298 if (patch_coord.contain_pos(r)) {
299 vec_acc.push_back(r);
300 hpart_acc.push_back(part_hpart_insert[i]);
301 u_acc.push_back(u);
302 }
303 }
304
305 if (vec_acc.size() == 0) {
306 return;
307 }
308
309 log += shambase::format(
310 "\n rank = {} patch id={}, add N={} particles, coords = {} {}",
312 p.id_patch,
313 vec_acc.size(),
314 patch_coord.lower,
315 patch_coord.upper);
316
317 PatchDataLayer tmp(sched.get_layout_ptr_old());
318 tmp.resize(vec_acc.size());
319 tmp.fields_raz();
320
321 {
322 u32 len = vec_acc.size();
324 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("xyz"));
325 sycl::buffer<Tvec> buf(vec_acc.data(), len);
326 f.override(buf, len);
327 }
328
329 {
330 u32 len = vec_acc.size();
332 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("hpart"));
333 sycl::buffer<Tscal> buf(hpart_acc.data(), len);
334 f.override(buf, len);
335 }
336
337 {
338 u32 len = u_acc.size();
340 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("uint"));
341 sycl::buffer<Tscal> buf(u_acc.data(), len);
342 f.override(buf, len);
343 }
344
345 pdat.insert_elements(tmp);
346
347 sched.check_patchdata_locality_correctness();
348
349 std::string log_gathered = "";
350 shamalgs::collective::gather_str(log, log_gathered);
351
352 if (shamcomm::world_rank() == 0) {
353 logger::info_ln("Model", "Push particles : ", log_gathered);
354 }
355 log = "";
356
357 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
358 .update_load_balancing();
359
360 post_insert_data<Tvec>(sched);
361 });
362}
363
364template<class Tvec, template<class> class SPHKernel>
366 std::vector<Tvec> &part_pos_insert,
367 std::vector<Tscal> &part_hpart_insert,
368 std::vector<Tscal> &part_u_insert,
369 std::vector<Tvec> &part_B_on_rho_insert,
370 std::vector<Tscal> &part_psi_on_ch_insert) {
371 StackEntry stack_loc{};
372
373 using namespace shamrock::patch;
374
375 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
376
377 std::string log = "";
378
379 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
380 PatchCoordTransform<Tvec> ptransf = sched.get_sim_box().get_patch_transform<Tvec>();
381
382 shammath::CoordRange<Tvec> patch_coord = ptransf.to_obj_coord(p);
383
384 std::vector<Tvec> vec_acc;
385 std::vector<Tscal> hpart_acc;
386 std::vector<Tscal> u_acc;
387 std::vector<Tvec> B_on_rho_acc;
388 std::vector<Tscal> psi_on_ch_acc;
389 for (u32 i = 0; i < part_pos_insert.size(); i++) {
390 Tvec r = part_pos_insert[i];
391 Tscal u = part_u_insert[i];
392 if (patch_coord.contain_pos(r)) {
393 vec_acc.push_back(r);
394 hpart_acc.push_back(part_hpart_insert[i]);
395 u_acc.push_back(u);
396 B_on_rho_acc.push_back(part_B_on_rho_insert[i]);
397 psi_on_ch_acc.push_back(part_psi_on_ch_insert[i]);
398 }
399 }
400
401 if (vec_acc.size() == 0) {
402 return;
403 }
404
405 log += shambase::format(
406 "\n rank = {} patch id={}, add N={} particles, coords = {} {}",
408 p.id_patch,
409 vec_acc.size(),
410 patch_coord.lower,
411 patch_coord.upper);
412
413 PatchDataLayer tmp(sched.get_layout_ptr_old());
414 tmp.resize(vec_acc.size());
415 tmp.fields_raz();
416
417 {
418 u32 len = vec_acc.size();
420 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("xyz"));
421 sycl::buffer<Tvec> buf(vec_acc.data(), len);
422 f.override(buf, len);
423 }
424
425 {
426 u32 len = vec_acc.size();
428 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("hpart"));
429 sycl::buffer<Tscal> buf(hpart_acc.data(), len);
430 f.override(buf, len);
431 }
432
433 {
434 u32 len = u_acc.size();
436 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("uint"));
437 sycl::buffer<Tscal> buf(u_acc.data(), len);
438 f.override(buf, len);
439 }
440
441 {
442 u32 len = vec_acc.size();
444 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("B/rho"));
445 sycl::buffer<Tvec> buf(B_on_rho_acc.data(), len);
446 f.override(buf, len);
447 }
448
449 {
450 u32 len = vec_acc.size();
452 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("psi/ch"));
453 sycl::buffer<Tscal> buf(psi_on_ch_acc.data(), len);
454 f.override(buf, len);
455 }
456
457 pdat.insert_elements(tmp);
458
459 sched.check_patchdata_locality_correctness();
460
461 std::string log_gathered = "";
462 shamalgs::collective::gather_str(log, log_gathered);
463
464 if (shamcomm::world_rank() == 0) {
465 logger::info_ln("Model", "Push particles MHD : ", log_gathered);
466 }
467 log = "";
468
469 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
470 .update_load_balancing();
471
472 post_insert_data<Tvec>(sched);
473 });
474}
475
476template<class Tvec, template<class> class SPHKernel>
478 Tscal dr, std::pair<Tvec, Tvec> _box) -> std::pair<Tvec, Tvec> {
479 StackEntry stack_loc{};
480
481 using Lattice = shammath::LatticeHCP<Tvec>;
482 using LatticeIter = typename shammath::LatticeHCP<Tvec>::Iterator;
483
485 auto [idxs_min, idxs_max] = Lattice::get_box_index_bounds(dr, box.lower, box.upper);
486
487 auto [idxs_min_per, idxs_max_per] = Lattice::nearest_periodic_box_indices(idxs_min, idxs_max);
488
489 shammath::CoordRange<Tvec> ret = Lattice::get_periodic_box(dr, idxs_min_per, idxs_max_per);
490
491 return {ret.lower, ret.upper};
492}
493
494template<class Tvec, template<class> class SPHKernel>
496 Tscal dr, std::pair<Tvec, Tvec> _box) {
497 shambase::Timer time_setup;
498 time_setup.start();
499
500 StackEntry stack_loc{};
501
503
504 using namespace shamrock::patch;
505
506 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
507
508 using Lattice = shammath::LatticeHCP<Tvec>;
509 using LatticeIter = typename shammath::LatticeHCP<Tvec>::IteratorDiscontinuous;
510
511 auto [idxs_min, idxs_max] = Lattice::get_box_index_bounds(dr, box.lower, box.upper);
512
513 LatticeIter gen = LatticeIter(dr, idxs_min, idxs_max);
514
515 u64 acc_count = 0;
516
517 std::string log = "";
518 while (!gen.is_done()) {
519
520 // loc maximum count of insert part
521 u64 loc_sum_ins_cnt = 0;
522 // sum_node( loc_sum_ins_cnt )
523 u64 max_loc_sum_ins_cnt = 0;
524
525 do {
526 std::vector<Tvec> to_ins = gen.next_n(sched.crit_patch_split * 2);
527 acc_count += to_ins.size();
528
529 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
530 PatchCoordTransform<Tvec> ptransf = sched.get_sim_box().get_patch_transform<Tvec>();
531
532 shammath::CoordRange<Tvec> patch_coord = ptransf.to_obj_coord(p);
533
534 std::vector<Tvec> vec_acc;
535 for (Tvec r : to_ins) {
536 if (patch_coord.contain_pos(r)) {
537 vec_acc.push_back(r);
538 }
539 }
540
541 // update max insert_count
542 loc_sum_ins_cnt += vec_acc.size();
543
544 if (vec_acc.size() == 0) {
545 return;
546 }
547
548 log += shambase::format(
549 "\n rank = {} patch id={}, add N={} particles, coords = {} {}",
551 p.id_patch,
552 vec_acc.size(),
553 patch_coord.lower,
554 patch_coord.upper);
555
556 // reserve space to avoid allocating during copy
557 pdat.reserve(vec_acc.size());
558
559 PatchDataLayer tmp(sched.get_layout_ptr_old());
560 tmp.resize(vec_acc.size());
561 tmp.fields_raz();
562
563 {
564 u32 len = vec_acc.size();
566 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("xyz"));
567 // sycl::buffer<Tvec> buf(vec_acc.data(), len);
568 f.override(vec_acc, len);
569 }
570
571 {
573 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("hpart"));
574 f.override(dr);
575 }
576
577 pdat.insert_elements(tmp);
578 });
579
580 max_loc_sum_ins_cnt = shamalgs::collective::allreduce_max(loc_sum_ins_cnt);
581
582 if (shamcomm::world_rank() == 0) {
583 logger::info_ln(
584 "Model",
585 "--> insertion loop : max loc insert count = ",
586 max_loc_sum_ins_cnt,
587 "sum =",
588 acc_count);
589 }
590 } while (!gen.is_done() && max_loc_sum_ins_cnt < sched.crit_patch_split * 8);
591
592 sched.check_patchdata_locality_correctness();
593
594 // if(logger::details::loglevel >= shamcomm::logs::log_info){
595 // std::string log_gathered = "";
596 // shamalgs::collective::gather_str(log, log_gathered);
597 //
598 // if (shamcomm::world_rank() == 0) {
599 // shamlog_debug_ln("Model", "Push particles : ", log_gathered);
600 // }
601 // }
602 log = "";
603
604 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
605 .update_load_balancing();
606 post_insert_data<Tvec>(sched);
607 }
608
609 if (true) {
610 modules::ParticleReordering<Tvec, u32, SPHKernel>(ctx, solver.solver_config, solver.storage)
611 .reorder_particles();
612 }
613
614 time_setup.end();
615 if (shamcomm::world_rank() == 0) {
616 logger::info_ln("Model", "add_cube_hcp took :", time_setup.elasped_sec(), "s");
617 }
618}
619
620template<class Tvec, template<class> class SPHKernel>
622 Tscal dr, std::pair<Tvec, Tvec> _box) {
623 shambase::Timer time_setup;
624 time_setup.start();
625 StackEntry stack_loc{};
626
628 using namespace shamrock::patch;
629
630 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
631
632 using Lattice = shammath::LatticeHCP<Tvec>;
633 using LatticeIter = typename shammath::LatticeHCP<Tvec>::IteratorDiscontinuous;
634
635 auto [idxs_min, idxs_max] = Lattice::get_box_index_bounds(dr, box.lower, box.upper);
636
637 LatticeIter gen = LatticeIter(dr, idxs_min, idxs_max);
638
639 shamrock::DataInserterUtility inserter(sched);
640
641 auto push_current_data = [&](std::vector<Tvec> pos_data) {
642 PatchDataLayer tmp(sched.get_layout_ptr_old());
643 tmp.resize(pos_data.size());
644 tmp.fields_raz();
645
646 {
647 u32 len = pos_data.size();
649 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("xyz"));
650 // sycl::buffer<Tvec> buf(pos_data.data(), len);
651 f.override(pos_data, len);
652 }
653
654 {
656 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("hpart"));
657 f.override(dr);
658 }
659
660 inserter.push_patch_data<Tvec>(tmp, "xyz", sched.crit_patch_split * 8, [&]() {
661 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(
662 ctx, solver.solver_config, solver.storage)
663 .update_load_balancing();
664 });
665 pos_data.clear();
666 };
667
668 u32 insert_step = sched.crit_patch_split * 8;
669
670 auto [bmin, bmax] = sched.patch_data.sim_box.get_bounding_box<Tvec>();
671
672 auto has_pdat = [&]() {
673 bool ret = false;
674 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
675 ret = true;
676 });
677 return ret;
678 };
679
680 // Every MPI rank should be synchroneous on gen state
681 while (!gen.is_done()) {
682
683 u64 loc_gen_count = (has_pdat()) ? insert_step : 0;
684
685 auto gen_info = shamalgs::collective::fetch_view(loc_gen_count);
686
687 u64 skip_start = gen_info.head_offset;
688 u64 gen_cnt = loc_gen_count;
689 u64 skip_end = gen_info.total_byte_count - loc_gen_count - gen_info.head_offset;
690
691 shamlog_debug_ln(
692 "Gen",
693 "generate : ",
694 skip_start,
695 gen_cnt,
696 skip_end,
697 "total",
698 skip_start + gen_cnt + skip_end);
699 gen.skip(skip_start);
700 auto tmp = gen.next_n(gen_cnt);
701 gen.skip(skip_end);
702
703 std::vector<Tvec> pos_data;
704 for (Tvec r : tmp) {
705 if (Patch::is_in_patch_converted(r, bmin, bmax)) {
706 pos_data.push_back(r);
707 }
708 }
709
710 push_current_data(pos_data);
711
712 shamlog_debug_ln("Gen", "gen.is_done()", gen.is_done());
713 }
714
715 time_setup.end();
716 if (shamcomm::world_rank() == 0) {
717 logger::info_ln("Model", "add_cube_hcp took :", time_setup.elasped_sec(), "s");
718 }
719}
720
721template<class Tvec>
723 public:
724 using Tscal = shambase::VecComponent<Tvec>;
726
728 bool done = false;
729 Tvec center;
730 Tscal central_mass;
731 u64 Npart;
732 Tscal r_in;
733 Tscal r_out;
734 Tscal disc_mass;
735 Tscal p;
736 Tscal H_r_in;
737 Tscal q;
738 Tscal G;
739
740 u64 current_index;
741
742 std::mt19937 eng;
743
744 std::function<Tscal(Tscal)> sigma_profile;
745 std::function<Tscal(Tscal)> cs_profile;
746 std::function<Tscal(Tscal)> rot_profile;
747 std::function<Tscal(Tscal)> vel_full_corr;
748
749 public:
751 Tvec center,
752 Tscal central_mass,
753 u64 Npart,
754 Tscal r_in,
755 Tscal r_out,
756 Tscal disc_mass,
757 Tscal p,
758 Tscal H_r_in,
759 Tscal q,
760 Tscal G,
761 std::mt19937 eng,
762 std::function<Tscal(Tscal)> sigma_profile,
763 std::function<Tscal(Tscal)> cs_profile,
764 std::function<Tscal(Tscal)> rot_profile)
765 : current_index(0), Npart(Npart), center(center), central_mass(central_mass),
766 r_in(r_in), r_out(r_out), disc_mass(disc_mass), p(p), H_r_in(H_r_in), q(q), G(G),
767 eng(eng), sigma_profile(sigma_profile), cs_profile(cs_profile),
768 rot_profile(rot_profile) {
769
770 if (Npart == 0) {
771 done = true;
772 }
773 }
774
775 inline bool is_done() { return done; } // just to make sure the result is not tempered with
776
778
779 constexpr Tscal _2pi = 2 * shambase::constants::pi<Tscal>;
780
781 auto f_func = [&](Tscal r) {
782 return r * sigma_profile(r);
783 };
784
785 Tscal fmax = f_func(r_out);
786
787 auto find_r = [&]() {
788 while (true) {
789 Tscal u2 = shamalgs::primitives::mock_value<Tscal>(eng, 0, fmax);
790 Tscal r = shamalgs::primitives::mock_value<Tscal>(eng, r_in, r_out);
791 if (u2 < f_func(r)) {
792 return r;
793 }
794 }
795 };
796
797 auto theta = shamalgs::primitives::mock_value<Tscal>(eng, 0, _2pi);
798 auto Gauss = shamalgs::random::mock_gaussian<Tscal>(eng);
799 Tscal aspin = 2.;
800
801 Tscal r = find_r();
802
803 Tscal vk = rot_profile(r);
804 Tscal cs = cs_profile(r);
805 Tscal sigma = sigma_profile(r);
806
807 Tscal Omega_Kep = sycl::sqrt(G * central_mass / (r * r * r));
808
809 // Tscal H_r = cs/vk;
810 // Tscal H = H_r * r;
811 Tscal H = sycl::sqrt(2.) * 3. * cs
812 / Omega_Kep; // factor taken from phantom, to fasten thermalizing
813
814 Tscal z = H * Gauss;
815
816 auto pos = sycl::vec<Tscal, 3>{r * sycl::cos(theta), z, r * sycl::sin(theta)};
817
818 auto etheta = sycl::vec<Tscal, 3>{-pos.z(), 0, pos.x()};
819 etheta /= sycl::length(etheta);
820
821 auto vel = vk * etheta;
822
823 // Tscal rho = (sigma / (H * shambase::constants::pi2_sqrt<Tscal>))*
824 // sycl::exp(- z*z / (2*H*H));
825
826 Tscal fs = 1. - sycl::sqrt(r_in / r);
827 Tscal rho = (sigma * fs) * sycl::exp(-z * z / (2 * H * H));
828
829 Out out{pos, vel, cs, rho};
830
831 // increase counter + check if finished
832 current_index++;
833 if (current_index == Npart) {
834 done = true;
835 }
836
837 return out;
838 }
839
840 inline std::vector<Out> next_n(u32 nmax) {
841 std::vector<Out> ret{};
842 for (u32 i = 0; i < nmax; i++) {
843 if (done) {
844 break;
845 }
846
847 ret.push_back(next());
848 }
849 return ret;
850 }
851 };
852};
853
854template<class Tvec, template<class> class SPHKernel>
856 Tvec center,
857 Tscal central_mass,
858 u32 Npart,
859 Tscal r_in,
860 Tscal r_out,
861 Tscal disc_mass,
862 Tscal p,
863 Tscal H_r_in,
864 Tscal q,
865 std::mt19937 eng) {
866
867 Tscal eos_gamma;
868 using Config = SolverConfig;
869 using SolverConfigEOS = typename Config::EOSConfig;
870 using SolverEOS_Adiabatic = typename SolverConfigEOS::Adiabatic;
871 if (SolverEOS_Adiabatic *eos_config
872 = std::get_if<SolverEOS_Adiabatic>(&solver.solver_config.eos_config.config)) {
873
874 eos_gamma = eos_config->gamma;
875
876 } else {
877 // dirty hack for disc setup in locally isothermal
878 eos_gamma = 2;
879 // shambase::throw_unimplemented();
880 }
881
882 auto sigma_profile = [=](Tscal r) {
883 // we setup with an adimensional mass since it is monte carlo
884 constexpr Tscal sigma_0 = 1;
885 return sigma_0 * sycl::pow(r / r_in, -p);
886 };
887
888 auto cs_law = [=](Tscal r) {
889 return sycl::pow(r / r_in, -q);
890 };
891
892 auto kep_profile = [&](Tscal r) {
893 Tscal G = solver.solver_config.get_constant_G();
894 return sycl::sqrt(G * central_mass / r);
895 };
896
897 auto rot_profile = [&](Tscal r) -> Tscal {
898 // carefull: needs r in cylindrical
899 Tscal G = solver.solver_config.get_constant_G();
900 Tscal c = solver.solver_config.get_constant_c();
901 Tscal aspin = 2.;
902 Tscal term = G * central_mass / r;
903 Tscal term_fs = 1. - sycl::sqrt(r_in / r);
904 Tscal term_pr
905 = -sycl::pown(cs_law(r), 2) * (1.5 + p + q); // NO CORRECTION from fs term, bad response
906 Tscal term_bh = 0.; //- (2. * aspin / sycl::pow(c, 3)) * sycl::pow(G * central_mass / r, 2);
907 Tscal det = sycl::pown(term_bh, 2) + 4. * (term + term_pr);
908 Tscal Rg = G * central_mass / sycl::pown(c, 2);
909 Tscal vkep = sqrt(G * central_mass / r);
910
911 Tscal vphi = 0.5 * (term_bh + sycl::sqrt(det));
912
913 return vphi;
914 };
915
916 auto cs_profile = [&](Tscal r) {
917 Tscal cs_in = (H_r_in * r_in / r) * kep_profile(r_in); // H_r_in*rot_profile(r_in);
918 return cs_law(r) * cs_in;
919 };
920
921 auto get_hfact = []() -> Tscal {
922 return Kernel::hfactd;
923 };
924
925 auto int_rho_h = [&](Tscal h) -> Tscal {
926 return shamrock::sph::rho_h(solver.solver_config.gpart_mass, h, Kernel::hfactd);
927 };
928
929 Tscal part_mass = disc_mass / Npart;
930
931 shambase::Timer time_setup;
932 time_setup.start();
933
934 StackEntry stack_loc{};
935
936 using namespace shamrock::patch;
937
938 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
939
941 using DIter = typename BigDiscUtils<Tvec>::DiscIterator;
942
943 Tscal G = solver.solver_config.get_constant_G();
944 DIter gen = DIter(
945 center,
946 central_mass,
947 Npart,
948 r_in,
949 r_out,
950 disc_mass,
951 p,
952 H_r_in,
953 q,
954 G,
955 eng,
956 sigma_profile,
957 cs_profile,
958 rot_profile);
959
960 u64 acc_count = 0;
961
962 std::string log = "";
963 while (!gen.is_done()) {
964
965 // loc maximum count of insert part
966 u64 loc_sum_ins_cnt = 0;
967 // sum_node( loc_sum_ins_cnt )
968 u64 max_loc_sum_ins_cnt = 0;
969
970 do {
971 std::vector<Out> to_ins = gen.next_n(sched.crit_patch_split * 2);
972 acc_count += to_ins.size();
973
974 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
975 PatchCoordTransform<Tvec> ptransf = sched.get_sim_box().get_patch_transform<Tvec>();
976
977 shammath::CoordRange<Tvec> patch_coord = ptransf.to_obj_coord(p);
978
979 std::vector<Out> part_list;
980 for (Out r : to_ins) {
981 if (patch_coord.contain_pos(r.pos)) {
982 // add all part to insert in a vector
983 part_list.push_back(r);
984 }
985 }
986
987 // update max insert_count
988 loc_sum_ins_cnt += part_list.size();
989
990 if (part_list.size() == 0) {
991 return;
992 }
993
994 log += shambase::format(
995 "\n rank = {} patch id={}, add N={} particles, coords = {} {}",
997 p.id_patch,
998 part_list.size(),
999 patch_coord.lower,
1000 patch_coord.upper);
1001
1002 // extract the pos from part_list
1003 std::vector<Tvec> vec_pos;
1004 std::vector<Tvec> vec_vel;
1005 std::vector<Tscal> vec_u;
1006 std::vector<Tscal> vec_h;
1007 std::vector<Tscal> vec_cs;
1008
1009 for (Out o : part_list) {
1010 vec_pos.push_back(o.pos);
1011 vec_vel.push_back(o.velocity);
1012 vec_u.push_back(o.cs * o.cs / (/*solver.eos_gamma * */ (eos_gamma - 1)));
1013 vec_h.push_back(shamrock::sph::h_rho(part_mass, o.rho * 0.1, Kernel::hfactd));
1014 vec_cs.push_back(o.cs);
1015 }
1016
1017 // reserve space to avoid allocating during copy
1018 pdat.reserve(vec_pos.size());
1019
1020 PatchDataLayer tmp(sched.get_layout_ptr_old());
1021 tmp.resize(vec_pos.size());
1022 tmp.fields_raz();
1023
1024 {
1025 u32 len = vec_pos.size();
1027 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("xyz"));
1028 sycl::buffer<Tvec> buf(vec_pos.data(), len);
1029 f.override(buf, len);
1030 }
1031
1032 {
1033 u32 len = vec_pos.size();
1035 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("hpart"));
1036 sycl::buffer<Tscal> buf(vec_h.data(), len);
1037 f.override(buf, len);
1038 }
1039
1040 {
1041 u32 len = vec_pos.size();
1043 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("uint"));
1044 sycl::buffer<Tscal> buf(vec_u.data(), len);
1045 f.override(buf, len);
1046 }
1047
1048 if (solver.solver_config.is_eos_locally_isothermal()) {
1049 u32 len = vec_pos.size();
1051 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("soundspeed"));
1052 sycl::buffer<Tscal> buf(vec_cs.data(), len);
1053 f.override(buf, len);
1054 }
1055
1056 {
1057 u32 len = vec_pos.size();
1059 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("vxyz"));
1060 sycl::buffer<Tvec> buf(vec_vel.data(), len);
1061 f.override(buf, len);
1062 }
1063
1064 pdat.insert_elements(tmp);
1065 });
1066
1067 max_loc_sum_ins_cnt = shamalgs::collective::allreduce_max(loc_sum_ins_cnt);
1068
1069 if (shamcomm::world_rank() == 0) {
1070 logger::info_ln(
1071 "Model",
1072 "--> insertion loop : max loc insert count = ",
1073 max_loc_sum_ins_cnt,
1074 "sum =",
1075 acc_count);
1076 }
1077 } while (!gen.is_done() && max_loc_sum_ins_cnt < sched.crit_patch_split * 8);
1078
1079 sched.check_patchdata_locality_correctness();
1080
1081 // if(logger::details::loglevel >= shamcomm::logs::log_info){
1082 // std::string log_gathered = "";
1083 // shamalgs::collective::gather_str(log, log_gathered);
1084 //
1085 // if (shamcomm::world_rank() == 0) {
1086 // shamlog_debug_ln("Model", "Push particles : ", log_gathered);
1087 // }
1088 // }
1089 log = "";
1090
1091 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
1092 .update_load_balancing();
1093 post_insert_data<Tvec>(sched);
1094 }
1095
1096 if (true) {
1097 modules::ParticleReordering<Tvec, u32, SPHKernel>(ctx, solver.solver_config, solver.storage)
1098 .reorder_particles();
1099 }
1100
1101 time_setup.end();
1102 if (shamcomm::world_rank() == 0) {
1103 logger::info_ln("Model", "add_big_disc took :", time_setup.elasped_sec(), "s");
1104 }
1105}
1106
1107template<class Tvec, template<class> class SPHKernel>
1109 Tscal dr, std::pair<Tvec, Tvec> _box) {
1110 StackEntry stack_loc{};
1111
1112 shammath::CoordRange<Tvec> box = _box;
1113
1114 using namespace shamrock::patch;
1115
1116 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
1117
1118 std::string log = "";
1119
1120 auto make_sliced = [&]() {
1121 std::vector<Tvec> vec_lst;
1122 generic::setup::generators::add_particles_fcc(
1123 dr,
1124 {box.lower, box.upper},
1125 [&](Tvec r) {
1126 return box.contain_pos(r);
1127 },
1128 [&](Tvec r, Tscal h) {
1129 vec_lst.push_back(r);
1130 });
1131
1132 std::vector<std::vector<Tvec>> sliced_buf;
1133
1134 u32 sz_buf = sched.crit_patch_split * 4;
1135
1136 std::vector<Tvec> cur_buf;
1137 for (u32 i = 0; i < vec_lst.size(); i++) {
1138 cur_buf.push_back(vec_lst[i]);
1139
1140 if (cur_buf.size() > sz_buf) {
1141 sliced_buf.push_back(std::exchange(cur_buf, std::vector<Tvec>{}));
1142 }
1143 }
1144
1145 if (cur_buf.size() > 0) {
1146 sliced_buf.push_back(std::exchange(cur_buf, std::vector<Tvec>{}));
1147 }
1148
1149 return sliced_buf;
1150 };
1151
1152 std::vector<std::vector<Tvec>> sliced_buf = make_sliced();
1153
1154 for (std::vector<Tvec> to_ins : sliced_buf) {
1155
1156 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
1157 PatchCoordTransform<Tvec> ptransf = sched.get_sim_box().get_patch_transform<Tvec>();
1158
1159 shammath::CoordRange<Tvec> patch_coord = ptransf.to_obj_coord(p);
1160
1161 std::vector<Tvec> vec_acc;
1162 for (Tvec r : to_ins) {
1163 if (patch_coord.contain_pos(r)) {
1164 vec_acc.push_back(r);
1165 }
1166 }
1167
1168 if (vec_acc.size() == 0) {
1169 return;
1170 }
1171
1172 log += shambase::format(
1173 "\n rank = {} patch id={}, add N={} particles, coords = {} {}",
1175 p.id_patch,
1176 vec_acc.size(),
1177 patch_coord.lower,
1178 patch_coord.upper);
1179
1180 PatchDataLayer tmp(sched.get_layout_ptr_old());
1181 tmp.resize(vec_acc.size());
1182 tmp.fields_raz();
1183
1184 {
1185 u32 len = vec_acc.size();
1187 = tmp.get_field<Tvec>(sched.pdl_old().get_field_idx<Tvec>("xyz"));
1188 sycl::buffer<Tvec> buf(vec_acc.data(), len);
1189 f.override(buf, len);
1190 }
1191
1192 {
1194 = tmp.get_field<Tscal>(sched.pdl_old().get_field_idx<Tscal>("hpart"));
1195 f.override(dr);
1196 }
1197
1198 pdat.insert_elements(tmp);
1199 });
1200
1201 sched.check_patchdata_locality_correctness();
1202
1203 std::string log_gathered = "";
1204 shamalgs::collective::gather_str(log, log_gathered);
1205
1206 if (shamcomm::world_rank() == 0) {
1207 logger::info_ln("Model", "Push particles : ", log_gathered);
1208 }
1209 log = "";
1210
1211 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
1212 .update_load_balancing();
1213 post_insert_data<Tvec>(sched);
1214 }
1215}
1216
1217template<class Tvec, template<class> class SPHKernel>
1219 PhantomDump &phdump, bool bypass_error) -> SolverConfig {
1220 StackEntry stack_loc{};
1221 SolverConfig conf{};
1222
1223 auto massoftype = phdump.read_header_floats<Tscal>("massoftype");
1224
1225 conf.gpart_mass = massoftype[0];
1226 conf.cfl_config.cfl_cour = phdump.read_header_float<Tscal>("C_cour");
1227 conf.cfl_config.cfl_force = phdump.read_header_float<Tscal>("C_force");
1228
1229 conf.eos_config = get_shamrock_eosconfig<Tvec>(phdump, bypass_error);
1230 conf.artif_viscosity = get_shamrock_avconfig<Tvec>(phdump);
1231
1232 conf.set_units(get_shamrock_units<Tscal>(phdump));
1233
1234 conf.boundary_config = get_shamrock_boundary_config<Tvec>(phdump);
1235
1236 return conf;
1237}
1238
1239template<class Tvec, template<class> class SPHKernel>
1241 PhantomDump &phdump, Tscal hpart_fact_load) {
1242 StackEntry stack_loc{};
1243
1244 bool has_coord_in_header = true;
1245
1246 Tscal xmin, xmax, ymin, ymax, zmin, zmax;
1247 has_coord_in_header = phdump.has_header_entry("xmin");
1248
1249 std::string log = "";
1250
1251 std::vector<Tvec> xyz, vxyz;
1252 std::vector<Tscal> h, u, alpha;
1253
1254 {
1255 std::vector<Tscal> x, y, z, vx, vy, vz;
1256
1257 phdump.blocks[0].fill_vec("x", x);
1258 phdump.blocks[0].fill_vec("y", y);
1259 phdump.blocks[0].fill_vec("z", z);
1260
1261 if (has_coord_in_header) {
1262 xmin = phdump.read_header_float<f64>("xmin");
1263 xmax = phdump.read_header_float<f64>("xmax");
1264 ymin = phdump.read_header_float<f64>("ymin");
1265 ymax = phdump.read_header_float<f64>("ymax");
1266 zmin = phdump.read_header_float<f64>("zmin");
1267 zmax = phdump.read_header_float<f64>("zmax");
1268
1269 resize_simulation_box({{xmin, ymin, zmin}, {xmax, ymax, zmax}});
1270 } else {
1271 Tscal box_tolerance = 1.2;
1272
1273 xmin = *std::min_element(x.begin(), x.end());
1274 xmax = *std::max_element(x.begin(), x.end());
1275 ymin = *std::min_element(y.begin(), y.end());
1276 ymax = *std::max_element(y.begin(), y.end());
1277 zmin = *std::min_element(z.begin(), z.end());
1278 zmax = *std::max_element(z.begin(), z.end());
1279
1280 Tvec bm = {xmin, ymin, zmin};
1281 Tvec bM = {xmax, ymax, zmax};
1282
1283 Tvec center = (bm + bM) * 0.5;
1284
1285 Tvec d = (bM - bm) * 0.5;
1286
1287 // expand the box
1288 d *= box_tolerance;
1289
1290 resize_simulation_box({center - d, center + d});
1291 }
1292
1293 phdump.blocks[0].fill_vec("h", h);
1294
1295 phdump.blocks[0].fill_vec("vx", vx);
1296 phdump.blocks[0].fill_vec("vy", vy);
1297 phdump.blocks[0].fill_vec("vz", vz);
1298
1299 phdump.blocks[0].fill_vec("u", u);
1300 phdump.blocks[0].fill_vec("alpha", alpha);
1301
1302 for (u32 i = 0; i < x.size(); i++) {
1303 xyz.push_back({x[i], y[i], z[i]});
1304 }
1305 for (u32 i = 0; i < vx.size(); i++) {
1306 vxyz.push_back({vx[i], vy[i], vz[i]});
1307 }
1308 }
1309
1310 // Load time infos
1311 f64 time_phdump = phdump.read_header_float<f64>("time");
1312 solver.solver_config.set_time(time_phdump);
1313
1314 using namespace shamrock::patch;
1315
1316 PatchScheduler &sched = shambase::get_check_ref(ctx.sched);
1317
1318 u32 sz_buf = sched.crit_patch_split * 4;
1319
1320 u32 Ntot = xyz.size();
1321
1322 std::vector<u64> insert_ranges;
1323 insert_ranges.push_back(0);
1324 for (u64 i = sz_buf; i < Ntot; i += sz_buf) {
1325 insert_ranges.push_back(i);
1326 }
1327 insert_ranges.push_back(Ntot);
1328
1329 for (u64 krange = 0; krange < insert_ranges.size() - 1; krange++) {
1330 u64 start_id = insert_ranges[krange];
1331 u64 end_id = insert_ranges[krange + 1];
1332
1333 u64 Nloc = end_id - start_id;
1334
1335 sched.for_each_local_patchdata([&](const Patch &p, PatchDataLayer &pdat) {
1336 PatchCoordTransform<Tvec> ptransf = sched.get_sim_box().get_patch_transform<Tvec>();
1337
1338 shammath::CoordRange<Tvec> patch_coord = ptransf.to_obj_coord(p);
1339
1340 std::vector<u64> sel_index;
1341 for (u64 i = start_id; i < end_id; i++) {
1342 Tvec r = xyz[i];
1343 Tscal h_ = h[i];
1344 if (patch_coord.contain_pos(r) && (h_ >= 0)) {
1345 sel_index.push_back(i);
1346 }
1347 }
1348
1349 if (sel_index.size() == 0) {
1350 return;
1351 }
1352
1353 log += shambase::format(
1354 "\n rank = {} patch id={}, add N={} particles, coords = {} {}",
1356 p.id_patch,
1357 sel_index.size(),
1358 patch_coord.lower,
1359 patch_coord.upper);
1360
1361 std::vector<Tvec> ins_xyz, ins_vxyz;
1362 std::vector<Tscal> ins_h, ins_u, ins_alpha;
1363 for (u64 i : sel_index) {
1364 ins_xyz.push_back(xyz[i]);
1365 }
1366 for (u64 i : sel_index) {
1367 ins_vxyz.push_back(vxyz[i]);
1368 }
1369 for (u64 i : sel_index) {
1370 ins_h.push_back(h[i] * hpart_fact_load);
1371 }
1372 if (u.size() > 0) {
1373 for (u64 i : sel_index) {
1374 ins_u.push_back(u[i]);
1375 }
1376 }
1377 if (alpha.size() > 0) {
1378 for (u64 i : sel_index) {
1379 ins_alpha.push_back(alpha[i]);
1380 }
1381 }
1382
1383 PatchDataLayer ptmp(sched.get_layout_ptr_old());
1384 ptmp.resize(sel_index.size());
1385 ptmp.fields_raz();
1386
1387 ptmp.override_patch_field("xyz", ins_xyz);
1388 ptmp.override_patch_field("vxyz", ins_vxyz);
1389 ptmp.override_patch_field("hpart", ins_h);
1390
1391 if (ins_alpha.size() > 0) {
1392 ptmp.override_patch_field("alpha_AV", ins_alpha);
1393 }
1394
1395 if (ins_u.size() > 0) {
1396 ptmp.override_patch_field("uint", ins_u);
1397 }
1398
1399 pdat.insert_elements(ptmp);
1400 });
1401
1402 sched.check_patchdata_locality_correctness();
1403
1404 std::string log_gathered = "";
1405 shamalgs::collective::gather_str(log, log_gathered);
1406
1407 if (shamcomm::world_rank() == 0) {
1408 logger::info_ln("Model", "Push particles : ", log_gathered);
1409 }
1410 log = "";
1411
1412 modules::ComputeLoadBalanceValue<Tvec, SPHKernel>(ctx, solver.solver_config, solver.storage)
1413 .update_load_balancing();
1414
1415 post_insert_data<Tvec>(sched);
1416
1417 // add sinks
1418
1419 PhantomDumpBlock &sink_block = phdump.blocks[1];
1420 {
1421 std::vector<Tscal> xsink, ysink, zsink;
1422 std::vector<Tscal> vxsink, vysink, vzsink;
1423 std::vector<Tscal> mass;
1424 std::vector<Tscal> Racc;
1425
1426 sink_block.fill_vec("x", xsink);
1427 sink_block.fill_vec("y", ysink);
1428 sink_block.fill_vec("z", zsink);
1429 sink_block.fill_vec("vx", vxsink);
1430 sink_block.fill_vec("vy", vysink);
1431 sink_block.fill_vec("vz", vzsink);
1432 sink_block.fill_vec("m", mass);
1433 sink_block.fill_vec("h", Racc);
1434
1435 for (u32 i = 0; i < xsink.size(); i++) {
1436 add_sink(
1437 mass[i],
1438 {xsink[i], ysink[i], zsink[i]},
1439 {vxsink[i], vysink[i], vzsink[i]},
1440 Racc[i]);
1441 }
1442 }
1443 }
1444}
1445
1446template<class Tvec, template<class> class SPHKernel>
1448 PhantomDumpBlock &block, shamrock::patch::PatchDataLayer &pdat) {
1449
1450 std::vector<Tvec> xyz = pdat.fetch_data<Tvec>("xyz");
1451
1452 u64 xid = block.get_ref_fort_real("x");
1453 u64 yid = block.get_ref_fort_real("y");
1454 u64 zid = block.get_ref_fort_real("z");
1455
1456 for (auto vec : xyz) {
1457 block.blocks_fort_real[xid].vals.push_back(vec.x());
1458 block.blocks_fort_real[yid].vals.push_back(vec.y());
1459 block.blocks_fort_real[zid].vals.push_back(vec.z());
1460 }
1461
1462 std::vector<Tscal> h = pdat.fetch_data<Tscal>("hpart");
1463 u64 hid = block.get_ref_f32("h");
1464 for (auto h_ : h) {
1465 block.blocks_f32[hid].vals.push_back(h_);
1466 }
1467
1468 if (solver.solver_config.has_field_alphaAV()) {
1469 std::vector<Tscal> alpha = pdat.fetch_data<Tscal>("alpha_AV");
1470 u64 aid = block.get_ref_f32("alpha");
1471 for (auto alp_ : alpha) {
1472 block.blocks_f32[aid].vals.push_back(alp_);
1473 }
1474 }
1475
1476 if (solver.solver_config.has_field_divv()) {
1477 std::vector<Tscal> vecdivv = pdat.fetch_data<Tscal>("divv");
1478 u64 divvid = block.get_ref_f32("divv");
1479 for (auto d_ : vecdivv) {
1480 block.blocks_f32[divvid].vals.push_back(d_);
1481 }
1482 }
1483
1484 std::vector<Tvec> vxyz = pdat.fetch_data<Tvec>("vxyz");
1485
1486 u64 vxid = block.get_ref_fort_real("vx");
1487 u64 vyid = block.get_ref_fort_real("vy");
1488 u64 vzid = block.get_ref_fort_real("vz");
1489
1490 for (auto vec : vxyz) {
1491 block.blocks_fort_real[vxid].vals.push_back(vec.x());
1492 block.blocks_fort_real[vyid].vals.push_back(vec.y());
1493 block.blocks_fort_real[vzid].vals.push_back(vec.z());
1494 }
1495
1496 std::vector<Tscal> u = pdat.fetch_data<Tscal>("uint");
1497 u64 uid = block.get_ref_fort_real("u");
1498 for (auto u_ : u) {
1499 block.blocks_fort_real[uid].vals.push_back(u_);
1500 }
1501
1502 block.tot_count = block.blocks_fort_real[xid].vals.size();
1503}
1504
1505template<class Tvec, template<class> class SPHKernel>
1507 StackEntry stack_loc{};
1508
1509 PhantomDump dump;
1510
1511 bool bypass_error_check = false;
1512
1513 auto get_sink_count = [&]() -> int {
1514 if (solver.storage.sinks.is_empty()) {
1515 return 0;
1516 } else {
1517 return int(solver.storage.sinks.get().size());
1518 }
1519 };
1520
1521 dump.override_magic_number();
1522 dump.iversion = 1;
1523 dump.fileid = shambase::format("{:100s}", "FT:Phantom Shamrock writer");
1524
1525 u32 Ntot = get_total_part_count();
1526 dump.table_header_fort_int.add("nparttot", Ntot);
1527 dump.table_header_fort_int.add("ntypes", 8);
1528 dump.table_header_fort_int.add("npartoftype", Ntot);
1529 dump.table_header_fort_int.add("npartoftype", 0);
1530 dump.table_header_fort_int.add("npartoftype", 0);
1531 dump.table_header_fort_int.add("npartoftype", 0);
1532 dump.table_header_fort_int.add("npartoftype", 0);
1533 dump.table_header_fort_int.add("npartoftype", 0);
1534 dump.table_header_fort_int.add("npartoftype", 0);
1535 dump.table_header_fort_int.add("npartoftype", 0);
1536
1537 dump.table_header_i64.add("nparttot", Ntot);
1538 dump.table_header_i64.add("ntypes", 8);
1539 dump.table_header_i64.add("npartoftype", Ntot);
1540 dump.table_header_i64.add("npartoftype", 0);
1541 dump.table_header_i64.add("npartoftype", 0);
1542 dump.table_header_i64.add("npartoftype", 0);
1543 dump.table_header_i64.add("npartoftype", 0);
1544 dump.table_header_i64.add("npartoftype", 0);
1545 dump.table_header_i64.add("npartoftype", 0);
1546 dump.table_header_i64.add("npartoftype", 0);
1547
1548 dump.table_header_fort_int.add("nblocks", 1);
1549 dump.table_header_fort_int.add("nptmass", get_sink_count());
1550 dump.table_header_fort_int.add("ndustlarge", 0);
1551 dump.table_header_fort_int.add("ndustsmall", 0);
1552 dump.table_header_fort_int.add("idust", 7);
1553 dump.table_header_fort_int.add("idtmax_n", 1);
1554 dump.table_header_fort_int.add("idtmax_frac", 0);
1555 dump.table_header_fort_int.add("idumpfile", 0);
1556 dump.table_header_fort_int.add("majorv", 2023);
1557 dump.table_header_fort_int.add("minorv", 0);
1558 dump.table_header_fort_int.add("microv", 0);
1559 dump.table_header_fort_int.add("isink", 0);
1560
1561 dump.table_header_i32.add("iexternalforce", 0);
1562
1563 write_shamrock_eos_in_phantom_dump(solver.solver_config.eos_config, dump, bypass_error_check);
1564
1565 dump.table_header_fort_real.add("time", solver.solver_config.get_time());
1566 dump.table_header_fort_real.add("dtmax", solver.solver_config.get_dt_sph());
1567
1568 dump.table_header_fort_real.add("rhozero", 0);
1569 dump.table_header_fort_real.add("hfact", Kernel::hfactd);
1570 dump.table_header_fort_real.add("tolh", 0.0001);
1571 dump.table_header_fort_real.add("C_cour", solver.solver_config.cfl_config.cfl_cour);
1572 dump.table_header_fort_real.add("C_force", solver.solver_config.cfl_config.cfl_force);
1573 dump.table_header_fort_real.add("alpha", 0);
1574 dump.table_header_fort_real.add("alphau", 1);
1575 dump.table_header_fort_real.add("alphaB", 1);
1576
1577 dump.table_header_fort_real.add("massoftype", solver.solver_config.gpart_mass);
1578 dump.table_header_fort_real.add("massoftype", 0);
1579 dump.table_header_fort_real.add("massoftype", 0);
1580 dump.table_header_fort_real.add("massoftype", 0);
1581 dump.table_header_fort_real.add("massoftype", 0);
1582 dump.table_header_fort_real.add("massoftype", 0);
1583 dump.table_header_fort_real.add("massoftype", 0);
1584 dump.table_header_fort_real.add("massoftype", 0);
1585
1586 dump.table_header_fort_real.add("Bextx", 0);
1587 dump.table_header_fort_real.add("Bexty", 0);
1588 dump.table_header_fort_real.add("Bextz", 0);
1589 dump.table_header_fort_real.add("dum", 0);
1590
1591 PatchScheduler &sched = shambase::get_check_ref(solver.context.sched);
1592
1593 auto box_size = sched.get_box_volume<Tvec>();
1594
1595 write_shamrock_boundaries_in_phantom_dump(
1596 solver.solver_config.boundary_config, box_size, dump, bypass_error_check);
1597
1598 dump.table_header_fort_real.add("get_conserv", -1);
1599 dump.table_header_fort_real.add("etot_in", 0.59762);
1600 dump.table_header_fort_real.add("angtot_in", 0.0189694);
1601 dump.table_header_fort_real.add("totmom_in", 0.0306284);
1602
1603 write_shamrock_units_in_phantom_dump(solver.solver_config.unit_sys, dump, bypass_error_check);
1604
1605 PhantomDumpBlock block_part;
1606
1607 {
1608 NamedStackEntry stack_loc{"gather data"};
1609 std::vector<std::unique_ptr<shamrock::patch::PatchDataLayer>> gathered
1610 = ctx.allgather_data();
1611
1612 for (auto &dat : gathered) {
1613 add_pdat_to_phantom_block(block_part, shambase::get_check_ref(dat));
1614 }
1615 }
1616
1617 dump.blocks.push_back(std::move(block_part));
1618
1619 if (!solver.storage.sinks.is_empty()) {
1620
1621 auto &sinks = solver.storage.sinks.get();
1622 // add sinks to block 1
1623 PhantomDumpBlock sink_block;
1624
1625 u64 xid = sink_block.get_ref_fort_real("x");
1626 u64 yid = sink_block.get_ref_fort_real("y");
1627 u64 zid = sink_block.get_ref_fort_real("z");
1628 u64 mid = sink_block.get_ref_fort_real("m");
1629 u64 hid = sink_block.get_ref_fort_real("h");
1630 u64 vxid = sink_block.get_ref_fort_real("vx");
1631 u64 vyid = sink_block.get_ref_fort_real("vy");
1632 u64 vzid = sink_block.get_ref_fort_real("vz");
1633
1634 for (SinkParticle<Tvec> s : sinks) {
1635 sink_block.blocks_fort_real[xid].vals.push_back(s.pos.x());
1636 sink_block.blocks_fort_real[yid].vals.push_back(s.pos.y());
1637 sink_block.blocks_fort_real[zid].vals.push_back(s.pos.z());
1638 sink_block.blocks_fort_real[mid].vals.push_back(s.mass);
1639 sink_block.blocks_fort_real[hid].vals.push_back(s.accretion_radius);
1640 sink_block.blocks_fort_real[vxid].vals.push_back(s.velocity.x());
1641 sink_block.blocks_fort_real[vyid].vals.push_back(s.velocity.y());
1642 sink_block.blocks_fort_real[vzid].vals.push_back(s.velocity.z());
1643 }
1644
1645 sink_block.tot_count = sinks.size();
1646
1647 dump.blocks.push_back(std::move(sink_block));
1648 }
1649
1650 return dump;
1651}
1652
1653using namespace shammath;
1654
1658
constexpr const char * vxyz
3-velocity field
constexpr const char * xyz
Position field (3D coordinates)
Header file describing a Node Instance.
MPI scheduler.
double f64
Alias for double.
std::uint32_t u32
32 bit unsigned integer
std::uint64_t u64
64 bit unsigned integer
The MPI scheduler.
SchedulerPatchData patch_data
handle the data of the patches of the scheduler
u64 crit_patch_split
splitting limit (if load value > crit_patch_split => patch split)
PatchTree patch_tree
handle the tree structure of the patches
void scheduler_step(bool do_split_merge, bool do_load_balancing)
scheduler step
SchedulerPatchList patch_list
handle the list of the patches of the scheduler
std::unordered_set< u64 > owned_patch_id
(owned_patch_id = patch_list.build_local())
void add_root_patch()
add patch to the scheduler
std::unordered_set< u64 > build_local()
select owned patches owned by the node to rebuild local
void build_local_idx_map()
recompute id_patch_to_local_idx
Class Timer measures the time elapsed since the timer was started.
Definition time.hpp:96
void end()
Stops the timer and stores the elapsed time in nanoseconds.
Definition time.hpp:111
f64 elasped_sec() const
Converts the stored nanosecond time to a floating point representation in seconds.
Definition time.hpp:123
void start()
Starts the timer.
Definition time.hpp:106
Iterator utility to generate the lattice.
Iterator utility to generate the lattice.
utility for generating HCP crystal lattices
Vector class based on std::array storage and mdspan.
Definition matrix.hpp:96
The shamrock SPH model.
Definition Model.hpp:55
void init()
Initialise the model and all the related data structures (patch scheduler in particular)
Definition Model.cpp:56
Class to insert data in the PatchScheduler.
Utility class used to move the objects between patches.
u32 get_field_idx(const std::string &field_name) const
Get the field id if matching name & type.
PatchDataLayer container class, the layout is described in patchdata_layout.
std::vector< T > fetch_data(std::string key)
Fetch data of a patchdata field into a std::vector.
std::tuple< T, T > get_bounding_box() const
Get the stored bounding box of the domain.
Definition SimBox.hpp:247
PatchCoordTransform< T > get_patch_transform() const
Get a PatchCoordTransform object that describes the conversion between patch coordinates and domain c...
Definition SimBox.hpp:285
shamrock::patch::SimulationBoxInfo sim_box
simulation box geometry info
shambase::DistributedData< PatchData > owned_data
map container for patchdata owned by the current node (layout : id_patch,data)
This header file contains utility functions related to exception handling in the code.
std::vector< int > vector_allgatherv(const std::vector< T > &send_vec, const MPI_Datatype &send_type, std::vector< T > &recv_vec, const MPI_Datatype &recv_type, const MPI_Comm comm)
allgatherv on vector with size query (size querying variant of vector_allgatherv_ks) //TODO add fault...
Definition exchanges.hpp:98
void gather_str(const std::string &send_vec, std::string &recv_vec)
Gathers a string from all nodes and store the result in a std::string.
void throw_with_loc(std::string message, SourceLocation loc=SourceLocation{})
Throw an exception and append the source location to it.
T & get_check_ref(const std::unique_ptr< T > &ptr, SourceLocation loc=SourceLocation())
Takes a std::unique_ptr and returns a reference to the object it holds. It throws a std::runtime_erro...
Definition memory.hpp:110
i32 world_rank()
Gives the rank of the current process in the MPI communicator.
Definition worldInfo.cpp:40
namespace for math utility
Definition AABB.hpp:26
void write_shamrock_units_in_phantom_dump(std::optional< shamunits::UnitSystem< Tscal > > &units, PhantomDump &dump, bool bypass_error)
Write shamrock units config into the phantom dump.
void write_shamrock_eos_in_phantom_dump(EOSConfig< Tvec > &cfg, PhantomDump &dump, bool bypass_error)
Write the eos config to th phantom dump header.
constexpr u32 u32_max
u32 max value
sph kernels
This file contains the definition for the stacktrace related functionality.
Class representing a Phantom dump file.
Patch object that contain generic patch information.
Definition Patch.hpp:33