Shamrock 2025.10.0
Astrophysical Code
Loading...
Searching...
No Matches
Solver.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
17#include "shambase/assert.hpp"
20#include "shambase/memory.hpp"
22#include "shambase/string.hpp"
23#include "shambase/time.hpp"
31#include "shambackends/math.hpp"
32#include "shamcomm/logs.hpp"
34#include "shamcomm/wrapper.hpp"
72#include "shamphys/mhd.hpp"
103#include <memory>
104#include <stdexcept>
105#include <vector>
106
107template<class Tvec, template<class> class Kern>
109
110 shamrock::patch::PatchDataLayerLayout &pdl = scheduler().pdl_old();
111 bool has_B_field = solver_config.has_field_B_on_rho();
112 bool has_psi_field = solver_config.has_field_psi_on_ch();
113 bool has_epsilon_field = solver_config.dust_config.has_epsilon_field();
114 bool has_deltav_field = solver_config.dust_config.has_deltav_field();
115
116 using namespace shamrock::solvergraph;
117
118 SolverGraph &solver_graph = storage.solver_graph;
119
120 solver_graph.register_edge(
121 "scheduler_patchdata", PatchDataLayerRefs("patchdatas", "\\mathbb{U}_{\\rm patch}"));
122 solver_graph.register_edge("part_counts", Indexes<u32>("Npart_patch", "N_{\\rm part}_p"));
123
124 solver_graph.register_edge("dt", IDataEdge<Tscal>("dt", "dt"));
125 solver_graph.register_edge("dt_half", IDataEdge<Tscal>("dt_half", "\\frac{dt}{2}"));
126 solver_graph.register_edge("gpart_mass", ScalarEdge<Tscal>("m", "m"));
127
128 solver_graph.register_edge("xyz", FieldRefs<Tvec>("xyz", "\\mathbf{r}"));
129 solver_graph.register_edge("vxyz", FieldRefs<Tvec>("vxyz", "\\mathbf{v}"));
130 solver_graph.register_edge("axyz", FieldRefs<Tvec>("axyz", "\\mathbf{a}"));
131 solver_graph.register_edge("uint", FieldRefs<Tscal>("uint", "u_{\\rm int}"));
132 solver_graph.register_edge("duint", FieldRefs<Tscal>("duint", "du_{\\rm int}"));
133 solver_graph.register_edge("hpart", FieldRefs<Tscal>("hpart", "h_{\\rm part}"));
134
135 if (has_B_field) {
136 solver_graph.register_edge("B/rho", FieldRefs<Tvec>("B/rho", "B_{\\rho}"));
137 solver_graph.register_edge("dB/rho", FieldRefs<Tvec>("dB/rho", "dB_{\\rho}"));
138 }
139 if (has_psi_field) {
140 solver_graph.register_edge("psi/ch", FieldRefs<Tscal>("psi/ch", "\\psi_{\\rm ch}"));
141 solver_graph.register_edge("dpsi/ch", FieldRefs<Tscal>("dpsi/ch", "d\\psi_{\\rm ch}"));
142 }
143 if (has_epsilon_field) {
144 solver_graph.register_edge("epsilon", FieldRefs<Tscal>("epsilon", "\\epsilon"));
145 solver_graph.register_edge("dtepsilon", FieldRefs<Tscal>("dtepsilon", "d\\epsilon"));
146 }
147 if (has_deltav_field) {
148 solver_graph.register_edge("deltav", FieldRefs<Tvec>("deltav", "\\Delta v"));
149 solver_graph.register_edge("dtdeltav", FieldRefs<Tvec>("dtdeltav", "d\\Delta v"));
150 }
151
152 {
153 auto set_gpart_mass = solver_graph.register_node(
154 "set_gpart_mass", NodeSetEdge<ScalarEdge<Tscal>>([&](ScalarEdge<Tscal> &gpart_mass) {
155 gpart_mass.value = solver_config.gpart_mass;
156 }));
157 shambase::get_check_ref(set_gpart_mass)
158 .set_edges(solver_graph.get_edge_ptr<ScalarEdge<Tscal>>("gpart_mass"));
159 }
160
162 // attach fields to scheduler
164 {
165 std::vector<std::shared_ptr<shamrock::solvergraph::INode>> attach_field_sequence;
166
167 {
168 auto set_scheduler_patchdata = solver_graph.register_node(
169 "set_scheduler_patchdata",
170 NodeSetEdge<PatchDataLayerRefs>([&](PatchDataLayerRefs &scheduler_patchdata) {
171 scheduler_patchdata.free_alloc();
172 scheduler().for_each_patchdata_nonempty(
173 [&](const shamrock::patch::Patch &p,
175 scheduler_patchdata.patchdatas.add_obj(p.id_patch, std::ref(pdat));
176 });
177 }));
178 shambase::get_check_ref(set_scheduler_patchdata)
179 .set_edges(solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"));
180 attach_field_sequence.push_back(set_scheduler_patchdata);
181 }
182
183 {
184 auto attach_part_counts
185 = solver_graph.register_node("attach_part_counts", GetObjCntFromLayer{});
186 shambase::get_check_ref(attach_part_counts)
187 .set_edges(
188 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
189 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"));
190 attach_field_sequence.push_back(attach_part_counts);
191 }
192
193 {
194 auto attach_xyz
195 = solver_graph.register_node("attach_xyz", GetFieldRefFromLayer<Tvec>(pdl, "xyz"));
196 shambase::get_check_ref(attach_xyz)
197 .set_edges(
198 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
199 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("xyz"));
200 attach_field_sequence.push_back(attach_xyz);
201 }
202
203 {
204 auto attach_vxyz = solver_graph.register_node(
205 "attach_vxyz", GetFieldRefFromLayer<Tvec>(pdl, "vxyz"));
206 shambase::get_check_ref(attach_vxyz)
207 .set_edges(
208 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
209 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("vxyz"));
210 attach_field_sequence.push_back(attach_vxyz);
211 }
212
213 {
214 auto attach_axyz = solver_graph.register_node(
215 "attach_axyz", GetFieldRefFromLayer<Tvec>(pdl, "axyz"));
216 shambase::get_check_ref(attach_axyz)
217 .set_edges(
218 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
219 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("axyz"));
220 attach_field_sequence.push_back(attach_axyz);
221 }
222
223 {
224 auto attach_uint = solver_graph.register_node(
225 "attach_uint", GetFieldRefFromLayer<Tscal>(pdl, "uint"));
226 shambase::get_check_ref(attach_uint)
227 .set_edges(
228 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
229 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("uint"));
230 attach_field_sequence.push_back(attach_uint);
231 }
232
233 {
234 auto attach_duint = solver_graph.register_node(
235 "attach_duint", GetFieldRefFromLayer<Tscal>(pdl, "duint"));
236 shambase::get_check_ref(attach_duint)
237 .set_edges(
238 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
239 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("duint"));
240 attach_field_sequence.push_back(attach_duint);
241 }
242
243 {
244 auto attach_hpart = solver_graph.register_node(
245 "attach_hpart", GetFieldRefFromLayer<Tscal>(pdl, "hpart"));
246 shambase::get_check_ref(attach_hpart)
247 .set_edges(
248 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
249 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("hpart"));
250 attach_field_sequence.push_back(attach_hpart);
251 }
252
253 if (has_B_field) {
254 auto attach_B_on_rho = solver_graph.register_node(
255 "attach_B_on_rho", GetFieldRefFromLayer<Tvec>(pdl, "B/rho"));
256 shambase::get_check_ref(attach_B_on_rho)
257 .set_edges(
258 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
259 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("B/rho"));
260 attach_field_sequence.push_back(attach_B_on_rho);
261 }
262
263 if (has_B_field) {
264 auto attach_dB_on_rho = solver_graph.register_node(
265 "attach_dB_on_rho", GetFieldRefFromLayer<Tvec>(pdl, "dB/rho"));
266 shambase::get_check_ref(attach_dB_on_rho)
267 .set_edges(
268 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
269 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("dB/rho"));
270 attach_field_sequence.push_back(attach_dB_on_rho);
271 }
272
273 if (has_psi_field) {
274 auto attach_psi_on_ch = solver_graph.register_node(
275 "attach_psi_on_ch", GetFieldRefFromLayer<Tscal>(pdl, "psi/ch"));
276 shambase::get_check_ref(attach_psi_on_ch)
277 .set_edges(
278 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
279 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("psi/ch"));
280 attach_field_sequence.push_back(attach_psi_on_ch);
281 }
282
283 if (has_psi_field) {
284 auto attach_dpsi_on_ch = solver_graph.register_node(
285 "attach_dpsi_on_ch", GetFieldRefFromLayer<Tscal>(pdl, "dpsi/ch"));
286 shambase::get_check_ref(attach_dpsi_on_ch)
287 .set_edges(
288 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
289 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("dpsi/ch"));
290 attach_field_sequence.push_back(attach_dpsi_on_ch);
291 }
292
293 if (has_epsilon_field) {
294 auto attach_epsilon = solver_graph.register_node(
295 "attach_epsilon", GetFieldRefFromLayer<Tscal>(pdl, "epsilon"));
296 shambase::get_check_ref(attach_epsilon)
297 .set_edges(
298 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
299 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("epsilon"));
300 attach_field_sequence.push_back(attach_epsilon);
301 }
302
303 if (has_epsilon_field) {
304 auto attach_dtepsilon = solver_graph.register_node(
305 "attach_dtepsilon", GetFieldRefFromLayer<Tscal>(pdl, "dtepsilon"));
306 shambase::get_check_ref(attach_dtepsilon)
307 .set_edges(
308 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
309 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("dtepsilon"));
310 attach_field_sequence.push_back(attach_dtepsilon);
311 }
312
313 if (has_deltav_field) {
314 auto attach_deltav = solver_graph.register_node(
315 "attach_deltav", GetFieldRefFromLayer<Tvec>(pdl, "deltav"));
316 shambase::get_check_ref(attach_deltav)
317 .set_edges(
318 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
319 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("deltav"));
320 attach_field_sequence.push_back(attach_deltav);
321 }
322
323 if (has_deltav_field) {
324 auto attach_dtdeltav = solver_graph.register_node(
325 "attach_dtdeltav", GetFieldRefFromLayer<Tvec>(pdl, "dtdeltav"));
326 shambase::get_check_ref(attach_dtdeltav)
327 .set_edges(
328 solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata"),
329 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("dtdeltav"));
330 attach_field_sequence.push_back(attach_dtdeltav);
331 }
332
333 solver_graph.register_node(
334 "attach fields to scheduler",
335 OperationSequence("attach fields", std::move(attach_field_sequence)));
336 }
337
339 // leapfrog predictor
341
342 {
343
344 auto make_half_step_sequence = [&](std::string prefix) {
345 std::vector<std::shared_ptr<shamrock::solvergraph::INode>> half_step_sequence;
346
347 {
348 auto half_step_vxyz = solver_graph.register_node(
350 shambase::get_check_ref(half_step_vxyz)
351 .set_edges(
352 solver_graph.get_edge_ptr<IDataEdge<Tscal>>("dt_half"),
353 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("axyz"),
354 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"),
355 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("vxyz"));
356 half_step_sequence.push_back(half_step_vxyz);
357 }
358
359 {
360 auto half_step_uint = solver_graph.register_node(
362 shambase::get_check_ref(half_step_uint)
363 .set_edges(
364 solver_graph.get_edge_ptr<IDataEdge<Tscal>>("dt_half"),
365 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("duint"),
366 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"),
367 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("uint"));
368 half_step_sequence.push_back(half_step_uint);
369 }
370
371 if (has_B_field) {
372 auto half_step_B_on_rho = solver_graph.register_node(
374 shambase::get_check_ref(half_step_B_on_rho)
375 .set_edges(
376 solver_graph.get_edge_ptr<IDataEdge<Tscal>>("dt_half"),
377 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("dB/rho"),
378 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"),
379 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("B/rho"));
380 half_step_sequence.push_back(half_step_B_on_rho);
381 }
382
383 if (has_psi_field) {
384 auto half_step_psi_on_ch = solver_graph.register_node(
385 prefix + "_psi_on_ch", shammodels::common::modules::ForwardEuler<Tscal>{});
386 shambase::get_check_ref(half_step_psi_on_ch)
387 .set_edges(
388 solver_graph.get_edge_ptr<IDataEdge<Tscal>>("dt_half"),
389 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("dpsi/ch"),
390 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"),
391 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("psi/ch"));
392 half_step_sequence.push_back(half_step_psi_on_ch);
393 }
394
395 if (has_epsilon_field) {
396 auto half_step_epsilon = solver_graph.register_node(
398 shambase::get_check_ref(half_step_epsilon)
399 .set_edges(
400 solver_graph.get_edge_ptr<IDataEdge<Tscal>>("dt_half"),
401 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("dtepsilon"),
402 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"),
403 solver_graph.get_edge_ptr<FieldRefs<Tscal>>("epsilon"));
404 half_step_sequence.push_back(half_step_epsilon);
405 }
406
407 if (has_deltav_field) {
408 auto half_step_deltav = solver_graph.register_node(
410 shambase::get_check_ref(half_step_deltav)
411 .set_edges(
412 solver_graph.get_edge_ptr<IDataEdge<Tscal>>("dt_half"),
413 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("dtdeltav"),
414 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"),
415 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("deltav"));
416 half_step_sequence.push_back(half_step_deltav);
417 }
418
419 return OperationSequence("half step", std::move(half_step_sequence));
420 };
421
422 solver_graph.register_node("half_step1", make_half_step_sequence("half_step1"));
423 solver_graph.register_node("half_step2", make_half_step_sequence("half_step2"));
424
425 {
426 auto full_step_xyz = solver_graph.register_node(
428 shambase::get_check_ref(full_step_xyz)
429 .set_edges(
430 solver_graph.get_edge_ptr<IDataEdge<Tscal>>("dt"),
431 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("vxyz"),
432 solver_graph.get_edge_ptr<Indexes<u32>>("part_counts"),
433 solver_graph.get_edge_ptr<FieldRefs<Tvec>>("xyz"));
434 }
435
436 {
437 auto leapfrog_predictor = solver_graph.register_node(
438 "leapfrog predictor",
440 "leapfrog predictor",
441 {
442 solver_graph.get_node_ptr_base("half_step1"),
443 solver_graph.get_node_ptr_base("full_step_xyz"),
444 solver_graph.get_node_ptr_base("half_step2"),
445 }));
446 }
447 }
448
450 // Part killing step
452 bool do_part_killing_step = solver_config.particle_killing.kill_list.size() > 0;
453
454 if (do_part_killing_step) {
455
456 auto patchdatas = solver_graph.get_edge_ptr<PatchDataLayerRefs>("scheduler_patchdata");
457 auto xyz_edge = solver_graph.get_edge_ptr<FieldRefs<Tvec>>("xyz");
458
459 auto part_to_remove = solver_graph.register_edge(
460 "part_to_remove", DistributedBuffers<u32>("part_to_remove", "part_to_remove"));
461
462 std::vector<std::shared_ptr<shamrock::solvergraph::INode>> part_kill_sequence{};
463
464 {
465
466 auto empty_part_to_remove
467 = solver_graph.register_node("empty_part_to_remove", NodeFreeAlloc{});
468 shambase::get_check_ref(empty_part_to_remove).set_edges(part_to_remove);
469 part_kill_sequence.push_back(empty_part_to_remove);
470 }
471
472 using kill_t = typename ParticleKillingConfig<Tvec>::kill_t;
473 using kill_sphere = typename ParticleKillingConfig<Tvec>::Sphere;
474
475 // selectors
476 for (kill_t &kill_obj : solver_config.particle_killing.kill_list) {
477 if (kill_sphere *kill_info = std::get_if<kill_sphere>(&kill_obj)) {
478
480 kill_info->center, kill_info->radius);
481 node_selector.set_edges(xyz_edge, part_to_remove);
482
483 part_kill_sequence.push_back(
484 std::make_shared<decltype(node_selector)>(std::move(node_selector)));
485 }
486 }
487
488 { // killing
489 modules::KillParticles node_killer{};
490 node_killer.set_edges(part_to_remove, patchdatas);
491
492 part_kill_sequence.push_back(
493 std::make_shared<decltype(node_killer)>(std::move(node_killer)));
494 }
495
496 solver_graph.register_node(
497 "part killing step",
498 OperationSequence("part killing step", std::move(part_kill_sequence)));
499 }
500
501 {
502 std::vector<std::shared_ptr<shamrock::solvergraph::INode>> seq{};
503
504 seq.push_back(solver_graph.get_node_ptr_base("set_gpart_mass"));
505 seq.push_back(solver_graph.get_node_ptr_base("attach fields to scheduler"));
506 seq.push_back(solver_graph.get_node_ptr_base("leapfrog predictor"));
507 if (do_part_killing_step) {
508 seq.push_back(solver_graph.get_node_ptr_base("part killing step"));
509 }
510
511 storage.solver_sequence = solver_graph.register_node(
512 "time_step", OperationSequence("time step", std::move(seq)));
513 }
514
515 storage.part_counts
516 = std::make_shared<shamrock::solvergraph::Indexes<u32>>("part_counts", "N_{\\rm part}");
517
518 storage.part_counts_with_ghost = std::make_shared<shamrock::solvergraph::Indexes<u32>>(
519 "part_counts_with_ghost", "N_{\\rm part, with ghost}");
520
521 storage.patch_rank_owner = std::make_shared<shamrock::solvergraph::RankGetter>(
522 [&](u64 patch_id) -> u32 {
523 return scheduler().get_patch_rank_owner(patch_id);
524 },
525 "patch_rank_owner",
526 "rank");
527
528 // merged ghost spans
529 storage.positions_with_ghosts
530 = std::make_shared<shamrock::solvergraph::FieldRefs<Tvec>>("part_pos", "\\mathbf{r}");
531 storage.hpart_with_ghosts
532 = std::make_shared<shamrock::solvergraph::FieldRefs<Tscal>>("h_part", "h");
533
534 storage.neigh_cache
535 = std::make_shared<shammodels::sph::solvergraph::NeighCache>("neigh_cache", "neigh");
536
537 storage.omega = std::make_shared<shamrock::solvergraph::Field<Tscal>>(1, "omega", "\\Omega");
538
539 if (solver_config.has_field_alphaAV()) {
540 storage.alpha_av_updated = std::make_shared<shamrock::solvergraph::Field<Tscal>>(
541 1, "alpha_av_updated", "\\alpha_{\\rm AV}");
542 }
543
544 storage.pressure = std::make_shared<shamrock::solvergraph::Field<Tscal>>(1, "pressure", "P");
545 storage.soundspeed
546 = std::make_shared<shamrock::solvergraph::Field<Tscal>>(1, "soundspeed", "c_s");
547
548 storage.exchange_gz_alpha
549 = std::make_shared<shamrock::solvergraph::ExchangeGhostField<Tscal>>();
550 storage.exchange_gz_node
551 = std::make_shared<shamrock::solvergraph::ExchangeGhostLayer>(storage.ghost_layout);
552 storage.exchange_gz_positions
553 = std::make_shared<shamrock::solvergraph::ExchangeGhostLayer>(storage.xyzh_ghost_layout);
554}
555
556template<class Tvec, template<class> class Kern>
558 std::string filename, bool add_patch_world_id) {
559
560 modules::VTKDump(context, solver_config).do_dump(filename, add_patch_world_id);
561}
562
564// Debug interface dump
566
567namespace shammodels::sph {
568
569 template<class Tvec>
571 using Tscal = shambase::VecComponent<Tvec>;
572
573 u64 nobj;
574 f64 gpart_mass;
575
576 sycl::buffer<Tvec> &buf_xyz;
577 sycl::buffer<Tscal> &buf_hpart;
578 sycl::buffer<Tvec> &buf_vxyz;
579 };
580
581 template<class Tvec>
582 void fill_blocks(PhantomDumpBlock &block, Debug_ph_dump<Tvec> &info) {
583
584 using Tscal = shambase::VecComponent<Tvec>;
585 std::vector<Tvec> xyz = shamalgs::memory::buf_to_vec(info.buf_xyz, info.nobj);
586
587 u64 xid = block.get_ref_fort_real("x");
588 u64 yid = block.get_ref_fort_real("y");
589 u64 zid = block.get_ref_fort_real("z");
590
591 for (auto vec : xyz) {
592 block.blocks_fort_real[xid].vals.push_back(vec.x());
593 block.blocks_fort_real[yid].vals.push_back(vec.y());
594 block.blocks_fort_real[zid].vals.push_back(vec.z());
595 }
596
597 std::vector<Tscal> h = shamalgs::memory::buf_to_vec(info.buf_hpart, info.nobj);
598 u64 hid = block.get_ref_f32("h");
599 for (auto h_ : h) {
600 block.blocks_f32[hid].vals.push_back(h_);
601 }
602
603 std::vector<Tvec> vxyz = shamalgs::memory::buf_to_vec(info.buf_vxyz, info.nobj);
604
605 u64 vxid = block.get_ref_fort_real("vx");
606 u64 vyid = block.get_ref_fort_real("vy");
607 u64 vzid = block.get_ref_fort_real("vz");
608
609 for (auto vec : vxyz) {
610 block.blocks_fort_real[vxid].vals.push_back(vec.x());
611 block.blocks_fort_real[vyid].vals.push_back(vec.y());
612 block.blocks_fort_real[vzid].vals.push_back(vec.z());
613 }
614
615 block.tot_count = block.blocks_fort_real[xid].vals.size();
616 }
617
618 template<class Tvec>
619 shammodels::sph::PhantomDump make_interface_debug_phantom_dump(Debug_ph_dump<Tvec> info) {
620
621 using Tscal = shambase::VecComponent<Tvec>;
622 PhantomDump dump;
623
625 dump.iversion = 1;
626 dump.fileid = shambase::format("{:100s}", "FT:Phantom Shamrock writer");
627
628 u32 Ntot = info.nobj;
629 dump.table_header_fort_int.add("nparttot", Ntot);
630 dump.table_header_fort_int.add("ntypes", 8);
631 dump.table_header_fort_int.add("npartoftype", Ntot);
632 dump.table_header_fort_int.add("npartoftype", 0);
633 dump.table_header_fort_int.add("npartoftype", 0);
634 dump.table_header_fort_int.add("npartoftype", 0);
635 dump.table_header_fort_int.add("npartoftype", 0);
636 dump.table_header_fort_int.add("npartoftype", 0);
637 dump.table_header_fort_int.add("npartoftype", 0);
638 dump.table_header_fort_int.add("npartoftype", 0);
639
640 dump.table_header_i64.add("nparttot", Ntot);
641 dump.table_header_i64.add("ntypes", 8);
642 dump.table_header_i64.add("npartoftype", Ntot);
643 dump.table_header_i64.add("npartoftype", 0);
644 dump.table_header_i64.add("npartoftype", 0);
645 dump.table_header_i64.add("npartoftype", 0);
646 dump.table_header_i64.add("npartoftype", 0);
647 dump.table_header_i64.add("npartoftype", 0);
648 dump.table_header_i64.add("npartoftype", 0);
649 dump.table_header_i64.add("npartoftype", 0);
650
651 dump.table_header_fort_int.add("nblocks", 1);
652 dump.table_header_fort_int.add("nptmass", 0);
653 dump.table_header_fort_int.add("ndustlarge", 0);
654 dump.table_header_fort_int.add("ndustsmall", 0);
655 dump.table_header_fort_int.add("idust", 7);
656 dump.table_header_fort_int.add("idtmax_n", 1);
657 dump.table_header_fort_int.add("idtmax_frac", 0);
658 dump.table_header_fort_int.add("idumpfile", 0);
659 dump.table_header_fort_int.add("majorv", 2023);
660 dump.table_header_fort_int.add("minorv", 0);
661 dump.table_header_fort_int.add("microv", 0);
662 dump.table_header_fort_int.add("isink", 0);
663
664 dump.table_header_i32.add("iexternalforce", 0);
665 dump.table_header_i32.add("ieos", 2);
666 dump.table_header_fort_real.add("gamma", 1.66667);
667 dump.table_header_fort_real.add("RK2", 0);
668 dump.table_header_fort_real.add("polyk2", 0);
669 dump.table_header_fort_real.add("qfacdisc", 0.75);
670 dump.table_header_fort_real.add("qfacdisc2", 0.75);
671
672 dump.table_header_fort_real.add("time", 0);
673 dump.table_header_fort_real.add("dtmax", 0.1);
674
675 dump.table_header_fort_real.add("rhozero", 0);
676 dump.table_header_fort_real.add("hfact", 1.2);
677 dump.table_header_fort_real.add("tolh", 0.0001);
678 dump.table_header_fort_real.add("C_cour", 0);
679 dump.table_header_fort_real.add("C_force", 0);
680 dump.table_header_fort_real.add("alpha", 0);
681 dump.table_header_fort_real.add("alphau", 1);
682 dump.table_header_fort_real.add("alphaB", 1);
683
684 dump.table_header_fort_real.add("massoftype", info.gpart_mass);
685 dump.table_header_fort_real.add("massoftype", 0);
686 dump.table_header_fort_real.add("massoftype", 0);
687 dump.table_header_fort_real.add("massoftype", 0);
688 dump.table_header_fort_real.add("massoftype", 0);
689 dump.table_header_fort_real.add("massoftype", 0);
690 dump.table_header_fort_real.add("massoftype", 0);
691 dump.table_header_fort_real.add("massoftype", 0);
692
693 dump.table_header_fort_real.add("Bextx", 0);
694 dump.table_header_fort_real.add("Bexty", 0);
695 dump.table_header_fort_real.add("Bextz", 0);
696 dump.table_header_fort_real.add("dum", 0);
697
698 dump.table_header_fort_real.add("get_conserv", -1);
699 dump.table_header_fort_real.add("etot_in", 0.59762);
700 dump.table_header_fort_real.add("angtot_in", 0.0189694);
701 dump.table_header_fort_real.add("totmom_in", 0.0306284);
702
703 dump.table_header_f64.add("udist", 1);
704 dump.table_header_f64.add("umass", 1);
705 dump.table_header_f64.add("utime", 1);
706 dump.table_header_f64.add("umagfd", 3.54491);
707
708 PhantomDumpBlock block_part;
709
710 fill_blocks(block_part, info);
711
712 dump.blocks.push_back(std::move(block_part));
713
714 return dump;
715 }
716
717} // namespace shammodels::sph
718
719template<class Tvec, template<class> class Kern>
721 StackEntry stack_loc{};
722
724 _sptree.attach_buf();
725 storage.serial_patch_tree.set(std::move(_sptree));
726}
727
733template<class Tvec, template<class> class Kern>
735
736 StackEntry stack_loc{};
737
738 shamlog_debug_ln("SphSolver", "apply position boundary");
739
740 PatchScheduler &sched = scheduler();
741
742 shamrock::SchedulerUtility integrators(sched);
744
745 auto &pdl = sched.pdl_old();
746
747 const u32 ixyz = pdl.get_field_idx<Tvec>("xyz");
748 const u32 ivxyz = pdl.get_field_idx<Tvec>("vxyz");
749 auto [bmin, bmax] = sched.get_box_volume<Tvec>();
750
751 using SolverConfigBC = typename Config::BCConfig;
752 using SolverBCFree = typename SolverConfigBC::Free;
753 using SolverBCPeriodic = typename SolverConfigBC::Periodic;
754 using SolverBCShearingPeriodic = typename SolverConfigBC::ShearingPeriodic;
755 if (SolverBCFree *c = std::get_if<SolverBCFree>(&solver_config.boundary_config.config)) {
756 if (shamcomm::world_rank() == 0) {
757 logger::info_ln("PositionUpdated", "free boundaries skipping geometry update");
758 }
759 } else if (
760 SolverBCPeriodic *c
761 = std::get_if<SolverBCPeriodic>(&solver_config.boundary_config.config)) {
762 integrators.fields_apply_periodicity(ixyz, std::pair{bmin, bmax});
763 } else if (
764 SolverBCShearingPeriodic *c
765 = std::get_if<SolverBCShearingPeriodic>(&solver_config.boundary_config.config)) {
766 integrators.fields_apply_shearing_periodicity(
767 ixyz,
768 ivxyz,
769 std::pair{bmin, bmax},
770 c->shear_base,
771 c->shear_dir,
772 c->shear_speed * time_val,
773 c->shear_speed);
774 }
775
776 reatrib.reatribute_patch_objects(storage.serial_patch_tree.get(), "xyz");
777}
778
779template<class Tvec, template<class> class Kern>
781
782 StackEntry stack_loc{};
783
784 using SPHUtils = sph::SPHUtilities<Tvec, Kernel>;
785 SPHUtils sph_utils(scheduler());
786
787 storage.ghost_patch_cache.set(sph_utils.build_interf_cache(
788 storage.ghost_handler.get(),
789 storage.serial_patch_tree.get(),
790 solver_config.htol_up_coarse_cycle));
791
792 // storage.ghost_handler.get().gen_debug_patch_ghost(storage.ghost_patch_cache.get());
793}
794
795template<class Tvec, template<class> class Kern>
797 StackEntry stack_loc{};
798 storage.ghost_patch_cache.reset();
799}
800
801template<class Tvec, template<class> class Kern>
803
804 StackEntry stack_loc{};
805
806 storage.merged_xyzh.set(storage.ghost_handler.get().build_comm_merge_positions(
807 storage.ghost_patch_cache.get(),
808 storage.exchange_gz_positions,
809 solver_config.show_ghost_zone_graph));
810
811 { // set element counts
812 shambase::get_check_ref(storage.part_counts).indexes
813 = storage.merged_xyzh.get().template map<u32>(
814 [&](u64 id, shamrock::patch::PatchDataLayer &mpdat) {
815 return scheduler().patch_data.get_pdat(id).get_obj_cnt();
816 });
817 }
818
819 { // set element counts
820 shambase::get_check_ref(storage.part_counts_with_ghost).indexes
821 = storage.merged_xyzh.get().template map<u32>(
822 [&](u64 id, shamrock::patch::PatchDataLayer &mpdat) {
823 return mpdat.get_obj_cnt();
824 });
825 }
826
827 { // Attach spans to block coords
828 shambase::get_check_ref(storage.positions_with_ghosts)
829 .set_refs(storage.merged_xyzh.get()
830 .template map<std::reference_wrapper<PatchDataField<Tvec>>>(
831 [&](u64 id, shamrock::patch::PatchDataLayer &mpdat) {
832 return std::ref(mpdat.get_field<Tvec>(0));
833 }));
834
835 shambase::get_check_ref(storage.hpart_with_ghosts)
836 .set_refs(storage.merged_xyzh.get()
837 .template map<std::reference_wrapper<PatchDataField<Tscal>>>(
838 [&](u64 id, shamrock::patch::PatchDataLayer &mpdat) {
839 return std::ref(mpdat.get_field<Tscal>(1));
840 }));
841 }
842}
843
844template<class Tvec, template<class> class Kern>
848
849template<class Tvec, template<class> class Kern>
851 StackEntry stack_loc{};
852 storage.merged_pos_trees.reset();
853}
854
855template<class Tvec, template<class> class Kern>
857
858 StackEntry stack_loc{};
859 using namespace shamrock::patch;
860 PatchDataLayerLayout &pdl = scheduler().pdl_old();
861 const u32 ixyz = pdl.get_field_idx<Tvec>("xyz");
862 const u32 ivxyz = pdl.get_field_idx<Tvec>("vxyz");
863 const u32 iaxyz = pdl.get_field_idx<Tvec>("axyz");
864 const u32 iuint = pdl.get_field_idx<Tscal>("uint");
865 const u32 iduint = pdl.get_field_idx<Tscal>("duint");
866
867 bool has_B_field = solver_config.has_field_B_on_rho();
868 bool has_psi_field = solver_config.has_field_psi_on_ch();
869 bool has_epsilon_field = solver_config.dust_config.has_epsilon_field();
870 bool has_deltav_field = solver_config.dust_config.has_deltav_field();
871
872 const u32 iB_on_rho = (has_B_field) ? pdl.get_field_idx<Tvec>("B/rho") : 0;
873 const u32 idB_on_rho = (has_B_field) ? pdl.get_field_idx<Tvec>("dB/rho") : 0;
874 const u32 ipsi_on_ch = (has_psi_field) ? pdl.get_field_idx<Tscal>("psi/ch") : 0;
875 const u32 idpsi_on_ch = (has_psi_field) ? pdl.get_field_idx<Tscal>("dpsi/ch") : 0;
876
877 const u32 iepsilon = (has_epsilon_field) ? pdl.get_field_idx<Tscal>("epsilon") : 0;
878 const u32 idtepsilon = (has_epsilon_field) ? pdl.get_field_idx<Tscal>("dtepsilon") : 0;
879 const u32 ideltav = (has_deltav_field) ? pdl.get_field_idx<Tvec>("deltav") : 0;
880 const u32 idtdeltav = (has_deltav_field) ? pdl.get_field_idx<Tvec>("dtdeltav") : 0;
881
882 shamrock::SchedulerUtility utility(scheduler());
883
884 // forward euler step f dt/2
885 shamlog_debug_ln("sph::BasicGas", "forward euler step f dt/2");
886 utility.fields_forward_euler<Tvec>(ivxyz, iaxyz, dt / 2);
887 utility.fields_forward_euler<Tscal>(iuint, iduint, dt / 2);
888
889 if (has_B_field) {
890 utility.fields_forward_euler<Tvec>(iB_on_rho, idB_on_rho, dt / 2);
891 }
892 if (has_psi_field) {
893 utility.fields_forward_euler<Tscal>(ipsi_on_ch, idpsi_on_ch, dt / 2);
894 }
895
896 // forward euler step positions dt
897 shamlog_debug_ln("sph::BasicGas", "forward euler step positions dt");
898 utility.fields_forward_euler<Tvec>(ixyz, ivxyz, dt);
899
900 // forward euler step f dt/2
901 shamlog_debug_ln("sph::BasicGas", "forward euler step f dt/2");
902 utility.fields_forward_euler<Tvec>(ivxyz, iaxyz, dt / 2);
903 utility.fields_forward_euler<Tscal>(iuint, iduint, dt / 2);
904
905 if (has_B_field) {
906 utility.fields_forward_euler<Tvec>(iB_on_rho, idB_on_rho, dt / 2);
907 }
908 if (has_psi_field) {
909 utility.fields_forward_euler<Tscal>(ipsi_on_ch, idpsi_on_ch, dt / 2);
910 }
911 if (has_epsilon_field) {
912 utility.fields_forward_euler<Tscal>(iepsilon, idtepsilon, dt / 2);
913 }
914 if (has_deltav_field) {
915 utility.fields_forward_euler<Tvec>(ideltav, idtdeltav, dt / 2);
916 }
917}
918
919template<class Tvec, template<class> class Kern>
921 StackEntry stack_loc{};
922
923 using namespace shamrock;
924 using namespace shamrock::patch;
925
927 using SPHUtils = sph::SPHUtilities<Tvec, Kernel>;
928
929 SPHUtils sph_utils(scheduler());
930 shamrock::SchedulerUtility utility(scheduler());
931
932 PatchDataLayerLayout &pdl = scheduler().pdl_old();
933 const u32 ihpart = pdl.get_field_idx<Tscal>("hpart");
934
935 ComputeField<Tscal> _epsilon_h, _h_old;
936
937 auto should_set_omega_mask = std::make_shared<shamrock::solvergraph::Field<u32>>(
938 1, "should_set_omega_mask", "should_set_omega_mask");
939
940 u32 hstep_cnt = 0;
941 u32 hstep_max = solver_config.h_max_subcycles_count;
942 for (; hstep_cnt < hstep_max; hstep_cnt++) {
943
944 gen_ghost_handler(time_val + dt);
945 build_ghost_cache();
946 merge_position_ghost();
947 build_merged_pos_trees();
948 compute_presteps_rint();
949 start_neighbors_cache();
950
951 _epsilon_h = utility.make_compute_field<Tscal>("epsilon_h", 1, Tscal(100));
952 _h_old = utility.save_field<Tscal>(ihpart, "h_old");
953
954 Tscal max_eps_h;
955
956 if (solver_config.gpart_mass == 0) {
958 "invalid gpart_mass {}, this configuration can not converge.\n"
959 "Please set it using either model.set_particle_mass(pmass) or "
960 "cfg.set_particle_mass(pmass)",
961 solver_config.gpart_mass));
962 }
963
964 // sizes
965 std::shared_ptr<shamrock::solvergraph::Indexes<u32>> sizes
966 = std::make_shared<shamrock::solvergraph::Indexes<u32>>("", "");
967 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
968 sizes->indexes.add_obj(p.id_patch, pdat.get_obj_cnt());
969 });
970
971 // neigh cache
972 auto &neigh_cache = storage.neigh_cache;
973
974 // positions
975 auto &pos_merged = storage.positions_with_ghosts;
976
977 // old smoothing length field
978 std::shared_ptr<shamrock::solvergraph::FieldRefs<Tscal>> hold
979 = std::make_shared<shamrock::solvergraph::FieldRefs<Tscal>>("", "");
981 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
982 auto &field = _h_old.get_field(p.id_patch);
983 hold_refs.add_obj(p.id_patch, std::ref(field));
984 });
985 hold->set_refs(hold_refs);
986
987 // new smoothing length field
988 std::shared_ptr<shamrock::solvergraph::FieldRefs<Tscal>> hnew
989 = std::make_shared<shamrock::solvergraph::FieldRefs<Tscal>>("", "");
991 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
992 auto &field = pdat.get_field<Tscal>(ihpart);
993 hnew_refs.add_obj(p.id_patch, std::ref(field));
994 });
995 hnew->set_refs(hnew_refs);
996
997 // epsilon field
998 std::shared_ptr<shamrock::solvergraph::FieldRefs<Tscal>> eps_h
999 = std::make_shared<shamrock::solvergraph::FieldRefs<Tscal>>("", "");
1001 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
1002 auto &field = _epsilon_h.get_field(p.id_patch);
1003 eps_h_refs.add_obj(p.id_patch, std::ref(field));
1004 });
1005 eps_h->set_refs(eps_h_refs);
1006
1007 std::shared_ptr<shamrock::solvergraph::INode> smth_h_iter_ptr;
1008
1009 using h_conf_density_based = typename SmoothingLengthConfig::DensityBased;
1010 using h_conf_neigh_lim = typename SmoothingLengthConfig::DensityBasedNeighLim;
1011
1012 if (h_conf_density_based *conf
1013 = std::get_if<h_conf_density_based>(&solver_config.smoothing_length_config.config)) {
1014 std::shared_ptr<shammodels::sph::modules::IterateSmoothingLengthDensity<Tvec, Kernel>>
1015 smth_h_iter = std::make_shared<
1017 solver_config.gpart_mass,
1018 solver_config.htol_up_coarse_cycle,
1019 solver_config.htol_up_fine_cycle);
1020 smth_h_iter->set_edges(sizes, neigh_cache, pos_merged, hold, hnew, eps_h);
1021 smth_h_iter_ptr = smth_h_iter;
1022 } else if (
1023 h_conf_neigh_lim *conf
1024 = std::get_if<h_conf_neigh_lim>(&solver_config.smoothing_length_config.config)) {
1025 std::shared_ptr<
1027 smth_h_iter_neigh_lim = std::make_shared<
1029 solver_config.gpart_mass,
1030 solver_config.htol_up_coarse_cycle,
1031 solver_config.htol_up_fine_cycle,
1032 conf->max_neigh_count);
1033 smth_h_iter_neigh_lim->set_edges(
1034 sizes, neigh_cache, pos_merged, hold, hnew, eps_h, should_set_omega_mask);
1035 smth_h_iter_ptr = smth_h_iter_neigh_lim;
1036 } else {
1037 shambase::throw_with_loc<std::runtime_error>("Invalid smoothing length configuration");
1038 }
1039 // iterate smoothing length
1040
1041 std::shared_ptr<shamrock::solvergraph::ScalarEdge<bool>> is_converged
1042 = std::make_shared<shamrock::solvergraph::ScalarEdge<bool>>("", "");
1043
1045 smth_h_iter_ptr, solver_config.epsilon_h, solver_config.h_iter_per_subcycles, false);
1046 loop_smth_h_iter.set_edges(eps_h, is_converged);
1047
1048 loop_smth_h_iter.evaluate();
1049
1050 if (!is_converged->value) {
1051
1052 Tscal largest_h = 0;
1053
1054 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
1055 largest_h = sham::max(largest_h, pdat.get_field<Tscal>(ihpart).compute_max());
1056 });
1057 Tscal global_largest_h = shamalgs::collective::allreduce_max(largest_h);
1058
1059 std::string add_info = "";
1060 u64 cnt_unconverged = 0;
1061 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
1062 auto res
1063 = _epsilon_h.get_field(p.id_patch).get_ids_buf_where([](auto access, u32 id) {
1064 return access[id] == -1;
1065 });
1066
1067 if (hstep_cnt == hstep_max - 1) {
1068 if (std::get<0>(res)) {
1069 add_info += "\n patch " + std::to_string(p.id_patch) + " ";
1070 add_info += "errored parts : \n";
1071 sycl::buffer<u32> &idx_err = *std::get<0>(res);
1072
1073 sham::DeviceBuffer<Tvec> &xyz = pdat.get_field_buf_ref<Tvec>(0);
1074 sham::DeviceBuffer<Tscal> &hpart = pdat.get_field_buf_ref<Tscal>(ihpart);
1075
1076 auto pos = xyz.copy_to_stdvec();
1077 auto h = hpart.copy_to_stdvec();
1078
1079 {
1080 sycl::host_accessor acc{idx_err};
1081 for (u32 i = 0; i < idx_err.size(); i++) {
1082 add_info += shambase::format(
1083 "{} - pos : {}, hpart : {}\n", acc[i], pos[acc[i]], h[acc[i]]);
1084 }
1085 }
1086 }
1087 }
1088
1089 cnt_unconverged += std::get<1>(res);
1090 });
1091
1092 u64 global_cnt_unconverged = shamalgs::collective::allreduce_sum(cnt_unconverged);
1093
1094 if (shamcomm::world_rank() == 0) {
1095 logger::warn_ln(
1096 "Smoothinglength",
1097 "smoothing length is not converged, rerunning the iterator ...\n largest h "
1098 "=",
1099 global_largest_h,
1100 "unconverged cnt =",
1101 global_cnt_unconverged,
1102 add_info);
1103 }
1104
1105 reset_ghost_handler();
1106 clear_ghost_cache();
1107
1108 shambase::get_check_ref(storage.part_counts).free_alloc();
1109 shambase::get_check_ref(storage.part_counts_with_ghost).free_alloc();
1110 shambase::get_check_ref(storage.positions_with_ghosts).free_alloc();
1111 shambase::get_check_ref(storage.hpart_with_ghosts).free_alloc();
1112
1113 storage.merged_xyzh.reset();
1114
1115 clear_merged_pos_trees();
1116 reset_presteps_rint();
1117 reset_neighbors_cache();
1118
1119 // scheduler().for_each_patchdata_nonempty([&](Patch cur_p, PatchData &pdat) {
1120 // pdat.synchronize_buf();
1121 // });
1122
1123 continue;
1124 }
1125
1126 // The hpart is not valid anymore in ghost zones since we iterated it's value
1127 shambase::get_check_ref(storage.hpart_with_ghosts).free_alloc();
1128
1129 _epsilon_h.reset();
1130 _h_old.reset();
1131 break;
1132 }
1133
1134 if (hstep_cnt == hstep_max) {
1135 logger::err_ln("SPH", "the h iterator is not converged after", hstep_cnt, "iterations");
1136 }
1137
1138 std::shared_ptr<shamrock::solvergraph::FieldRefs<Tscal>> hnew_edge
1139 = std::make_shared<shamrock::solvergraph::FieldRefs<Tscal>>("", "");
1141 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
1142 auto &field = pdat.get_field<Tscal>(ihpart);
1143 hnew_refs.add_obj(p.id_patch, std::ref(field));
1144 });
1145 hnew_edge->set_refs(hnew_refs);
1146
1147 modules::NodeComputeOmega<Tvec, Kern> compute_omega{solver_config.gpart_mass};
1148 compute_omega.set_edges(
1149 storage.part_counts,
1150 storage.neigh_cache,
1151 storage.positions_with_ghosts,
1152 hnew_edge,
1153 storage.omega);
1154 compute_omega.evaluate();
1155
1156 if (solver_config.smoothing_length_config.is_density_based_neigh_lim()) {
1157 // if the h limiter is triggered, omega does not hold it's sense of dh/dr anymore
1158 // so we set it to 1, this effectively is equivalent of disabling the energy correction
1159 // term corresponding to dh/dr
1160 modules::SetWhenMask<Tscal> set_omega_mask{1};
1161 set_omega_mask.set_edges(storage.part_counts, should_set_omega_mask, storage.omega);
1162 set_omega_mask.evaluate();
1163 }
1164}
1165
1166template<class Tvec, template<class> class Kern>
1168
1169 storage.ghost_layout = std::make_shared<shamrock::patch::PatchDataLayerLayout>();
1170
1172 = shambase::get_check_ref(storage.ghost_layout);
1173
1174 solver_config.set_ghost_layout(ghost_layout);
1175
1176 storage.xyzh_ghost_layout = std::make_shared<shamrock::patch::PatchDataLayerLayout>();
1177 storage.xyzh_ghost_layout->template add_field<Tvec>("xyz", 1);
1178 storage.xyzh_ghost_layout->template add_field<Tscal>("hpart", 1);
1179}
1180
1181template<class Tvec, template<class> class Kern>
1183
1184 StackEntry stack_loc{};
1185
1186 auto &xyzh_merged = storage.merged_xyzh.get();
1187 auto dev_sched = shamsys::instance::get_compute_scheduler_ptr();
1188
1189 storage.rtree_rint_field.set(
1190 storage.merged_pos_trees.get().template map<shamtree::KarrasRadixTreeField<Tscal>>(
1191 [&](u64 id, RTree &rtree) -> shamtree::KarrasRadixTreeField<Tscal> {
1192 shamrock::patch::PatchDataLayer &tmp = xyzh_merged.get(id);
1193 auto &buf = tmp.get_field_buf_ref<Tscal>(1);
1194 auto buf_int = shamtree::new_empty_karras_radix_tree_field<Tscal>();
1195
1196 auto ret = shamtree::compute_tree_field_max_field<Tscal>(
1197 rtree.structure,
1198 rtree.reduced_morton_set.get_leaf_cell_iterator(),
1199 std::move(buf_int),
1200 buf);
1201
1202 // the old tree used to increase the size of the hmax of the tree nodes by the
1203 // tolerance so we do it also with the new tree, maybe we should move that somewhere
1204 // else.
1205 sham::kernel_call(
1206 dev_sched->get_queue(),
1207 sham::MultiRef{},
1208 sham::MultiRef{ret.buf_field},
1209 ret.buf_field.get_size(),
1210 [htol = solver_config.htol_up_coarse_cycle](u32 i, Tscal *h_tree) {
1211 h_tree[i] *= htol;
1212 });
1213
1214 return std::move(ret);
1215 }));
1216}
1217
1218template<class Tvec, template<class> class Kern>
1220 storage.rtree_rint_field.reset();
1221}
1222
1223template<class Tvec, template<class> class Kern>
1225 if (solver_config.use_two_stage_search) {
1227 context, solver_config, storage)
1228 .start_neighbors_cache_2stages();
1229 } else {
1231 context, solver_config, storage)
1232 .start_neighbors_cache();
1233 }
1234
1235 if (solver_config.show_neigh_stats) {
1236 auto &pos_merged = storage.positions_with_ghosts;
1237 auto &neigh_cache = storage.neigh_cache;
1238 auto &hpart_with_ghosts = storage.hpart_with_ghosts;
1239 auto &part_counts = storage.part_counts;
1240
1241 modules::ComputeNeighStats<Tvec> compute_neigh_stats(Kernel::Rkern);
1242
1243 compute_neigh_stats.set_edges(part_counts, neigh_cache, pos_merged, hpart_with_ghosts);
1244 compute_neigh_stats.evaluate();
1245 }
1246}
1247
1248template<class Tvec, template<class> class Kern>
1250 // storage.neighbors_cache.reset();
1251}
1252
1253template<class Tvec, template<class> class Kern>
1255
1256 StackEntry stack_loc{};
1257
1258 shambase::Timer timer_interf;
1259 timer_interf.start();
1260
1261 using namespace shamrock;
1262 using namespace shamrock::patch;
1263
1264 bool has_alphaAV_field = solver_config.has_field_alphaAV();
1265 bool has_soundspeed_field = solver_config.ghost_has_soundspeed();
1266
1267 bool has_B_field = solver_config.has_field_B_on_rho();
1268 bool has_psi_field = solver_config.has_field_psi_on_ch();
1269 bool has_curlB_field = solver_config.has_field_curlB();
1270 bool has_epsilon_field = solver_config.dust_config.has_epsilon_field();
1271 bool has_deltav_field = solver_config.dust_config.has_deltav_field();
1272
1273 PatchDataLayerLayout &pdl = scheduler().pdl_old();
1274 const u32 ixyz = pdl.get_field_idx<Tvec>("xyz");
1275 const u32 ivxyz = pdl.get_field_idx<Tvec>("vxyz");
1276 const u32 iaxyz = pdl.get_field_idx<Tvec>("axyz");
1277 const u32 iuint = pdl.get_field_idx<Tscal>("uint");
1278 const u32 iduint = pdl.get_field_idx<Tscal>("duint");
1279 const u32 ihpart = pdl.get_field_idx<Tscal>("hpart");
1280
1281 const u32 ialpha_AV = (has_alphaAV_field) ? pdl.get_field_idx<Tscal>("alpha_AV") : 0;
1282 const u32 isoundspeed = (has_soundspeed_field) ? pdl.get_field_idx<Tscal>("soundspeed") : 0;
1283
1284 const u32 iB_on_rho = (has_B_field) ? pdl.get_field_idx<Tvec>("B/rho") : 0;
1285 const u32 idB_on_rho = (has_B_field) ? pdl.get_field_idx<Tvec>("dB/rho") : 0;
1286 const u32 ipsi_on_ch = (has_psi_field) ? pdl.get_field_idx<Tscal>("psi/ch") : 0;
1287 const u32 idpsi_on_ch = (has_psi_field) ? pdl.get_field_idx<Tscal>("dpsi/ch") : 0;
1288 const u32 icurlB = (has_curlB_field) ? pdl.get_field_idx<Tvec>("curlB") : 0;
1289
1290 bool do_MHD_debug = solver_config.do_MHD_debug();
1291 const u32 imag_pressure = (do_MHD_debug) ? pdl.get_field_idx<Tvec>("mag_pressure") : -1;
1292 const u32 imag_tension = (do_MHD_debug) ? pdl.get_field_idx<Tvec>("mag_tension") : -1;
1293 const u32 igas_pressure = (do_MHD_debug) ? pdl.get_field_idx<Tvec>("gas_pressure") : -1;
1294 const u32 itensile_corr = (do_MHD_debug) ? pdl.get_field_idx<Tvec>("tensile_corr") : -1;
1295 const u32 ipsi_propag = (do_MHD_debug) ? pdl.get_field_idx<Tscal>("psi_propag") : -1;
1296 const u32 ipsi_diff = (do_MHD_debug) ? pdl.get_field_idx<Tscal>("psi_diff") : -1;
1297 const u32 ipsi_cons = (do_MHD_debug) ? pdl.get_field_idx<Tscal>("psi_cons") : -1;
1298 const u32 iu_mhd = (do_MHD_debug) ? pdl.get_field_idx<Tscal>("u_mhd") : -1;
1299
1300 const u32 iepsilon = (has_epsilon_field) ? pdl.get_field_idx<Tscal>("epsilon") : 0;
1301 const u32 ideltav = (has_deltav_field) ? pdl.get_field_idx<Tvec>("deltav") : 0;
1302
1303 auto &ghost_layout_ptr = storage.ghost_layout;
1304 shamrock::patch::PatchDataLayerLayout &ghost_layout = shambase::get_check_ref(ghost_layout_ptr);
1305 u32 ihpart_interf = ghost_layout.get_field_idx<Tscal>("hpart");
1306 u32 iuint_interf = ghost_layout.get_field_idx<Tscal>("uint");
1307 u32 ivxyz_interf = ghost_layout.get_field_idx<Tvec>("vxyz");
1308 u32 iomega_interf = ghost_layout.get_field_idx<Tscal>("omega");
1309
1310 const u32 iaxyz_interf
1311 = (solver_config.has_axyz_in_ghost()) ? ghost_layout.get_field_idx<Tvec>("axyz") : 0;
1312
1313 const u32 isoundspeed_interf
1314 = (has_soundspeed_field) ? ghost_layout.get_field_idx<Tscal>("soundspeed") : 0;
1315
1316 const u32 iB_interf = (has_B_field) ? ghost_layout.get_field_idx<Tvec>("B/rho") : 0;
1317 const u32 ipsi_interf = (has_psi_field) ? ghost_layout.get_field_idx<Tscal>("psi/ch") : 0;
1318 const u32 icurlB_interf = (has_curlB_field) ? ghost_layout.get_field_idx<Tvec>("curlB") : 0;
1319
1320 const u32 iepsilon_interf
1321 = (has_epsilon_field) ? ghost_layout.get_field_idx<Tscal>("epsilon") : 0;
1322 const u32 ideltav_interf = (has_deltav_field) ? ghost_layout.get_field_idx<Tvec>("deltav") : 0;
1323
1324 using InterfaceBuildInfos = typename sph::BasicSPHGhostHandler<Tvec>::InterfaceBuildInfos;
1325
1326 sph::BasicSPHGhostHandler<Tvec> &ghost_handle = storage.ghost_handler.get();
1328
1329 auto pdat_interf = ghost_handle.template build_interface_native<PatchDataLayer>(
1330 storage.ghost_patch_cache.get(),
1331 [&](u64 sender, u64, InterfaceBuildInfos binfo, sham::DeviceBuffer<u32> &buf_idx, u32 cnt) {
1332 PatchDataLayer pdat(ghost_layout_ptr);
1333
1334 pdat.reserve(cnt);
1335
1336 return pdat;
1337 });
1338
1339 ghost_handle.template modify_interface_native<PatchDataLayer>(
1340 storage.ghost_patch_cache.get(),
1341 pdat_interf,
1342 [&](u64 sender,
1343 u64,
1344 InterfaceBuildInfos binfo,
1345 sham::DeviceBuffer<u32> &buf_idx,
1346 u32 cnt,
1347 PatchDataLayer &pdat) {
1348 PatchDataLayer &sender_patch = scheduler().patch_data.get_pdat(sender);
1349 PatchDataField<Tscal> &sender_omega = omega.get(sender);
1350
1351 sender_patch.get_field<Tscal>(ihpart).append_subset_to(
1352 buf_idx, cnt, pdat.get_field<Tscal>(ihpart_interf));
1353 sender_patch.get_field<Tscal>(iuint).append_subset_to(
1354 buf_idx, cnt, pdat.get_field<Tscal>(iuint_interf));
1355
1356 if (solver_config.has_axyz_in_ghost()) {
1357 sender_patch.get_field<Tvec>(iaxyz).append_subset_to(
1358 buf_idx, cnt, pdat.get_field<Tvec>(iaxyz_interf));
1359 }
1360
1361 sender_patch.get_field<Tvec>(ivxyz).append_subset_to(
1362 buf_idx, cnt, pdat.get_field<Tvec>(ivxyz_interf));
1363
1364 sender_omega.append_subset_to(buf_idx, cnt, pdat.get_field<Tscal>(iomega_interf));
1365
1366 if (has_soundspeed_field) {
1367 sender_patch.get_field<Tscal>(isoundspeed)
1368 .append_subset_to(buf_idx, cnt, pdat.get_field<Tscal>(isoundspeed_interf));
1369 }
1370
1371 if (has_B_field) {
1372 sender_patch.get_field<Tvec>(iB_on_rho).append_subset_to(
1373 buf_idx, cnt, pdat.get_field<Tvec>(iB_interf));
1374 }
1375
1376 if (has_psi_field) {
1377 sender_patch.get_field<Tscal>(ipsi_on_ch)
1378 .append_subset_to(buf_idx, cnt, pdat.get_field<Tscal>(ipsi_interf));
1379 }
1380
1381 if (has_curlB_field) {
1382 sender_patch.get_field<Tvec>(icurlB).append_subset_to(
1383 buf_idx, cnt, pdat.get_field<Tvec>(icurlB_interf));
1384 }
1385
1386 if (has_epsilon_field) {
1387 sender_patch.get_field<Tscal>(iepsilon).append_subset_to(
1388 buf_idx, cnt, pdat.get_field<Tscal>(iepsilon_interf));
1389 }
1390
1391 if (has_deltav_field) {
1392 sender_patch.get_field<Tvec>(ideltav).append_subset_to(
1393 buf_idx, cnt, pdat.get_field<Tvec>(ideltav_interf));
1394 }
1395 });
1396
1397 ghost_handle.template modify_interface_native<PatchDataLayer>(
1398 storage.ghost_patch_cache.get(),
1399 pdat_interf,
1400 [&](u64 sender,
1401 u64,
1402 InterfaceBuildInfos binfo,
1403 sham::DeviceBuffer<u32> &buf_idx,
1404 u32 cnt,
1405 PatchDataLayer &pdat) {
1406 if (sycl::length(binfo.offset_speed) > 0) {
1407 pdat.get_field<Tvec>(ivxyz_interf).apply_offset(binfo.offset_speed);
1408 }
1409 });
1410
1411 shambase::DistributedDataShared<PatchDataLayer> interf_pdat = ghost_handle.communicate_pdat(
1412 ghost_layout_ptr,
1413 std::move(pdat_interf),
1414 storage.exchange_gz_node,
1415 solver_config.show_ghost_zone_graph);
1416
1417 std::map<u64, u64> sz_interf_map;
1418 interf_pdat.for_each([&](u64 s, u64 r, PatchDataLayer &pdat_interf) {
1419 sz_interf_map[r] += pdat_interf.get_obj_cnt();
1420 });
1421
1422 storage.merged_patchdata_ghost.set(
1423 ghost_handle.template merge_native<PatchDataLayer, PatchDataLayer>(
1424 std::move(interf_pdat),
1426 PatchDataLayer pdat_new(ghost_layout_ptr);
1427
1428 u32 or_elem = pdat.get_obj_cnt();
1429 pdat_new.reserve(or_elem + sz_interf_map[p.id_patch]);
1430 u32 total_elements = or_elem;
1431
1432 PatchDataField<Tscal> &cur_omega = omega.get(p.id_patch);
1433
1434 pdat_new.get_field<Tscal>(ihpart_interf).insert(pdat.get_field<Tscal>(ihpart));
1435 pdat_new.get_field<Tscal>(iuint_interf).insert(pdat.get_field<Tscal>(iuint));
1436 pdat_new.get_field<Tvec>(ivxyz_interf).insert(pdat.get_field<Tvec>(ivxyz));
1437
1438 if (solver_config.has_axyz_in_ghost()) {
1439 pdat_new.get_field<Tvec>(iaxyz_interf).insert(pdat.get_field<Tvec>(iaxyz));
1440 }
1441
1442 pdat_new.get_field<Tscal>(iomega_interf).insert(cur_omega);
1443
1444 if (has_soundspeed_field) {
1445 pdat_new.get_field<Tscal>(isoundspeed_interf)
1446 .insert(pdat.get_field<Tscal>(isoundspeed));
1447 }
1448
1449 if (has_B_field) {
1450 pdat_new.get_field<Tvec>(iB_interf).insert(pdat.get_field<Tvec>(iB_on_rho));
1451 }
1452
1453 if (has_psi_field) {
1454 pdat_new.get_field<Tscal>(ipsi_interf)
1455 .insert(pdat.get_field<Tscal>(ipsi_on_ch));
1456 }
1457
1458 if (has_curlB_field) {
1459 pdat_new.get_field<Tvec>(icurlB_interf).insert(pdat.get_field<Tvec>(icurlB));
1460 }
1461
1462 if (has_epsilon_field) {
1463 pdat_new.get_field<Tscal>(iepsilon_interf)
1464 .insert(pdat.get_field<Tscal>(iepsilon));
1465 }
1466
1467 if (has_deltav_field) {
1468 pdat_new.get_field<Tvec>(ideltav_interf).insert(pdat.get_field<Tvec>(ideltav));
1469 }
1470
1471 pdat_new.check_field_obj_cnt_match();
1472
1473 return pdat_new;
1474 },
1475 [](PatchDataLayer &pdat, PatchDataLayer &pdat_interf) {
1476 pdat.insert_elements(pdat_interf);
1477 }));
1478
1479 timer_interf.end();
1480 storage.timings_details.interface += timer_interf.elasped_sec();
1481}
1482
1483template<class Tvec, template<class> class Kern>
1485 storage.merged_patchdata_ghost.reset();
1486}
1487
1489// start artificial viscosity section //////////////////////////////////////////////////////////////
1491
1492template<class Tvec, template<class> class Kern>
1494
1495 sph::modules::UpdateViscosity<Tvec, Kern>(context, solver_config, storage)
1496 .update_artificial_viscosity(dt);
1497}
1498
1500// end artificial viscosity section ////////////////////////////////////////////////////////////////
1502
1503template<class Tvec, template<class> class Kern>
1508
1509template<class Tvec, template<class> class Kern>
1511 shambase::get_check_ref(storage.pressure).free_alloc();
1512 shambase::get_check_ref(storage.soundspeed).free_alloc();
1513}
1514
1515template<class Tvec, template<class> class Kern>
1517
1518 StackEntry stack_loc{};
1519
1520 using namespace shamrock;
1521 using namespace shamrock::patch;
1522 shamrock::SchedulerUtility utility(scheduler());
1523 PatchDataLayerLayout &pdl = scheduler().pdl_old();
1524
1525 bool has_B_field = solver_config.has_field_B_on_rho();
1526 bool has_psi_field = solver_config.has_field_psi_on_ch();
1527 bool has_epsilon_field = solver_config.dust_config.has_epsilon_field();
1528 bool has_deltav_field = solver_config.dust_config.has_deltav_field();
1529
1530 const u32 iduint = pdl.get_field_idx<Tscal>("duint");
1531 const u32 iaxyz = pdl.get_field_idx<Tvec>("axyz");
1532 const u32 idB_on_rho = (has_B_field) ? pdl.get_field_idx<Tvec>("dB/rho") : 0;
1533 const u32 idpsi_on_ch = (has_psi_field) ? pdl.get_field_idx<Tscal>("dpsi/ch") : 0;
1534
1535 shamlog_debug_ln("sph::BasicGas", "save old fields");
1536 storage.old_axyz.set(utility.save_field<Tvec>(iaxyz, "axyz_old"));
1537 storage.old_duint.set(utility.save_field<Tscal>(iduint, "duint_old"));
1538
1539 if (has_B_field) {
1540 storage.old_dB_on_rho.set(utility.save_field<Tvec>(idB_on_rho, "dB/rho_old"));
1541 }
1542 if (has_psi_field) {
1543 storage.old_dpsi_on_ch.set(utility.save_field<Tscal>(idpsi_on_ch, "dpsi/ch_old"));
1544 }
1545 if (has_epsilon_field) {
1546 storage.old_dtepsilon.set(
1547 utility.save_field<Tscal>(pdl.get_field_idx<Tscal>("dtepsilon"), "dtepsilon_old"));
1548 }
1549 if (has_deltav_field) {
1550 storage.old_dtdeltav.set(
1551 utility.save_field<Tvec>(pdl.get_field_idx<Tvec>("dtdeltav"), "dtdeltav_old"));
1552 }
1553}
1554
1555template<class Tvec, template<class> class Kern>
1557
1558 modules::UpdateDerivs<Tvec, Kern> derivs(context, solver_config, storage);
1559 derivs.update_derivs();
1560
1561 modules::ExternalForces<Tvec, Kern> ext_forces(context, solver_config, storage);
1562 ext_forces.add_ext_forces();
1563}
1564
1565template<class Tvec, template<class> class Kern>
1567 return false;
1568}
1569
1570template<class Tvec, template<class> class Kern>
1572 modules::ComputeLoadBalanceValue<Tvec, Kern>(context, solver_config, storage)
1573 .update_load_balancing();
1574 scheduler().scheduler_step(false, false);
1575}
1576
1577template<class Tvec, template<class> class Kern>
1579
1580 // has to be first since there is a barrier that may mess the other timers
1581 shamsys::SystemMetrics system_metrics_start = shamsys::get_system_metrics();
1582
1583 sham::MemPerfInfos mem_perf_infos_start = sham::details::get_mem_perf_info();
1584 f64 mpi_timer_start = shamcomm::mpi::get_timer("total");
1585
1586 for (auto &callbacks : timestep_callbacks) {
1587 if (callbacks.step_begin_callback) {
1588 shambase::get_check_ref(callbacks.step_begin_callback)();
1589 }
1590 }
1591
1592 Tscal t_current = solver_config.get_time();
1593 Tscal dt = solver_config.get_dt_sph();
1594
1595 StackEntry stack_loc{};
1596
1597 if (shamcomm::world_rank() == 0) {
1598 shamcomm::logs::raw_ln(
1599 shambase::format("---------------- t = {}, dt = {} ----------------", t_current, dt));
1600 }
1601
1602 shambase::Timer tstep;
1603 tstep.start();
1604
1605 // if(shamcomm::world_rank() == 0) std::cout << scheduler().dump_status() << std::endl;
1606 modules::ComputeLoadBalanceValue<Tvec, Kern>(context, solver_config, storage)
1607 .update_load_balancing();
1608 scheduler().scheduler_step(true, true);
1609 modules::ComputeLoadBalanceValue<Tvec, Kern>(context, solver_config, storage)
1610 .update_load_balancing();
1611 // if(shamcomm::world_rank() == 0) std::cout << scheduler().dump_status() << std::endl;
1612 scheduler().scheduler_step(false, false);
1613 // if(shamcomm::world_rank() == 0) std::cout << scheduler().dump_status() << std::endl;
1614
1616
1617 using namespace shamrock;
1618 using namespace shamrock::patch;
1619
1620 bool has_B_field = solver_config.has_field_B_on_rho();
1621 bool has_psi_field = solver_config.has_field_psi_on_ch();
1622 bool has_epsilon_field = solver_config.dust_config.has_epsilon_field();
1623 bool has_deltav_field = solver_config.dust_config.has_deltav_field();
1624
1625 PatchDataLayerLayout &pdl = scheduler().pdl_old();
1626
1627 const u32 ixyz = pdl.get_field_idx<Tvec>("xyz");
1628 const u32 ivxyz = pdl.get_field_idx<Tvec>("vxyz");
1629 const u32 iaxyz = pdl.get_field_idx<Tvec>("axyz");
1630 const u32 iuint = pdl.get_field_idx<Tscal>("uint");
1631 const u32 iduint = pdl.get_field_idx<Tscal>("duint");
1632 const u32 ihpart = pdl.get_field_idx<Tscal>("hpart");
1633 const u32 iB_on_rho = (has_B_field) ? pdl.get_field_idx<Tvec>("B/rho") : 0;
1634 const u32 idB_on_rho = (has_B_field) ? pdl.get_field_idx<Tvec>("dB/rho") : 0;
1635 const u32 ipsi_on_ch = (has_psi_field) ? pdl.get_field_idx<Tscal>("psi/ch") : 0;
1636 const u32 idpsi_on_ch = (has_psi_field) ? pdl.get_field_idx<Tscal>("dpsi/ch") : 0;
1637 const u32 iepsilon = (has_epsilon_field) ? pdl.get_field_idx<Tscal>("epsilon") : 0;
1638 const u32 idtepsilon = (has_epsilon_field) ? pdl.get_field_idx<Tscal>("dtepsilon") : 0;
1639 const u32 ideltav = (has_deltav_field) ? pdl.get_field_idx<Tvec>("deltav") : 0;
1640 const u32 idtdeltav = (has_deltav_field) ? pdl.get_field_idx<Tvec>("dtdeltav") : 0;
1641
1642 shamrock::SchedulerUtility utility(scheduler());
1643
1644 modules::SinkParticlesUpdate<Tvec, Kern> sink_update(context, solver_config, storage);
1645 modules::ExternalForces<Tvec, Kern> ext_forces(context, solver_config, storage);
1646
1647 sink_update.accrete_particles(dt);
1648 ext_forces.point_mass_accrete_particles();
1649
1650 sink_update.predictor_step(dt);
1651
1652 {
1653 // beginning of SolverGraph migration
1654
1655 using namespace shamrock::solvergraph;
1656
1657 SolverGraph &solver_graph = storage.solver_graph;
1658
1659 // change the graph inputs
1660 {
1661 solver_graph.get_edge_ref<IDataEdge<Tscal>>("dt").data = dt;
1662 solver_graph.get_edge_ref<IDataEdge<Tscal>>("dt_half").data = dt / 2.0;
1663 }
1664
1666 // Solver evaluation
1668
1669 shambase::get_check_ref(storage.solver_sequence).evaluate();
1670 }
1671
1672 sink_update.compute_ext_forces();
1673
1674 ext_forces.compute_ext_forces_indep_v();
1675
1676 gen_serial_patch_tree();
1677
1678 apply_position_boundary(t_current + dt);
1679
1680 u64 Npart_all = scheduler().get_total_obj_count();
1681
1682 if (solver_config.enable_particle_reordering
1683 && solve_logs.step_count % solver_config.particle_reordering_step_freq == 0) {
1684 logger::info_ln("SPH", "Reordering particles at step ", solve_logs.step_count);
1685 modules::ParticleReordering<Tvec, u_morton, Kern>(context, solver_config, storage)
1687 }
1688
1689 {
1690 // update part counts and spans since particles have been moved and thus
1691 // new patch can be non-empty/empty
1692 using namespace shamrock::solvergraph;
1693 SolverGraph &solver_graph = storage.solver_graph;
1694 solver_graph.get_node_ref_base("attach fields to scheduler").evaluate();
1695 }
1696
1697 sph_prestep(t_current, dt);
1698
1700
1701 // Here we will add self grav to the external forces indep of vel (this will be moved into a
1702 // sperate module later)
1703 if (solver_config.self_grav_config.is_sg_on()) {
1704
1706
1709 constant_G.data = solver_config.get_constant_G();
1710 });
1711
1712 set_constant_G.set_edges(constant_G);
1713
1715
1717 [&](shamrock::solvergraph::FieldRefs<Tvec> &field_xyz_edge) {
1719 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
1720 auto &field = pdat.get_field<Tvec>(ixyz);
1721 field_xyz_refs.add_obj(p.id_patch, std::ref(field));
1722 });
1723 field_xyz_edge.set_refs(field_xyz_refs);
1724 });
1725 set_field_xyz.set_edges(field_xyz);
1726
1727 const u32 iaxyz_ext = pdl.get_field_idx<Tvec>("axyz_ext");
1728
1729 auto field_axyz_ext = shamrock::solvergraph::FieldRefs<Tvec>::make_shared("", "");
1730
1732 set_field_axyz_ext([&](shamrock::solvergraph::FieldRefs<Tvec> &field_axyz_ext_edge) {
1733 shamrock::solvergraph::DDPatchDataFieldRef<Tvec> field_axyz_ext_refs = {};
1734 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
1735 auto &field = pdat.get_field<Tvec>(iaxyz_ext);
1736 field_axyz_ext_refs.add_obj(p.id_patch, std::ref(field));
1737 });
1738 field_axyz_ext_edge.set_refs(field_axyz_ext_refs);
1739 });
1740 set_field_axyz_ext.set_edges(field_axyz_ext);
1741
1743
1746 sizes.indexes = {};
1747 scheduler().for_each_patchdata_nonempty([&](const Patch p, PatchDataLayer &pdat) {
1748 sizes.indexes.add_obj(p.id_patch, pdat.get_obj_cnt());
1749 });
1750 });
1751 set_sizes.set_edges(sizes);
1752
1754
1757 gpart_mass.data = solver_config.gpart_mass;
1758 });
1759
1760 set_gpart_mass.set_edges(gpart_mass);
1761
1762 set_gpart_mass.evaluate();
1763 set_constant_G.evaluate();
1764 set_field_xyz.evaluate();
1765 set_field_axyz_ext.evaluate();
1766 set_sizes.evaluate();
1767
1768 Tscal eps_grav = shambase::get_check_ref(
1769 std::get_if<SelfGravConfig::SofteningPlummer>(
1770 &solver_config.self_grav_config.softening_mode))
1771 .epsilon;
1772
1773 if (solver_config.self_grav_config.is_none()) {
1774 // do nothing
1775 } else if (solver_config.self_grav_config.is_direct()) {
1776
1778 std::get_if<SelfGravConfig::Direct>(&solver_config.self_grav_config.config));
1779
1780 modules::SGDirectPlummer<Tvec> self_gravity_direct_node(
1781 eps_grav, direct_config.reference_mode);
1782 self_gravity_direct_node.set_edges(
1783 sizes, gpart_mass, constant_G, field_xyz, field_axyz_ext);
1784 self_gravity_direct_node.evaluate();
1785
1786 } else if (solver_config.self_grav_config.is_mm()) {
1787
1789 std::get_if<SelfGravConfig::MM>(&solver_config.self_grav_config.config));
1790
1791 auto run_sg_mm = [&](auto mm_order_tag) {
1792 constexpr u32 order = decltype(mm_order_tag)::value;
1793 modules::SGMMPlummer<Tvec, order> self_gravity_mm_node(
1794 eps_grav, mm_config.opening_angle, mm_config.reduction_level);
1795 self_gravity_mm_node.set_edges(
1796 sizes, gpart_mass, constant_G, field_xyz, field_axyz_ext);
1797 self_gravity_mm_node.evaluate();
1798 };
1799
1800 switch (mm_config.order) {
1801 case 1 : run_sg_mm(std::integral_constant<u32, 1>{}); break;
1802 case 2 : run_sg_mm(std::integral_constant<u32, 2>{}); break;
1803 case 3 : run_sg_mm(std::integral_constant<u32, 3>{}); break;
1804 case 4 : run_sg_mm(std::integral_constant<u32, 4>{}); break;
1805 case 5 : run_sg_mm(std::integral_constant<u32, 5>{}); break;
1807 }
1808
1809 } else if (solver_config.self_grav_config.is_fmm()) {
1810
1812 std::get_if<SelfGravConfig::FMM>(&solver_config.self_grav_config.config));
1813
1814 auto run_sg_fmm = [&](auto fmm_order_tag) {
1815 constexpr u32 order = decltype(fmm_order_tag)::value;
1816 modules::SGFMMPlummer<Tvec, order> self_gravity_mm_node(
1817 eps_grav, fmm_config.opening_angle, fmm_config.reduction_level);
1818 self_gravity_mm_node.set_edges(
1819 sizes, gpart_mass, constant_G, field_xyz, field_axyz_ext);
1820 self_gravity_mm_node.evaluate();
1821 };
1822
1823 switch (fmm_config.order) {
1824 case 1 : run_sg_fmm(std::integral_constant<u32, 1>{}); break;
1825 case 2 : run_sg_fmm(std::integral_constant<u32, 2>{}); break;
1826 case 3 : run_sg_fmm(std::integral_constant<u32, 3>{}); break;
1827 case 4 : run_sg_fmm(std::integral_constant<u32, 4>{}); break;
1828 case 5 : run_sg_fmm(std::integral_constant<u32, 5>{}); break;
1830 }
1831
1832 } else if (solver_config.self_grav_config.is_sfmm()) {
1833
1835 std::get_if<SelfGravConfig::SFMM>(&solver_config.self_grav_config.config));
1836
1837 auto run_sg_sfmm = [&](auto sfmm_order_tag) {
1838 constexpr u32 order = decltype(sfmm_order_tag)::value;
1839 modules::SGSFMMPlummer<Tvec, order> self_gravity_mm_node(
1840 eps_grav,
1841 sfmm_config.opening_angle,
1842 sfmm_config.leaf_lowering,
1843 sfmm_config.reduction_level);
1844 self_gravity_mm_node.set_edges(
1845 sizes, gpart_mass, constant_G, field_xyz, field_axyz_ext);
1846 self_gravity_mm_node.evaluate();
1847 };
1848
1849 switch (sfmm_config.order) {
1850 case 1 : run_sg_sfmm(std::integral_constant<u32, 1>{}); break;
1851 case 2 : run_sg_sfmm(std::integral_constant<u32, 2>{}); break;
1852 case 3 : run_sg_sfmm(std::integral_constant<u32, 3>{}); break;
1853 case 4 : run_sg_sfmm(std::integral_constant<u32, 4>{}); break;
1854 case 5 : run_sg_sfmm(std::integral_constant<u32, 5>{}); break;
1856 }
1857
1858 } else {
1860 "Self gravity config not supported, current state is : \n"
1861 + nlohmann::json{solver_config.self_grav_config}.dump(4));
1862 }
1863 }
1864
1865 sph::BasicSPHGhostHandler<Tvec> &ghost_handle = storage.ghost_handler.get();
1866 auto &merged_xyzh = storage.merged_xyzh.get();
1867 shambase::DistributedData<RTree> &trees = storage.merged_pos_trees.get();
1868 // ComputeField<Tscal> &omega = storage.omega.get();
1869
1871 = shambase::get_check_ref(storage.ghost_layout.get());
1872 u32 ihpart_interf = ghost_layout.get_field_idx<Tscal>("hpart");
1873 u32 iuint_interf = ghost_layout.get_field_idx<Tscal>("uint");
1874 u32 ivxyz_interf = ghost_layout.get_field_idx<Tvec>("vxyz");
1875 u32 iomega_interf = ghost_layout.get_field_idx<Tscal>("omega");
1876 u32 iB_on_rho_interf = (has_B_field) ? ghost_layout.get_field_idx<Tvec>("B/rho") : 0;
1877 u32 ipsi_on_rho_interf = (has_psi_field) ? ghost_layout.get_field_idx<Tscal>("psi/ch") : 0;
1878
1879 using RTreeField = RadixTreeField<Tscal>;
1881
1882 Tscal next_cfl = 0;
1883
1884 u32 corrector_iter_cnt = 0;
1885 bool need_rerun_corrector = false;
1886 do {
1887
1888 reset_merge_ghosts_fields();
1889 reset_eos_fields();
1890
1891 if (corrector_iter_cnt == 50) {
1893 "the corrector has made over 50 loops, either their is a bug, either you are using "
1894 "a dt that is too large");
1895 }
1896
1897 // communicate fields
1898 communicate_merge_ghosts_fields();
1899
1900 if (solver_config.has_field_alphaAV()) {
1901
1902 std::shared_ptr<shamrock::solvergraph::PatchDataLayerRefs> patchdatas
1903 = std::make_shared<shamrock::solvergraph::PatchDataLayerRefs>(
1904 "patchdata_layer_ref", "patchdata_layer_ref");
1905
1906 auto node_set_edge = scheduler().get_node_set_edge_patchdata_layer_refs();
1907 node_set_edge->set_edges(patchdatas);
1908 node_set_edge->evaluate();
1909
1911 scheduler().get_layout_ptr_old(), "alpha_AV");
1912 node_copy.set_edges(patchdatas, storage.alpha_av_updated);
1913 node_copy.evaluate();
1914 }
1915
1916 if (solver_config.has_field_dtdivv()) {
1917
1918 if (solver_config.combined_dtdiv_divcurlv_compute) {
1919 if (solver_config.has_field_dtdivv()) {
1920 sph::modules::DiffOperatorDtDivv<Tvec, Kern>(context, solver_config, storage)
1921 .update_dtdivv(true);
1922 }
1923 } else {
1924
1925 if (solver_config.has_field_divv()) {
1926 sph::modules::DiffOperators<Tvec, Kern>(context, solver_config, storage)
1927 .update_divv();
1928 }
1929
1930 if (solver_config.has_field_curlv()) {
1931 sph::modules::DiffOperators<Tvec, Kern>(context, solver_config, storage)
1932 .update_curlv();
1933 }
1934
1935 if (solver_config.has_field_dtdivv()) {
1936 sph::modules::DiffOperatorDtDivv<Tvec, Kern>(context, solver_config, storage)
1937 .update_dtdivv(false);
1938 }
1939 }
1940
1941 } else {
1942 if (solver_config.has_field_divv()) {
1943 sph::modules::DiffOperators<Tvec, Kern>(context, solver_config, storage)
1944 .update_divv();
1945 }
1946
1947 if (solver_config.has_field_curlv()) {
1948 sph::modules::DiffOperators<Tvec, Kern>(context, solver_config, storage)
1949 .update_curlv();
1950 }
1951 }
1952
1953 // if (solver_config.has_field_divB()) {
1954 // sph::modules::DiffOperatorsB<Tvec, Kern>(context, solver_config, storage)
1955 // .update_divB();
1956 // }
1957
1958 // if (solver_config.has_field_curlB()) {
1959 // sph::modules::DiffOperatorsB<Tvec, Kern>(context, solver_config, storage)
1960 // .update_curlB();
1961 // }
1962 update_artificial_viscosity(dt);
1963
1964 if (solver_config.has_field_alphaAV()) {
1965
1967 = shambase::get_check_ref(storage.alpha_av_updated);
1968
1969 using InterfaceBuildInfos =
1971
1972 shambase::Timer time_interf;
1973 time_interf.start();
1974
1975 auto field_interf = ghost_handle.template build_interface_native<PatchDataField<Tscal>>(
1976 storage.ghost_patch_cache.get(),
1977 [&](u64 sender,
1978 u64 /*receiver*/,
1979 InterfaceBuildInfos binfo,
1980 sham::DeviceBuffer<u32> &buf_idx,
1981 u32 cnt) -> PatchDataField<Tscal> {
1982 PatchDataField<Tscal> &sender_field = comp_field_send.get_field(sender);
1983
1984 return sender_field.make_new_from_subset(buf_idx, cnt);
1985 });
1986
1988 = ghost_handle.communicate_pdatfield(
1989 std::move(field_interf), 1, storage.exchange_gz_alpha);
1990
1992 = ghost_handle.template merge_native<PatchDataField<Tscal>, PatchDataField<Tscal>>(
1993 std::move(interf_pdat),
1995 PatchDataField<Tscal> &receiver_field
1996 = comp_field_send.get_field(p.id_patch);
1997 return receiver_field.duplicate();
1998 },
1999 [](PatchDataField<Tscal> &mpdat, PatchDataField<Tscal> &pdat_interf) {
2000 mpdat.insert(pdat_interf);
2001 });
2002
2003 time_interf.end();
2004 storage.timings_details.interface += time_interf.elasped_sec();
2005
2006 storage.alpha_av_ghost.set(std::move(merged_field));
2007 }
2008
2009 // compute pressure
2010 compute_eos_fields();
2011
2012 constexpr bool debug_interfaces = false;
2013 if constexpr (debug_interfaces) {
2014
2015 if (solver_config.do_debug_dump) {
2016
2018 = storage.merged_patchdata_ghost.get();
2019
2020 scheduler().for_each_patchdata_nonempty([&](Patch cur_p, PatchDataLayer &pdat) {
2021 MergedPatchData &merged_patch = mpdat.get(cur_p.id_patch);
2022 PatchDataLayer &mpdat = merged_patch.pdat;
2023
2024 sycl::buffer<Tvec> &buf_xyz = shambase::get_check_ref(
2025 merged_xyzh.get(cur_p.id_patch).field_pos.get_buf());
2026 sycl::buffer<Tvec> &buf_vxyz = mpdat.get_field_buf_ref<Tvec>(ivxyz_interf);
2027 sycl::buffer<Tscal> &buf_hpart = mpdat.get_field_buf_ref<Tscal>(ihpart_interf);
2028
2029 u32 total_elements = shambase::get_check_ref(storage.part_counts_with_ghost)
2030 .indexes.get(cur_p.id_patch);
2031 SHAM_ASSERT(merged_patch.total_elements == total_elements);
2032
2034 total_elements,
2035 solver_config.gpart_mass,
2036
2037 buf_xyz,
2038 buf_hpart,
2039 buf_vxyz};
2040
2041 make_interface_debug_phantom_dump(info).gen_file().write_to_file(
2042 solver_config.debug_dump_filename);
2043 logger::raw_ln("writing : ", solver_config.debug_dump_filename);
2044 });
2045 }
2046 }
2047
2048 // compute force
2049 shamlog_debug_ln("sph::BasicGas", "compute force");
2050
2051 // save old acceleration
2052 prepare_corrector();
2053
2054 update_derivs();
2055
2056 bool has_luminosity = solver_config.compute_luminosity;
2057
2058 if (has_luminosity) {
2059 const u32 iluminosity = pdl.get_field_idx<Tscal>("luminosity");
2060
2061 shambase::get_check_ref(storage.hpart_with_ghosts)
2062 .set_refs(storage.merged_xyzh.get()
2063 .template map<std::reference_wrapper<PatchDataField<Tscal>>>(
2064 [&](u64 id, shamrock::patch::PatchDataLayer &mpdat) {
2065 return std::ref(mpdat.get_field<Tscal>(
2066 1)); // hpart is at index 1 in merged_xyzh
2067 }));
2068
2069 auto uint_with_ghost = shamrock::solvergraph::FieldRefs<Tscal>::make_shared("", "");
2070
2071 shambase::get_check_ref(storage.hpart_with_ghosts)
2072 .set_refs(storage.merged_xyzh.get()
2073 .template map<std::reference_wrapper<PatchDataField<Tscal>>>(
2074 [&](u64 id, shamrock::patch::PatchDataLayer &mpdat) {
2075 return std::ref(mpdat.get_field<Tscal>(1));
2076 }));
2077
2079 set_uint_with_ghost_refs(
2080 [&](shamrock::solvergraph::FieldRefs<Tscal> &field_uint_with_ghost_edge) {
2082 = storage.merged_patchdata_ghost.get();
2083
2084 shamrock::solvergraph::DDPatchDataFieldRef<Tscal> field_uint_with_ghost_refs
2085 = {};
2086
2087 scheduler().for_each_patchdata_nonempty(
2088 [&](const Patch p, PatchDataLayer &pdat) {
2089 PatchDataLayer &mpdat = mpdats.get(p.id_patch);
2090
2091 auto &field = mpdat.get_field<Tscal>(iuint_interf);
2092 field_uint_with_ghost_refs.add_obj(p.id_patch, std::ref(field));
2093 });
2094
2095 field_uint_with_ghost_edge.set_refs(field_uint_with_ghost_refs);
2096 });
2097
2098 set_uint_with_ghost_refs.set_edges(uint_with_ghost);
2099
2101
2103 set_luminosity_refs(
2104 [&](shamrock::solvergraph::FieldRefs<Tscal> &field_luminosity_edge) {
2106 = storage.merged_patchdata_ghost.get();
2107
2109 = {};
2110
2111 scheduler().for_each_patchdata_nonempty(
2112 [&](const Patch p, PatchDataLayer &pdat) {
2113 auto &field = pdat.get_field<Tscal>(iluminosity);
2114 field_luminosity_refs.add_obj(p.id_patch, std::ref(field));
2115 });
2116 field_luminosity_edge.set_refs(field_luminosity_refs);
2117 });
2118
2119 set_luminosity_refs.set_edges(luminosity);
2120
2121 set_uint_with_ghost_refs.evaluate();
2122 set_luminosity_refs.evaluate();
2123
2124 Tscal alpha_u = solver_config.artif_viscosity.get_alpha_u().value();
2125
2127 solver_config.gpart_mass, alpha_u};
2128
2129 compute_luminosity.set_edges(
2130 storage.part_counts,
2131 storage.neigh_cache,
2132 storage.positions_with_ghosts,
2133 storage.hpart_with_ghosts,
2134 storage.omega,
2135 uint_with_ghost,
2136 storage.pressure,
2137 luminosity);
2138
2139 compute_luminosity.evaluate();
2140 }
2141
2142 modules::ConservativeCheck<Tvec, Kern> cv_check(context, solver_config, storage);
2143 cv_check.check_conservation();
2144
2145 ComputeField<Tscal> vepsilon_v_sq
2146 = utility.make_compute_field<Tscal>("vmean epsilon_v^2", 1);
2147 ComputeField<Tscal> uepsilon_u_sq
2148 = utility.make_compute_field<Tscal>("umean epsilon_u^2", 1);
2149
2150 // corrector
2151 shamlog_debug_ln("sph::BasicGas", "leapfrog corrector");
2152 utility.fields_leapfrog_corrector<Tvec>(
2153 ivxyz, iaxyz, storage.old_axyz.get(), vepsilon_v_sq, dt / 2);
2154 utility.fields_leapfrog_corrector<Tscal>(
2155 iuint, iduint, storage.old_duint.get(), uepsilon_u_sq, dt / 2);
2156
2157 if (solver_config.has_field_B_on_rho()) {
2158 ComputeField<Tscal> BOR_epsilon_BOR_sq
2159 = utility.make_compute_field<Tscal>("B/rho epsilon_B/rho^2", 1);
2160 utility.fields_leapfrog_corrector<Tvec>(
2161 iB_on_rho, idB_on_rho, storage.old_dB_on_rho.get(), BOR_epsilon_BOR_sq, dt / 2);
2162 }
2163 if (solver_config.has_field_B_on_rho()) {
2164 ComputeField<Tscal> POC_epsilon_POC_sq
2165 = utility.make_compute_field<Tscal>("psi/ch epsilon_psi/ch^2", 1);
2166 utility.fields_leapfrog_corrector<Tscal>(
2167 ipsi_on_ch, idpsi_on_ch, storage.old_dpsi_on_ch.get(), POC_epsilon_POC_sq, dt / 2);
2168 }
2169
2170 if (solver_config.dust_config.has_epsilon_field()) {
2171 ComputeField<Tscal> epsilon_epsilon_sq
2172 = utility.make_compute_field<Tscal>("epsilon epsilon^2", 1);
2173 utility.fields_leapfrog_corrector<Tscal>(
2174 iepsilon, idtepsilon, storage.old_dtepsilon.get(), epsilon_epsilon_sq, dt / 2);
2175 }
2176
2177 if (solver_config.dust_config.has_deltav_field()) {
2178 ComputeField<Tscal> epsilon_deltav_sq
2179 = utility.make_compute_field<Tscal>("deltav deltav^2", 1);
2180 utility.fields_leapfrog_corrector<Tvec>(
2181 ideltav, idtdeltav, storage.old_dtdeltav.get(), epsilon_deltav_sq, dt / 2);
2182 }
2183
2184 storage.old_axyz.reset();
2185 storage.old_duint.reset();
2186 if (solver_config.has_field_B_on_rho()) {
2187 storage.old_dB_on_rho.reset();
2188 }
2189 if (solver_config.has_field_B_on_rho()) {
2190 storage.old_dpsi_on_ch.reset();
2191 }
2192
2193 if (solver_config.dust_config.has_epsilon_field()) {
2194 storage.old_dtepsilon.reset();
2195 }
2196
2197 if (solver_config.dust_config.has_deltav_field()) {
2198 storage.old_dtdeltav.reset();
2199 }
2200
2201 Tscal rank_veps_v = sycl::sqrt(vepsilon_v_sq.compute_rank_max());
2203 // compute means //////////////////////////
2205
2206 Tscal sum_vsq = utility.compute_rank_dot_sum<Tvec>(ivxyz);
2207
2208 Tscal vmean_sq = shamalgs::collective::allreduce_sum(sum_vsq) / Tscal(Npart_all);
2209
2210 Tscal vmean = sycl::sqrt(vmean_sq);
2211
2212 Tscal rank_eps_v = rank_veps_v / vmean;
2213
2214 if (vmean <= 0) {
2215 rank_eps_v = 0;
2216 }
2217
2218 Tscal eps_v = shamalgs::collective::allreduce_max(rank_eps_v);
2219
2220 shamlog_debug_ln("BasicGas", "epsilon v :", eps_v);
2221
2222 if (eps_v > 1e-2) {
2223 if (shamcomm::world_rank() == 0) {
2224 logger::warn_ln(
2225 "BasicGasSPH",
2226 shambase::format(
2227 "the corrector tolerance are broken the step will "
2228 "be re rerunned\n eps_v = {}",
2229 eps_v));
2230 }
2231 need_rerun_corrector = true;
2232 solver_config.time_state.cfl_multiplier /= 2;
2233
2234 // logger::info_ln("rerun corrector ...");
2235 } else {
2236 need_rerun_corrector = false;
2237 }
2238
2239 if (!need_rerun_corrector) {
2240
2241 sink_update.corrector_step(dt);
2242
2243 // write back alpha av field
2244 if (solver_config.has_field_alphaAV()) {
2245
2246 const u32 ialpha_AV = pdl.get_field_idx<Tscal>("alpha_AV");
2247 shamrock::solvergraph::Field<Tscal> &alpha_av_updated
2248 = shambase::get_check_ref(storage.alpha_av_updated);
2249
2250 scheduler().for_each_patchdata_nonempty([&](Patch cur_p, PatchDataLayer &pdat) {
2251 sham::DeviceBuffer<Tscal> &buf_alpha_av
2252 = pdat.get_field<Tscal>(ialpha_AV).get_buf();
2253 sham::DeviceBuffer<Tscal> &buf_alpha_av_updated
2254 = alpha_av_updated.get_field(cur_p.id_patch).get_buf();
2255
2256 auto &q = shamsys::instance::get_compute_scheduler().get_queue();
2257 sham::EventList depends_list;
2258
2259 auto alpha_av = buf_alpha_av.get_write_access(depends_list);
2260 auto alpha_av_updated = buf_alpha_av_updated.get_read_access(depends_list);
2261
2262 auto e = q.submit(depends_list, [&](sycl::handler &cgh) {
2263 shambase::parallel_for(
2264 cgh, pdat.get_obj_cnt(), "write back alpha_av", [=](i32 id_a) {
2265 alpha_av[id_a] = alpha_av_updated[id_a];
2266 });
2267 });
2268
2269 buf_alpha_av.complete_event_state(e);
2270 buf_alpha_av_updated.complete_event_state(e);
2271 });
2272 }
2273
2274 shamlog_debug_ln("BasicGas", "computing next CFL");
2275
2276 ComputeField<Tscal> vsig_max_dt = utility.make_compute_field<Tscal>("vsig_a", 1);
2277 std::unique_ptr<ComputeField<Tscal>> vclean_dt;
2278 if (has_psi_field) {
2279 vclean_dt = std::make_unique<ComputeField<Tscal>>(
2280 utility.make_compute_field<Tscal>("vclean_a", 1));
2281 }
2282
2284 = storage.merged_patchdata_ghost.get();
2285
2286 scheduler().for_each_patchdata_nonempty([&](Patch cur_p, PatchDataLayer &pdat) {
2287 PatchDataLayer &mpdat = mpdats.get(cur_p.id_patch);
2288
2290 = merged_xyzh.get(cur_p.id_patch).template get_field_buf_ref<Tvec>(0);
2291 sham::DeviceBuffer<Tvec> &buf_vxyz = mpdat.get_field_buf_ref<Tvec>(ivxyz_interf);
2292 sham::DeviceBuffer<Tscal> &buf_hpart
2293 = mpdat.get_field_buf_ref<Tscal>(ihpart_interf);
2294 sham::DeviceBuffer<Tscal> &buf_uint = mpdat.get_field_buf_ref<Tscal>(iuint_interf);
2295 sham::DeviceBuffer<Tscal> &buf_pressure
2296 = shambase::get_check_ref(storage.pressure).get_field(cur_p.id_patch).get_buf();
2297 sham::DeviceBuffer<Tscal> &cs_buf = shambase::get_check_ref(storage.soundspeed)
2298 .get_field(cur_p.id_patch)
2299 .get_buf();
2300
2301 sham::DeviceBuffer<Tscal> &vsig_buf = vsig_max_dt.get_buf_check(cur_p.id_patch);
2302
2303 sycl::range range_npart{pdat.get_obj_cnt()};
2304
2305 tree::ObjectCache &pcache
2306 = shambase::get_check_ref(storage.neigh_cache).get_cache(cur_p.id_patch);
2307
2309
2310 {
2311
2312 auto &q = shamsys::instance::get_compute_scheduler().get_queue();
2313 sham::EventList depends_list;
2314
2315 auto xyz = buf_xyz.get_read_access(depends_list);
2316 auto vxyz = buf_vxyz.get_read_access(depends_list);
2317 auto hpart = buf_hpart.get_read_access(depends_list);
2318 auto u = buf_uint.get_read_access(depends_list);
2319 auto pressure = buf_pressure.get_read_access(depends_list);
2320 auto cs = cs_buf.get_read_access(depends_list);
2321 auto vsig = vsig_buf.get_write_access(depends_list);
2322 auto particle_looper_ptrs = pcache.get_read_access(depends_list);
2323
2324 NamedStackEntry tmppp{"compute vsig"};
2325 auto e = q.submit(depends_list, [&](sycl::handler &cgh) {
2326 const Tscal pmass = solver_config.gpart_mass;
2327 const Tscal alpha_u = 1.0;
2328 const Tscal alpha_AV = 1.0;
2329 const Tscal beta_AV = 2.0;
2330
2331 tree::ObjectCacheIterator particle_looper(particle_looper_ptrs);
2332
2333 constexpr Tscal Rker2 = Kernel::Rkern * Kernel::Rkern;
2334
2335 shambase::parallel_for(
2336 cgh, pdat.get_obj_cnt(), "compute vsig", [=](i32 id_a) {
2337 using namespace shamrock::sph;
2338
2339 Tvec sum_axyz = {0, 0, 0};
2340 Tscal sum_du_a = 0;
2341 Tscal h_a = hpart[id_a];
2342
2343 Tvec xyz_a = xyz[id_a];
2344 Tvec vxyz_a = vxyz[id_a];
2345
2346 Tscal rho_a = rho_h(pmass, h_a, Kernel::hfactd);
2347 Tscal rho_a_sq = rho_a * rho_a;
2348 Tscal rho_a_inv = 1. / rho_a;
2349
2350 Tscal P_a = pressure[id_a];
2351
2352 const Tscal u_a = u[id_a];
2353
2354 Tscal cs_a = cs[id_a];
2355
2356 Tscal vsig_max = 0;
2357
2358 particle_looper.for_each_object(id_a, [&](u32 id_b) {
2359 // compute only omega_a
2360 Tvec dr = xyz_a - xyz[id_b];
2361 Tscal rab2 = sycl::dot(dr, dr);
2362 Tscal h_b = hpart[id_b];
2363
2364 if (rab2 > h_a * h_a * Rker2 && rab2 > h_b * h_b * Rker2) {
2365 return;
2366 }
2367
2368 Tscal rab = sycl::sqrt(rab2);
2369 Tvec vxyz_b = vxyz[id_b];
2370 Tvec v_ab = vxyz_a - vxyz_b;
2371 const Tscal u_b = u[id_b];
2372
2373 Tvec r_ab_unit = dr / rab;
2374
2375 if (rab < 1e-9) {
2376 r_ab_unit = {0, 0, 0};
2377 }
2378
2379 Tscal rho_b = rho_h(pmass, h_b, Kernel::hfactd);
2380 Tscal P_b = pressure[id_b];
2381 Tscal cs_b = cs[id_b];
2382 Tscal v_ab_r_ab = sycl::dot(v_ab, r_ab_unit);
2383 Tscal abs_v_ab_r_ab = sycl::fabs(v_ab_r_ab);
2384
2386 // internal energy update
2387 // scalar : f32 | vector : f32_3
2388 const Tscal alpha_a = alpha_AV;
2389 const Tscal alpha_b = alpha_AV;
2390
2391 Tscal vsig_a = alpha_a * cs_a + beta_AV * abs_v_ab_r_ab;
2392
2393 vsig_max = sycl::fmax(vsig_max, vsig_a);
2394 });
2395
2396 vsig[id_a] = vsig_max;
2397 });
2398 });
2399
2400 if (has_psi_field) {
2401 NamedStackEntry tmppp{"compute vclean"};
2402 Tscal const mu_0 = solver_config.get_constant_mu_0();
2403 sham::DeviceBuffer<Tscal> &vclean_buf
2404 = shambase::get_check_ref(vclean_dt).get_buf_check(cur_p.id_patch);
2405
2406 Tvec *B_on_rho = mpdat.get_field_buf_ref<Tvec>(iB_on_rho_interf)
2407 .get_write_access(depends_list);
2408
2409 auto vclean = vclean_buf.get_write_access(depends_list);
2410
2411 auto e = q.submit(depends_list, [&](sycl::handler &cgh) {
2412 const Tscal pmass = solver_config.gpart_mass;
2413
2414 tree::ObjectCacheIterator particle_looper(particle_looper_ptrs);
2415
2416 constexpr Tscal Rker2 = Kernel::Rkern * Kernel::Rkern;
2417
2418 shambase::parallel_for(
2419 cgh, pdat.get_obj_cnt(), "compute vclean", [=](i32 id_a) {
2420 using namespace shamrock::sph;
2421
2422 Tscal h_a = hpart[id_a];
2423 Tscal rho_a = rho_h(pmass, h_a, Kernel::hfactd);
2424 const Tscal u_a = u[id_a];
2425 Tscal cs_a = cs[id_a];
2426 Tvec B_a = B_on_rho[id_a] * rho_a;
2427
2428 Tscal vclean_a = shamphys::MHD_physics<Tvec, Tscal>::v_shock(
2429 cs_a, B_a, rho_a, mu_0);
2430
2431 vclean[id_a] = vclean_a;
2432 });
2433 });
2434 mpdat.get_field_buf_ref<Tvec>(iB_on_rho_interf).complete_event_state(e);
2435 vclean_buf.complete_event_state(e);
2436 };
2437
2438 buf_xyz.complete_event_state(e);
2439 buf_vxyz.complete_event_state(e);
2440 buf_hpart.complete_event_state(e);
2441 buf_uint.complete_event_state(e);
2442 buf_pressure.complete_event_state(e);
2443 cs_buf.complete_event_state(e);
2444 vsig_buf.complete_event_state(e);
2445
2446 sham::EventList resulting_events;
2447 resulting_events.add_event(e);
2448 pcache.complete_event_state(resulting_events);
2449 }
2450 });
2451
2452 ComputeField<Tscal> cfl_dt = utility.make_compute_field<Tscal>("cfl_dt", 1);
2453
2454 scheduler().for_each_patchdata_nonempty([&](Patch cur_p, PatchDataLayer &pdat) {
2455 PatchDataLayer &mpdat = mpdats.get(cur_p.id_patch);
2456
2457 sham::DeviceBuffer<Tvec> &buf_axyz = pdat.get_field<Tvec>(iaxyz).get_buf();
2458 sham::DeviceBuffer<Tscal> &buf_hpart
2459 = mpdat.get_field<Tscal>(ihpart_interf).get_buf();
2460 sham::DeviceBuffer<Tscal> &vsig_buf = vsig_max_dt.get_buf_check(cur_p.id_patch);
2461 sham::DeviceBuffer<Tscal> &cfl_dt_buf = cfl_dt.get_buf_check(cur_p.id_patch);
2462
2463 auto &q = shamsys::instance::get_compute_scheduler().get_queue();
2464 sham::EventList depends_list;
2465
2466 auto hpart = buf_hpart.get_read_access(depends_list);
2467 auto a = buf_axyz.get_read_access(depends_list);
2468 auto vsig = vsig_buf.get_read_access(depends_list);
2469 auto cfl_dt = cfl_dt_buf.get_write_access(depends_list);
2470
2471 auto e = q.submit(depends_list, [&](sycl::handler &cgh) {
2472 Tscal C_cour = solver_config.cfl_config.cfl_cour
2473 * solver_config.time_state.cfl_multiplier;
2474 Tscal C_force = solver_config.cfl_config.cfl_force
2475 * solver_config.time_state.cfl_multiplier;
2476
2477 cgh.parallel_for(sycl::range<1>{pdat.get_obj_cnt()}, [=](sycl::item<1> item) {
2478 Tscal h_a = hpart[item];
2479 Tscal vsig_a = vsig[item];
2480 Tscal abs_a_a = sycl::length(a[item]);
2481
2482 Tscal dt_c = C_cour * h_a / vsig_a;
2483 Tscal dt_f = C_force * sycl::sqrt(h_a / abs_a_a);
2484
2485 cfl_dt[item] = sycl::min(dt_c, dt_f);
2486 });
2487 });
2488
2489 if (has_psi_field) {
2490 sham::DeviceBuffer<Tscal> &vclean_buf
2491 = shambase::get_check_ref(vclean_dt).get_buf_check(cur_p.id_patch);
2492 auto vclean = vclean_buf.get_read_access(depends_list);
2493 auto e = q.submit(depends_list, [&](sycl::handler &cgh) {
2494 Tscal C_cour = solver_config.cfl_config.cfl_cour
2495 * solver_config.time_state.cfl_multiplier;
2496
2497 cgh.parallel_for(
2498 sycl::range<1>{pdat.get_obj_cnt()}, [=](sycl::item<1> item) {
2499 Tscal h_a = hpart[item];
2500 Tscal vclean_a = vclean[item];
2501
2502 Tscal dt_divB_cleaning = C_cour * h_a / vclean_a;
2503
2504 cfl_dt[item] = sycl::min(cfl_dt[item], dt_divB_cleaning);
2505 });
2506 });
2507 vclean_buf.complete_event_state(e);
2508 };
2509
2510 buf_hpart.complete_event_state(e);
2511 buf_axyz.complete_event_state(e);
2512 vsig_buf.complete_event_state(e);
2513 cfl_dt_buf.complete_event_state(e);
2514 });
2515
2516 Tscal rank_dt = cfl_dt.compute_rank_min();
2517
2518 if (solver_config.should_save_dt_to_fields()) {
2519
2520 const u32 idt_part = pdl.get_field_idx<Tscal>("dt_part");
2521
2522 scheduler().for_each_patchdata_nonempty([&](Patch cur_p, PatchDataLayer &pdat) {
2523 sham::DeviceBuffer<Tscal> &buf_dt_part
2524 = pdat.get_field_buf_ref<Tscal>(idt_part);
2525 sham::DeviceBuffer<Tscal> &buf_dt = cfl_dt.get_buf_check(cur_p.id_patch);
2526
2527 buf_dt_part.copy_from(buf_dt, pdat.get_obj_cnt());
2528 });
2529 }
2530
2531 Tscal sink_sink_cfl = shambase::get_infty<Tscal>();
2532 if (!storage.sinks.is_empty()) {
2533 // sink sink CFL
2534
2535 Tscal G = solver_config.get_constant_G();
2536
2537 Tscal C_force
2538 = solver_config.cfl_config.cfl_force * solver_config.time_state.cfl_multiplier;
2539 Tscal eta_phi = solver_config.cfl_config.eta_sink;
2540
2541 std::vector<SinkParticle<Tvec>> &sink_parts = storage.sinks.get();
2542
2543 for (u32 i = 0; i < sink_parts.size(); i++) {
2544 SinkParticle<Tvec> &s_i = sink_parts[i];
2545 Tscal sink_sink_cfl_i = shambase::get_infty<Tscal>();
2546
2547 Tvec f_i = s_i.ext_acceleration;
2548
2549 Tscal grad_phi_i_sq = sham::dot(f_i, f_i); // m^2.s^-4
2550
2551 if (grad_phi_i_sq == 0) {
2552 continue;
2553 }
2554
2555 for (u32 j = 0; j < sink_parts.size(); j++) {
2556 SinkParticle<Tvec> &s_j = sink_parts[j];
2557
2558 if (i == j) {
2559 continue;
2560 }
2561
2562 Tvec rij = s_i.pos - s_j.pos;
2563 Tscal rij_scal = sycl::length(rij);
2564
2565 Tscal phi_ij = G * s_j.mass / rij_scal; // J / kg = m^2.s^-2
2566 Tscal term_ij = sham::abs(phi_ij) / grad_phi_i_sq; // s^2
2567 Tscal dt_ij = C_force * eta_phi * sycl::sqrt(term_ij); // s
2568
2569 sink_sink_cfl_i = sham::min(sink_sink_cfl_i, dt_ij);
2570 }
2571
2572 sink_sink_cfl = sham::min(sink_sink_cfl, sink_sink_cfl_i);
2573 }
2574
2575 sink_sink_cfl = shamalgs::collective::allreduce_min(sink_sink_cfl);
2576 }
2577
2578 shamlog_debug_ln("BasigGas", "rank", shamcomm::world_rank(), "found cfl dt =", rank_dt);
2579
2580 Tscal hydro_cfl = shamalgs::collective::allreduce_min(rank_dt);
2581
2582 if (shamcomm::world_rank() == 0) {
2583 shamlog_info_ln("SPH", "CFL hydro =", hydro_cfl, "sink sink =", sink_sink_cfl);
2584 }
2585
2586 next_cfl = sham::min(hydro_cfl, sink_sink_cfl);
2587
2588 if (shamcomm::world_rank() == 0) {
2589 logger::info_ln(
2590 "sph::Model",
2591 "cfl dt =",
2592 next_cfl,
2593 "cfl multiplier :",
2594 solver_config.time_state.cfl_multiplier);
2595 }
2596
2597 // this should not be needed idealy, but we need the pressure on the ghosts and
2598 // we don't want to communicate it as it can be recomputed from the other fields
2599 // hence we copy the soundspeed at the end of the step to a field in the patchdata
2600 if (solver_config.has_field_soundspeed()) {
2601
2602 const u32 isoundspeed = pdl.get_field_idx<Tscal>("soundspeed");
2603
2604 scheduler().for_each_patchdata_nonempty([&](Patch cur_p, PatchDataLayer &pdat) {
2605 sham::DeviceBuffer<Tscal> &buf_cs = pdat.get_field_buf_ref<Tscal>(isoundspeed);
2606 sham::DeviceBuffer<Tscal> &buf_cs_in
2607 = shambase::get_check_ref(storage.soundspeed)
2608 .get_field(cur_p.id_patch)
2609 .get_buf();
2610
2611 sycl::range range_npart{pdat.get_obj_cnt()};
2612
2614
2615 auto &q = shamsys::instance::get_compute_scheduler().get_queue();
2616 sham::EventList depends_list;
2617
2618 auto cs_in = buf_cs_in.get_read_access(depends_list);
2619 auto cs = buf_cs.get_write_access(depends_list);
2620
2621 auto e = q.submit(depends_list, [&](sycl::handler &cgh) {
2622 const Tscal pmass = solver_config.gpart_mass;
2623
2624 cgh.parallel_for(
2625 sycl::range<1>{pdat.get_obj_cnt()}, [=](sycl::item<1> item) {
2626 cs[item] = cs_in[item];
2627 });
2628 });
2629
2630 buf_cs_in.complete_event_state(e);
2631 buf_cs.complete_event_state(e);
2632 });
2633 }
2634
2635 } // if (!need_rerun_corrector) {
2636
2637 corrector_iter_cnt++;
2638
2639 if (solver_config.has_field_alphaAV()) {
2640 storage.alpha_av_ghost.reset();
2641 }
2642
2643 } while (need_rerun_corrector);
2644
2645 reset_merge_ghosts_fields();
2646 reset_eos_fields();
2647
2648 // if delta too big jump to compute force
2649
2650 tstep.end();
2651
2652 for (auto it = timestep_callbacks.rbegin(); it != timestep_callbacks.rend(); ++it) {
2653 if (it->step_end_callback) {
2654 shambase::get_check_ref(it->step_end_callback)();
2655 }
2656 }
2657
2658 f64 delta_mpi_timer = shamcomm::mpi::get_timer("total") - mpi_timer_start;
2659 sham::MemPerfInfos mem_perf_infos_end = sham::details::get_mem_perf_info();
2660
2662 shamsys::SystemMetrics system_metrics_end = shamsys::get_system_metrics();
2663 shamsys::SystemMetrics system_metrics_delta = system_metrics_end - system_metrics_start;
2664
2665 f64 t_dev_alloc
2666 = (mem_perf_infos_end.time_alloc_device - mem_perf_infos_start.time_alloc_device)
2667 + (mem_perf_infos_end.time_free_device - mem_perf_infos_start.time_free_device);
2668 f64 t_host_alloc = (mem_perf_infos_end.time_alloc_host - mem_perf_infos_start.time_alloc_host)
2669 + (mem_perf_infos_end.time_free_host - mem_perf_infos_start.time_free_host);
2670
2671 u64 rank_count = scheduler().get_rank_count();
2672 f64 rate = f64(rank_count) / tstep.elasped_sec();
2673
2674 u64 npatch = scheduler().patch_list.local.size();
2675
2676 // logger::info_ln("SPHSolver", "process rate : ", rate, "particle.s-1");
2677
2678 std::string log_step = report_perf_timestep(
2679 rate,
2680 rank_count,
2681 npatch,
2682 tstep.elasped_sec(),
2683 delta_mpi_timer,
2684 t_dev_alloc,
2685 t_host_alloc,
2686 mem_perf_infos_end.max_allocated_byte_device,
2687 mem_perf_infos_end.max_allocated_byte_host,
2688 system_metrics_delta,
2689 shamsys::has_reporter());
2690
2691 if (shamcomm::world_rank() == 0) {
2692 logger::info_ln("sph::Model", log_step);
2693 logger::info_ln(
2694 "sph::Model", "estimated rate :", dt * (3600 / tstep.elasped_sec()), "(tsim/hr)");
2695 }
2696
2697 solve_logs.register_log(
2698 {t_current, // f64 solver_t;
2699 dt, // f64 solver_dt;
2700 shamcomm::world_rank(), // i32 world_rank;
2701 rank_count, // u64 rank_count;
2702 rate, // f64 rate;
2703 tstep.elasped_sec(), // f64 elasped_sec;
2705 system_metrics_delta});
2706
2707 storage.timings_details.reset();
2708
2709 reset_serial_patch_tree();
2710 reset_ghost_handler();
2711
2712 shambase::get_check_ref(storage.part_counts).free_alloc();
2713 shambase::get_check_ref(storage.part_counts_with_ghost).free_alloc();
2714 shambase::get_check_ref(storage.positions_with_ghosts).free_alloc();
2715 shambase::get_check_ref(storage.hpart_with_ghosts).free_alloc();
2716 storage.merged_xyzh.reset();
2717 shambase::get_check_ref(storage.omega).free_alloc();
2718 clear_merged_pos_trees();
2719 clear_ghost_cache();
2720 reset_presteps_rint();
2721 reset_neighbors_cache();
2722
2723 shambase::get_check_ref(storage.neigh_cache).free_alloc();
2724
2725 solver_config.set_next_dt(next_cfl);
2726 solver_config.set_time(t_current + dt);
2727
2728 auto get_next_cfl_mult = [&]() {
2729 Tscal cfl_m = solver_config.time_state.cfl_multiplier;
2730 Tscal stiff = solver_config.cfl_config.cfl_multiplier_stiffness;
2731
2732 return (cfl_m * stiff + 1.) / (stiff + 1.);
2733 };
2734
2735 solver_config.time_state.cfl_multiplier = get_next_cfl_mult();
2736
2737 TimestepLog log;
2738 log.rank = shamcomm::world_rank();
2739 log.rate = rate;
2740 log.npart = rank_count;
2741 log.tcompute = tstep.elasped_sec();
2742
2743 return log;
2744}
2745
2746using namespace shammath;
2747
2751
A module to compute and display statistics on neighbor counts for SPH particles.
Defines the CopyPatchDataFieldFromLayer class for copying fields between patch data layers.
Defines the DistributedBuffers class for managing distributed device buffers in a solver graph.
Implements a forward Euler integration step as a solver graph node.
Defines the GetFieldRefFromLayer class for extracting field references from patch data layers.
Defines the GetObjCntFromLayer class for extracting object counts from patch data layers.
Declares the GetParticlesOutsideSphere module for removing particles.
Declares the IterateSmoothingLengthDensityNeighLim module for iterating smoothing length based on the...
Declares the IterateSmoothingLengthDensity module for iterating smoothing length based on the SPH den...
Declares the KillParticles module for removing particles.
Declares the LoopSmoothingLengthIter module for looping over the smoothing length iteration until con...
Field variant object to instanciate a variant on the patch types.
Header file describing a Node Instance.
Node that applies a custom function to modify connected edges.
Defines the PatchDataLayerRefs class for managing distributed references to patch data layers.
MPI scheduler.
Header file for the patch struct and related function.
Declare a class to register and retrieve nodes and edges from a unique container.
double f64
Alias for double.
std::uint32_t u32
32 bit unsigned integer
std::uint64_t u64
64 bit unsigned integer
std::int32_t i32
32 bit integer
Shamrock assertion utility.
#define SHAM_ASSERT(x)
Shorthand for SHAM_ASSERT_NAMED without a message.
Definition assert.hpp:67
The MPI scheduler.
A buffer allocated in USM (Unified Shared Memory)
void complete_event_state(sycl::event e) const
Complete the event state of the buffer.
T * get_write_access(sham::EventList &depends_list, SourceLocation src_loc=SourceLocation{})
Get a read-write pointer to the buffer's data.
void copy_from(const DeviceBuffer< T, new_target > &other, size_t copy_size)
Copies the content of another buffer to this one.
const T * get_read_access(sham::EventList &depends_list, SourceLocation src_loc=SourceLocation{}) const
Get a read-only pointer to the buffer's data.
DeviceQueue & get_queue(u32 id=0)
Get a reference to a DeviceQueue.
Class to manage a list of SYCL events.
Definition EventList.hpp:31
void add_event(sycl::event e)
Add an event to the list of events.
Definition EventList.hpp:87
Container for objects shared between two distributed data elements.
void for_each(std::function< void(u64, u64, T &)> &&f)
Apply a function to all stored objects.
Represents a collection of objects distributed across patches identified by a u64 id.
iterator add_obj(u64 id, T &&obj)
Adds a new object to the collection.
T & get(u64 id)
Returns a reference to an object in the collection.
void write_to_file(std::string fname)
Write the Fortran formatted file to disk.
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
Vector class based on std::array storage and mdspan.
Definition matrix.hpp:96
handle basic utilities dealing with SPH
The shamrock SPH model.
Definition Solver.hpp:61
void reset_presteps_rint()
Resets tree radius interval field.
Definition Solver.cpp:1219
void reset_merge_ghosts_fields()
Resets merged ghost field data.
Definition Solver.cpp:1484
void update_sync_load_values()
Updates load balancing values and synchronizes patch ownership.
Definition Solver.cpp:1571
bool apply_corrector(Tscal dt, u64 Npart_all)
Definition Solver.cpp:1566
void update_derivs()
Updates time derivatives and applies external forces.
Definition Solver.cpp:1556
void merge_position_ghost()
Merges ghost particle positions from neighboring patches.
Definition Solver.cpp:802
void reset_eos_fields()
Frees memory allocated for EOS fields.
Definition Solver.cpp:1510
void do_predictor_leapfrog(Tscal dt)
Performs predictor step for leapfrog integration.
Definition Solver.cpp:856
void prepare_corrector()
Saves old derivative fields for predictor-corrector integration.
Definition Solver.cpp:1516
void build_ghost_cache()
Builds ghost particle interface cache for inter-patch communication.
Definition Solver.cpp:780
void update_artificial_viscosity(Tscal dt)
Updates artificial viscosity coefficients for shock capturing.
Definition Solver.cpp:1493
TimestepLog evolve_once()
Performs one complete SPH timestep evolution.
Definition Solver.cpp:1578
void vtk_do_dump(std::string filename, bool add_patch_world_id)
Writes VTK dump file for visualization.
Definition Solver.cpp:557
void build_merged_pos_trees()
Builds spatial BVH trees for merged positions including ghosts.
Definition Solver.cpp:845
void clear_merged_pos_trees()
Clears merged position trees to free memory.
Definition Solver.cpp:850
void init_solver_graph()
Initializes the solver graph for computation pipeline.
Definition Solver.cpp:108
void sph_prestep(Tscal time_val, Tscal dt)
Performs pre-step operations for SPH timestep.
Definition Solver.cpp:920
void compute_presteps_rint()
Computes maximum smoothing length in tree nodes for neighbor search.
Definition Solver.cpp:1182
void compute_eos_fields()
Computes equation of state fields (pressure, sound speed)
Definition Solver.cpp:1504
void apply_position_boundary(Tscal time_val)
Applies position-based boundary conditions.
Definition Solver.cpp:734
void reset_neighbors_cache()
Resets neighbor cache.
Definition Solver.cpp:1249
void communicate_merge_ghosts_fields()
Communicates and merges ghost particle fields across processes.
Definition Solver.cpp:1254
void clear_ghost_cache()
Clears ghost particle cache to free memory.
Definition Solver.cpp:796
void init_ghost_layout()
Initializes data layout for ghost particle fields.
Definition Solver.cpp:1167
void start_neighbors_cache()
Builds neighbor particle cache for SPH calculations.
Definition Solver.cpp:1224
Module for constructing spatial tree structures for SPH neighbor searches.
void build_merged_pos_trees()
Builds compressed leaf BVH trees for merged particle positions including ghosts.
Module for computing equation of state quantities.
void compute_eos()
Computes pressure and sound speed from equation of state.
Module for checking conservation of physical quantities.
void check_conservation()
Verifies conservation of mass, momentum, and energy.
void add_ext_forces()
add external forces to the particle acceleration, note that forces dependant on velocity shlould be a...
void compute_ext_forces_indep_v()
is ran once per timestep, it computes the forces that are independant of velocity
Module for reordering particles to improve cache locality.
void reorder_particles()
Reorders particles by Morton code for improved memory access patterns.
Module for writing VTK format output files.
Definition VTKDump.hpp:33
void do_dump(std::string filename, bool add_patch_world_id)
Writes particle data to VTK file for visualization.
Definition VTKDump.cpp:37
Utility class used to move the objects between patches.
void reatribute_patch_objects(SerialPatchTree< T > &sptree, std::string position_field)
Reattribute objects based on a given position field.
ComputeField< T > make_compute_field(std::string new_name, u32 nvar)
create a compute field and init it to zeros
ComputeField< T > save_field(u32 field_idx, std::string new_name)
save a field in patchdata to a compute field
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.
Interface for a solver graph edge representing a field as spans.
PatchDataField< T > & get_field(u64 id) const
Get the underlying PatchDataField at the given id.
void evaluate()
Evaluate the node.
Definition INode.hpp:109
A node that simply frees the allocation of the connected node.
A node that applies a custom function to modify connected edges.
void set_edges(std::shared_ptr< IEdge > to_set)
Set the edges of the node.
virtual void free_alloc() override
Free allocated memory.
A graph container for managing solver nodes and edges with type-safe access.
std::shared_ptr< INode > & get_node_ptr_base(const std::string &name)
Retrieve a node by name as a shared pointer to the base interface.
T & get_edge_ref(const std::string &name)
Get a typed reference to an edge by name.
std::shared_ptr< T > get_edge_ptr(const std::string &name)
Get a typed shared pointer to an edge by name.
std::shared_ptr< T > register_edge(const std::string &name, T &&edge)
Register an edge with automatic type deduction and shared pointer creation.
std::shared_ptr< T > register_node(const std::string &name, T &&node)
Register a node with automatic type deduction and shared pointer creation.
INode & get_node_ref_base(const std::string &name)
Get a reference to a node by name through the base interface.
A Compressed Leaf Bounding Volume Hierarchy (CLBVH) for neighborhood queries.
A data structure representing a Karras Radix Tree Field.
This header file contains utility functions related to exception handling in the code.
MPI string gather / allgather helpers (declarations; implementations in shamalgs/src/collective/gathe...
This file contains the declaration of the memory handling and its methods.
std::vector< T > buf_to_vec(sycl::buffer< T > &buf, u32 len)
Convert a sycl::buffer to a std::vector
Definition memory.cpp:34
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
void throw_unimplemented(SourceLocation loc=SourceLocation{})
Throw a std::runtime_error saying that the function is unimplemented.
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
namespace for the sph model
namespace for the main framework
Definition __init__.py:1
void info(std::string module_name, Types... var2)
Prints a log message with multiple arguments.
Definition logs.hpp:133
file containing formulas for sph forces
sph kernels
f64 get_wtime()
Returns the current wall clock time in seconds.
Structure to store the performance informations about memory allocation and deallocation.
f64 time_alloc_host
Time spent allocating memory on the host.
size_t max_allocated_byte_host
max bytes allocated on the host
f64 time_free_device
Time spent deallocating memory on the device.
size_t max_allocated_byte_device
max bytes allocated on the device
f64 time_alloc_device
Time spent allocating memory on the device.
f64 time_free_host
Time spent deallocating memory on the host.
A class that references multiple buffers or similar objects.
A class to represent a single block of data in a Phantom dump.
u64 get_ref_f32(std::string s)
Gets the index of a block of type f32 with the given name.
u64 get_ref_fort_real(std::string s)
Gets the index of a block of type fort_real with the given name.
i64 tot_count
The total number of values in the block.
std::vector< PhantomDumpBlockArray< fort_real > > blocks_fort_real
The blocks of values of type fort_real.
std::vector< PhantomDumpBlockArray< f32 > > blocks_f32
The blocks of values of type f32.
Class representing a Phantom dump file.
void override_magic_number()
Overrides the magic numbers used in the PhantomDump struct.
shambase::FortranIOFile gen_file()
Generates a Phantom dump file from the current state of the object.
BCConfig< Tvec > BCConfig
Configuration of the boundary conditions.
Patch object that contain generic patch information.
Definition Patch.hpp:33
u64 id_patch
unique key that identify the patch
Definition Patch.hpp:86
Functions related to the MPI communicator.