00001 #include <stdio.h>
00002 #include <stdlib.h>
00003 #include <string.h>
00004 #include <math.h>
00005 #include <time.h>
00006 #include <mpi.h>
00007
00008 #include "allvars.h"
00009 #include "proto.h"
00010
00011
00026 static int last;
00027
00028
00029
00031 #define NTAB 1000
00032
00033 static float tabfac, shortrange_table[NTAB], shortrange_table_potential[NTAB];
00034
00036 static int first_flag = 0;
00037
00038
00039
00040
00041 #ifdef PERIODIC
00042
00043 #define NEAREST(x) (((x)>boxhalf)?((x)-boxsize):(((x)<-boxhalf)?((x)+boxsize):(x)))
00044
00045 #define EN 64
00046
00049 static FLOAT fcorrx[EN + 1][EN + 1][EN + 1];
00050 static FLOAT fcorry[EN + 1][EN + 1][EN + 1];
00051 static FLOAT fcorrz[EN + 1][EN + 1][EN + 1];
00052 static FLOAT potcorr[EN + 1][EN + 1][EN + 1];
00053 static double fac_intp;
00054 #endif
00055
00056
00057
00061 int force_treebuild(int npart)
00062 {
00063 Numnodestree = force_treebuild_single(npart);
00064
00065 force_update_pseudoparticles();
00066
00067 force_flag_localnodes();
00068
00069 TimeOfLastTreeConstruction = All.Time;
00070
00071 return Numnodestree;
00072 }
00073
00074
00075
00093 int force_treebuild_single(int npart)
00094 {
00095 int i, j, subnode = 0, parent, numnodes;
00096 int nfree, th, nn, no;
00097 struct NODE *nfreep;
00098 double lenhalf, epsilon;
00099 peanokey key;
00100
00101
00102
00103 nfree = All.MaxPart;
00104 nfreep = &Nodes[nfree];
00105
00106 nfreep->len = DomainLen;
00107 for(j = 0; j < 3; j++)
00108 nfreep->center[j] = DomainCenter[j];
00109 for(j = 0; j < 8; j++)
00110 nfreep->u.suns[j] = -1;
00111
00112
00113 numnodes = 1;
00114 nfreep++;
00115 nfree++;
00116
00117
00118
00119
00120
00121
00122
00123 force_create_empty_nodes(All.MaxPart, 0, 1, 0, 0, 0, &numnodes, &nfree);
00124
00125
00126
00127
00128
00129
00130
00131 nfreep = &Nodes[nfree];
00132 parent = -1;
00133
00134
00135
00136 for(i = 0; i < npart; i++)
00137 {
00138
00139
00140
00141
00142 epsilon = All.ForceSoftening[P[i].Type];
00143
00144 key = peano_hilbert_key((P[i].Pos[0] - DomainCorner[0]) * DomainFac,
00145 (P[i].Pos[1] - DomainCorner[1]) * DomainFac,
00146 (P[i].Pos[2] - DomainCorner[2]) * DomainFac, BITS_PER_DIMENSION);
00147
00148 no = 0;
00149 while(TopNodes[no].Daughter >= 0)
00150 no = TopNodes[no].Daughter + (key - TopNodes[no].StartKey) / (TopNodes[no].Size / 8);
00151
00152 no = TopNodes[no].Leaf;
00153 th = DomainNodeIndex[no];
00154
00155 while(1)
00156 {
00157 if(th >= All.MaxPart)
00158 {
00159 subnode = 0;
00160 if(P[i].Pos[0] > Nodes[th].center[0])
00161 subnode += 1;
00162 if(P[i].Pos[1] > Nodes[th].center[1])
00163 subnode += 2;
00164 if(P[i].Pos[2] > Nodes[th].center[2])
00165 subnode += 4;
00166
00167 nn = Nodes[th].u.suns[subnode];
00168
00169 if(nn >= 0)
00170 {
00171 parent = th;
00172 th = nn;
00173 }
00174 else
00175 {
00176
00177
00178
00179 Nodes[th].u.suns[subnode] = i;
00180 break;
00181 }
00182 }
00183 else
00184 {
00185
00186
00187
00188 Nodes[parent].u.suns[subnode] = nfree;
00189
00190 nfreep->len = 0.5 * Nodes[parent].len;
00191 lenhalf = 0.25 * Nodes[parent].len;
00192
00193 if(subnode & 1)
00194 nfreep->center[0] = Nodes[parent].center[0] + lenhalf;
00195 else
00196 nfreep->center[0] = Nodes[parent].center[0] - lenhalf;
00197
00198 if(subnode & 2)
00199 nfreep->center[1] = Nodes[parent].center[1] + lenhalf;
00200 else
00201 nfreep->center[1] = Nodes[parent].center[1] - lenhalf;
00202
00203 if(subnode & 4)
00204 nfreep->center[2] = Nodes[parent].center[2] + lenhalf;
00205 else
00206 nfreep->center[2] = Nodes[parent].center[2] - lenhalf;
00207
00208 nfreep->u.suns[0] = -1;
00209 nfreep->u.suns[1] = -1;
00210 nfreep->u.suns[2] = -1;
00211 nfreep->u.suns[3] = -1;
00212 nfreep->u.suns[4] = -1;
00213 nfreep->u.suns[5] = -1;
00214 nfreep->u.suns[6] = -1;
00215 nfreep->u.suns[7] = -1;
00216
00217
00218 subnode = 0;
00219 if(P[th].Pos[0] > nfreep->center[0])
00220 subnode += 1;
00221 if(P[th].Pos[1] > nfreep->center[1])
00222 subnode += 2;
00223 if(P[th].Pos[2] > nfreep->center[2])
00224 subnode += 4;
00225 #ifndef NOTREERND
00226 if(nfreep->len < 1.0e-3 * epsilon)
00227 {
00228
00229
00230
00231
00232
00233 subnode = (int) (8.0 * get_random_number((0xffff & P[i].ID) + P[i].GravCost));
00234 P[i].GravCost += 1;
00235 if(subnode >= 8)
00236 subnode = 7;
00237 }
00238 #endif
00239 nfreep->u.suns[subnode] = th;
00240
00241 th = nfree;
00242
00243
00244
00245 numnodes++;
00246 nfree++;
00247 nfreep++;
00248
00249 if((numnodes) >= MaxNodes)
00250 {
00251 printf("task %d: maximum number %d of tree-nodes reached.\n", ThisTask, MaxNodes);
00252 printf("for particle %d\n", i);
00253 dump_particles();
00254 endrun(1);
00255 }
00256 }
00257 }
00258 }
00259
00260
00261
00262 force_insert_pseudo_particles();
00263
00264
00265
00266 last = -1;
00267
00268 force_update_node_recursive(All.MaxPart, -1, -1);
00269
00270 if(last >= All.MaxPart)
00271 {
00272 if(last >= All.MaxPart + MaxNodes)
00273 Nextnode[last - MaxNodes] = -1;
00274 else
00275 Nodes[last].u.d.nextnode = -1;
00276 }
00277 else
00278 Nextnode[last] = -1;
00279
00280 return numnodes;
00281 }
00282
00283
00284
00292 void force_create_empty_nodes(int no, int topnode, int bits, int x, int y, int z, int *nodecount,
00293 int *nextfree)
00294 {
00295 int i, j, k, n, sub, count;
00296
00297 if(TopNodes[topnode].Daughter >= 0)
00298 {
00299 for(i = 0; i < 2; i++)
00300 for(j = 0; j < 2; j++)
00301 for(k = 0; k < 2; k++)
00302 {
00303 sub = 7 & peano_hilbert_key((x << 1) + i, (y << 1) + j, (z << 1) + k, bits);
00304
00305 count = i + 2 * j + 4 * k;
00306
00307 Nodes[no].u.suns[count] = *nextfree;
00308
00309
00310 Nodes[*nextfree].len = 0.5 * Nodes[no].len;
00311 Nodes[*nextfree].center[0] = Nodes[no].center[0] + (2 * i - 1) * 0.25 * Nodes[no].len;
00312 Nodes[*nextfree].center[1] = Nodes[no].center[1] + (2 * j - 1) * 0.25 * Nodes[no].len;
00313 Nodes[*nextfree].center[2] = Nodes[no].center[2] + (2 * k - 1) * 0.25 * Nodes[no].len;
00314
00315 for(n = 0; n < 8; n++)
00316 Nodes[*nextfree].u.suns[n] = -1;
00317
00318 if(TopNodes[TopNodes[topnode].Daughter + sub].Daughter == -1)
00319 DomainNodeIndex[TopNodes[TopNodes[topnode].Daughter + sub].Leaf] = *nextfree;
00320
00321 *nextfree = *nextfree + 1;
00322 *nodecount = *nodecount + 1;
00323
00324 if((*nodecount) >= MaxNodes)
00325 {
00326 printf("task %d: maximum number %d of tree-nodes reached.\n", ThisTask, MaxNodes);
00327 printf("in create empty nodes\n");
00328 dump_particles();
00329 endrun(11);
00330 }
00331
00332 force_create_empty_nodes(*nextfree - 1, TopNodes[topnode].Daughter + sub,
00333 bits + 1, 2 * x + i, 2 * y + j, 2 * z + k, nodecount, nextfree);
00334 }
00335 }
00336 }
00337
00338
00339
00346 void force_insert_pseudo_particles(void)
00347 {
00348 int i, index, subnode, nn, th;
00349
00350 for(i = 0; i < NTopleaves; i++)
00351 {
00352 index = DomainNodeIndex[i];
00353
00354 DomainMoment[i].mass = 0;
00355 DomainMoment[i].s[0] = Nodes[index].center[0];
00356 DomainMoment[i].s[1] = Nodes[index].center[1];
00357 DomainMoment[i].s[2] = Nodes[index].center[2];
00358 }
00359
00360 for(i = 0; i < NTopleaves; i++)
00361 {
00362 if(i < DomainMyStart || i > DomainMyLast)
00363 {
00364 th = All.MaxPart;
00365
00366 while(1)
00367 {
00368 if(th >= All.MaxPart)
00369 {
00370 if(th >= All.MaxPart + MaxNodes)
00371 endrun(888);
00372
00373 subnode = 0;
00374 if(DomainMoment[i].s[0] > Nodes[th].center[0])
00375 subnode += 1;
00376 if(DomainMoment[i].s[1] > Nodes[th].center[1])
00377 subnode += 2;
00378 if(DomainMoment[i].s[2] > Nodes[th].center[2])
00379 subnode += 4;
00380
00381 nn = Nodes[th].u.suns[subnode];
00382
00383 if(nn >= 0)
00384 {
00385 th = nn;
00386 }
00387 else
00388 {
00389
00390
00391
00392 Nodes[th].u.suns[subnode] = All.MaxPart + MaxNodes + i;
00393
00394 break;
00395 }
00396 }
00397 else
00398 {
00399 endrun(889);
00400 }
00401 }
00402 }
00403 }
00404 }
00405
00406
00422 void force_update_node_recursive(int no, int sib, int father)
00423 {
00424 int j, jj, p, pp, nextsib, suns[8];
00425 FLOAT hmax;
00426
00427 #ifdef UNEQUALSOFTENINGS
00428 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00429 int maxsofttype, diffsoftflag;
00430 #else
00431 FLOAT maxsoft;
00432 #endif
00433 #endif
00434 struct particle_data *pa;
00435 double s[3], vs[3], mass;
00436
00437 if(no >= All.MaxPart && no < All.MaxPart + MaxNodes)
00438 {
00439 for(j = 0; j < 8; j++)
00440 suns[j] = Nodes[no].u.suns[j];
00441
00442 if(last >= 0)
00443 {
00444 if(last >= All.MaxPart)
00445 {
00446 if(last >= All.MaxPart + MaxNodes)
00447 Nextnode[last - MaxNodes] = no;
00448 else
00449 Nodes[last].u.d.nextnode = no;
00450 }
00451 else
00452 Nextnode[last] = no;
00453 }
00454
00455 last = no;
00456
00457 mass = 0;
00458 s[0] = 0;
00459 s[1] = 0;
00460 s[2] = 0;
00461 vs[0] = 0;
00462 vs[1] = 0;
00463 vs[2] = 0;
00464 hmax = 0;
00465 #ifdef UNEQUALSOFTENINGS
00466 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00467 maxsofttype = 7;
00468 diffsoftflag = 0;
00469 #else
00470 maxsoft = 0;
00471 #endif
00472 #endif
00473
00474 for(j = 0; j < 8; j++)
00475 {
00476 if((p = suns[j]) >= 0)
00477 {
00478
00479 for(jj = j + 1; jj < 8; jj++)
00480 if((pp = suns[jj]) >= 0)
00481 break;
00482
00483 if(jj < 8)
00484 nextsib = pp;
00485 else
00486 nextsib = sib;
00487
00488 force_update_node_recursive(p, nextsib, no);
00489
00490
00491 if(p >= All.MaxPart)
00492 {
00493 if(p >= All.MaxPart + MaxNodes)
00494 {
00495
00496
00497
00498
00499 }
00500 else
00501 {
00502 mass += Nodes[p].u.d.mass;
00503 s[0] += Nodes[p].u.d.mass * Nodes[p].u.d.s[0];
00504 s[1] += Nodes[p].u.d.mass * Nodes[p].u.d.s[1];
00505 s[2] += Nodes[p].u.d.mass * Nodes[p].u.d.s[2];
00506 vs[0] += Nodes[p].u.d.mass * Extnodes[p].vs[0];
00507 vs[1] += Nodes[p].u.d.mass * Extnodes[p].vs[1];
00508 vs[2] += Nodes[p].u.d.mass * Extnodes[p].vs[2];
00509
00510 if(Extnodes[p].hmax > hmax)
00511 hmax = Extnodes[p].hmax;
00512
00513 #ifdef UNEQUALSOFTENINGS
00514 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00515 diffsoftflag |= (Nodes[p].u.d.bitflags >> 5) & 1;
00516
00517 if(maxsofttype == 7)
00518 {
00519 maxsofttype = (Nodes[p].u.d.bitflags >> 2) & 7;
00520 }
00521 else
00522 {
00523 if(((Nodes[p].u.d.bitflags >> 2) & 7) != 7)
00524 {
00525 if(All.ForceSoftening[((Nodes[p].u.d.bitflags >> 2) & 7)] >
00526 All.ForceSoftening[maxsofttype])
00527 {
00528 maxsofttype = ((Nodes[p].u.d.bitflags >> 2) & 7);
00529 diffsoftflag = 1;
00530 }
00531 else
00532 {
00533 if(All.ForceSoftening[((Nodes[p].u.d.bitflags >> 2) & 7)] <
00534 All.ForceSoftening[maxsofttype])
00535 diffsoftflag = 1;
00536 }
00537 }
00538 }
00539 #else
00540 if(Nodes[p].maxsoft > maxsoft)
00541 maxsoft = Nodes[p].maxsoft;
00542 #endif
00543 #endif
00544 }
00545 }
00546 else
00547 {
00548 pa = &P[p];
00549
00550 mass += pa->Mass;
00551 s[0] += pa->Mass * pa->Pos[0];
00552 s[1] += pa->Mass * pa->Pos[1];
00553 s[2] += pa->Mass * pa->Pos[2];
00554 vs[0] += pa->Mass * pa->Vel[0];
00555 vs[1] += pa->Mass * pa->Vel[1];
00556 vs[2] += pa->Mass * pa->Vel[2];
00557
00558 #ifdef UNEQUALSOFTENINGS
00559 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00560 if(maxsofttype == 7)
00561 {
00562 maxsofttype = pa->Type;
00563 }
00564 else
00565 {
00566 if(All.ForceSoftening[pa->Type] > All.ForceSoftening[maxsofttype])
00567 {
00568 maxsofttype = pa->Type;
00569 diffsoftflag = 1;
00570 }
00571 else
00572 {
00573 if(All.ForceSoftening[pa->Type] < All.ForceSoftening[maxsofttype])
00574 diffsoftflag = 1;
00575 }
00576 }
00577 #else
00578 if(pa->Type == 0)
00579 {
00580 if(SphP[p].Hsml > maxsoft)
00581 maxsoft = SphP[p].Hsml;
00582 }
00583 else
00584 {
00585 if(All.ForceSoftening[pa->Type] > maxsoft)
00586 maxsoft = All.ForceSoftening[pa->Type];
00587 }
00588 #endif
00589 #endif
00590 if(pa->Type == 0)
00591 if(SphP[p].Hsml > hmax)
00592 hmax = SphP[p].Hsml;
00593 }
00594 }
00595 }
00596
00597
00598 if(mass)
00599 {
00600 s[0] /= mass;
00601 s[1] /= mass;
00602 s[2] /= mass;
00603 vs[0] /= mass;
00604 vs[1] /= mass;
00605 vs[2] /= mass;
00606 }
00607 else
00608 {
00609 s[0] = Nodes[no].center[0];
00610 s[1] = Nodes[no].center[1];
00611 s[2] = Nodes[no].center[2];
00612 }
00613
00614 Nodes[no].u.d.s[0] = s[0];
00615 Nodes[no].u.d.s[1] = s[1];
00616 Nodes[no].u.d.s[2] = s[2];
00617 Nodes[no].u.d.mass = mass;
00618
00619
00620 #ifdef UNEQUALSOFTENINGS
00621 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00622 Nodes[no].u.d.bitflags = 4 * maxsofttype + 32 * diffsoftflag;
00623 #else
00624 Nodes[no].u.d.bitflags = 0;
00625 Nodes[no].maxsoft = maxsoft;
00626 #endif
00627 #else
00628 Nodes[no].u.d.bitflags = 0;
00629 #endif
00630
00631
00632 Extnodes[no].vs[0] = vs[0];
00633 Extnodes[no].vs[1] = vs[1];
00634 Extnodes[no].vs[2] = vs[2];
00635 Extnodes[no].hmax = hmax;
00636
00637 Nodes[no].u.d.sibling = sib;
00638 Nodes[no].u.d.father = father;
00639 }
00640 else
00641 {
00642 if(last >= 0)
00643 {
00644 if(last >= All.MaxPart)
00645 {
00646 if(last >= All.MaxPart + MaxNodes)
00647 Nextnode[last - MaxNodes] = no;
00648 else
00649 Nodes[last].u.d.nextnode = no;
00650 }
00651 else
00652 Nextnode[last] = no;
00653 }
00654
00655 last = no;
00656
00657 if(no < All.MaxPart)
00658 Father[no] = father;
00659 }
00660
00661 }
00662
00663
00664
00671 void force_update_pseudoparticles(void)
00672 {
00673 force_exchange_pseudodata();
00674
00675 force_treeupdate_pseudos();
00676 }
00677
00678
00679
00684 void force_exchange_pseudodata(void)
00685 {
00686 int i, no;
00687 MPI_Status status;
00688 int level, sendTask, recvTask;
00689
00690 for(i = DomainMyStart; i <= DomainMyLast; i++)
00691 {
00692 no = DomainNodeIndex[i];
00693
00694
00695 DomainMoment[i].s[0] = Nodes[no].u.d.s[0];
00696 DomainMoment[i].s[1] = Nodes[no].u.d.s[1];
00697 DomainMoment[i].s[2] = Nodes[no].u.d.s[2];
00698 DomainMoment[i].vs[0] = Extnodes[no].vs[0];
00699 DomainMoment[i].vs[1] = Extnodes[no].vs[1];
00700 DomainMoment[i].vs[2] = Extnodes[no].vs[2];
00701 DomainMoment[i].mass = Nodes[no].u.d.mass;
00702 #ifdef UNEQUALSOFTENINGS
00703 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00704 DomainMoment[i].bitflags = Nodes[no].u.d.bitflags;
00705 #else
00706 DomainMoment[i].maxsoft = Nodes[no].maxsoft;
00707 #endif
00708 #endif
00709 }
00710
00711
00712
00713 for(level = 1; level < (1 << PTask); level++)
00714 {
00715 sendTask = ThisTask;
00716 recvTask = ThisTask ^ level;
00717
00718 if(recvTask < NTask)
00719 MPI_Sendrecv(&DomainMoment[DomainStartList[sendTask]],
00720 (DomainEndList[sendTask] - DomainStartList[sendTask] + 1) * sizeof(struct DomainNODE),
00721 MPI_BYTE, recvTask, TAG_DMOM,
00722 &DomainMoment[DomainStartList[recvTask]],
00723 (DomainEndList[recvTask] - DomainStartList[recvTask] + 1) * sizeof(struct DomainNODE),
00724 MPI_BYTE, recvTask, TAG_DMOM, MPI_COMM_WORLD, &status);
00725 }
00726
00727 }
00728
00732 void force_treeupdate_pseudos(void)
00733 {
00734 int i, k, no;
00735 FLOAT sold[3], vsold[3], snew[3], vsnew[3], massold, massnew, mm;
00736
00737 #ifdef UNEQUALSOFTENINGS
00738 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00739 int maxsofttype, diffsoftflag;
00740 #else
00741 FLOAT maxsoft;
00742 #endif
00743 #endif
00744
00745 for(i = 0; i < NTopleaves; i++)
00746 if(i < DomainMyStart || i > DomainMyLast)
00747 {
00748 no = DomainNodeIndex[i];
00749
00750 for(k = 0; k < 3; k++)
00751 {
00752 sold[k] = Nodes[no].u.d.s[k];
00753 vsold[k] = Extnodes[no].vs[k];
00754 }
00755 massold = Nodes[no].u.d.mass;
00756
00757 for(k = 0; k < 3; k++)
00758 {
00759 snew[k] = DomainMoment[i].s[k];
00760 vsnew[k] = DomainMoment[i].vs[k];
00761 }
00762 massnew = DomainMoment[i].mass;
00763
00764
00765 #ifdef UNEQUALSOFTENINGS
00766 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00767 maxsofttype = (DomainMoment[i].bitflags >> 2) & 7;
00768 diffsoftflag = (DomainMoment[i].bitflags >> 5) & 1;
00769 #else
00770 maxsoft = DomainMoment[i].maxsoft;
00771 #endif
00772 #endif
00773 do
00774 {
00775 mm = Nodes[no].u.d.mass + massnew - massold;
00776 for(k = 0; k < 3; k++)
00777 {
00778 if(mm > 0)
00779 {
00780 Nodes[no].u.d.s[k] =
00781 (Nodes[no].u.d.mass * Nodes[no].u.d.s[k] + massnew * snew[k] - massold * sold[k]) / mm;
00782 Extnodes[no].vs[k] =
00783 (Nodes[no].u.d.mass * Extnodes[no].vs[k] + massnew * vsnew[k] -
00784 massold * vsold[k]) / mm;
00785 }
00786 }
00787 Nodes[no].u.d.mass = mm;
00788
00789
00790 #ifdef UNEQUALSOFTENINGS
00791 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
00792 diffsoftflag |= (Nodes[no].u.d.bitflags >> 5) & 1;
00793
00794 if(maxsofttype == 7)
00795 maxsofttype = (Nodes[no].u.d.bitflags >> 2) & 7;
00796 else
00797 {
00798 if(((Nodes[no].u.d.bitflags >> 2) & 7) != 7)
00799 {
00800 if(All.ForceSoftening[((Nodes[no].u.d.bitflags >> 2) & 7)] >
00801 All.ForceSoftening[maxsofttype])
00802 {
00803 maxsofttype = ((Nodes[no].u.d.bitflags >> 2) & 7);
00804 diffsoftflag = 1;
00805 }
00806 else
00807 {
00808 if(All.ForceSoftening[((Nodes[no].u.d.bitflags >> 2) & 7)] <
00809 All.ForceSoftening[maxsofttype])
00810 diffsoftflag = 1;
00811 }
00812 }
00813 }
00814
00815 Nodes[no].u.d.bitflags = (Nodes[no].u.d.bitflags & 3) + 4 * maxsofttype + 32 * diffsoftflag;
00816 #else
00817 if(Nodes[no].maxsoft < maxsoft)
00818 Nodes[no].maxsoft = maxsoft;
00819 maxsoft = Nodes[no].maxsoft;
00820 #endif
00821 #endif
00822 no = Nodes[no].u.d.father;
00823
00824 }
00825 while(no >= 0);
00826 }
00827 }
00828
00829
00830
00834 void force_flag_localnodes(void)
00835 {
00836 int no, i;
00837
00838
00839
00840 for(i = 0; i < NTopleaves; i++)
00841 {
00842 no = DomainNodeIndex[i];
00843
00844 while(no >= 0)
00845 {
00846 if((Nodes[no].u.d.bitflags & 1))
00847 break;
00848
00849 Nodes[no].u.d.bitflags |= 1;
00850
00851 no = Nodes[no].u.d.father;
00852 }
00853 }
00854
00855
00856
00857 for(i = DomainMyStart; i <= DomainMyLast; i++)
00858 {
00859
00860
00861
00862 {
00863 no = DomainNodeIndex[i];
00864
00865 while(no >= 0)
00866 {
00867 if((Nodes[no].u.d.bitflags & 2))
00868 break;
00869
00870 Nodes[no].u.d.bitflags |= 2;
00871
00872 no = Nodes[no].u.d.father;
00873 }
00874 }
00875 }
00876 }
00877
00878
00879
00885 void force_update_len(void)
00886 {
00887 int i, no;
00888 MPI_Status status;
00889 int level, sendTask, recvTask;
00890
00891 force_update_node_len_local();
00892
00893
00894 for(i = DomainMyStart; i <= DomainMyLast; i++)
00895 {
00896 no = DomainNodeIndex[i];
00897
00898 DomainTreeNodeLen[i] = Nodes[no].len;
00899 }
00900
00901 for(level = 1; level < (1 << PTask); level++)
00902 {
00903 sendTask = ThisTask;
00904 recvTask = ThisTask ^ level;
00905
00906 if(recvTask < NTask)
00907 MPI_Sendrecv(&DomainTreeNodeLen[DomainStartList[sendTask]],
00908 (DomainEndList[sendTask] - DomainStartList[sendTask] + 1) * sizeof(FLOAT),
00909 MPI_BYTE, recvTask, TAG_NODELEN,
00910 &DomainTreeNodeLen[DomainStartList[recvTask]],
00911 (DomainEndList[recvTask] - DomainStartList[recvTask] + 1) * sizeof(FLOAT),
00912 MPI_BYTE, recvTask, TAG_NODELEN, MPI_COMM_WORLD, &status);
00913 }
00914
00915
00916 force_update_node_len_toptree();
00917 }
00918
00919
00923 void force_update_node_len_local(void)
00924 {
00925 int i, p, k, no;
00926 FLOAT dist, distmax;
00927
00928 for(i = 0; i < NumPart; i++)
00929 {
00930 no = Father[i];
00931
00932 for(k = 0, distmax = 0; k < 3; k++)
00933 {
00934 dist = P[i].Pos[k] - Nodes[no].center[k];
00935 if(dist < 0)
00936 dist = -dist;
00937 if(dist > distmax)
00938 distmax = dist;
00939 }
00940
00941 if(distmax + distmax > Nodes[no].len)
00942 {
00943 Nodes[no].len = distmax + distmax;
00944 p = Nodes[no].u.d.father;
00945
00946 while(p >= 0)
00947 {
00948 distmax = Nodes[p].center[0] - Nodes[no].center[0];
00949 if(distmax < 0)
00950 distmax = -distmax;
00951 distmax = distmax + distmax + Nodes[no].len;
00952
00953 if(0.999999 * distmax > Nodes[p].len)
00954 {
00955 Nodes[p].len = distmax;
00956 no = p;
00957 p = Nodes[p].u.d.father;
00958 }
00959 else
00960 break;
00961 }
00962 }
00963 }
00964 }
00965
00966
00970 void force_update_node_len_toptree(void)
00971 {
00972 int i, no, p;
00973 FLOAT distmax;
00974
00975 for(i = 0; i < NTopleaves; i++)
00976 if(i < DomainMyStart || i > DomainMyLast)
00977 {
00978 no = DomainNodeIndex[i];
00979
00980 if(Nodes[no].len < DomainTreeNodeLen[i])
00981 Nodes[no].len = DomainTreeNodeLen[i];
00982
00983 p = Nodes[no].u.d.father;
00984
00985 while(p >= 0)
00986 {
00987 distmax = Nodes[p].center[0] - Nodes[no].center[0];
00988 if(distmax < 0)
00989 distmax = -distmax;
00990 distmax = distmax + distmax + Nodes[no].len;
00991
00992 if(0.999999 * distmax > Nodes[p].len)
00993 {
00994 Nodes[p].len = distmax;
00995 no = p;
00996 p = Nodes[p].u.d.father;
00997 }
00998 else
00999 break;
01000 }
01001 }
01002 }
01003
01004
01005
01006
01014 void force_update_hmax(void)
01015 {
01016 int i, no;
01017 MPI_Status status;
01018 int level, sendTask, recvTask;
01019
01020 force_update_node_hmax_local();
01021
01022 for(i = DomainMyStart; i <= DomainMyLast; i++)
01023 {
01024 no = DomainNodeIndex[i];
01025
01026 DomainHmax[i] = Extnodes[no].hmax;
01027 }
01028
01029
01030
01031 for(level = 1; level < (1 << PTask); level++)
01032 {
01033 sendTask = ThisTask;
01034 recvTask = ThisTask ^ level;
01035
01036 if(recvTask < NTask)
01037 MPI_Sendrecv(&DomainHmax[DomainStartList[sendTask]],
01038 (DomainEndList[sendTask] - DomainStartList[sendTask] + 1) * sizeof(FLOAT),
01039 MPI_BYTE, recvTask, TAG_HMAX,
01040 &DomainHmax[DomainStartList[recvTask]],
01041 (DomainEndList[recvTask] - DomainStartList[recvTask] + 1) * sizeof(FLOAT),
01042 MPI_BYTE, recvTask, TAG_HMAX, MPI_COMM_WORLD, &status);
01043 }
01044
01045
01046 force_update_node_hmax_toptree();
01047 }
01048
01051 void force_update_node_hmax_local(void)
01052 {
01053 int i, p, no;
01054
01055 for(i = 0; i < N_gas; i++)
01056 {
01057
01058 no = Father[i];
01059
01060 if(SphP[i].Hsml > Extnodes[no].hmax)
01061 {
01062
01063 Extnodes[no].hmax = SphP[i].Hsml;
01064 p = Nodes[no].u.d.father;
01065
01066 while(p >= 0)
01067 {
01068 if(Extnodes[no].hmax > Extnodes[p].hmax)
01069 {
01070 Extnodes[p].hmax = Extnodes[no].hmax;
01071 no = p;
01072 p = Nodes[p].u.d.father;
01073 }
01074 else
01075 break;
01076 }
01077 }
01078
01079 }
01080 }
01081
01082
01083
01084
01087 void force_update_node_hmax_toptree(void)
01088 {
01089
01090 int i, no, p;
01091
01092
01093 for(i = 0; i < NTopleaves; i++)
01094 if(i < DomainMyStart || i > DomainMyLast)
01095 {
01096 no = DomainNodeIndex[i];
01097
01098 if(Extnodes[no].hmax < DomainHmax[i])
01099 Extnodes[no].hmax = DomainHmax[i];
01100
01101 p = Nodes[no].u.d.father;
01102
01103 while(p >= 0)
01104 {
01105 if(Extnodes[no].hmax > Extnodes[p].hmax)
01106 {
01107 Extnodes[p].hmax = Extnodes[no].hmax;
01108 no = p;
01109 p = Nodes[p].u.d.father;
01110 }
01111 else
01112 break;
01113 }
01114 }
01115 }
01116
01117
01118
01119
01120
01126 int force_treeevaluate(int target, int mode, double *ewaldcountsum)
01127 {
01128 struct NODE *nop = 0;
01129 int no, ninteractions, ptype;
01130 double r2, dx, dy, dz, mass, r, fac, u, h, h_inv, h3_inv;
01131 double acc_x, acc_y, acc_z, pos_x, pos_y, pos_z, aold;
01132 #if defined(UNEQUALSOFTENINGS) && !defined(ADAPTIVE_GRAVSOFT_FORGAS)
01133 int maxsofttype;
01134 #endif
01135 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01136 double soft = 0;
01137 #endif
01138 #ifdef PERIODIC
01139 double boxsize, boxhalf;
01140
01141 boxsize = All.BoxSize;
01142 boxhalf = 0.5 * All.BoxSize;
01143 #endif
01144
01145 acc_x = 0;
01146 acc_y = 0;
01147 acc_z = 0;
01148 ninteractions = 0;
01149
01150 if(mode == 0)
01151 {
01152 pos_x = P[target].Pos[0];
01153 pos_y = P[target].Pos[1];
01154 pos_z = P[target].Pos[2];
01155 ptype = P[target].Type;
01156 aold = All.ErrTolForceAcc * P[target].OldAcc;
01157 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01158 if(ptype == 0)
01159 soft = SphP[target].Hsml;
01160 #endif
01161 }
01162 else
01163 {
01164 pos_x = GravDataGet[target].u.Pos[0];
01165 pos_y = GravDataGet[target].u.Pos[1];
01166 pos_z = GravDataGet[target].u.Pos[2];
01167 #ifdef UNEQUALSOFTENINGS
01168 ptype = GravDataGet[target].Type;
01169 #else
01170 ptype = P[0].Type;
01171 #endif
01172 aold = All.ErrTolForceAcc * GravDataGet[target].w.OldAcc;
01173 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01174 if(ptype == 0)
01175 soft = GravDataGet[target].Soft;
01176 #endif
01177 }
01178
01179
01180
01181 #ifndef UNEQUALSOFTENINGS
01182 h = All.ForceSoftening[ptype];
01183 h_inv = 1.0 / h;
01184 h3_inv = h_inv * h_inv * h_inv;
01185 #endif
01186 no = All.MaxPart;
01187
01188 while(no >= 0)
01189 {
01190 if(no < All.MaxPart)
01191 {
01192
01193
01194
01195 dx = P[no].Pos[0] - pos_x;
01196 dy = P[no].Pos[1] - pos_y;
01197 dz = P[no].Pos[2] - pos_z;
01198
01199 mass = P[no].Mass;
01200 }
01201 else
01202 {
01203 if(no >= All.MaxPart + MaxNodes)
01204 {
01205 if(mode == 0)
01206 {
01207 Exportflag[DomainTask[no - (All.MaxPart + MaxNodes)]] = 1;
01208 }
01209 no = Nextnode[no - MaxNodes];
01210 continue;
01211 }
01212 nop = &Nodes[no];
01213 dx = nop->u.d.s[0] - pos_x;
01214 dy = nop->u.d.s[1] - pos_y;
01215 dz = nop->u.d.s[2] - pos_z;
01216
01217 mass = nop->u.d.mass;
01218 }
01219 #ifdef PERIODIC
01220 dx = NEAREST(dx);
01221 dy = NEAREST(dy);
01222 dz = NEAREST(dz);
01223 #endif
01224 r2 = dx * dx + dy * dy + dz * dz;
01225
01226 if(no < All.MaxPart)
01227 {
01228 #ifdef UNEQUALSOFTENINGS
01229 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01230 if(ptype == 0)
01231 h = soft;
01232 else
01233 h = All.ForceSoftening[ptype];
01234
01235 if(P[no].Type == 0)
01236 {
01237 if(h < SphP[no].Hsml)
01238 h = SphP[no].Hsml;
01239 }
01240 else
01241 {
01242 if(h < All.ForceSoftening[P[no].Type])
01243 h = All.ForceSoftening[P[no].Type];
01244 }
01245 #else
01246 h = All.ForceSoftening[ptype];
01247 if(h < All.ForceSoftening[P[no].Type])
01248 h = All.ForceSoftening[P[no].Type];
01249 #endif
01250 #endif
01251 no = Nextnode[no];
01252 }
01253 else
01254 {
01255 if(mode == 1)
01256 {
01257 if((nop->u.d.bitflags & 3) == 1)
01258
01259
01260
01261 {
01262 no = nop->u.d.sibling;
01263 continue;
01264 }
01265 }
01266
01267
01268 if(All.ErrTolTheta)
01269 {
01270 if(nop->len * nop->len > r2 * All.ErrTolTheta * All.ErrTolTheta)
01271 {
01272
01273 no = nop->u.d.nextnode;
01274 continue;
01275 }
01276 }
01277 else
01278 {
01279 if(mass * nop->len * nop->len > r2 * r2 * aold)
01280 {
01281
01282 no = nop->u.d.nextnode;
01283 continue;
01284 }
01285
01286
01287
01288 if(fabs(nop->center[0] - pos_x) < 0.60 * nop->len)
01289 {
01290 if(fabs(nop->center[1] - pos_y) < 0.60 * nop->len)
01291 {
01292 if(fabs(nop->center[2] - pos_z) < 0.60 * nop->len)
01293 {
01294 no = nop->u.d.nextnode;
01295 continue;
01296 }
01297 }
01298 }
01299 }
01300
01301 #ifdef UNEQUALSOFTENINGS
01302 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
01303 h = All.ForceSoftening[ptype];
01304 maxsofttype = (nop->u.d.bitflags >> 2) & 7;
01305 if(maxsofttype == 7)
01306 {
01307 if(mass > 0)
01308 endrun(986);
01309 no = nop->u.d.nextnode;
01310 continue;
01311 }
01312 else
01313 {
01314 if(h < All.ForceSoftening[maxsofttype])
01315 {
01316 h = All.ForceSoftening[maxsofttype];
01317 if(r2 < h * h)
01318 {
01319 if(((nop->u.d.bitflags >> 5) & 1))
01320 {
01321 no = nop->u.d.nextnode;
01322 continue;
01323 }
01324 }
01325 }
01326 }
01327 #else
01328 if(ptype == 0)
01329 h = soft;
01330 else
01331 h = All.ForceSoftening[ptype];
01332
01333 if(h < nop->maxsoft)
01334 {
01335 h = nop->maxsoft;
01336 if(r2 < h * h)
01337 {
01338 no = nop->u.d.nextnode;
01339 continue;
01340 }
01341 }
01342 #endif
01343 #endif
01344
01345 no = nop->u.d.sibling;
01346
01347 if(mode == 1)
01348 {
01349 if(((nop->u.d.bitflags) & 1))
01350 continue;
01351 }
01352 }
01353
01354 r = sqrt(r2);
01355
01356 if(r >= h)
01357 fac = mass / (r2 * r);
01358 else
01359 {
01360 #ifdef UNEQUALSOFTENINGS
01361 h_inv = 1.0 / h;
01362 h3_inv = h_inv * h_inv * h_inv;
01363 #endif
01364 u = r * h_inv;
01365 if(u < 0.5)
01366 fac = mass * h3_inv * (10.666666666667 + u * u * (32.0 * u - 38.4));
01367 else
01368 fac =
01369 mass * h3_inv * (21.333333333333 - 48.0 * u +
01370 38.4 * u * u - 10.666666666667 * u * u * u - 0.066666666667 / (u * u * u));
01371 }
01372
01373 acc_x += dx * fac;
01374 acc_y += dy * fac;
01375 acc_z += dz * fac;
01376
01377 ninteractions++;
01378 }
01379
01380
01381
01382 if(mode == 0)
01383 {
01384 P[target].GravAccel[0] = acc_x;
01385 P[target].GravAccel[1] = acc_y;
01386 P[target].GravAccel[2] = acc_z;
01387 P[target].GravCost = ninteractions;
01388 }
01389 else
01390 {
01391 GravDataResult[target].u.Acc[0] = acc_x;
01392 GravDataResult[target].u.Acc[1] = acc_y;
01393 GravDataResult[target].u.Acc[2] = acc_z;
01394 GravDataResult[target].w.Ninteractions = ninteractions;
01395 }
01396
01397 #ifdef PERIODIC
01398 *ewaldcountsum += force_treeevaluate_ewald_correction(target, mode, pos_x, pos_y, pos_z, aold);
01399 #endif
01400
01401 return ninteractions;
01402 }
01403
01404
01405
01406
01407
01408
01409 #ifdef PMGRID
01410
01420 int force_treeevaluate_shortrange(int target, int mode)
01421 {
01422 struct NODE *nop = 0;
01423 int no, ptype, ninteractions, tabindex;
01424 double r2, dx, dy, dz, mass, r, fac, u, h, h_inv, h3_inv;
01425 double acc_x, acc_y, acc_z, pos_x, pos_y, pos_z, aold;
01426 double eff_dist;
01427 double rcut, asmth, asmthfac, rcut2, dist;
01428 #if defined(UNEQUALSOFTENINGS) && !defined(ADAPTIVE_GRAVSOFT_FORGAS)
01429 int maxsofttype;
01430 #endif
01431 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01432 double soft = 0;
01433 #endif
01434 #ifdef PERIODIC
01435 double boxsize, boxhalf;
01436
01437 boxsize = All.BoxSize;
01438 boxhalf = 0.5 * All.BoxSize;
01439 #endif
01440
01441
01442 acc_x = 0;
01443 acc_y = 0;
01444 acc_z = 0;
01445 ninteractions = 0;
01446
01447 if(mode == 0)
01448 {
01449 pos_x = P[target].Pos[0];
01450 pos_y = P[target].Pos[1];
01451 pos_z = P[target].Pos[2];
01452 ptype = P[target].Type;
01453 aold = All.ErrTolForceAcc * P[target].OldAcc;
01454 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01455 if(ptype == 0)
01456 soft = SphP[target].Hsml;
01457 #endif
01458 }
01459 else
01460 {
01461 pos_x = GravDataGet[target].u.Pos[0];
01462 pos_y = GravDataGet[target].u.Pos[1];
01463 pos_z = GravDataGet[target].u.Pos[2];
01464 #ifdef UNEQUALSOFTENINGS
01465 ptype = GravDataGet[target].Type;
01466 #else
01467 ptype = P[0].Type;
01468 #endif
01469 aold = All.ErrTolForceAcc * GravDataGet[target].w.OldAcc;
01470 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01471 if(ptype == 0)
01472 soft = GravDataGet[target].Soft;
01473 #endif
01474 }
01475
01476 rcut = All.Rcut[0];
01477 asmth = All.Asmth[0];
01478 #ifdef PLACEHIGHRESREGION
01479 if(((1 << ptype) & (PLACEHIGHRESREGION)))
01480 {
01481 rcut = All.Rcut[1];
01482 asmth = All.Asmth[1];
01483 }
01484 #endif
01485 rcut2 = rcut * rcut;
01486
01487 asmthfac = 0.5 / asmth * (NTAB / 3.0);
01488
01489 #ifndef UNEQUALSOFTENINGS
01490 h = All.ForceSoftening[ptype];
01491 h_inv = 1.0 / h;
01492 h3_inv = h_inv * h_inv * h_inv;
01493 #endif
01494 no = All.MaxPart;
01495
01496 while(no >= 0)
01497 {
01498 if(no < All.MaxPart)
01499 {
01500
01501 dx = P[no].Pos[0] - pos_x;
01502 dy = P[no].Pos[1] - pos_y;
01503 dz = P[no].Pos[2] - pos_z;
01504 #ifdef PERIODIC
01505 dx = NEAREST(dx);
01506 dy = NEAREST(dy);
01507 dz = NEAREST(dz);
01508 #endif
01509 r2 = dx * dx + dy * dy + dz * dz;
01510
01511 mass = P[no].Mass;
01512 #ifdef UNEQUALSOFTENINGS
01513 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
01514 if(ptype == 0)
01515 h = soft;
01516 else
01517 h = All.ForceSoftening[ptype];
01518
01519 if(P[no].Type == 0)
01520 {
01521 if(h < SphP[no].Hsml)
01522 h = SphP[no].Hsml;
01523 }
01524 else
01525 {
01526 if(h < All.ForceSoftening[P[no].Type])
01527 h = All.ForceSoftening[P[no].Type];
01528 }
01529 #else
01530 h = All.ForceSoftening[ptype];
01531 if(h < All.ForceSoftening[P[no].Type])
01532 h = All.ForceSoftening[P[no].Type];
01533 #endif
01534 #endif
01535 no = Nextnode[no];
01536 }
01537 else
01538 {
01539 if(no >= All.MaxPart + MaxNodes)
01540 {
01541 if(mode == 0)
01542 {
01543 Exportflag[DomainTask[no - (All.MaxPart + MaxNodes)]] = 1;
01544 }
01545 no = Nextnode[no - MaxNodes];
01546 continue;
01547 }
01548
01549 nop = &Nodes[no];
01550
01551 if(mode == 1)
01552 {
01553 if((nop->u.d.bitflags & 3) == 1)
01554
01555
01556
01557
01558 {
01559 no = nop->u.d.sibling;
01560 continue;
01561 }
01562 }
01563
01564 mass = nop->u.d.mass;
01565
01566 dx = nop->u.d.s[0] - pos_x;
01567 dy = nop->u.d.s[1] - pos_y;
01568 dz = nop->u.d.s[2] - pos_z;
01569 #ifdef PERIODIC
01570 dx = NEAREST(dx);
01571 dy = NEAREST(dy);
01572 dz = NEAREST(dz);
01573 #endif
01574 r2 = dx * dx + dy * dy + dz * dz;
01575
01576 if(r2 > rcut2)
01577 {
01578
01579 eff_dist = rcut + 0.5 * nop->len;
01580 #ifdef PERIODIC
01581 dist = NEAREST(nop->center[0] - pos_x);
01582 #else
01583 dist = nop->center[0] - pos_x;
01584 #endif
01585 if(dist < -eff_dist || dist > eff_dist)
01586 {
01587 no = nop->u.d.sibling;
01588 continue;
01589 }
01590 #ifdef PERIODIC
01591 dist = NEAREST(nop->center[1] - pos_y);
01592 #else
01593 dist = nop->center[1] - pos_y;
01594 #endif
01595 if(dist < -eff_dist || dist > eff_dist)
01596 {
01597 no = nop->u.d.sibling;
01598 continue;
01599 }
01600 #ifdef PERIODIC
01601 dist = NEAREST(nop->center[2] - pos_z);
01602 #else
01603 dist = nop->center[2] - pos_z;
01604 #endif
01605 if(dist < -eff_dist || dist > eff_dist)
01606 {
01607 no = nop->u.d.sibling;
01608 continue;
01609 }
01610 }
01611
01612
01613 if(All.ErrTolTheta)
01614 {
01615 if(nop->len * nop->len > r2 * All.ErrTolTheta * All.ErrTolTheta)
01616 {
01617
01618 no = nop->u.d.nextnode;
01619 continue;
01620 }
01621 }
01622 else
01623 {
01624 if(mass * nop->len * nop->len > r2 * r2 * aold)
01625 {
01626
01627 no = nop->u.d.nextnode;
01628 continue;
01629 }
01630
01631
01632
01633 if(fabs(nop->center[0] - pos_x) < 0.60 * nop->len)
01634 {
01635 if(fabs(nop->center[1] - pos_y) < 0.60 * nop->len)
01636 {
01637 if(fabs(nop->center[2] - pos_z) < 0.60 * nop->len)
01638 {
01639 no = nop->u.d.nextnode;
01640 continue;
01641 }
01642 }
01643 }
01644 }
01645
01646 #ifdef UNEQUALSOFTENINGS
01647 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
01648 h = All.ForceSoftening[ptype];
01649 maxsofttype = (nop->u.d.bitflags >> 2) & 7;
01650 if(maxsofttype == 7)
01651 {
01652 if(mass > 0)
01653 endrun(987);
01654 no = nop->u.d.nextnode;
01655 continue;
01656 }
01657 else
01658 {
01659 if(h < All.ForceSoftening[maxsofttype])
01660 {
01661 h = All.ForceSoftening[maxsofttype];
01662 if(r2 < h * h)
01663 {
01664 if(((nop->u.d.bitflags >> 5) & 1))
01665 {
01666 no = nop->u.d.nextnode;
01667
01668 continue;
01669 }
01670 }
01671 }
01672 }
01673 #else
01674 if(ptype == 0)
01675 h = soft;
01676 else
01677 h = All.ForceSoftening[ptype];
01678
01679 if(h < nop->maxsoft)
01680 {
01681 h = nop->maxsoft;
01682 if(r2 < h * h)
01683 {
01684 no = nop->u.d.nextnode;
01685 continue;
01686 }
01687 }
01688 #endif
01689 #endif
01690 no = nop->u.d.sibling;
01691
01692 if(mode == 1)
01693 {
01694 if(((nop->u.d.bitflags) & 1))
01695 continue;
01696 }
01697 }
01698
01699 r = sqrt(r2);
01700
01701 if(r >= h)
01702 fac = mass / (r2 * r);
01703 else
01704 {
01705 #ifdef UNEQUALSOFTENINGS
01706 h_inv = 1.0 / h;
01707 h3_inv = h_inv * h_inv * h_inv;
01708 #endif
01709 u = r * h_inv;
01710 if(u < 0.5)
01711 fac = mass * h3_inv * (10.666666666667 + u * u * (32.0 * u - 38.4));
01712 else
01713 fac =
01714 mass * h3_inv * (21.333333333333 - 48.0 * u +
01715 38.4 * u * u - 10.666666666667 * u * u * u - 0.066666666667 / (u * u * u));
01716 }
01717
01718 tabindex = (int) (asmthfac * r);
01719
01720 if(tabindex < NTAB)
01721 {
01722 fac *= shortrange_table[tabindex];
01723
01724 acc_x += dx * fac;
01725 acc_y += dy * fac;
01726 acc_z += dz * fac;
01727
01728 ninteractions++;
01729 }
01730 }
01731
01732
01733
01734 if(mode == 0)
01735 {
01736 P[target].GravAccel[0] = acc_x;
01737 P[target].GravAccel[1] = acc_y;
01738 P[target].GravAccel[2] = acc_z;
01739 P[target].GravCost = ninteractions;
01740 }
01741 else
01742 {
01743 GravDataResult[target].u.Acc[0] = acc_x;
01744 GravDataResult[target].u.Acc[1] = acc_y;
01745 GravDataResult[target].u.Acc[2] = acc_z;
01746 GravDataResult[target].w.Ninteractions = ninteractions;
01747 }
01748
01749 return ninteractions;
01750 }
01751
01752 #endif
01753
01754
01755
01756 #ifdef PERIODIC
01757
01775 int force_treeevaluate_ewald_correction(int target, int mode, double pos_x, double pos_y, double pos_z,
01776 double aold)
01777 {
01778 struct NODE *nop = 0;
01779 int no, cost;
01780 double dx, dy, dz, mass, r2;
01781 int signx, signy, signz;
01782 int i, j, k, openflag;
01783 double u, v, w;
01784 double f1, f2, f3, f4, f5, f6, f7, f8;
01785 double acc_x, acc_y, acc_z;
01786 double boxsize, boxhalf;
01787
01788 boxsize = All.BoxSize;
01789 boxhalf = 0.5 * All.BoxSize;
01790
01791 acc_x = 0;
01792 acc_y = 0;
01793 acc_z = 0;
01794 cost = 0;
01795
01796 no = All.MaxPart;
01797
01798 while(no >= 0)
01799 {
01800 if(no < All.MaxPart)
01801 {
01802
01803
01804
01805 dx = P[no].Pos[0] - pos_x;
01806 dy = P[no].Pos[1] - pos_y;
01807 dz = P[no].Pos[2] - pos_z;
01808 mass = P[no].Mass;
01809 }
01810 else
01811 {
01812 if(no >= All.MaxPart + MaxNodes)
01813 {
01814 if(mode == 0)
01815 {
01816 Exportflag[DomainTask[no - (All.MaxPart + MaxNodes)]] = 1;
01817 }
01818
01819 no = Nextnode[no - MaxNodes];
01820 continue;
01821 }
01822
01823 nop = &Nodes[no];
01824 dx = nop->u.d.s[0] - pos_x;
01825 dy = nop->u.d.s[1] - pos_y;
01826 dz = nop->u.d.s[2] - pos_z;
01827 mass = nop->u.d.mass;
01828 }
01829
01830 dx = NEAREST(dx);
01831 dy = NEAREST(dy);
01832 dz = NEAREST(dz);
01833
01834 if(no < All.MaxPart)
01835 no = Nextnode[no];
01836 else
01837 {
01838 openflag = 0;
01839
01840 r2 = dx * dx + dy * dy + dz * dz;
01841
01842 if(All.ErrTolTheta)
01843 {
01844 if(nop->len * nop->len > r2 * All.ErrTolTheta * All.ErrTolTheta)
01845 {
01846 openflag = 1;
01847 }
01848 }
01849 else
01850 {
01851 if(mass * nop->len * nop->len > r2 * r2 * aold)
01852 {
01853 openflag = 1;
01854 }
01855 else
01856 {
01857 if(fabs(nop->center[0] - pos_x) < 0.60 * nop->len)
01858 {
01859 if(fabs(nop->center[1] - pos_y) < 0.60 * nop->len)
01860 {
01861 if(fabs(nop->center[2] - pos_z) < 0.60 * nop->len)
01862 {
01863 openflag = 1;
01864 }
01865 }
01866 }
01867 }
01868 }
01869
01870 if(openflag)
01871 {
01872
01873
01874 u = nop->center[0] - pos_x;
01875 if(u > boxhalf)
01876 u -= boxsize;
01877 if(u < -boxhalf)
01878 u += boxsize;
01879
01880 if(fabs(u) > 0.5 * (boxsize - nop->len))
01881 {
01882 no = nop->u.d.nextnode;
01883 continue;
01884 }
01885
01886 u = nop->center[1] - pos_y;
01887 if(u > boxhalf)
01888 u -= boxsize;
01889 if(u < -boxhalf)
01890 u += boxsize;
01891
01892 if(fabs(u) > 0.5 * (boxsize - nop->len))
01893 {
01894 no = nop->u.d.nextnode;
01895 continue;
01896 }
01897
01898 u = nop->center[2] - pos_z;
01899 if(u > boxhalf)
01900 u -= boxsize;
01901 if(u < -boxhalf)
01902 u += boxsize;
01903
01904 if(fabs(u) > 0.5 * (boxsize - nop->len))
01905 {
01906 no = nop->u.d.nextnode;
01907 continue;
01908 }
01909
01910
01911
01912
01913 if(nop->len > 0.20 * boxsize)
01914 {
01915
01916 no = nop->u.d.nextnode;
01917 continue;
01918 }
01919 }
01920
01921 no = nop->u.d.sibling;
01922
01923 if(mode == 1)
01924 {
01925 if((nop->u.d.bitflags & 1))
01926 continue;
01927 }
01928 }
01929
01930
01931
01932 if(dx < 0)
01933 {
01934 dx = -dx;
01935 signx = +1;
01936 }
01937 else
01938 signx = -1;
01939
01940 if(dy < 0)
01941 {
01942 dy = -dy;
01943 signy = +1;
01944 }
01945 else
01946 signy = -1;
01947
01948 if(dz < 0)
01949 {
01950 dz = -dz;
01951 signz = +1;
01952 }
01953 else
01954 signz = -1;
01955
01956 u = dx * fac_intp;
01957 i = (int) u;
01958 if(i >= EN)
01959 i = EN - 1;
01960 u -= i;
01961 v = dy * fac_intp;
01962 j = (int) v;
01963 if(j >= EN)
01964 j = EN - 1;
01965 v -= j;
01966 w = dz * fac_intp;
01967 k = (int) w;
01968 if(k >= EN)
01969 k = EN - 1;
01970 w -= k;
01971
01972
01973
01974 f1 = (1 - u) * (1 - v) * (1 - w);
01975 f2 = (1 - u) * (1 - v) * (w);
01976 f3 = (1 - u) * (v) * (1 - w);
01977 f4 = (1 - u) * (v) * (w);
01978 f5 = (u) * (1 - v) * (1 - w);
01979 f6 = (u) * (1 - v) * (w);
01980 f7 = (u) * (v) * (1 - w);
01981 f8 = (u) * (v) * (w);
01982
01983 acc_x += mass * signx * (fcorrx[i][j][k] * f1 +
01984 fcorrx[i][j][k + 1] * f2 +
01985 fcorrx[i][j + 1][k] * f3 +
01986 fcorrx[i][j + 1][k + 1] * f4 +
01987 fcorrx[i + 1][j][k] * f5 +
01988 fcorrx[i + 1][j][k + 1] * f6 +
01989 fcorrx[i + 1][j + 1][k] * f7 + fcorrx[i + 1][j + 1][k + 1] * f8);
01990
01991 acc_y += mass * signy * (fcorry[i][j][k] * f1 +
01992 fcorry[i][j][k + 1] * f2 +
01993 fcorry[i][j + 1][k] * f3 +
01994 fcorry[i][j + 1][k + 1] * f4 +
01995 fcorry[i + 1][j][k] * f5 +
01996 fcorry[i + 1][j][k + 1] * f6 +
01997 fcorry[i + 1][j + 1][k] * f7 + fcorry[i + 1][j + 1][k + 1] * f8);
01998
01999 acc_z += mass * signz * (fcorrz[i][j][k] * f1 +
02000 fcorrz[i][j][k + 1] * f2 +
02001 fcorrz[i][j + 1][k] * f3 +
02002 fcorrz[i][j + 1][k + 1] * f4 +
02003 fcorrz[i + 1][j][k] * f5 +
02004 fcorrz[i + 1][j][k + 1] * f6 +
02005 fcorrz[i + 1][j + 1][k] * f7 + fcorrz[i + 1][j + 1][k + 1] * f8);
02006 cost++;
02007 }
02008
02009
02010
02011
02012 if(mode == 0)
02013 {
02014 P[target].GravAccel[0] += acc_x;
02015 P[target].GravAccel[1] += acc_y;
02016 P[target].GravAccel[2] += acc_z;
02017 P[target].GravCost += cost;
02018 }
02019 else
02020 {
02021 GravDataResult[target].u.Acc[0] += acc_x;
02022 GravDataResult[target].u.Acc[1] += acc_y;
02023 GravDataResult[target].u.Acc[2] += acc_z;
02024 GravDataResult[target].w.Ninteractions += cost;
02025 }
02026
02027 return cost;
02028 }
02029
02030 #endif
02031
02032
02033
02034
02035
02036
02041 void force_treeevaluate_potential(int target, int mode)
02042 {
02043 struct NODE *nop = 0;
02044 int no, ptype;
02045 double r2, dx, dy, dz, mass, r, u, h, h_inv, wp;
02046 double pot, pos_x, pos_y, pos_z, aold;
02047 #if defined(UNEQUALSOFTENINGS) && !defined(ADAPTIVE_GRAVSOFT_FORGAS)
02048 int maxsofttype;
02049 #endif
02050 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02051 double soft = 0;
02052 #endif
02053 #ifdef PERIODIC
02054 double boxsize, boxhalf;
02055
02056 boxsize = All.BoxSize;
02057 boxhalf = 0.5 * All.BoxSize;
02058 #endif
02059
02060 pot = 0;
02061
02062 if(mode == 0)
02063 {
02064 pos_x = P[target].Pos[0];
02065 pos_y = P[target].Pos[1];
02066 pos_z = P[target].Pos[2];
02067 ptype = P[target].Type;
02068 aold = All.ErrTolForceAcc * P[target].OldAcc;
02069 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02070 if(ptype == 0)
02071 soft = SphP[target].Hsml;
02072 #endif
02073 }
02074 else
02075 {
02076 pos_x = GravDataGet[target].u.Pos[0];
02077 pos_y = GravDataGet[target].u.Pos[1];
02078 pos_z = GravDataGet[target].u.Pos[2];
02079 #ifdef UNEQUALSOFTENINGS
02080 ptype = GravDataGet[target].Type;
02081 #else
02082 ptype = P[0].Type;
02083 #endif
02084 aold = All.ErrTolForceAcc * GravDataGet[target].w.OldAcc;
02085 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02086 if(ptype == 0)
02087 soft = GravDataGet[target].Soft;
02088 #endif
02089 }
02090
02091
02092 #ifndef UNEQUALSOFTENINGS
02093 h = All.ForceSoftening[ptype];
02094 h_inv = 1.0 / h;
02095 #endif
02096 no = All.MaxPart;
02097
02098 while(no >= 0)
02099 {
02100 if(no < All.MaxPart)
02101 {
02102
02103
02104
02105 dx = P[no].Pos[0] - pos_x;
02106 dy = P[no].Pos[1] - pos_y;
02107 dz = P[no].Pos[2] - pos_z;
02108 mass = P[no].Mass;
02109 }
02110 else
02111 {
02112 if(no >= All.MaxPart + MaxNodes)
02113 {
02114 if(mode == 0)
02115 {
02116 Exportflag[DomainTask[no - (All.MaxPart + MaxNodes)]] = 1;
02117 }
02118 no = Nextnode[no - MaxNodes];
02119 continue;
02120 }
02121
02122 nop = &Nodes[no];
02123 dx = nop->u.d.s[0] - pos_x;
02124 dy = nop->u.d.s[1] - pos_y;
02125 dz = nop->u.d.s[2] - pos_z;
02126 mass = nop->u.d.mass;
02127 }
02128
02129 #ifdef PERIODIC
02130 dx = NEAREST(dx);
02131 dy = NEAREST(dy);
02132 dz = NEAREST(dz);
02133 #endif
02134 r2 = dx * dx + dy * dy + dz * dz;
02135
02136 if(no < All.MaxPart)
02137 {
02138 #ifdef UNEQUALSOFTENINGS
02139 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02140 if(ptype == 0)
02141 h = soft;
02142 else
02143 h = All.ForceSoftening[ptype];
02144
02145 if(P[no].Type == 0)
02146 {
02147 if(h < SphP[no].Hsml)
02148 h = SphP[no].Hsml;
02149 }
02150 else
02151 {
02152 if(h < All.ForceSoftening[P[no].Type])
02153 h = All.ForceSoftening[P[no].Type];
02154 }
02155 #else
02156 h = All.ForceSoftening[ptype];
02157 if(h < All.ForceSoftening[P[no].Type])
02158 h = All.ForceSoftening[P[no].Type];
02159 #endif
02160 #endif
02161 no = Nextnode[no];
02162 }
02163 else
02164 {
02165 if(mode == 1)
02166 {
02167 if((nop->u.d.bitflags & 3) == 1)
02168
02169
02170
02171
02172 {
02173 no = nop->u.d.sibling;
02174 continue;
02175 }
02176 }
02177
02178 if(All.ErrTolTheta)
02179 {
02180 if(nop->len * nop->len > r2 * All.ErrTolTheta * All.ErrTolTheta)
02181 {
02182
02183 no = nop->u.d.nextnode;
02184 continue;
02185 }
02186 }
02187 else
02188 {
02189 if(mass * nop->len * nop->len > r2 * r2 * aold)
02190 {
02191
02192 no = nop->u.d.nextnode;
02193 continue;
02194 }
02195
02196 if(fabs(nop->center[0] - pos_x) < 0.60 * nop->len)
02197 {
02198 if(fabs(nop->center[1] - pos_y) < 0.60 * nop->len)
02199 {
02200 if(fabs(nop->center[2] - pos_z) < 0.60 * nop->len)
02201 {
02202 no = nop->u.d.nextnode;
02203 continue;
02204 }
02205 }
02206 }
02207 }
02208 #ifdef UNEQUALSOFTENINGS
02209 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
02210 h = All.ForceSoftening[ptype];
02211 maxsofttype = (nop->u.d.bitflags >> 2) & 7;
02212 if(maxsofttype == 7)
02213 {
02214 if(mass > 0)
02215 endrun(988);
02216 no = nop->u.d.nextnode;
02217 continue;
02218 }
02219 else
02220 {
02221 if(h < All.ForceSoftening[maxsofttype])
02222 {
02223 h = All.ForceSoftening[maxsofttype];
02224 if(r2 < h * h)
02225 {
02226 if(((nop->u.d.bitflags >> 5) & 1))
02227 {
02228 no = nop->u.d.nextnode;
02229 continue;
02230 }
02231 }
02232 }
02233 }
02234 #else
02235 if(ptype == 0)
02236 h = soft;
02237 else
02238 h = All.ForceSoftening[ptype];
02239
02240 if(h < nop->maxsoft)
02241 {
02242 h = nop->maxsoft;
02243 if(r2 < h * h)
02244 {
02245 no = nop->u.d.nextnode;
02246 continue;
02247 }
02248 }
02249 #endif
02250 #endif
02251
02252 no = nop->u.d.sibling;
02253
02254 if(mode == 1)
02255 {
02256 if(((nop->u.d.bitflags) & 1))
02257 continue;
02258 }
02259 }
02260
02261 r = sqrt(r2);
02262
02263 if(r >= h)
02264 pot -= mass / r;
02265 else
02266 {
02267 #ifdef UNEQUALSOFTENINGS
02268 h_inv = 1.0 / h;
02269 #endif
02270 u = r * h_inv;
02271
02272 if(u < 0.5)
02273 wp = -2.8 + u * u * (5.333333333333 + u * u * (6.4 * u - 9.6));
02274 else
02275 wp =
02276 -3.2 + 0.066666666667 / u + u * u * (10.666666666667 +
02277 u * (-16.0 + u * (9.6 - 2.133333333333 * u)));
02278
02279 pot += mass * h_inv * wp;
02280 }
02281 #ifdef PERIODIC
02282 pot += mass * ewald_pot_corr(dx, dy, dz);
02283 #endif
02284 }
02285
02286
02287
02288 if(mode == 0)
02289 P[target].Potential = pot;
02290 else
02291 GravDataResult[target].u.Potential = pot;
02292 }
02293
02294
02295
02296
02297 #ifdef PMGRID
02298
02302 void force_treeevaluate_potential_shortrange(int target, int mode)
02303 {
02304 struct NODE *nop = 0;
02305 int no, ptype, tabindex;
02306 double r2, dx, dy, dz, mass, r, u, h, h_inv, wp;
02307 double pot, pos_x, pos_y, pos_z, aold;
02308 double eff_dist, fac, rcut, asmth, asmthfac;
02309 double dxx, dyy, dzz;
02310 #if defined(UNEQUALSOFTENINGS) && !defined(ADAPTIVE_GRAVSOFT_FORGAS)
02311 int maxsofttype;
02312 #endif
02313 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02314 double soft = 0;
02315 #endif
02316
02317 #ifdef PERIODIC
02318 double boxsize, boxhalf;
02319
02320 boxsize = All.BoxSize;
02321 boxhalf = 0.5 * All.BoxSize;
02322 #endif
02323
02324 pot = 0;
02325
02326 if(mode == 0)
02327 {
02328 pos_x = P[target].Pos[0];
02329 pos_y = P[target].Pos[1];
02330 pos_z = P[target].Pos[2];
02331 ptype = P[target].Type;
02332 aold = All.ErrTolForceAcc * P[target].OldAcc;
02333 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02334 if(ptype == 0)
02335 soft = SphP[target].Hsml;
02336 #endif
02337 }
02338 else
02339 {
02340 pos_x = GravDataGet[target].u.Pos[0];
02341 pos_y = GravDataGet[target].u.Pos[1];
02342 pos_z = GravDataGet[target].u.Pos[2];
02343 #ifdef UNEQUALSOFTENINGS
02344 ptype = GravDataGet[target].Type;
02345 #else
02346 ptype = P[0].Type;
02347 #endif
02348 aold = All.ErrTolForceAcc * GravDataGet[target].w.OldAcc;
02349 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02350 if(ptype == 0)
02351 soft = GravDataGet[target].Soft;
02352 #endif
02353 }
02354
02355
02356 rcut = All.Rcut[0];
02357 asmth = All.Asmth[0];
02358 #ifdef PLACEHIGHRESREGION
02359 if(((1 << ptype) & (PLACEHIGHRESREGION)))
02360 {
02361 rcut = All.Rcut[1];
02362 asmth = All.Asmth[1];
02363 }
02364 #endif
02365 asmthfac = 0.5 / asmth * (NTAB / 3.0);
02366
02367 #ifndef UNEQUALSOFTENINGS
02368 h = All.ForceSoftening[ptype];
02369 h_inv = 1.0 / h;
02370 #endif
02371
02372 no = All.MaxPart;
02373
02374 while(no >= 0)
02375 {
02376 if(no < All.MaxPart)
02377 {
02378
02379
02380
02381 dx = P[no].Pos[0] - pos_x;
02382 dy = P[no].Pos[1] - pos_y;
02383 dz = P[no].Pos[2] - pos_z;
02384 mass = P[no].Mass;
02385 }
02386 else
02387 {
02388 if(no >= All.MaxPart + MaxNodes)
02389 {
02390 if(mode == 0)
02391 {
02392 Exportflag[DomainTask[no - (All.MaxPart + MaxNodes)]] = 1;
02393 }
02394 no = Nextnode[no - MaxNodes];
02395 continue;
02396 }
02397
02398 nop = &Nodes[no];
02399 dx = nop->u.d.s[0] - pos_x;
02400 dy = nop->u.d.s[1] - pos_y;
02401 dz = nop->u.d.s[2] - pos_z;
02402 mass = nop->u.d.mass;
02403 }
02404
02405 #ifdef PERIODIC
02406 dx = NEAREST(dx);
02407 dy = NEAREST(dy);
02408 dz = NEAREST(dz);
02409 #endif
02410 r2 = dx * dx + dy * dy + dz * dz;
02411
02412 if(no < All.MaxPart)
02413 {
02414 #ifdef UNEQUALSOFTENINGS
02415 #ifdef ADAPTIVE_GRAVSOFT_FORGAS
02416 if(ptype == 0)
02417 h = soft;
02418 else
02419 h = All.ForceSoftening[ptype];
02420
02421 if(P[no].Type == 0)
02422 {
02423 if(h < SphP[no].Hsml)
02424 h = SphP[no].Hsml;
02425 }
02426 else
02427 {
02428 if(h < All.ForceSoftening[P[no].Type])
02429 h = All.ForceSoftening[P[no].Type];
02430 }
02431 #else
02432 h = All.ForceSoftening[ptype];
02433 if(h < All.ForceSoftening[P[no].Type])
02434 h = All.ForceSoftening[P[no].Type];
02435 #endif
02436 #endif
02437 no = Nextnode[no];
02438 }
02439 else
02440 {
02441
02442 if(no >= All.MaxPart + MaxNodes)
02443 {
02444 if(mode == 0)
02445 {
02446 Exportflag[DomainTask[no - (All.MaxPart + MaxNodes)]] = 1;
02447 }
02448 no = Nextnode[no - MaxNodes];
02449 continue;
02450 }
02451
02452 if(mode == 1)
02453 {
02454 if((nop->u.d.bitflags & 3) == 1)
02455 {
02456 no = nop->u.d.sibling;
02457 continue;
02458 }
02459 }
02460
02461 eff_dist = rcut + 0.5 * nop->len;
02462
02463 dxx = nop->center[0] - pos_x;
02464 dyy = nop->center[1] - pos_y;
02465 dzz = nop->center[2] - pos_z;
02466 #ifdef PERIODIC
02467 dxx = NEAREST(dxx);
02468 dyy = NEAREST(dyy);
02469 dzz = NEAREST(dzz);
02470 #endif
02471 if(dxx < -eff_dist || dxx > eff_dist)
02472 {
02473 no = nop->u.d.sibling;
02474 continue;
02475 }
02476
02477 if(dyy < -eff_dist || dyy > eff_dist)
02478 {
02479 no = nop->u.d.sibling;
02480 continue;
02481 }
02482
02483 if(dzz < -eff_dist || dzz > eff_dist)
02484 {
02485 no = nop->u.d.sibling;
02486 continue;
02487 }
02488
02489 if(All.ErrTolTheta)
02490 {
02491 if(nop->len * nop->len > r2 * All.ErrTolTheta * All.ErrTolTheta)
02492 {
02493
02494 no = nop->u.d.nextnode;
02495 continue;
02496 }
02497 }
02498 else
02499 {
02500 if(mass * nop->len * nop->len > r2 * r2 * aold)
02501 {
02502
02503 no = nop->u.d.nextnode;
02504 continue;
02505 }
02506
02507 if(fabs(nop->center[0] - pos_x) < 0.60 * nop->len)
02508 {
02509 if(fabs(nop->center[1] - pos_y) < 0.60 * nop->len)
02510 {
02511 if(fabs(nop->center[2] - pos_z) < 0.60 * nop->len)
02512 {
02513 no = nop->u.d.nextnode;
02514 continue;
02515 }
02516 }
02517 }
02518 }
02519
02520 #ifdef UNEQUALSOFTENINGS
02521 #ifndef ADAPTIVE_GRAVSOFT_FORGAS
02522 h = All.ForceSoftening[ptype];
02523 maxsofttype = (nop->u.d.bitflags >> 2) & 7;
02524 if(maxsofttype == 7)
02525 {
02526 if(mass > 0)
02527 endrun(989);
02528 no = nop->u.d.nextnode;
02529 continue;
02530 }
02531 else
02532 {
02533 if(h < All.ForceSoftening[maxsofttype])
02534 {
02535 h = All.ForceSoftening[maxsofttype];
02536 if(r2 < h * h)
02537 {
02538
02539
02540
02541 if(((nop->u.d.bitflags >> 5) & 1))
02542 {
02543 no = nop->u.d.nextnode;
02544 continue;
02545 }
02546 }
02547 }
02548 }
02549 #else
02550 if(ptype == 0)
02551 h = soft;
02552 else
02553 h = All.ForceSoftening[ptype];
02554
02555 if(h < nop->maxsoft)
02556 {
02557 h = nop->maxsoft;
02558 if(r2 < h * h)
02559 {
02560 no = nop->u.d.nextnode;
02561 continue;
02562 }
02563 }
02564 #endif
02565 #endif
02566 no = nop->u.d.sibling;
02567
02568 if(mode == 1)
02569 {
02570 if(((nop->u.d.bitflags) & 1))
02571 continue;
02572 }
02573 }
02574
02575 r = sqrt(r2);
02576
02577 tabindex = (int) (r * asmthfac);
02578
02579 if(tabindex < NTAB)
02580 {
02581 fac = shortrange_table_potential[tabindex];
02582
02583 if(r >= h)
02584 pot -= fac * mass / r;
02585 else
02586 {
02587 #ifdef UNEQUALSOFTENINGS
02588 h_inv = 1.0 / h;
02589 #endif
02590 u = r * h_inv;
02591
02592 if(u < 0.5)
02593 wp = -2.8 + u * u * (5.333333333333 + u * u * (6.4 * u - 9.6));
02594 else
02595 wp =
02596 -3.2 + 0.066666666667 / u + u * u * (10.666666666667 +
02597 u * (-16.0 + u * (9.6 - 2.133333333333 * u)));
02598 pot += fac * mass * h_inv * wp;
02599 }
02600 }
02601 }
02602
02603
02604
02605 if(mode == 0)
02606 P[target].Potential = pot;
02607 else
02608 GravDataResult[target].u.Potential = pot;
02609 }
02610
02611 #endif
02612
02613
02614
02620 void force_treeallocate(int maxnodes, int maxpart)
02621 {
02622 int i;
02623 size_t bytes;
02624 double allbytes = 0;
02625 double u;
02626
02627 MaxNodes = maxnodes;
02628
02629 if(!(Nodes_base = malloc(bytes = (MaxNodes + 1) * sizeof(struct NODE))))
02630 {
02631 printf("failed to allocate memory for %d tree-nodes (%g MB).\n", MaxNodes, bytes / (1024.0 * 1024.0));
02632 endrun(3);
02633 }
02634 allbytes += bytes;
02635
02636 if(!(Extnodes_base = malloc(bytes = (MaxNodes + 1) * sizeof(struct extNODE))))
02637 {
02638 printf("failed to allocate memory for %d tree-extnodes (%g MB).\n", MaxNodes,
02639 bytes / (1024.0 * 1024.0));
02640 endrun(3);
02641 }
02642 allbytes += bytes;
02643
02644 Nodes = Nodes_base - All.MaxPart;
02645 Extnodes = Extnodes_base - All.MaxPart;
02646
02647 if(!(Nextnode = malloc(bytes = (maxpart + MAXTOPNODES) * sizeof(int))))
02648 {
02649 printf("Failed to allocate %d spaces for 'Nextnode' array (%g MB)\n", maxpart + MAXTOPNODES,
02650 bytes / (1024.0 * 1024.0));
02651 exit(0);
02652 }
02653 allbytes += bytes;
02654
02655 if(!(Father = malloc(bytes = (maxpart) * sizeof(int))))
02656 {
02657 printf("Failed to allocate %d spaces for 'Father' array (%g MB)\n", maxpart, bytes / (1024.0 * 1024.0));
02658 exit(0);
02659 }
02660 allbytes += bytes;
02661
02662 if(first_flag == 0)
02663 {
02664 first_flag = 1;
02665
02666 if(ThisTask == 0)
02667 printf("\nAllocated %g MByte for BH-tree. %d\n\n", allbytes / (1024.0 * 1024.0),
02668 sizeof(struct NODE) + sizeof(struct extNODE));
02669
02670 tabfac = NTAB / 3.0;
02671
02672 for(i = 0; i < NTAB; i++)
02673 {
02674 u = 3.0 / NTAB * (i + 0.5);
02675 shortrange_table[i] = erfc(u) + 2.0 * u / sqrt(M_PI) * exp(-u * u);
02676 shortrange_table_potential[i] = erfc(u);
02677 }
02678 }
02679 }
02680
02681
02685 void force_treefree(void)
02686 {
02687 free(Father);
02688 free(Nextnode);
02689 free(Extnodes_base);
02690 free(Nodes_base);
02691 }
02692
02693
02694
02695
02701 #ifdef FORCETEST
02702 int force_treeevaluate_direct(int target, int mode)
02703 {
02704 double epsilon;
02705 double h, h_inv, dx, dy, dz, r, r2, u, r_inv, fac;
02706 int i, ptype;
02707 double pos_x, pos_y, pos_z;
02708 double acc_x, acc_y, acc_z;
02709
02710 #ifdef PERIODIC
02711 double fcorr[3];
02712 #endif
02713 #ifdef PERIODIC
02714 double boxsize, boxhalf;
02715
02716 boxsize = All.BoxSize;
02717 boxhalf = 0.5 * All.BoxSize;
02718 #endif
02719
02720 acc_x = 0;
02721 acc_y = 0;
02722 acc_z = 0;
02723
02724 if(mode == 0)
02725 {
02726 pos_x = P[target].Pos[0];
02727 pos_y = P[target].Pos[1];
02728 pos_z = P[target].Pos[2];
02729 ptype = P[target].Type;
02730 }
02731 else
02732 {
02733 pos_x = GravDataGet[target].u.Pos[0];
02734 pos_y = GravDataGet[target].u.Pos[1];
02735 pos_z = GravDataGet[target].u.Pos[2];
02736 #ifdef UNEQUALSOFTENINGS
02737 ptype = GravDataGet[target].Type;
02738 #else
02739 ptype = P[0].Type;
02740 #endif
02741 }
02742
02743 for(i = 0; i < NumPart; i++)
02744 {
02745 epsilon = dmax(All.ForceSoftening[P[i].Type], All.ForceSoftening[ptype]);
02746
02747 h = epsilon;
02748 h_inv = 1 / h;
02749
02750 dx = P[i].Pos[0] - pos_x;
02751 dy = P[i].Pos[1] - pos_y;
02752 dz = P[i].Pos[2] - pos_z;
02753
02754 #ifdef PERIODIC
02755 while(dx > boxhalf)
02756 dx -= boxsize;
02757 while(dy > boxhalf)
02758 dy -= boxsize;
02759 while(dz > boxhalf)
02760 dz -= boxsize;
02761 while(dx < -boxhalf)
02762 dx += boxsize;
02763 while(dy < -boxhalf)
02764 dy += boxsize;
02765 while(dz < -boxhalf)
02766 dz += boxsize;
02767 #endif
02768 r2 = dx * dx + dy * dy + dz * dz;
02769
02770 r = sqrt(r2);
02771
02772 u = r * h_inv;
02773
02774 if(u >= 1)
02775 {
02776 r_inv = 1 / r;
02777
02778 fac = P[i].Mass * r_inv * r_inv * r_inv;
02779 }
02780 else
02781 {
02782 if(u < 0.5)
02783 fac = P[i].Mass * h_inv * h_inv * h_inv * (10.666666666667 + u * u * (32.0 * u - 38.4));
02784 else
02785 fac =
02786 P[i].Mass * h_inv * h_inv * h_inv * (21.333333333333 -
02787 48.0 * u + 38.4 * u * u -
02788 10.666666666667 * u * u *
02789 u - 0.066666666667 / (u * u * u));
02790 }
02791
02792 acc_x += dx * fac;
02793 acc_y += dy * fac;
02794 acc_z += dz * fac;
02795
02796 #ifdef PERIODIC
02797 if(u > 1.0e-5)
02798 {
02799 ewald_corr(dx, dy, dz, fcorr);
02800
02801 acc_x += P[i].Mass * fcorr[0];
02802 acc_y += P[i].Mass * fcorr[1];
02803 acc_z += P[i].Mass * fcorr[2];
02804 }
02805 #endif
02806 }
02807
02808
02809 if(mode == 0)
02810 {
02811 P[target].GravAccelDirect[0] = acc_x;
02812 P[target].GravAccelDirect[1] = acc_y;
02813 P[target].GravAccelDirect[2] = acc_z;
02814 }
02815 else
02816 {
02817 GravDataResult[target].u.Acc[0] = acc_x;
02818 GravDataResult[target].u.Acc[1] = acc_y;
02819 GravDataResult[target].u.Acc[2] = acc_z;
02820 }
02821
02822
02823 return NumPart;
02824 }
02825 #endif
02826
02827
02833 void dump_particles(void)
02834 {
02835 FILE *fd;
02836 char buffer[200];
02837 int i;
02838
02839 sprintf(buffer, "particles%d.dat", ThisTask);
02840 fd = fopen(buffer, "w");
02841 my_fwrite(&NumPart, 1, sizeof(int), fd);
02842
02843 for(i = 0; i < NumPart; i++)
02844 my_fwrite(&P[i].Pos[0], 3, sizeof(FLOAT), fd);
02845
02846 for(i = 0; i < NumPart; i++)
02847 my_fwrite(&P[i].Vel[0], 3, sizeof(FLOAT), fd);
02848
02849 for(i = 0; i < NumPart; i++)
02850 my_fwrite(&P[i].ID, 1, sizeof(int), fd);
02851
02852 fclose(fd);
02853 }
02854
02855
02856
02857 #ifdef PERIODIC
02858
02872 void ewald_init(void)
02873 {
02874 int i, j, k, beg, len, size, n, task, count;
02875 double x[3], force[3];
02876 char buf[200];
02877 FILE *fd;
02878
02879 if(ThisTask == 0)
02880 {
02881 printf("initialize Ewald correction...\n");
02882 fflush(stdout);
02883 }
02884
02885 #ifdef DOUBLEPRECISION
02886 sprintf(buf, "ewald_spc_table_%d_dbl.dat", EN);
02887 #else
02888 sprintf(buf, "ewald_spc_table_%d.dat", EN);
02889 #endif
02890
02891 if((fd = fopen(buf, "r")))
02892 {
02893 if(ThisTask == 0)
02894 {
02895 printf("\nreading Ewald tables from file `%s'\n", buf);
02896 fflush(stdout);
02897 }
02898
02899 my_fread(&fcorrx[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02900 my_fread(&fcorry[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02901 my_fread(&fcorrz[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02902 my_fread(&potcorr[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02903 fclose(fd);
02904 }
02905 else
02906 {
02907 if(ThisTask == 0)
02908 {
02909 printf("\nNo Ewald tables in file `%s' found.\nRecomputing them...\n", buf);
02910 fflush(stdout);
02911 }
02912
02913
02914
02915 size = (EN + 1) * (EN + 1) * (EN + 1) / NTask;
02916
02917
02918 beg = ThisTask * size;
02919 len = size;
02920 if(ThisTask == (NTask - 1))
02921 len = (EN + 1) * (EN + 1) * (EN + 1) - beg;
02922
02923 for(i = 0, count = 0; i <= EN; i++)
02924 for(j = 0; j <= EN; j++)
02925 for(k = 0; k <= EN; k++)
02926 {
02927 n = (i * (EN + 1) + j) * (EN + 1) + k;
02928 if(n >= beg && n < (beg + len))
02929 {
02930 if(ThisTask == 0)
02931 {
02932 if((count % (len / 20)) == 0)
02933 {
02934 printf("%4.1f percent done\n", count / (len / 100.0));
02935 fflush(stdout);
02936 }
02937 }
02938
02939 x[0] = 0.5 * ((double) i) / EN;
02940 x[1] = 0.5 * ((double) j) / EN;
02941 x[2] = 0.5 * ((double) k) / EN;
02942
02943 ewald_force(i, j, k, x, force);
02944
02945 fcorrx[i][j][k] = force[0];
02946 fcorry[i][j][k] = force[1];
02947 fcorrz[i][j][k] = force[2];
02948
02949 if(i + j + k == 0)
02950 potcorr[i][j][k] = 2.8372975;
02951 else
02952 potcorr[i][j][k] = ewald_psi(x);
02953
02954 count++;
02955 }
02956 }
02957
02958 for(task = 0; task < NTask; task++)
02959 {
02960 beg = task * size;
02961 len = size;
02962 if(task == (NTask - 1))
02963 len = (EN + 1) * (EN + 1) * (EN + 1) - beg;
02964
02965 #ifdef DOUBLEPRECISION
02966 MPI_Bcast(&fcorrx[0][0][beg], len, MPI_DOUBLE, task, MPI_COMM_WORLD);
02967 MPI_Bcast(&fcorry[0][0][beg], len, MPI_DOUBLE, task, MPI_COMM_WORLD);
02968 MPI_Bcast(&fcorrz[0][0][beg], len, MPI_DOUBLE, task, MPI_COMM_WORLD);
02969 MPI_Bcast(&potcorr[0][0][beg], len, MPI_DOUBLE, task, MPI_COMM_WORLD);
02970 #else
02971 MPI_Bcast(&fcorrx[0][0][beg], len, MPI_FLOAT, task, MPI_COMM_WORLD);
02972 MPI_Bcast(&fcorry[0][0][beg], len, MPI_FLOAT, task, MPI_COMM_WORLD);
02973 MPI_Bcast(&fcorrz[0][0][beg], len, MPI_FLOAT, task, MPI_COMM_WORLD);
02974 MPI_Bcast(&potcorr[0][0][beg], len, MPI_FLOAT, task, MPI_COMM_WORLD);
02975 #endif
02976 }
02977
02978 if(ThisTask == 0)
02979 {
02980 printf("\nwriting Ewald tables to file `%s'\n", buf);
02981 fflush(stdout);
02982
02983 if((fd = fopen(buf, "w")))
02984 {
02985 my_fwrite(&fcorrx[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02986 my_fwrite(&fcorry[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02987 my_fwrite(&fcorrz[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02988 my_fwrite(&potcorr[0][0][0], sizeof(FLOAT), (EN + 1) * (EN + 1) * (EN + 1), fd);
02989 fclose(fd);
02990 }
02991 }
02992 }
02993
02994 fac_intp = 2 * EN / All.BoxSize;
02995
02996 for(i = 0; i <= EN; i++)
02997 for(j = 0; j <= EN; j++)
02998 for(k = 0; k <= EN; k++)
02999 {
03000 potcorr[i][j][k] /= All.BoxSize;
03001 fcorrx[i][j][k] /= All.BoxSize * All.BoxSize;
03002 fcorry[i][j][k] /= All.BoxSize * All.BoxSize;
03003 fcorrz[i][j][k] /= All.BoxSize * All.BoxSize;
03004 }
03005
03006 if(ThisTask == 0)
03007 {
03008 printf("initialization of periodic boundaries finished.\n");
03009 fflush(stdout);
03010 }
03011 }
03012
03013
03020 #ifdef FORCETEST
03021 void ewald_corr(double dx, double dy, double dz, double *fper)
03022 {
03023 int signx, signy, signz;
03024 int i, j, k;
03025 double u, v, w;
03026 double f1, f2, f3, f4, f5, f6, f7, f8;
03027
03028 if(dx < 0)
03029 {
03030 dx = -dx;
03031 signx = +1;
03032 }
03033 else
03034 signx = -1;
03035
03036 if(dy < 0)
03037 {
03038 dy = -dy;
03039 signy = +1;
03040 }
03041 else
03042 signy = -1;
03043
03044 if(dz < 0)
03045 {
03046 dz = -dz;
03047 signz = +1;
03048 }
03049 else
03050 signz = -1;
03051
03052 u = dx * fac_intp;
03053 i = (int) u;
03054 if(i >= EN)
03055 i = EN - 1;
03056 u -= i;
03057 v = dy * fac_intp;
03058 j = (int) v;
03059 if(j >= EN)
03060 j = EN - 1;
03061 v -= j;
03062 w = dz * fac_intp;
03063 k = (int) w;
03064 if(k >= EN)
03065 k = EN - 1;
03066 w -= k;
03067
03068 f1 = (1 - u) * (1 - v) * (1 - w);
03069 f2 = (1 - u) * (1 - v) * (w);
03070 f3 = (1 - u) * (v) * (1 - w);
03071 f4 = (1 - u) * (v) * (w);
03072 f5 = (u) * (1 - v) * (1 - w);
03073 f6 = (u) * (1 - v) * (w);
03074 f7 = (u) * (v) * (1 - w);
03075 f8 = (u) * (v) * (w);
03076
03077 fper[0] = signx * (fcorrx[i][j][k] * f1 +
03078 fcorrx[i][j][k + 1] * f2 +
03079 fcorrx[i][j + 1][k] * f3 +
03080 fcorrx[i][j + 1][k + 1] * f4 +
03081 fcorrx[i + 1][j][k] * f5 +
03082 fcorrx[i + 1][j][k + 1] * f6 +
03083 fcorrx[i + 1][j + 1][k] * f7 + fcorrx[i + 1][j + 1][k + 1] * f8);
03084
03085 fper[1] = signy * (fcorry[i][j][k] * f1 +
03086 fcorry[i][j][k + 1] * f2 +
03087 fcorry[i][j + 1][k] * f3 +
03088 fcorry[i][j + 1][k + 1] * f4 +
03089 fcorry[i + 1][j][k] * f5 +
03090 fcorry[i + 1][j][k + 1] * f6 +
03091 fcorry[i + 1][j + 1][k] * f7 + fcorry[i + 1][j + 1][k + 1] * f8);
03092
03093 fper[2] = signz * (fcorrz[i][j][k] * f1 +
03094 fcorrz[i][j][k + 1] * f2 +
03095 fcorrz[i][j + 1][k] * f3 +
03096 fcorrz[i][j + 1][k + 1] * f4 +
03097 fcorrz[i + 1][j][k] * f5 +
03098 fcorrz[i + 1][j][k + 1] * f6 +
03099 fcorrz[i + 1][j + 1][k] * f7 + fcorrz[i + 1][j + 1][k + 1] * f8);
03100 }
03101 #endif
03102
03103
03110 double ewald_pot_corr(double dx, double dy, double dz)
03111 {
03112 int i, j, k;
03113 double u, v, w;
03114 double f1, f2, f3, f4, f5, f6, f7, f8;
03115
03116 if(dx < 0)
03117 dx = -dx;
03118
03119 if(dy < 0)
03120 dy = -dy;
03121
03122 if(dz < 0)
03123 dz = -dz;
03124
03125 u = dx * fac_intp;
03126 i = (int) u;
03127 if(i >= EN)
03128 i = EN - 1;
03129 u -= i;
03130 v = dy * fac_intp;
03131 j = (int) v;
03132 if(j >= EN)
03133 j = EN - 1;
03134 v -= j;
03135 w = dz * fac_intp;
03136 k = (int) w;
03137 if(k >= EN)
03138 k = EN - 1;
03139 w -= k;
03140
03141 f1 = (1 - u) * (1 - v) * (1 - w);
03142 f2 = (1 - u) * (1 - v) * (w);
03143 f3 = (1 - u) * (v) * (1 - w);
03144 f4 = (1 - u) * (v) * (w);
03145 f5 = (u) * (1 - v) * (1 - w);
03146 f6 = (u) * (1 - v) * (w);
03147 f7 = (u) * (v) * (1 - w);
03148 f8 = (u) * (v) * (w);
03149
03150 return potcorr[i][j][k] * f1 +
03151 potcorr[i][j][k + 1] * f2 +
03152 potcorr[i][j + 1][k] * f3 +
03153 potcorr[i][j + 1][k + 1] * f4 +
03154 potcorr[i + 1][j][k] * f5 +
03155 potcorr[i + 1][j][k + 1] * f6 + potcorr[i + 1][j + 1][k] * f7 + potcorr[i + 1][j + 1][k + 1] * f8;
03156 }
03157
03158
03159
03163 double ewald_psi(double x[3])
03164 {
03165 double alpha, psi;
03166 double r, sum1, sum2, hdotx;
03167 double dx[3];
03168 int i, n[3], h[3], h2;
03169
03170 alpha = 2.0;
03171
03172 for(n[0] = -4, sum1 = 0; n[0] <= 4; n[0]++)
03173 for(n[1] = -4; n[1] <= 4; n[1]++)
03174 for(n[2] = -4; n[2] <= 4; n[2]++)
03175 {
03176 for(i = 0; i < 3; i++)
03177 dx[i] = x[i] - n[i];
03178
03179 r = sqrt(dx[0] * dx[0] + dx[1] * dx[1] + dx[2] * dx[2]);
03180 sum1 += erfc(alpha * r) / r;
03181 }
03182
03183 for(h[0] = -4, sum2 = 0; h[0] <= 4; h[0]++)
03184 for(h[1] = -4; h[1] <= 4; h[1]++)
03185 for(h[2] = -4; h[2] <= 4; h[2]++)
03186 {
03187 hdotx = x[0] * h[0] + x[1] * h[1] + x[2] * h[2];
03188 h2 = h[0] * h[0] + h[1] * h[1] + h[2] * h[2];
03189 if(h2 > 0)
03190 sum2 += 1 / (M_PI * h2) * exp(-M_PI * M_PI * h2 / (alpha * alpha)) * cos(2 * M_PI * hdotx);
03191 }
03192
03193 r = sqrt(x[0] * x[0] + x[1] * x[1] + x[2] * x[2]);
03194
03195 psi = M_PI / (alpha * alpha) - sum1 - sum2 + 1 / r;
03196
03197 return psi;
03198 }
03199
03200
03204 void ewald_force(int iii, int jjj, int kkk, double x[3], double force[3])
03205 {
03206 double alpha, r2;
03207 double r, val, hdotx, dx[3];
03208 int i, h[3], n[3], h2;
03209
03210 alpha = 2.0;
03211
03212 for(i = 0; i < 3; i++)
03213 force[i] = 0;
03214
03215 if(iii == 0 && jjj == 0 && kkk == 0)
03216 return;
03217
03218 r2 = x[0] * x[0] + x[1] * x[1] + x[2] * x[2];
03219
03220 for(i = 0; i < 3; i++)
03221 force[i] += x[i] / (r2 * sqrt(r2));
03222
03223 for(n[0] = -4; n[0] <= 4; n[0]++)
03224 for(n[1] = -4; n[1] <= 4; n[1]++)
03225 for(n[2] = -4; n[2] <= 4; n[2]++)
03226 {
03227 for(i = 0; i < 3; i++)
03228 dx[i] = x[i] - n[i];
03229
03230 r = sqrt(dx[0] * dx[0] + dx[1] * dx[1] + dx[2] * dx[2]);
03231
03232 val = erfc(alpha * r) + 2 * alpha * r / sqrt(M_PI) * exp(-alpha * alpha * r * r);
03233
03234 for(i = 0; i < 3; i++)
03235 force[i] -= dx[i] / (r * r * r) * val;
03236 }
03237
03238 for(h[0] = -4; h[0] <= 4; h[0]++)
03239 for(h[1] = -4; h[1] <= 4; h[1]++)
03240 for(h[2] = -4; h[2] <= 4; h[2]++)
03241 {
03242 hdotx = x[0] * h[0] + x[1] * h[1] + x[2] * h[2];
03243 h2 = h[0] * h[0] + h[1] * h[1] + h[2] * h[2];
03244
03245 if(h2 > 0)
03246 {
03247 val = 2.0 / ((double) h2) * exp(-M_PI * M_PI * h2 / (alpha * alpha)) * sin(2 * M_PI * hdotx);
03248
03249 for(i = 0; i < 3; i++)
03250 force[i] -= h[i] * val;
03251 }
03252 }
03253 }
03254
03255 #endif