12#include "gadgetconfig.h"
24#include "../data/allvars.h"
25#include "../data/dtypes.h"
26#include "../data/intposconvert.h"
27#include "../data/mymalloc.h"
28#include "../domain/domain.h"
29#include "../fof/fof.h"
30#include "../fof/fof_io.h"
31#include "../gravtree/gravtree.h"
32#include "../logs/timer.h"
33#include "../main/simulation.h"
34#include "../mpi_utils/mpi_utils.h"
35#include "../ngbtree/ngbtree.h"
36#include "../sort/cxxsort.h"
37#include "../sort/parallel_sort.h"
38#include "../sort/peano.h"
39#include "../subfind/subfind.h"
40#include "../system/system.h"
41#include "../time_integration/timestep.h"
53template <
typename partset>
54void fof<partset>::fof_fof(
int num,
const char *grpcat_basename,
const char *grpcat_dirbasename,
double inner_distance)
58 mpi_printf(
"FOF: Begin to compute FoF group catalogue... (presently allocated=%g MB)\n",
Mem.
getAllocatedBytesInMB());
63 Tp->LinkL = fof_get_comoving_linking_length();
64 mpi_printf(
"FOF: Comoving linking length: %g\n", Tp->LinkL);
67#if defined(LIGHTCONE_PARTICLES_GROUPS)
68 Tp->DistanceOrigin = (
double *)
Mem.mymalloc(
"DistanceOrigin", Tp->NumPart *
sizeof(
double));
72 Tp->MinIDTask = (
int *)
Mem.mymalloc(
"MinIDTask", Tp->NumPart *
sizeof(
int));
73 Tp->Head = (
int *)
Mem.mymalloc(
"Head", Tp->NumPart *
sizeof(
int));
74 Tp->Next = (
int *)
Mem.mymalloc(
"Next", Tp->NumPart *
sizeof(
int));
75 Tp->Tail = (
int *)
Mem.mymalloc(
"Tail", Tp->NumPart *
sizeof(
int));
76 Tp->Len = (
int *)
Mem.mymalloc(
"Len", Tp->NumPart *
sizeof(
int));
80 for(
int i = 0; i < Tp->NumPart; i++)
82 Tp->Head[i] = Tp->Tail[i] = i;
85 Tp->MinID[i] = Tp->P[i].ID;
86 Tp->MinIDTask[i] = ThisTask;
88#if defined(LIGHTCONE_PARTICLES_GROUPS)
89 Tp->DistanceOrigin[i] = fof_distance_to_origin(i);
90 Tp->P[i].setFlagSaveDistance();
96 int npart = Tp->NumPart;
97 int *targetlist = NULL;
100 int *targetlist = (
int *)
Mem.mymalloc(
"targetlist", Tp->NumPart *
sizeof(
int));
101 for(
int i = 0; i < Tp->NumPart; i++)
102 if(is_type_primary_link_type(Tp->P[i].getType()))
103 targetlist[npart++] = i;
107 FoFNgbTree.treeallocate(Tp->NumPart, Tp, FoFDomain);
108 FoFNgbTree.treebuild(npart, targetlist);
111 double cputime = fof_find_groups();
112 mpi_printf(
"FOF: primary group finding took = %g sec\n", cputime);
115 cputime = fof_find_nearest_dmparticle();
116 mpi_printf(
"FOF: attaching gas and star particles to nearest dm particles took = %g sec\n", cputime);
119 FoFNgbTree.treefree();
121 Mem.myfree(targetlist);
124 Mem.myfree(Tp->Tail);
125 Mem.myfree(Tp->Next);
127#if defined(LIGHTCONE_PARTICLES_GROUPS)
128 double save_distance = 0;
129 if(inner_distance > 0)
130 save_distance = inner_distance + Tp->LinkL;
139 FOF_PList = (fof_particle_list *)
Mem.mymalloc_movable(&FOF_PList,
"FOF_PList", Tp->NumPart *
sizeof(fof_particle_list));
140 for(
int i = 0; i < Tp->NumPart; i++)
142 FOF_PList[i].MinID = Tp->MinID[Tp->Head[i]];
143 FOF_PList[i].MinIDTask = Tp->MinIDTask[Tp->Head[i]];
144 FOF_PList[i].Pindex = i;
146#if defined(LIGHTCONE_PARTICLES_GROUPS)
147 FOF_PList[i].DistanceOrigin = Tp->DistanceOrigin[Tp->Head[i]];
149 if(Tp->DistanceOrigin[Tp->Head[i]] < save_distance)
150 Tp->P[i].clearFlagSaveDistance();
155 Mem.myfree_movable(Tp->Head);
156 Mem.myfree_movable(Tp->MinIDTask);
157 Mem.myfree_movable(Tp->MinID);
158#if defined(LIGHTCONE_PARTICLES_GROUPS)
159 Mem.myfree_movable(Tp->DistanceOrigin);
166 FOF_GList = (fof_group_list *)
Mem.mymalloc_movable(&FOF_GList,
"FOF_GList",
sizeof(fof_group_list) * Tp->NumPart);
168 fof_compile_catalogue(inner_distance);
172 mpi_printf(
"FOF: compiling local group data and catalogue took = %g sec\n",
Logs.
timediff(t0, t1));
178 long long largestgroup = 0;
181 long long largestloc = 0;
183 for(
int i = 0; i < NgroupsExt; i++)
184 if(FOF_GList[i].Count > largestloc)
185 largestloc = FOF_GList[i].Count;
186 MPI_Allreduce(&largestloc, &largestgroup, 1, MPI_LONG_LONG, MPI_MAX, Communicator);
189 mpi_printf(
"FOF: Total number of FOF groups with at least %d particles: %lld\n",
FOF_GROUP_MIN_LEN, TotNgroups);
190 mpi_printf(
"FOF: Largest FOF group has %lld particles.\n", largestgroup);
191 mpi_printf(
"FOF: Total number of particles in FOF groups: %lld\n", TotNids);
198 MaxNgroups = NgroupsExt;
199 Group = (group_properties *)
Mem.mymalloc_movable(&Group,
"Group",
sizeof(group_properties) * MaxNgroups);
206 mycxxsort(FOF_GList, FOF_GList + NgroupsExt, fof_compare_FOF_GList_MinID);
209 long long count_nids = 0;
210 for(
int i = 0, start = 0; i < NgroupsExt; i++)
212 while(FOF_PList[start].MinID.get() < FOF_GList[i].MinID.get())
215 if(start > Tp->NumPart)
219 if(FOF_PList[start].MinID.get() != FOF_GList[i].MinID.get())
223 for(lenloc = 0; start + lenloc < Tp->NumPart;)
224 if(FOF_PList[start + lenloc].MinID.get() == FOF_GList[i].MinID.get())
229 Group[i].MinID = FOF_GList[i].MinID.get();
230 Group[i].MinIDTask = FOF_GList[i].MinIDTask;
231 Group[i].Len = FOF_GList[i].Count;
234 fof_compute_group_properties(i, start, lenloc);
237 count_nids += lenloc;
240 Mem.myfree_movable(FOF_GList);
245 sumup_longs(1, &count_nids, &totNids, Communicator);
246 if(totNids != TotNids)
247 Terminate(
"Task=%d Nids=%lld count_nids=%lld totNids=%lld TotNids=%lld\n", ThisTask, Nids, count_nids, totNids, TotNids);
250 fof_add_in_properties_of_group_segments();
253 mpi_printf(
"FOF: computation of group properties took = %g sec\n",
Logs.
timediff(t0, t1));
256 fof_assign_group_numbers();
258 Mem.myfree_movable(FOF_PList);
262 fof_finish_group_properties();
267 fof_assign_group_offset();
271 mpi_printf(
"FOF: Finished computing FoF groups. Complete work took %g sec (presently allocated=%g MB)\n",
Logs.
timediff(ta, tb),
280 sprintf(catname,
"%s_subhalo_tab", grpcat_basename);
282 subfind_find_subhalos(num, catname, grpcat_dirbasename);
297 sprintf(catname,
"%s_tab", grpcat_basename);
299 FoF_IO.fof_subfind_save_groups(num, catname, grpcat_dirbasename);
306 Mem.myfree_movable(Group);
316 fof_prepare_output_order();
318 mpi_printf(
"FOF: preparing output order of particles took %g sec\n",
Logs.
timediff(t0, t1));
327template <
typename partset>
328void fof<partset>::fof_prepare_output_order(
void)
331 for(
int i = 0; i <
NTYPES; i++)
334 for(
int i = 0; i < Tp->NumPart; i++)
337 Tp->PS[i].Type = Tp->P[i].getType();
339 ntype[Tp->P[i].getType()]++;
342#if defined(MERGERTREE) && defined(SUBFIND)
347 mycxxsort_parallel(Tp->PS, Tp->PS + Tp->NumPart, fof_compare_subfind_data_GroupNr_SubRankInNr_BndEgy, Communicator);
349 int *num_list = (
int *)
Mem.mymalloc(
"num_list", NTask *
sizeof(
int));
350 MPI_Allgather(&Tp->NumPart, 1, MPI_INT, num_list, 1, MPI_INT, Communicator);
352 int prev_non_empty_task = ThisTask - 1;
354 while(prev_non_empty_task >= 0)
355 if(num_list[prev_non_empty_task] == 0)
356 prev_non_empty_task--;
362 MPI_Allgather(&Tp->PS[Tp->NumPart > 0 ? Tp->NumPart - 1 : 0],
sizeof(
subfind_data), MPI_BYTE, aux_last_element,
sizeof(
subfind_data),
363 MPI_BYTE, Communicator);
365 int first_index = INT_MAX;
366 int first_task = NTask;
368 for(
int i = 0; i < Tp->NumPart; i++)
369 if(Tp->PS[i].GroupNr.get() <
HALONR_MAX && Tp->PS[i].SubRankInGr < INT_MAX)
375 int *first_list = (
int *)
Mem.mymalloc(
"first_list", NTask *
sizeof(
int));
376 MPI_Allgather(&first_index, 1, MPI_INT, first_list, 1, MPI_INT, Communicator);
377 for(
int n = 0; n < NTask; n++)
379 if(first_list[n] != INT_MAX)
381 first_index = first_list[n];
386 Mem.myfree(first_list);
391 for(
int i = 0; i < Tp->NumPart; i++)
393 if(Tp->PS[i].GroupNr.get() <
HALONR_MAX && Tp->PS[i].SubRankInGr < INT_MAX)
397 if(prev_non_empty_task >= 0)
399 if(Tp->PS[i].GroupNr.get() != aux_last_element[prev_non_empty_task].GroupNr.get() ||
400 Tp->PS[i].SubRankInGr !=
401 aux_last_element[prev_non_empty_task].SubRankInGr)
403 if(ThisTask != first_task || i != first_index)
411 else if(Tp->PS[i].GroupNr.get() != Tp->PS[i - 1].GroupNr.get() ||
412 Tp->PS[i].SubRankInGr != Tp->PS[i - 1].SubRankInGr)
414 if(ThisTask != first_task || i != first_index)
421 Tp->PS[i].SubhaloNr.set(subnr);
422 Tp->PS[i].RankInSubhalo.set(rank++);
424 if(subnr < 0 || subnr >= TotNsubhalos)
425 Terminate(
"i=%d NumPart=%d subnr=%lld PS[i].SubhaloNr.get()=%lld >= TotNsubhalos=%lld", i, Tp->NumPart, subnr,
426 (
long long)Tp->PS[i].SubhaloNr.get(), TotNsubhalos);
431 Tp->PS[i].RankInSubhalo.set(INT_MAX);
435 long long *subnr_list = (
long long *)
Mem.mymalloc(
"subnr_list", NTask *
sizeof(
long long));
436 MPI_Allgather(&subnr, 1, MPI_LONG_LONG, subnr_list, 1, MPI_LONG_LONG, Communicator);
438 long long subnr_prev = 0;
439 for(
int i = 0; i < ThisTask; i++)
440 subnr_prev += subnr_list[i];
442 for(
int i = 0; i < Tp->NumPart; i++)
443 if(Tp->PS[i].GroupNr.get() <
HALONR_MAX && Tp->PS[i].SubRankInGr < INT_MAX)
444 Tp->PS[i].SubhaloNr.set(Tp->PS[i].SubhaloNr.get() + subnr_prev);
447 long long *rank_list = (
long long *)
Mem.mymalloc(
"rank_list", NTask *
sizeof(
long long));
448 MPI_Allgather(&rank, 1, MPI_LONG_LONG, rank_list, 1, MPI_LONG_LONG, Communicator);
450 if(prev_non_empty_task >= 0)
452 long long rank = rank_list[prev_non_empty_task];
454 for(
int i = 0; i < Tp->NumPart; i++)
456 if(Tp->PS[i].GroupNr.get() <
HALONR_MAX && Tp->PS[i].SubRankInGr < INT_MAX)
460 if(prev_non_empty_task >= 0)
461 if(Tp->PS[i].GroupNr.get() != aux_last_element[prev_non_empty_task].GroupNr.get() ||
462 Tp->PS[i].SubRankInGr !=
463 aux_last_element[prev_non_empty_task].SubRankInGr)
466 else if(Tp->PS[i].GroupNr.get() != Tp->PS[i - 1].GroupNr.get() ||
467 Tp->PS[i].SubRankInGr != Tp->PS[i - 1].SubRankInGr)
470 Tp->PS[i].RankInSubhalo.set(rank++);
475 Mem.myfree(rank_list);
476 Mem.myfree(subnr_list);
477 Mem.myfree(aux_last_element);
478 Mem.myfree(num_list);
481 mycxxsort_parallel(Tp->PS, Tp->PS + Tp->NumPart, fof_compare_subfind_data_OriginTask_OriginIndex, Communicator);
485 for(
int i = 0; i < Tp->NumPart; i++)
487#if defined(RANDOMIZE_DOMAINCENTER) && defined(PERIODIC)
489 peano_hilbert_key(Tp->P[i].IntPos[0] - Tp->CurrentShiftVector[0], Tp->P[i].IntPos[1] - Tp->CurrentShiftVector[1],
499 Tp->PS[i].SubRankInGr = 0;
500 Tp->PS[i].v.DM_BindingEnergy = 0;
506 mycxxsort(Tp->PS, Tp->PS + Tp->NumPart, fof_compare_subfind_data_Type);
509 for(
int i = 0, off = 0; i <
NTYPES; i++)
511 mycxxsort_parallel(Tp->PS + off, Tp->PS + off + ntype[i], fof_compare_subfind_data_GroupNr_SubNr_Egy_Key, Communicator);
516 for(
int i = 0; i < Tp->NumPart; i++)
518 Tp->PS[i].TargetTask = ThisTask;
519 Tp->PS[i].TargetIndex = i;
523 mycxxsort_parallel(Tp->PS, Tp->PS + Tp->NumPart, fof_compare_subfind_data_OriginTask_OriginIndex, Communicator);
526 FoFDomain->particle_exchange_based_on_PS(Communicator);
530template <
typename partset>
531double fof<partset>::fof_get_comoving_linking_length(
void)
535 double mass = 0, masstot;
537 for(
int i = 0; i < Tp->NumPart; i++)
538 if(is_type_primary_link_type(Tp->P[i].getType()))
541 mass += Tp->P[i].getMass();
544 MPI_Allreduce(&mass, &masstot, 1, MPI_DOUBLE, MPI_SUM, Communicator);
550template <
typename partset>
551void fof<partset>::fof_compile_catalogue(
double inner_distance)
554 mycxxsort(FOF_PList, FOF_PList + Tp->NumPart, fof_compare_FOF_PList_MinID);
560 for(
int i = 0; i < Tp->NumPart; i++)
562 FOF_GList[i].MinID = FOF_PList[i].MinID;
563 FOF_GList[i].MinIDTask = FOF_PList[i].MinIDTask;
564 FOF_GList[i].Count = 1;
565#if defined(LIGHTCONE_PARTICLES_GROUPS)
566 FOF_GList[i].DistanceOrigin = FOF_PList[i].DistanceOrigin;
579 for(
int i = 1, start = 0; i < Tp->NumPart; i++)
581 if(FOF_GList[i].MinID.get() == FOF_GList[start].MinID.get())
583 if(FOF_GList[i].MinIDTask != FOF_GList[start].MinIDTask)
584 Terminate(
"FOF_GList[i].MinIDTask=%d != FOF_GList[start].MinIDTask=%d", FOF_GList[i].MinIDTask,
585 FOF_GList[start].MinIDTask);
587 FOF_GList[start].Count += FOF_GList[i].Count;
589#if defined(LIGHTCONE_PARTICLES_GROUPS)
590 if(FOF_GList[start].DistanceOrigin != FOF_GList[i].DistanceOrigin)
592 "start=%d i=%d Tp->NumPart=%d FOF_GList[start].DistanceOrigin=%g != FOF_GList[i].DistanceOrigin=%g MinID=%lld "
594 start, i, Tp->NumPart, FOF_GList[start].DistanceOrigin, FOF_GList[i].DistanceOrigin,
595 (
long long)FOF_GList[i].MinID.get(), FOF_GList[i].MinIDTask);
601 FOF_GList[start] = FOF_GList[i];
607 FOF_GList = (fof_group_list *)
Mem.myrealloc_movable(FOF_GList,
sizeof(fof_group_list) * NgroupsExt);
610 mycxxsort(FOF_GList, FOF_GList + NgroupsExt, fof_compare_FOF_GList_MinIDTask);
612 int *Send_count = (
int *)
Mem.mymalloc(
"Send_count",
sizeof(
int) * NTask);
613 int *Send_offset = (
int *)
Mem.mymalloc(
"Send_offset",
sizeof(
int) * NTask);
614 int *Recv_count = (
int *)
Mem.mymalloc(
"Recv_count",
sizeof(
int) * NTask);
615 int *Recv_offset = (
int *)
Mem.mymalloc(
"Recv_offset",
sizeof(
int) * NTask);
618 for(
int i = 0; i < NTask; i++)
620 for(
int i = 0; i < NgroupsExt; i++)
621 Send_count[FOF_GList[i].MinIDTask]++;
624 myMPI_Alltoall(Send_count, 1, MPI_INT, Recv_count, 1, MPI_INT, Communicator);
628 Recv_offset[0] = 0, Send_offset[0] = 0;
629 for(
int j = 0; j < NTask; j++)
633 nimport += Recv_count[j];
637 Send_offset[j] = Send_offset[j - 1] + Send_count[j - 1];
638 Recv_offset[j] = Recv_offset[j - 1] + Recv_count[j - 1];
643 fof_group_list *get_FOF_GList = (fof_group_list *)
Mem.mymalloc(
"get_FOF_GList", nimport *
sizeof(fof_group_list));
646 for(
int ngrp = 1; ngrp < (1 << PTask); ngrp++)
648 int recvTask = ThisTask ^ ngrp;
652 if(Send_count[recvTask] > 0 || Recv_count[recvTask] > 0)
655 myMPI_Sendrecv(&FOF_GList[Send_offset[recvTask]], Send_count[recvTask] *
sizeof(fof_group_list), MPI_BYTE, recvTask,
656 TAG_DENS_A, &get_FOF_GList[Recv_offset[recvTask]], Recv_count[recvTask] *
sizeof(fof_group_list), MPI_BYTE,
657 recvTask,
TAG_DENS_A, Communicator, MPI_STATUS_IGNORE);
664 for(
int i = 0; i < nimport; i++)
665 get_FOF_GList[i].MinIDTask = i;
668 mycxxsort(FOF_GList, FOF_GList + NgroupsExt, fof_compare_FOF_GList_MinID);
669 mycxxsort(get_FOF_GList, get_FOF_GList + nimport, fof_compare_FOF_GList_MinID);
676 for(
int i = 0; i < nimport; i++)
678 while(FOF_GList[start].MinID.get() < get_FOF_GList[i].MinID.get())
681 if(start >= NgroupsExt)
682 Terminate(
"start=%d >= NgroupsExt=%d", start, NgroupsExt);
685 if(FOF_GList[start].MinIDTask != ThisTask)
686 Terminate(
"FOF_GList[start].MinIDTask=%d != ThisTask=%d", FOF_GList[start].MinIDTask, ThisTask);
688 if(FOF_GList[start].MinID.get() != get_FOF_GList[i].MinID.get())
690 "FOF_GList[start].MinID != get_FOF_GList[i].MinID start=%d i=%d FOF_GList[start].MinID=%llu get_FOF_GList[i].MinID=%llu\n",
691 start, i, (
long long)FOF_GList[start].MinID.get(), (
long long)get_FOF_GList[i].MinID.get());
693 FOF_GList[start].Count += get_FOF_GList[i].Count;
698 for(
int i = 0; i < nimport; i++)
700 while(FOF_GList[start].MinID.get() < get_FOF_GList[i].MinID.get())
703 if(start >= NgroupsExt)
707 get_FOF_GList[i].Count = FOF_GList[start].Count;
714 mycxxsort(get_FOF_GList, get_FOF_GList + nimport, fof_compare_FOF_GList_MinIDTask);
715 mycxxsort(FOF_GList, FOF_GList + NgroupsExt, fof_compare_FOF_GList_MinIDTask);
718 for(
int i = 0; i < nimport; i++)
719 get_FOF_GList[i].MinIDTask = ThisTask;
722 for(
int ngrp = 1; ngrp < (1 << PTask); ngrp++)
724 int recvTask = ThisTask ^ ngrp;
728 if(Send_count[recvTask] > 0 || Recv_count[recvTask] > 0)
731 myMPI_Sendrecv(&get_FOF_GList[Recv_offset[recvTask]], Recv_count[recvTask] *
sizeof(fof_group_list), MPI_BYTE, recvTask,
732 TAG_DENS_A, &FOF_GList[Send_offset[recvTask]], Send_count[recvTask] *
sizeof(fof_group_list), MPI_BYTE,
733 recvTask,
TAG_DENS_A, Communicator, MPI_STATUS_IGNORE);
739 Mem.myfree(get_FOF_GList);
745#if defined(LIGHTCONE_PARTICLES_GROUPS)
746 double save_distance = 0;
747 if(inner_distance > 0)
748 save_distance = inner_distance + Tp->LinkL;
751 Ngroups = 0, Nids = 0;
752 for(
int i = 0; i < NgroupsExt; i++)
756 FOF_GList[i] = FOF_GList[NgroupsExt - 1];
760#if defined(LIGHTCONE_PARTICLES_GROUPS)
761 else if(FOF_GList[i].DistanceOrigin < save_distance)
763 FOF_GList[i] = FOF_GList[NgroupsExt - 1];
770 if(FOF_GList[i].MinIDTask == ThisTask)
773 Nids += FOF_GList[i].Count;
778 Mem.myfree(Recv_offset);
779 Mem.myfree(Recv_count);
780 Mem.myfree(Send_offset);
781 Mem.myfree(Send_count);
784 FOF_GList = (fof_group_list *)
Mem.myrealloc_movable(FOF_GList,
sizeof(fof_group_list) * NgroupsExt);
787template <
typename partset>
788void fof<partset>::fof_assign_group_numbers(
void)
790 mpi_printf(
"FOF: start assigning group numbers\n");
794 for(
int i = 0; i < NgroupsExt; i++)
795 Group[i].OriginTask = ThisTask;
799 mycxxsort_parallel(Group, Group + NgroupsExt, fof_compare_Group_Len_MinID_DiffOriginTaskMinIDTask, Communicator);
803 for(
int i = 0; i < NgroupsExt; i++)
805 if(Group[i].OriginTask == Group[i].MinIDTask)
808 Group[i].GroupNr = ngr - 1;
814 int *ngr_list = (
int *)
Mem.mymalloc(
"ngr_list",
sizeof(
int) * NTask);
815 MPI_Allgather(&ngr, 1, MPI_INT, ngr_list, 1, MPI_INT, Communicator);
817 long long ngr_sum = 0;
818 for(
int j = 0; j < ThisTask; j++)
819 ngr_sum += ngr_list[j];
821 Mem.myfree(ngr_list);
824 for(
int i = 0; i < NgroupsExt; i++)
825 Group[i].GroupNr += ngr_sum;
829 if(ngr_sum != TotNgroups)
830 Terminate(
"inconsistency ngr_sum=%lld\n", ngr_sum);
833 mycxxsort_parallel(Group, Group + NgroupsExt, fof_compare_Group_OriginTask_MinID, Communicator);
838 for(
int i = 0; i < Tp->NumPart; i++)
841 long long Nids_old = Nids;
847 for(
int i = 0, start = 0; i < NgroupsExt; i++)
849 while(FOF_PList[start].MinID.get() < Group[i].MinID)
852 if(start > Tp->NumPart)
856 if(FOF_PList[start].MinID.get() != Group[i].MinID)
857 Terminate(
"FOF_PList[start=%d].MinID=%lld != Group[i=%d].MinID=%lld", start, (
long long)FOF_PList[start].MinID.get(), i,
858 (
long long)Group[i].MinID);
861 for(lenloc = 0; start + lenloc < Tp->NumPart;)
862 if(FOF_PList[start + lenloc].MinID.get() == Group[i].MinID)
864 Tp->PS[FOF_PList[start + lenloc].Pindex].GroupNr.set(Group[i].GroupNr);
877 if(totNids != TotNids)
878 Terminate(
"Task=%d Nids=%lld Nids_old=%lld totNids=%lld TotNids=%lld\n", ThisTask, Nids, Nids_old, totNids, TotNids);
882 mpi_printf(
"FOF: Assigning of group numbers took = %g sec\n",
Logs.
timediff(t0, t1));
888template <
typename partset>
889void fof<partset>::fof_assign_group_offset(
void)
894 long long gtype_previous[
NTYPES];
896 for(
int i = 0; i <
NTYPES; i++)
899 for(
int i = 0; i < Ngroups; i++)
900 for(
int j = 0; j <
NTYPES; j++)
901 gtype_loc[j] += Group[i].LenType[j];
903 int *gtype_all = (
int *)
Mem.mymalloc(
"gtype_all",
NTYPES * NTask *
sizeof(
int));
904 MPI_Allgather(gtype_loc,
NTYPES, MPI_INT, gtype_all,
NTYPES, MPI_INT, Communicator);
906 for(
int i = 0; i <
NTYPES; i++)
907 gtype_previous[i] = 0;
909 for(
int i = 0; i < ThisTask; i++)
910 for(
int j = 0; j <
NTYPES; j++)
911 gtype_previous[j] += gtype_all[i *
NTYPES + j];
914 for(
int j = 0; j <
NTYPES; j++)
915 Group[0].OffsetType[j] = gtype_previous[j];
917 for(
int i = 1; i < Ngroups; i++)
918 for(
int j = 0; j <
NTYPES; j++)
919 Group[i].OffsetType[j] = Group[i - 1].OffsetType[j] + Group[i - 1].LenType[j];
921 Mem.myfree(gtype_all);
924template <
typename partset>
925void fof<partset>::fof_compute_group_properties(
int gr,
int start,
int len)
927 int start_index = FOF_PList[start].Pindex;
930 Group[gr].Ascale = 0;
934#if defined(SUBFIND_ORPHAN_TREATMENT)
935 Group[gr].LenPrevMostBnd = 0;
938 for(
int k = 0; k < 3; k++)
941 Group[gr].Vel[k] = 0;
942 Group[gr].FirstIntPos[k] = Tp->P[start_index].IntPos[k];
945 for(
int k = 0; k <
NTYPES; k++)
947 Group[gr].LenType[k] = 0;
948 Group[gr].MassType[k] = 0;
951 for(
int k = 0; k < len; k++)
953 int index = FOF_PList[start + k].Pindex;
955 Group[gr].Mass += Tp->P[index].getMass();
956 int type = Tp->P[index].getType();
958 Group[gr].Ascale += Tp->P[index].getMass() * Tp->P[index].getAscale();
960 Group[gr].LenType[type]++;
961 Group[gr].MassType[type] += Tp->P[index].getMass();
963#if defined(SUBFIND_ORPHAN_TREATMENT)
964 if(Tp->P[index].ID.is_previously_most_bound())
965 Group[gr].LenPrevMostBnd++;
969 if(Tp->P[index].getType() == 0)
970 Group[gr].Sfr += Tp->SphP[index].Sfr;
974 Tp->nearest_image_intpos_to_pos(Tp->P[index].IntPos, Tp->P[start_index].IntPos,
977 for(
int j = 0; j < 3; j++)
979 Group[gr].CM[j] += Tp->P[index].getMass() * xyz[j];
980 Group[gr].Vel[j] += Tp->P[index].getMass() * Tp->P[index].Vel[j];
985template <
typename partset>
986void fof<partset>::fof_add_in_properties_of_group_segments(
void)
988 int *Send_count = (
int *)
Mem.mymalloc(
"Send_count",
sizeof(
int) * NTask);
989 int *Send_offset = (
int *)
Mem.mymalloc(
"Send_offset",
sizeof(
int) * NTask);
990 int *Recv_count = (
int *)
Mem.mymalloc(
"Recv_count",
sizeof(
int) * NTask);
991 int *Recv_offset = (
int *)
Mem.mymalloc(
"Recv_offset",
sizeof(
int) * NTask);
994 mycxxsort(Group, Group + NgroupsExt, fof_compare_Group_MinIDTask);
997 for(
int i = 0; i < NTask; i++)
999 for(
int i = 0; i < NgroupsExt; i++)
1000 Send_count[Group[i].MinIDTask]++;
1002 myMPI_Alltoall(Send_count, 1, MPI_INT, Recv_count, 1, MPI_INT, Communicator);
1005 Recv_offset[0] = 0, Send_offset[0] = 0;
1007 for(
int j = 0; j < NTask; j++)
1011 nimport += Recv_count[j];
1015 Send_offset[j] = Send_offset[j - 1] + Send_count[j - 1];
1016 Recv_offset[j] = Recv_offset[j - 1] + Recv_count[j - 1];
1020 group_properties *get_Group = (group_properties *)
Mem.mymalloc(
"get_Group",
sizeof(group_properties) * nimport);
1022 for(
int ngrp = 1; ngrp < (1 << PTask); ngrp++)
1024 int recvTask = ThisTask ^ ngrp;
1026 if(recvTask < NTask)
1028 if(Send_count[recvTask] > 0 || Recv_count[recvTask] > 0)
1031 myMPI_Sendrecv(&Group[Send_offset[recvTask]], Send_count[recvTask] *
sizeof(group_properties), MPI_BYTE, recvTask,
1032 TAG_DENS_A, &get_Group[Recv_offset[recvTask]], Recv_count[recvTask] *
sizeof(group_properties), MPI_BYTE,
1033 recvTask,
TAG_DENS_A, Communicator, MPI_STATUS_IGNORE);
1039 mycxxsort(Group, Group + NgroupsExt, fof_compare_Group_MinID);
1040 mycxxsort(get_Group, get_Group + nimport, fof_compare_Group_MinID);
1044 for(
int i = 0; i < nimport; i++)
1046 while(Group[start].MinID < get_Group[i].MinID)
1049 if(start >= NgroupsExt)
1053 Group[start].Mass += get_Group[i].Mass;
1054 Group[start].Ascale += get_Group[i].Ascale;
1056 for(
int j = 0; j <
NTYPES; j++)
1058 Group[start].LenType[j] += get_Group[i].LenType[j];
1059 Group[start].MassType[j] += get_Group[i].MassType[j];
1061#if defined(SUBFIND_ORPHAN_TREATMENT)
1062 Group[start].LenPrevMostBnd += get_Group[i].LenPrevMostBnd;
1065 Group[start].Sfr += get_Group[i].Sfr;
1068 double tmpxyz[3] = {get_Group[i].CM[0] / get_Group[i].Mass, get_Group[i].CM[1] / get_Group[i].Mass,
1069 get_Group[i].CM[2] / get_Group[i].Mass};
1074 delta[0] += get_Group[i].FirstIntPos[0];
1075 delta[1] += get_Group[i].FirstIntPos[1];
1076 delta[2] += get_Group[i].FirstIntPos[2];
1078 Tp->constrain_intpos(delta);
1081 Tp->nearest_image_intpos_to_pos(delta, Group[start].FirstIntPos, xyz);
1083 for(
int j = 0; j < 3; j++)
1085 Group[start].CM[j] += get_Group[i].Mass * xyz[j];
1086 Group[start].Vel[j] += get_Group[i].Vel[j];
1090 Mem.myfree(get_Group);
1092 Mem.myfree(Recv_offset);
1093 Mem.myfree(Recv_count);
1094 Mem.myfree(Send_offset);
1095 Mem.myfree(Send_count);
1098template <
typename partset>
1099void fof<partset>::fof_finish_group_properties(
void)
1101 for(
int i = 0; i < NgroupsExt; i++)
1103 if(Group[i].MinIDTask == ThisTask)
1106 for(
int j = 0; j < 3; j++)
1108 Group[i].Vel[j] /= Group[i].Mass;
1109 cm[j] = Group[i].CM[j] / Group[i].Mass;
1112 Group[i].Ascale /= Group[i].Mass;
1117 delta[0] += Group[i].FirstIntPos[0];
1118 delta[1] += Group[i].FirstIntPos[1];
1119 delta[2] += Group[i].FirstIntPos[2];
1121 Tp->constrain_intpos(delta);
1123 fof_get_halo_position(delta, cm);
1125 Group[i].CM[0] = cm[0];
1126 Group[i].CM[1] = cm[1];
1127 Group[i].CM[2] = cm[2];
1132 for(
int j = 0; j < 3; j++)
1133 Group[i].Pos[j] = Group[i].CM[j];
1135 for(
int j = 0; j < 3; j++)
1136 Group[i].IntPos[j] = delta[j];
1140 int ngr = NgroupsExt;
1143 for(
int i = 0; i < ngr; i++)
1145 if(Group[i].MinIDTask != ThisTask)
1147 Group[i] = Group[ngr - 1];
1156 mycxxsort(Group, Group + Ngroups, fof_compare_Group_MinID);
1159#if defined(LIGHTCONE_PARTICLES_GROUPS)
1162double fof<simparticles>::fof_distance_to_origin(
int i)
1168double fof<lcparticles>::fof_distance_to_origin(
int i)
1176void fof<simparticles>::fof_get_halo_position(
MyIntPosType *intpos,
double *pos)
1178 Tp->intpos_to_pos(intpos, pos);
1181#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
1183void fof<lcparticles>::fof_get_halo_position(
MyIntPosType *intpos,
double *pos)
1187 Tp->nearest_image_intpos_to_pos(intpos, origin, pos);
1192#include "../data/simparticles.h"
1193template class fof<simparticles>;
1195#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
1196#include "../data/lcparticles.h"
1197template class fof<lcparticles>;
global_data_all_processes All
double timediff(double t0, double t1)
double getAllocatedBytesInMB(void)
#define FOF_GROUP_MIN_LEN
double mycxxsort(T *begin, T *end, Tcomp comp)
int32_t MySignedIntPosType
#define BITS_FOR_POSITIONS
#define TIMER_START(counter)
Starts the timer counter.
#define TIMER_STOP(counter)
Stops the timer counter.
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
int myMPI_Alltoall(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
int myMPI_Sendrecv(void *sendbuf, size_t sendcount, MPI_Datatype sendtype, int dest, int sendtag, void *recvbuf, size_t recvcount, MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm, MPI_Status *status)
void sumup_longs(int n, long long *src, long long *res, MPI_Comm comm)
expr pow(half base, half exp)
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)
peanokey peano_hilbert_key(MyIntPosType x, MyIntPosType y, MyIntPosType z, int bits)