GADGET-4
fmm.h
Go to the documentation of this file.
1/*******************************************************************************
2 * \copyright This file is part of the GADGET4 N-body/SPH code developed
3 * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4 * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5 *******************************************************************************/
6
12#ifndef FMM_H
13#define FMM_H
14
15#ifdef FMM
16
17#include "../data/symtensors.h"
18#include "../gravtree/gravtree.h"
19
20class fmm : public gravtree<simparticles>
21{
22 public:
23 void gravity_fmm(int timebin);
24
25 private:
26 long long interactioncountPP;
27 long long interactioncountPN;
28 long long interactioncountNN;
29 long long sum_NumForeignNodes;
30 long long sum_NumForeignPoints;
31 long long interactioncountEffective;
32
33 char *Topnode_depends_on_local_mass;
34
35 MyReal errTolTheta2;
36 MyReal errTolThetaMax2;
37 MyReal errTolForceAcc;
38
39#ifdef PRESERVE_SHMEM_BINARY_INVARIANCE
40 bool skip_actual_force_computation;
41#endif
42
43 struct fieldcoeff
44 {
45 MyReal phi;
46 vector<MyReal> dphi; /* first potential derivative at center of mass (i.e. minus the acceleration) */
47#if(MULTIPOLE_ORDER >= 2)
48 symtensor2<MyReal> d2phi; /* second derivatives */
49#endif
50#if(MULTIPOLE_ORDER >= 3)
51 symtensor3<MyReal> d3phi; /* third derivatives */
52#endif
53#if(MULTIPOLE_ORDER >= 4)
54 symtensor4<MyReal> d4phi; /* fourth derivatives */
55#endif
56#if(MULTIPOLE_ORDER >= 5)
57 symtensor5<MyReal> d5phi; /* fifth derivatives */
58#endif
59
60 MyReal interactions;
61 };
62
63 struct taylor_data
64 {
65 fieldcoeff coeff;
66 };
67 taylor_data *TaylorCoeff;
68
69 struct fmm_workstack_data
70 {
71 int Node1;
72 int Node2;
73 unsigned char ShmRank1;
74 unsigned char ShmRank2;
75 int MinTopLeafNode;
76 };
77
78 fmm_workstack_data *FMM_WorkStack;
79
80 static bool compare_fmm_workstack(const fmm_workstack_data &a, const fmm_workstack_data &b)
81 {
82 if(a.MinTopLeafNode < b.MinTopLeafNode)
83 return true;
84 if(a.MinTopLeafNode > b.MinTopLeafNode)
85 return false;
86
87 return a.Node1 < b.Node1;
88 }
89
90 inline void fmm_add_to_work_stack(int node1, int node2, unsigned char shmrank1, unsigned char shmrank2, int mintopleafnode)
91 {
92 if(NumOnWorkStack + NewOnWorkStack >= MaxOnWorkStack)
93 {
94 Terminate("we have run out of space: NumOnWorkStack=%d + NewOnWorkStack=%d >= MaxOnWorkStack=%d", NumOnWorkStack,
95 NewOnWorkStack, MaxOnWorkStack);
96 }
97
98 FMM_WorkStack[NumOnWorkStack + NewOnWorkStack].Node1 = node1;
99 FMM_WorkStack[NumOnWorkStack + NewOnWorkStack].Node2 = node2;
100 FMM_WorkStack[NumOnWorkStack + NewOnWorkStack].ShmRank1 = shmrank1;
101 FMM_WorkStack[NumOnWorkStack + NewOnWorkStack].ShmRank2 = shmrank2;
102 FMM_WorkStack[NumOnWorkStack + NewOnWorkStack].MinTopLeafNode = mintopleafnode;
103
104 NewOnWorkStack++;
105 }
106
107 inline bool fmm_depends_on_local_mass(int no, unsigned char shmrank)
108 {
109 if(no >= MaxPart && no < FirstNonTopLevelNode)
110 {
111 if(Topnode_depends_on_local_mass[no])
112 return true;
113 else
114 return false;
115 }
116 else if(no >= FirstNonTopLevelNode && no < MaxPart + MaxNodes)
117 {
118 if(shmrank == Shmem.Island_ThisTask)
119 return true;
120 else
121 return false;
122 }
123 else
124 return false;
125 }
126
127 void fmm_force_interact(int no_sink, int no_source, char type_sink, char type_source, unsigned char shmrank_sink,
128 unsigned char shmrank_source, int mintopleafnode, int committed);
129
130 void fmm_force_passdown(int no, unsigned char shmrank, taylor_data taylor_current);
131 void fmm_open_both(gravnode *noptr_sink, gravnode *noptr_source, int mintopleafnode, int committed);
132
133 void fmm_open_node(int no_particle, gravnode *nop, char type_particle, unsigned char shmrank_particle, int mintopleafnode,
134 int committed);
135 void fmm_particle_particle_interaction(int no_sink, int no_source, int type_sink, int type_source, unsigned char shmrank_sink,
136 unsigned char shmrank_source);
137
138 void fmm_particle_node_interaction(int no_sink, int no_source, int type_sink, int type_source, unsigned char shmrank_sink,
139 unsigned char shmrank_source, gravnode *noptr_source, vector<MyReal> &dxyz, MyReal &r2);
140
141 void fmm_node_node_interaction(int no_sink, int no_source, int type_sink, int type_source, unsigned char shmrank_sink,
142 unsigned char shmrank_source, gravnode *noptr_sink, gravnode *noptr_source, vector<MyReal> &dxyz,
143 MyReal &r2);
144
145 int fmm_evaluate_node_node_opening_criterion(gravnode *noptr_sink, gravnode *noptr_source, vector<MyReal> &dxyz, MyReal &r2);
146
147 int fmm_evaluate_particle_node_opening_criterion(int no_sink, char type_sink, unsigned char shmrank_sink, gravnode *nop_source,
148 vector<MyReal> &dxyz, MyReal &r2);
149
150 void fmm_determine_nodes_with_local_mass(int no, int sib);
151};
152
153#endif
154#endif
int Island_ThisTask
double MyReal
Definition: dtypes.h:82
#define Terminate(...)
Definition: macros.h:15
shmem Shmem
Definition: main.cc:45