Source

enzo-timing-woc / src / enzo / Grid.h

The trunk branch has multiple heads

   1
   2
   3
   4
   5
   6
   7
   8
   9
  10
  11
  12
  13
  14
  15
  16
  17
  18
  19
  20
  21
  22
  23
  24
  25
  26
  27
  28
  29
  30
  31
  32
  33
  34
  35
  36
  37
  38
  39
  40
  41
  42
  43
  44
  45
  46
  47
  48
  49
  50
  51
  52
  53
  54
  55
  56
  57
  58
  59
  60
  61
  62
  63
  64
  65
  66
  67
  68
  69
  70
  71
  72
  73
  74
  75
  76
  77
  78
  79
  80
  81
  82
  83
  84
  85
  86
  87
  88
  89
  90
  91
  92
  93
  94
  95
  96
  97
  98
  99
 100
 101
 102
 103
 104
 105
 106
 107
 108
 109
 110
 111
 112
 113
 114
 115
 116
 117
 118
 119
 120
 121
 122
 123
 124
 125
 126
 127
 128
 129
 130
 131
 132
 133
 134
 135
 136
 137
 138
 139
 140
 141
 142
 143
 144
 145
 146
 147
 148
 149
 150
 151
 152
 153
 154
 155
 156
 157
 158
 159
 160
 161
 162
 163
 164
 165
 166
 167
 168
 169
 170
 171
 172
 173
 174
 175
 176
 177
 178
 179
 180
 181
 182
 183
 184
 185
 186
 187
 188
 189
 190
 191
 192
 193
 194
 195
 196
 197
 198
 199
 200
 201
 202
 203
 204
 205
 206
 207
 208
 209
 210
 211
 212
 213
 214
 215
 216
 217
 218
 219
 220
 221
 222
 223
 224
 225
 226
 227
 228
 229
 230
 231
 232
 233
 234
 235
 236
 237
 238
 239
 240
 241
 242
 243
 244
 245
 246
 247
 248
 249
 250
 251
 252
 253
 254
 255
 256
 257
 258
 259
 260
 261
 262
 263
 264
 265
 266
 267
 268
 269
 270
 271
 272
 273
 274
 275
 276
 277
 278
 279
 280
 281
 282
 283
 284
 285
 286
 287
 288
 289
 290
 291
 292
 293
 294
 295
 296
 297
 298
 299
 300
 301
 302
 303
 304
 305
 306
 307
 308
 309
 310
 311
 312
 313
 314
 315
 316
 317
 318
 319
 320
 321
 322
 323
 324
 325
 326
 327
 328
 329
 330
 331
 332
 333
 334
 335
 336
 337
 338
 339
 340
 341
 342
 343
 344
 345
 346
 347
 348
 349
 350
 351
 352
 353
 354
 355
 356
 357
 358
 359
 360
 361
 362
 363
 364
 365
 366
 367
 368
 369
 370
 371
 372
 373
 374
 375
 376
 377
 378
 379
 380
 381
 382
 383
 384
 385
 386
 387
 388
 389
 390
 391
 392
 393
 394
 395
 396
 397
 398
 399
 400
 401
 402
 403
 404
 405
 406
 407
 408
 409
 410
 411
 412
 413
 414
 415
 416
 417
 418
 419
 420
 421
 422
 423
 424
 425
 426
 427
 428
 429
 430
 431
 432
 433
 434
 435
 436
 437
 438
 439
 440
 441
 442
 443
 444
 445
 446
 447
 448
 449
 450
 451
 452
 453
 454
 455
 456
 457
 458
 459
 460
 461
 462
 463
 464
 465
 466
 467
 468
 469
 470
 471
 472
 473
 474
 475
 476
 477
 478
 479
 480
 481
 482
 483
 484
 485
 486
 487
 488
 489
 490
 491
 492
 493
 494
 495
 496
 497
 498
 499
 500
 501
 502
 503
 504
 505
 506
 507
 508
 509
 510
 511
 512
 513
 514
 515
 516
 517
 518
 519
 520
 521
 522
 523
 524
 525
 526
 527
 528
 529
 530
 531
 532
 533
 534
 535
 536
 537
 538
 539
 540
 541
 542
 543
 544
 545
 546
 547
 548
 549
 550
 551
 552
 553
 554
 555
 556
 557
 558
 559
 560
 561
 562
 563
 564
 565
 566
 567
 568
 569
 570
 571
 572
 573
 574
 575
 576
 577
 578
 579
 580
 581
 582
 583
 584
 585
 586
 587
 588
 589
 590
 591
 592
 593
 594
 595
 596
 597
 598
 599
 600
 601
 602
 603
 604
 605
 606
 607
 608
 609
 610
 611
 612
 613
 614
 615
 616
 617
 618
 619
 620
 621
 622
 623
 624
 625
 626
 627
 628
 629
 630
 631
 632
 633
 634
 635
 636
 637
 638
 639
 640
 641
 642
 643
 644
 645
 646
 647
 648
 649
 650
 651
 652
 653
 654
 655
 656
 657
 658
 659
 660
 661
 662
 663
 664
 665
 666
 667
 668
 669
 670
 671
 672
 673
 674
 675
 676
 677
 678
 679
 680
 681
 682
 683
 684
 685
 686
 687
 688
 689
 690
 691
 692
 693
 694
 695
 696
 697
 698
 699
 700
 701
 702
 703
 704
 705
 706
 707
 708
 709
 710
 711
 712
 713
 714
 715
 716
 717
 718
 719
 720
 721
 722
 723
 724
 725
 726
 727
 728
 729
 730
 731
 732
 733
 734
 735
 736
 737
 738
 739
 740
 741
 742
 743
 744
 745
 746
 747
 748
 749
 750
 751
 752
 753
 754
 755
 756
 757
 758
 759
 760
 761
 762
 763
 764
 765
 766
 767
 768
 769
 770
 771
 772
 773
 774
 775
 776
 777
 778
 779
 780
 781
 782
 783
 784
 785
 786
 787
 788
 789
 790
 791
 792
 793
 794
 795
 796
 797
 798
 799
 800
 801
 802
 803
 804
 805
 806
 807
 808
 809
 810
 811
 812
 813
 814
 815
 816
 817
 818
 819
 820
 821
 822
 823
 824
 825
 826
 827
 828
 829
 830
 831
 832
 833
 834
 835
 836
 837
 838
 839
 840
 841
 842
 843
 844
 845
 846
 847
 848
 849
 850
 851
 852
 853
 854
 855
 856
 857
 858
 859
 860
 861
 862
 863
 864
 865
 866
 867
 868
 869
 870
 871
 872
 873
 874
 875
 876
 877
 878
 879
 880
 881
 882
 883
 884
 885
 886
 887
 888
 889
 890
 891
 892
 893
 894
 895
 896
 897
 898
 899
 900
 901
 902
 903
 904
 905
 906
 907
 908
 909
 910
 911
 912
 913
 914
 915
 916
 917
 918
 919
 920
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
/***********************************************************************
/
/  GRID CLASS
/
/  written by: Greg Bryan
/  date:       November, 1994
/  modified:   Many times by AK, DC, RH, JB, DR...
/
/  PURPOSE:
/
************************************************************************/

#ifndef GRID_DEFINED__
#define GRID_DEFINED__


#include "ProtoSubgrid.h"
#include "ListOfParticles.h"
#include "region.h"
#include "FastSiblingLocator.h"
#include "Star.h"

#ifdef FLUX_FIX
#include "TopGridData.h"
#endif

struct HierarchyEntry;

#include "EnzoArray.h"

#ifdef ANALYSIS_TOOLS
#   include "AnalyzeClusters.h"
#endif

#ifdef TRANSFER
#include "PhotonPackage.h"
#include "ListOfPhotonsToMove.h"
#endif /* TRANSFER */

//extern int CommunicationDirection;

extern int CommunicationDirection;
int FindField(int f, int farray[], int n);
struct LevelHierarchyEntry;

class grid
{
 private:
//
//  General grid class data
//
  int GridRank;                        // number of dimensions
  int GridDimension[MAX_DIMENSION];    // total dimensions of all grids
  int GridStartIndex[MAX_DIMENSION];   // starting index of the active region
                                       //   (zero based)
  int GridEndIndex[MAX_DIMENSION];     // stoping index of the active region
                                       //   (zero based)
  FLOAT GridLeftEdge[MAX_DIMENSION];   // starting pos (active problem space)
  FLOAT GridRightEdge[MAX_DIMENSION];  // ending pos (active problem space)
  float dtFixed;                       // current (fixed) timestep
  FLOAT Time;                          // current problem time
  FLOAT OldTime;                       // time corresponding to OldBaryonField
  int   SubgridsAreStatic;             // 
  int   ID;                            // Grid ID Number
//
//  Baryon grid data
//
  int    NumberOfBaryonFields;                        // active baryon fields
  float *BaryonField[MAX_NUMBER_OF_BARYON_FIELDS];    // pointers to arrays
  float *OldBaryonField[MAX_NUMBER_OF_BARYON_FIELDS]; // pointers to old arrays
  float *InterpolatedField[MAX_NUMBER_OF_BARYON_FIELDS]; // For RT and movies
  float *RandomForcingField[MAX_DIMENSION];           // pointers to arrays //AK
  int    FieldType[MAX_NUMBER_OF_BARYON_FIELDS];
  FLOAT *CellLeftEdge[MAX_DIMENSION];
  FLOAT *CellWidth[MAX_DIMENSION];
  fluxes *BoundaryFluxes;

  float  CourantSafetyNumber;                       // Hydro parameter
  int    PPMFlatteningParameter;                    // PPM parameter
  int    PPMDiffusionParameter;                     // PPM parameter
  int    PPMSteepeningParameter;                    // PPM parameter
//
//  Particle data
//
  int    NumberOfParticles;
  FLOAT *ParticlePosition[MAX_DIMENSION];  // pointers to position arrays
  float *ParticleVelocity[MAX_DIMENSION];  // pointers to velocity arrays
  float *ParticleAcceleration[MAX_DIMENSION+1];  // 
  float *ParticleMass;                     // pointer to mass array
  int   *ParticleNumber;                   // unique identifier
  int   *ParticleType;                     // type of particle
  float *ParticleAttribute[MAX_NUMBER_OF_PARTICLE_ATTRIBUTES];
//
//  Star particle data
//
  int NumberOfStars;
  Star *Stars;
//
//  Gravity data
// 
  float *PotentialField;
  float *AccelerationField[MAX_DIMENSION]; // cell cntr acceleration at n+1/2
  float *GravitatingMassField;
  FLOAT  GravitatingMassFieldLeftEdge[MAX_DIMENSION];
  int    GravitatingMassFieldDimension[MAX_DIMENSION];
  FLOAT  GravitatingMassFieldCellSize;     // all dimensions must be the same
  float *GravitatingMassFieldParticles;     // for particles only
  FLOAT  GravitatingMassFieldParticlesLeftEdge[MAX_DIMENSION];
  FLOAT  GravitatingMassFieldParticlesCellSize;
  int    GravitatingMassFieldParticlesDimension[MAX_DIMENSION];
  gravity_boundary_type GravityBoundaryType;
  float  PotentialSum;
//
//  Rebuild Hierarchy Temporaries
//
  bool  *FlaggingField;             // Boolean flagging field (for refinement)
  float *MassFlaggingField;         // Used by mass flagging criteria
  float *ParticleMassFlaggingField; // Used by particle mass flagging criteria
//
//  Parallel Information
//
  int ProcessorNumber;
//
// Movie Data Format
//
  int TimestepsSinceCreation; 	// Not really since creation anymore... 
  				// resets everytime the grid outputs
//
// Friends
//
  friend int ExternalBoundary::Prepare(grid *TopGrid);
  friend int ProtoSubgrid::CopyFlaggedZonesFromGrid(grid *Grid);
  friend class Star;

#ifdef TRANSFER
#include "PhotonGrid_Variables.h"
#endif

 public:

// -------------------------------------------------------------------------
//  Main hydro/AMR functions
//

/* Grid constructor (Set all data to null/default state). */

   grid();

/* Grid deconstructor (free up memory usage) */

   ~grid();

/* Read grid data from a file (returns: success/failure) */

   int ReadGrid(FILE *main_file_pointer, int GridID, 
		int ReadText = TRUE, int ReadData = TRUE);

/* Read grid data from a group file (returns: success/failure) */

   int Group_ReadGrid(FILE *main_file_pointer, int GridID, HDF5_hid_t file_id, 
		      int ReadText = TRUE, int ReadData = TRUE);

/* Get field or particle data based on name or integer 
   defined in typedefs.h. Details are in Grid_CreateFieldArray.C. */

   EnzoArrayBool *CreateFieldArrayBool(field_type field);
   EnzoArrayBool *CreateFieldArrayBool(char *field_name);

   EnzoArrayInt *CreateFieldArrayInt(field_type field);
   EnzoArrayInt *CreateFieldArrayInt(char *field_name);
  
   EnzoArrayFloat *CreateFieldArrayFloat(field_type field);
   EnzoArrayFloat *CreateFieldArrayFloat(char *field_name);
  
   EnzoArrayFLOAT *CreateFieldArrayFLOAT(field_type field);
   EnzoArrayFLOAT *CreateFieldArrayFLOAT(char *field_name);

/* Write unigrid cubes to a file (returns: success/failure) */

   int WriteCube(char *base_name, int grid_id, int TGdims[]);

/* Write grid data to a file (returns: success/failure) */

   int WriteGrid(FILE *main_file_pointer, char *base_name, int grid_id);

/* Write grid data to a group file (returns: success/failure) */

   int Group_WriteGrid(FILE *main_file_pointer, char *base_name, int grid_id, HDF5_hid_t file_id);

/* Write grid data to separate files (returns: success/failure) */

   int WriteGridX(FILE *main_file_pointer, char *base_name, int grid_id);

/* Write task memory map */

   int WriteMemoryMap(FILE *file_pointer, char *base_name, int grid_id);

/* Write grid-to-task map */

   int WriteTaskMap(FILE *file_pointer, char *base_name, int grid_id);

/* Write grid hierarchy only */

   int WriteStuff(FILE *main_file_pointer, char *base_name, int grid_id);

/* Interpolate to specified time and write unigrid cube data to a file
   (returns: success/failure). */

   int WriteCubeInterpolate(FLOAT WriteTime, char *base_name, int grid_id, int TGdims[]);

/* Interpolate to specified time and write grid data to a file
   (returns: success/failure). */

   int WriteGridInterpolate(FLOAT WriteTime, FILE *main_file_pointer, 
			    char *base_name, int grid_id);

/* Interpolate to specified time and write grid data to a group file
   (returns: success/failure). */

   int Group_WriteGridInterpolate(FLOAT WriteTime, FILE *main_file_pointer,
                            char *base_name, int grid_id, HDF5_hid_t file_id);

/* Compute the timestep constraint for this grid
    (for steps #3 and #4) */

   float ComputeTimeStep();

/* Set the timestep in this grid to the timestep in the argument
    (for step #3) */

   void SetTimeStep(float dt) {dtFixed = dt;};

/* Check timestep (dtFixed) against argument (return fail if dtFixed > dt).
    (for step #4) */

   int CheckTimeStep(float dt) {return ((dtFixed > dt) ? FAIL : SUCCESS);};

/* Return time, timestep */

   FLOAT ReturnTime() {return Time;};
   float ReturnTimeStep() {return dtFixed;};

  /* Return, set grid ID */

  void SetGridID(int id) { ID = id; };
  int GetGridID(void) { return ID; };
   
  /* Baryons: return field types. */

  int ReturnFieldType(int type[]) 
  {
    for (int i = 0; i < NumberOfBaryonFields; i++) type[i] = FieldType[i];
    return SUCCESS;
  };

/* Baryons: Interpolate (parental) grid in argument to current grid.
            (returns success or fail).
    (for step #16) */

   int InterpolateBoundaryFromParent(grid *ParentGrid);

/* Baryons: Copy current solution to Old solution (returns success/fail)
    (for step #16) */

   int CopyBaryonFieldToOldBaryonField();

/* Copy potential field to baryon potential for output purposes. */

   int CopyPotentialToBaryonField();

/* Baryons: Update boundary according to the external boundary values
    (for step #16) */

   int SetExternalBoundaryValues(ExternalBoundary *Exterior);

/* Baryons: solve hydro equations in this grid (returns: the fluxes of the
           subgrids in the argument).  Returns SUCCESS or FAIL.
    (for step #16) */

   int SolveHydroEquations(int CycleNumber, int NumberOfSubgrids, 
			   fluxes *SubgridFluxes[], int level);

/* Baryons: return pointer to the BoundaryFluxes of this grid */

   int ReturnFluxDims(fluxes &f, int RefinementFactors[]);

/* Baryons: prepare and clear the accumulated boundary fluxes for this grid.
    (for step #16) */

   void PrepareBoundaryFluxes();
   void ClearBoundaryFluxes();

/* Baryons: projected solution in current grid to the grid in the 
           argument which must have a lower resolution (i.e. downsample 
           the current grid to the appropriate level).
    (for step #18) */

   int ProjectSolutionToParentGrid(grid &ParentGrid);

/* Baryons: return boundary fluxes from this grid.  Downsample them to
           the refinement factors specified in the argument.
	   Returns FAIL or SUCCESS.
    (for step #19) */

   int GetProjectedBoundaryFluxes(grid *ParentGrid, fluxes &ProjectedFluxes);

/* Return the refinement factors as compared to the grid in the argument
   (integer version) (for step #19) */

   void ComputeRefinementFactors(grid *SubGrid, int RefinementFactors[]) {
     int dim;
     for (dim = 0; dim < GridRank; dim++) RefinementFactors[dim] = 
	 int( CellWidth[dim][0] / SubGrid->CellWidth[dim][0] + 0.5);
     for (dim = GridRank; dim < MAX_DIMENSION; dim++)
       RefinementFactors[dim] = 1;
   };

/* Return the refinement factors as compared to the grid in the argument
   (float version) (for step #19) */

   void ComputeRefinementFactorsFloat(grid *SubGrid, float Factors[]) {
     int dim;
     for (dim = 0; dim < GridRank; dim++) Factors[dim] = 
       (*CellWidth[dim]) / (*(SubGrid->CellWidth[dim]));;
     for (dim = GridRank; dim < MAX_DIMENSION; dim++)
       Factors[dim] = 1.0;
   };

/* Baryons: Search for redundant overlap between two sets of fluxes (other
            and refined).  If found, set the refined fluxes equal to the
	    initial fluxes so there will be no double corrections.(step #19) */

   void CorrectRedundantFluxes(fluxes *OtherFluxes, fluxes *InitialFluxes, 
                               fluxes *RefinedFluxes);

/* Baryons: correct for better flux estimates produced by subgrids
           (i.e given the initial flux estimates and the subgrid flux 
	   estimates, correct the grid to account for the subgrid 
	   flux estimates).  Returns SUCCESS or FAIL.
    (for step #19) */

#ifdef FLUX_FIX
   int CorrectForRefinedFluxes(fluxes *InitialFluxes, fluxes *RefinedFluxes,
			       fluxes *BoundaryFluxesThisTimeStep,
			       int SUBlingGrid,
			       TopGridData *MetaData);
#else
   int CorrectForRefinedFluxes(fluxes *InitialFluxes, fluxes *RefinedFluxes,
			       fluxes *BoundaryFluxesThisTimeStep);
#endif

/* Baryons: add the fluxes pointed to by the argument to the boundary fluxes
            of this grid (sort of for step #16).  Note that the two fluxes
	    must have the same size. */

   int AddToBoundaryFluxes(fluxes *BoundaryFluxesToBeAdded);

/* set new time (time += dt)
    (step #21) */

   void SetTimeNextTimestep() {Time += dtFixed;};

/* set time of this grid (used in setup) */

   void SetTime(FLOAT NewTime) {Time = NewTime;};

/* set hydro parameters (used in setup) */

   void SetHydroParameters(float co, int p1, int p2, int p3) 
     {
       CourantSafetyNumber    = co;
       PPMFlatteningParameter = p1;
       PPMDiffusionParameter  = p2;
       PPMSteepeningParameter = p3;
     }

/* Baryons: compute the pressure at the requested time. */

   int ComputePressure(FLOAT time, float *pressure);

/* Baryons: compute the pressure at the requested time using the dual energy
            formalism. */

   int ComputePressureDualEnergyFormalism(FLOAT time, float *pressure);

/* Baryons: compute the temperature. */

   int ComputeTemperatureField(float *temperature);

/* Baryons: compute X-ray emissivity in specified band. */

   int ComputeXrayEmissivity(float *temperature,
			     float *xray_emissivity, float keV1, float keV2,
			     char *XrayFileName);

/* Baryons: compute number density of ionized elements (just O7 and O8). */

   int ComputeElementalDensity(float *temperature, float *elemental_density,
			       int Type);

/* Baryons: compute the ratio of specific heats. */

   int ComputeGammaField(float *GammaField);

/* Baryons: compute the cooling time. */

   int ComputeCoolingTime(float *cooling_time);

/* Baryons & DualEnergyFormalism: Restore consistency between total and
                                  internal energy fields. */

   int RestoreEnergyConsistency(int Region);

/* Returns some grid info. */

   int ReturnGridInfo(int *Rank, int Dims[], FLOAT Left[], FLOAT Right[]);

/* Subtracts kinetic component from total energy. */

   int ConvertTotalEnergyToGasEnergy();

/* Sets the energy to provide Jean's level support (Zeus: returns coeff). */
   
   int SetMinimumSupport(float &MinimumSupportEnergyCoefficient);

/* Debugging support. */

   int DebugCheck(char *message = "Debug");

// -------------------------------------------------------------------------
// Functions used for analysis
//

/* Calculate the angular momentum of a grid (given center). */

   int CalculateAngularMomentum(FLOAT Center[], float AngularMomentum[],
				float MeanVelocity[], float DMVelocity[],
				FLOAT CenterOfMass[], FLOAT DMCofM[]);

/* Find and track density peaks. */

   int AnalyzeTrackPeaks(int level, int ReportLevel);

/* Project some of the fields to a plane. */

   int ProjectToPlane(FLOAT ProjectedFieldLeftEdge[], 
		      FLOAT ProjectedFieldRightEdge[],
		      int ProjectedFieldDims[], float *ProjectedField[], 
		      int ProjectionDimension, int ProjectionSmooth,
                      int NumberOfProjectedFields, int level,
		      int XrayUseLookupTable, float XrayLowerCutoffkeV,
		      float XrayUpperCutoffkeV, char *XrayFileName);

/* Set the fields to zero under the active region of the specified subgrid. */

   int ZeroSolutionUnderSubgrid(grid *Subgrid, int FieldsToZero, 
                                float Value = 1.0, int AllProcessors = FALSE);

/* Convert the grid data to particle data for output. */

   int OutputAsParticleData(FLOAT RegionLeftEdge[], FLOAT RegionRightEdge[],
                           ListOfParticles *ParticleList[NUM_PARTICLE_TYPES],
                           float BaseRadius);

/* Output star particles to a binary file */

   int OutputStarParticleInformation(FILE *StarFile);

/* Return some information about the grid. */

   int CollectGridInformation(int &GridMemory, float &GridVolume, 
                              int &NumberOfCells, float &AxialRatio,
                              int &CellsTotal, int &Particles);

/* Output grid information (for movie generation). */

   int OutputGridMovieData(FILE *Gridfptr, FILE *DMfptr, FILE *Starfptr,
			   FLOAT RegionLeftEdge[], FLOAT RegionRightEdge[],
			   FLOAT WriteOutTime, int NumberOfPoints[3],
			   int NumberOfValuesPerPoint[3],
			   char *PointValueNames[3][20], float BaseRadius);

/* Output movie data (sequential format) */

   int WriteNewMovieData(FLOAT RegionLeftEdge[], FLOAT RegionRightEdge[], 
			 FLOAT StopTime, int lastMovieStep, int &cycle);

   int ReturnMovieTimestep() { return TimestepsSinceCreation; };

/* Output tracer particle information (interpolated from baryon grid). */

   int TracerParticleOutputData(FILE *ptr, FLOAT WriteOutTime);

// -------------------------------------------------------------------------
// Functions for radiative cooling and multi-species rate equations
//

/* Solve the radiative cooling/heating equations  */

   int SolveRadiativeCooling();

/* Solve the rate equations. */

   int SolveRateEquations();

/* Solve the joint rate and radiative cooling/heating equations  */

   int SolveRateAndCoolEquations();

/* Compute densities of various species for RadiationFieldUpdate. */

   int RadiationComputeDensities(int level);

// -------------------------------------------------------------------------
// Functions for grid (re)generation.
//

/* Remove un-needed arrays before rebuilding. */

   void CleanUp();

/* Delete all the fields, but leave other grid data. */

   void DeleteAllFields();

/* Delete all the fields except for the particle data */

   void DeleteAllButParticles();

/* Sum particle mass flagging fields into ProcessorNumber if particles
   aren't local. */

   int SetParticleMassFlaggingField(int StartProc=0, int EndProc=0, int level=-1, 
				    int ParticleMassMethod=-1, int *SendProcs=NULL, 
				    int NumberOfSends=0);
   int CollectParticleMassFlaggingField(void);
   void ClearParticleMassFlaggingField(void);

/* Clear mass flagging field (gg #1) */

   void ClearMassFlaggingField();

/* Clear boolean flagging field (gg #0) */

   void ClearFlaggingField();

/* Set boolean flagging field */

   int SetFlaggingField(int &NumberOfFlaggedCells, int level);

/* Set flagging field from static regions */

   int SetFlaggingFieldStaticRegions(int level, int &NumberOfFlaggedCells);

/* Delete flagging field */

   void DeleteFlaggingField();

/* Particles: deposit particles living in this grid into the Mass Flagging
             field (gg #2) */

   void DepositParticlesToMassFlaggingField() {};

/* Particles: deposit particles to particle mass flagging field. */

   int DepositMustRefineParticles(int pmethod, int level);

/* baryons: add baryon density to mass flaggin field (so the mass flagging
            field contains the mass in the cell (not the density) 
            (gg #3) */

   int AddFieldMassToMassFlaggingField();

/* Flag all points that require refining  (and delete Mass Flagging Field).
     Returns the number of flagged cells.  Returns the number of flagged cells
     (gg #4) */

   int FlagCellsToBeRefinedByMass(int level, int method);

/* Flag all points that require refining by their slope.
     Returns the number of flagged cells.  Returns the number of flagged cells
     (gg #4) */

   int FlagCellsToBeRefinedBySlope();

/* Flag all points that require refinging by the presence of shocks.
     Returns the number of flagged cells.  Returns the number of flagged cells
     (gg #4) */

   int FlagCellsToBeRefinedByShocks();

/* Flag all points that require refining by the Jean's length criterion. */

   int FlagCellsToBeRefinedByJeansLength();

/* Flag all points that require refining by Shear. */

   int FlagCellsToBeRefinedByShear();

/* Flag all cells for which tcool < dx/sound_speed. */

   int FlagCellsToBeRefinedByCoolingTime();

/* Flag all cells which are near a must-refine particle. */

   int FlagCellsToBeRefinedByMustRefineParticles();

/* Flagging all cell adjacent to a previous flagged cell.  Also, remove all
   Flagged cells in the boundary zones and within one zone of the boundary. */

   int FlagBufferZones();

/* Identify new subgrids for this grid (and prove Fermat's last theorem too)
   (gg #5) */

   void IdentifyNewSubgrids(GridList &list);

/* Identify new subgrids for this grid (1x1x1 subgrids).
   (gg #5) */

   void IdentifyNewSubgridsSmall(GridList &list);

/* Coalesce neighbouring subgrids */

   // void CoalesceSubgrids(GridList &list);

/* Inherit properties (rank, baryon field types, etc.) from ParentGrid
   (gg # 5,6) */

   void InheritProperties(grid *ParentGrid);

/* set the grid dimensions, left, right edges and cell quantities based
   on arguments (gg #5,6) */

   void PrepareGrid(int Rank, int Dimensions[], 
		    FLOAT LeftEdge[], FLOAT RightEdge[], int NumParticles);

/* Allocates space for grids (dims and NumberOfBaryonFields must be set). */

   void AllocateGrids();

/* set the grid derived quantites (CellLeftEdge, CellWidth & BoundaryFluxes) */

   void PrepareGridDerivedQuantities();

/* baryons: interpolate field values from the Parent Grid (gg #6).
            Returns SUCCESS or FAIL. */

   int InterpolateFieldValues(grid *ParentGrid);

/* baryons: check for coincident zones between grids & copy if found.
            (correctly includes periodic boundary conditions). */

   int CheckForOverlap(grid *OtherGrid,
		       boundary_type LeftFaceBoundaryCondition[],
		       boundary_type RightFaceBoundaryCondition[],
		       int (grid::*CopyFunction)(grid *OtherGrid,
						 FLOAT EdgeOffset[]));

/* baryons: check for subgrids adjacent to external boundary with reflecting BCs. */

   int CheckForExternalReflections(
				   boundary_type LeftFaceBoundaryCondition[],
				   boundary_type RightFaceBoundaryCondition[]);

/* David Collins flux correction - July 2005 */
#ifdef FLUX_FIX
   int CheckForSharedFace(grid *OtherGrid,
			       boundary_type LeftFaceBoundaryCondition[],
			       boundary_type RightFaceBoundaryCondition[]);

   int CheckForSharedFaceHelper(grid *OtherGrid,
				     FLOAT EdgeOffset[MAX_DIMENSION]);
#endif

/* baryons: check for overlap between grids & return TRUE if it exists
            (correctly includes periodic boundary conditions). */

   int CheckForPossibleOverlap(grid *OtherGrid,
                       boundary_type LeftFaceBoundaryCondition[],
                       boundary_type RightFaceBoundaryCondition[]);
   int CheckForPossibleOverlapHelper(grid *OtherGrid,
                                        FLOAT EdgeOffset[MAX_DIMENSION]);

/* baryons: copy coincident zone from the (old) grid in the argument
            (gg #7).  Return SUCCESS or FAIL. */

   int CopyZonesFromGrid(grid *GridOnSameLevel, 
			 FLOAT EdgeOffset[MAX_DIMENSION]);

/* gravity: copy coincident potential field zones from grid in the argument
            (gg #7).  Return SUCCESS or FAIL. */

   int CopyPotentialField(grid *GridOnSameLevel, 
			  FLOAT EdgeOffset[MAX_DIMENSION]);

/* baryons: check for coincident zone from the (old) grid in the argument
            (gg #7).  Return SUCCESS or FAIL. */

   int CopyZonesFromGridCountOnly(grid *GridOnSameLevel, int &Overlap);

/* Returns whether or not the subgrids of this grid are static. */

   int AreSubgridsStatic() {return SubgridsAreStatic;};

/* Check the energy conservation. */

   int ComputeEnergy(float EnergySum[]);

/* These two routines add grids to the chaining mesh used in the
   FastSiblingLocator method and use the chaining mesh to find
   possible siblings. */

   int FastSiblingLocatorAddGrid(ChainingMeshStructure *mesh);

   int FastSiblingLocatorFindSiblings(ChainingMeshStructure *mesh,
                          SiblingGridList *list,
                          boundary_type LeftBoundaryCondition[],
                          boundary_type RightBoundaryCondition[]);

   /* hack: add density squared field to grid (used in ExtractSection). */

   void CreateDensitySquaredField() {
     int size = GridDimension[0]*GridDimension[1]*GridDimension[2];
     BaryonField[NumberOfBaryonFields] = new float[size];
     for (int i = 0; i < size; i++)
       BaryonField[NumberOfBaryonFields][i] = 
	 BaryonField[0][i]*BaryonField[0][i];
     FieldType[NumberOfBaryonFields++] = Density;
   };

// -------------------------------------------------------------------------
// Functions for use with gravity.
//

/* Set the gravity boundary type of a grid. */

   void SetGravityParameters(gravity_boundary_type Boundary) {
     GravityBoundaryType = Boundary;};
   gravity_boundary_type ReturnGravityBoundaryType() 
     {return GravityBoundaryType;};

/* Gravity: Initialize, the gravitating Mass Field
    (for steps #5, #6). */

   int InitializeGravitatingMassField(int RefinementFactor);

/* Gravity: Initialize, the particle component of the mass field. */

   int InitializeGravitatingMassFieldParticles(int RefinementFactor);

/* Gravity: allocate & clear the GravitatingMassField. */

   int ClearGravitatingMassField();

/* Gravity & baryons: Copy the parent density field to the extra boundary
      region of GravitatingMassField (if any). */

   int CopyParentToGravitatingFieldBoundary(grid *ParentGrid);

/* Gravity & Particles: allocate & clear the GravitatingMassFieldParticles. */

   int ClearGravitatingMassFieldParticles();

/* Baryons: add the baryon mass to the GravitatingMassField. */

   int AddBaryonsToGravitatingMassField();

/* Generic deposit particles/grids to grid (either GravitatingMassField or
   GravitatingMassFieldParticles depending on the value of DepositField). */

   int DepositPositions(FLOAT *Positions[], float *Mass, int Number, 
			int DepositField);

/* deposit particles/grids to grid (if they are on the grid). */

/* int DepositPositionsEdgeOff(float *Positions[], float *Mass, int Number);*/

/* Gravity: Difference potential to get acceleration field. */

   int ComputeAccelerationField(int DifferenceType, int level);

/* Gravity: Interpolate accelerations from other grid. */

   int InterpolateAccelerations(grid *FromGrid);

/* Gravity: Compute particle and grid accelerations. */

   int ComputeAccelerations(int level);

/* Particles: add overlapping ParticleMassField to Target's 
   GravitatingMassField. */

   int CopyOverlappingMassField(grid *TargetGrid, 
				FLOAT EdgeOffset[MAX_DIMENSION]);

/* Gravity: Allocate and make initial guess for PotentialField. */

   int PreparePotentialField(grid *ParentGrid);

/* Gravity: Allocate and make initial guess for PotentialField. */

   int SolveForPotential(int &Done, int level, FLOAT PotentialTime = -1);

/* Gravity: Prepare the Greens Function. */

   int PrepareGreensFunction();
   int PreparePeriodicGreensFunction(region *GreensRegion);

/* Gravity: Copy potential/density into/out of FFT regions. */

   int PrepareFFT(region *InitialRegion, int Field, int DomainDim[]);
   int FinishFFT(region *InitialRegion, int Field, int DomainDim[]);

/* Gravity: set the potential boundary for isolated BC's */

#ifdef ISOLATED_GRAVITY
   int SetIsolatedPotentialBoundary();
#endif /* ISOLATED_GRAVITY */

/* Gravity: Set the external acceleration fields. */

   int ComputeAccelerationFieldExternal();

/* Particles + Gravity: Clear ParticleAccleration. */

   int ClearParticleAccelerations();

/* Baryons + Gravity: Interpolate the AccelerationField in FromGrid to
             AccelerationFieldForCells at the GridPositions in this grid. */

   int InterpolateGridPositions(grid *FromGrid);

/* Particles + Gravity: Interpolate the AccelerationField in FromGrid to
             ParticleAcceleration at the ParticlePositions in this grid. */

   int InterpolateParticlePositions(grid *FromGrid, int DifferenceType);

/* Generic routine for interpolating particles/grid. */

   int InterpolatePositions(FLOAT *Positions[], int dim, float *Field, 
			    int Number);

/* Gravity: Delete GravitatingMassField. */

   void DeleteGravitatingMassField() {
     delete [] GravitatingMassField; 
     GravitatingMassField = NULL;
   };

/* Gravity: Delete AccelerationField. */

   void DeleteAccelerationField() {
     for (int dim = 0; dim < GridRank; dim++) {
       delete [] AccelerationField[dim];
       AccelerationField[dim] = NULL;
     }
   };

/* Gravity: Add fixed, external acceleration to baryons & particles. */

   int AddExternalAcceleration();

/* Gravity: deposit baryons into target GravitatingMassField. */

   int DepositBaryons(grid *TargetGrid, FLOAT DepositTime);

// -------------------------------------------------------------------------
// Functions for use with particles.
//

/* Particles: Deposit particles in the specified field (DepositField) of the
              TargetGrid at the given time. */

   int DepositParticlePositions(grid *TargetGrid, FLOAT DepositTime, 
				int DepositField);

   int DepositParticlePositionsLocal(FLOAT DepositTime, int DepositField);

/* Particles: add overlapping ParticleMassField to Target's 
   GravitatingMassField. */

   int AddOverlappingParticleMassField(grid *TargetGrid, 
				       FLOAT EdgeOffset[MAX_DIMENSION]);

/* Particles: Apply particle acceleration to velocity for particles in this 
              grid
    (for step #9) */

   int UpdateParticleVelocity(float TimeStep);

/* Particles: Update particle positions (push)
    (for step #13) */

   int UpdateParticlePosition(float TimeStep);

/* Particles: Move particles from TargetGrid to this grid. */

   int MoveAllParticles(int NumberOfGrids, grid* TargetGrids[]);

/* Particles: Move particles that lie within this grid from the TargetGrid
              to this grid. */

//   int MoveSubgridParticles(grid *TargetGrid);


   int MoveSubgridParticles(grid *TargetGrid,
                            int *Counter,
                            int *Number,
                            int *Type,
                            float *Mass,
                            FLOAT *Position[],
                            float *Velocity[],
                            float *Attribute[]);


/* Particles: same as above, but a version that is much more efficient. */

   int MoveSubgridParticlesFast(int NumberOfSubgrids, grid *ToGrids[],
				int AllLocal);

/* Particles: Clean up moved particles (from MoveSubgridParticles). */

   int CleanUpMovedParticles();

/* Particles: delete accleration fields. */

   void DeleteParticleAcceleration() {
     for (int dim = 0; dim < GridRank+ComputePotential; dim++) {
       delete [] ParticleAcceleration[dim];
       ParticleAcceleration[dim] = NULL;
     }
   };

/* Particles & Gravity: Delete GravitatingMassField. */

   void DeleteGravitatingMassFieldParticles() {
     delete [] GravitatingMassFieldParticles; 
     GravitatingMassFieldParticles = NULL;
     GravitatingMassFieldParticlesCellSize = FLOAT_UNDEFINED;
   };

/* Particles: return number of particles. */

   int ReturnNumberOfParticles() {return NumberOfParticles;};

/* Particles: set number of particles. */

   void SetNumberOfParticles(int num) {NumberOfParticles = num;};

/* Particles: delete particle fields and set null. */

   void DeleteParticles() {
     delete [] ParticleMass;
     delete [] ParticleNumber;
     delete [] ParticleType;
     ParticleMass = NULL;
     ParticleNumber = NULL;
     ParticleType = NULL;
     for (int dim = 0; dim < GridRank; dim++) {
       delete [] ParticlePosition[dim];
       delete [] ParticleVelocity[dim];
       ParticlePosition[dim] = NULL;
       ParticleVelocity[dim] = NULL;
     }
     for (int i = 0; i < NumberOfParticleAttributes; i++) {
       delete [] ParticleAttribute[i];
       ParticleAttribute[i] = NULL;
     }   
   };

/* Particles: allocate new particle fields. */

   void AllocateNewParticles(int NumberOfNewParticles) {
     ParticleMass = new float[NumberOfNewParticles];
     ParticleNumber = new int[NumberOfNewParticles];
     ParticleType = new int[NumberOfNewParticles];
     for (int dim = 0; dim < GridRank; dim++) {
       ParticlePosition[dim] = new FLOAT[NumberOfNewParticles];
       ParticleVelocity[dim] = new float[NumberOfNewParticles];
     }
     for (int i = 0; i < NumberOfParticleAttributes; i++)
       ParticleAttribute[i] = new float[NumberOfNewParticles];
   };

/* Particles: Copy pointers passed into into grid. */

   void SetParticlePointers(float *Mass, int *Number, int *Type,
                            FLOAT *Position[], 
			    float *Velocity[], float *Attribute[]) {
    ParticleMass   = Mass;
    ParticleNumber = Number;
    ParticleType   = Type;
    for (int dim = 0; dim < GridRank; dim++) {
      ParticlePosition[dim] = Position[dim];
      ParticleVelocity[dim] = Velocity[dim];
    }
    for (int i = 0; i < NumberOfParticleAttributes; i++)
      ParticleAttribute[i] = Attribute[i];
   };

/* Particles: Set new star particle index. */

   void SetNewParticleIndex(int &NumberCount, int BaseNumber) {
    for (int n = 0; n < NumberOfParticles; n++)
      if (ParticleNumber[n] == INT_UNDEFINED)
	ParticleNumber[n] = BaseNumber + NumberCount++;
   };

/* Particles: Add given number to particle index. */

   void AddToParticleNumber(int *Count) {
     if (MyProcessorNumber == ProcessorNumber)
       for (int n = 0; n < NumberOfParticles; n++)
	 ParticleNumber[n] += *Count;
     *Count += NumberOfParticles;
   }

  float ReturnTotalSinkMass() {
    float total = 0;
    double dx3 = CellWidth[0][0] * CellWidth[0][0] * CellWidth[0][0];
    if (MyProcessorNumber == ProcessorNumber)
      for (int n = 0; n < NumberOfParticles; n++)
	if (ParticleType[n] == PARTICLE_TYPE_MUST_REFINE)
	  total += ParticleMass[n] * dx3;
    return total;
  };

/* Particles: return particle type (1 - dm, 2 - star).   Note that this
   has now been superceded by a real particle type field. */

   int ReturnParticleType(int index) {
     if (NumberOfParticleAttributes > 0 && StarParticleCreation > 0)
       if (ParticleAttribute[0][index] > 0)
         return PARTICLE_TYPE_STAR;
     return PARTICLE_TYPE_DARK_MATTER;
   }

/* Particles: sort particle data in ascending order by number (id). */

void SortParticlesByNumber();

// -------------------------------------------------------------------------
// Communication functions
//

/* Set grid's processor number. */

  void SetProcessorNumber(int Proc) {
    ProcessorNumber = Proc;
  };

/* Return grid's processor number. */

  int ReturnProcessorNumber() {
    return ProcessorNumber;
  }

/* Send a region from a real grid to a 'fake' grid on another processor. */

  int CommunicationSendRegion(grid *ToGrid, int ToProcessor, int SendField, 
			     int NewOrOld, int RegionStart[], int RegionDim[]);

/* Send a region from a 'fake' grid to a real grid on another processor. */

  int CommunicationReceiveRegion(grid *ToGrid, int ToProcessor, 
				 int SendField, int NewOrOld, 
				 int RegionStart[], int RegionDim[],
				 int IncludeBoundary);

/* Move a grid from one processor to another. */

  int CommunicationMoveGrid(int ToProcessor, int MoveParticles = TRUE);

/* Send particles from one grid to another. */

  int CommunicationSendParticles(grid *ToGrid, int ToProcessor, 
				int FromStart, int FromNumber, int ToStart);

/* Transfer particle amount level 0 grids. */

#ifdef OPTIMIZED_CTP
  int CommunicationTransferParticles(grid* Grids[], int NumberOfGrids,
				     int ThisGridNum, int *&NumberToMove, 
				     int StartIndex, int EndIndex, 
				     particle_data *&List,
				     int *Layout, int *GridMap, 
				     int CopyDirection);
#else
  int CommunicationTransferParticles(grid* Grids[], int NumberOfGrids,
				     int ToGrid[6], int NumberToMove[6],
				     float_int *ParticleData[6], int CopyDirection);
#endif

  int CollectParticles(int GridNum, int* &NumberToMove, 
		       int &StartIndex, int &EndIndex, 
		       particle_data* &List, int CopyDirection);

  int TransferSubgridParticles(grid* Subgrids[], int NumberOfSubgrids, 
			       int* &NumberToMove, int StartIndex, 
			       int EndIndex, particle_data* &List, 
			       bool KeepLocal, bool ParticlesAreLocal,
			       int CopyDirection);


// -------------------------------------------------------------------------
// Helper functions (should be made private)
//

  /* This is a simple helper function which determines if the method
     should return immediately because of communication-mode reasons.
     Assumes that info is being sent from the "other-grid" processor to
     the "this-grid" processor (i.e. the one that holds this object). */

  int CommunicationMethodShouldExit(grid *OtherGrid) {

    /* Return if neither grid lives on this processor. */

    if (MyProcessorNumber != ProcessorNumber && 
        MyProcessorNumber != OtherGrid->ProcessorNumber)
      return SUCCESS;

    /* If the two grids are on the same processor then return if
       either in post-receive or receive modes to avoid duplicating method
       (i.e. action is only carried out if in send mode (or send-receive)). */

    if (ProcessorNumber == OtherGrid->ProcessorNumber)
      if (CommunicationDirection == COMMUNICATION_POST_RECEIVE ||
	  CommunicationDirection == COMMUNICATION_RECEIVE)
        return SUCCESS;

    /* If in send mode then exit if the send grid is not on this processor. */

    if (CommunicationDirection == COMMUNICATION_SEND &&
        MyProcessorNumber != OtherGrid->ProcessorNumber)
      return SUCCESS;

    /* If in either receive phase then exit if receive grid is not on this 
       processor. */

    if ((CommunicationDirection == COMMUNICATION_RECEIVE ||
         CommunicationDirection == COMMUNICATION_POST_RECEIVE) &&
        MyProcessorNumber != ProcessorNumber)
      return SUCCESS;

    return FAIL; /* i.e. method should not exit immediately. */
  }

/* Baryons: find certain commonly used variables from the list of fields. */

  int IdentifyPhysicalQuantities(int &DensNum, int &GENum,   int &Vel1Num, 
				 int &Vel2Num, int &Vel3Num, int &TENum);

/* Identify Multi-species fields. */

  int IdentifySpeciesFields(int &DeNum, int &HINum, int &HIINum, 
			    int &HeINum, int &HeIINum, int &HeIIINum,
			    int &HMNum, int &H2INum, int &H2IINum,
                            int &DINum, int &DIINum, int &HDINum);

/* Zeus Solver. */

  int ZeusSolver(float *gamma, int igamfield, int nhy, 
		 float dx[], float dy[], float dz[], 
		 int gravity, int NumberOfSubgrids, long_int GridGlobalStart[],
		 fluxes *SubgridFluxes[],
		 int NumberOfColours, int colnum[], int bottom,
		 float minsupecoef);

/* PPM Direct Euler Solver. */

  int PPMDirectEuler(int CycleNumber, int NumberOfSubgrids, 
                     fluxes *SubgridFluxes[], float *CellWidthTemp[], 
                     long_int GridGlobalStart[], int GravityOn,
		     int NumberOfColours, int colnum[]);

  int euler_sweep(int dim, int iter, int CycleNumber, int NumberOfSubgrids, 
                  fluxes *SubgridFluxes[], float *CellWidthTemp[],
                  long_int GridGlobalStart[], int GravityOn,
		  int NumberOfColours, int colnum[], float *temp, int tempsize);

// AccelerationHack

  int AccelerationHack;

#ifdef SAB
  //These should be moved later.
  //Used for boundary condition set of AccelerationField.
  int ActualNumberOfBaryonFields;
  int AttachAcceleration();
  int DetachAcceleration();

  int    ActualFieldType[MAX_NUMBER_OF_BARYON_FIELDS];
  float *ActualBaryonField[MAX_NUMBER_OF_BARYON_FIELDS];
  float *ActualOldBaryonField[MAX_NUMBER_OF_BARYON_FIELDS];
  float *OldAccelerationField[3];
#endif

// -------------------------------------------------------------------------
// Functions for Specific problems (usually problem generator functions).
//

/* Protostellar Collapse problem: initialize grid (returns SUCCESS or FAIL) */
  int ProtostellarCollapseInitializeGrid(float CoreDensity,
					 float CoreEnergy,
					 float CoreRadius,
					 float AngularVelocity);

/* ShockTube problem: initialize grid (returns SUCCESS or FAIL) */

  int ShockTubeInitializeGrid(int Direction, float Boundary, float Density[],
			      float Pressure[], float Velocity[]);

/* Initialize for a uniform grid (returns SUCCESS or FAIL) */

  int InitializeUniformGrid(float UniformDensity, float UniformTotalEnergy,
			    float UniformGasEnergy, float UniformVelocity[]);


/* Initialize a grid for the Double Mach reflection problem. */

  int DoubleMachInitializeGrid(float d0, float e0, float u0,float v0,float w0);


/* Initialize a grid for Implosion test problem */

  int ImplosionInitializeGrid(float ImplosionDiamondDensity,
			      float ImplosionDiamondTotalEnergy);


/* Initialize a grid for Sedov Explosion */

  int SedovBlastInitializeGrid(float SedovBlastInitialRadius,
                               float SedovBlastInnerTotalEnergy);

  int SedovBlastInitializeGrid3D(char * fname);

  /* Initialize a grid for a rotating cylinder collapse */
  int RotatingCylinderInitializeGrid(FLOAT RotatingCylinderRadius,
				     FLOAT RotatingCylinderCenterPosition[MAX_DIMENSION],
				     float RotatingCylinderLambda,
				     float RotatingCylinderOverdensity);

  /* Initialize a grid for the KH instability problem. */

  int KHInitializeGrid(float KHInnerDensity,
                       float KHInnerInternalEnergy,
                       float KHOuterInternalEnergy,
                       float KHPerturbationAmplitude,
                       float KHInnerVx, float KHOuterVx);

  /* Initialize a grid and set boundary for the 2D/3D Noh problem. */

  int NohInitializeGrid(float d0, float p0, float u0);
  int ComputeExternalNohBoundary();

/* Zeldovich Pancake: initial grid (returns SUCCESS or FAIL). */

  int ZeldovichPancakeInitializeGrid(int   ZeldovichPancakeDirection,
				     float ZeldovichPancakeCentralOffset,
				     float ZeldovichPancakeOmegaBaryonNow,
				     float ZeldovichPancakeOmegaCDMNow,
				     float ZeldovichPancakeCollapseRedshift,
				     float ZeldovichPancakeInitialTemperature);

/* 1D Pressureless Collapse: initialize grid. */

  int PressurelessCollapseInitializeGrid(int PressurelessCollapseDirection,
				    float PressurelessCollapseInitialDensity,
				      int PressurelessCollapseNumberOfCells);

/* Gravity Test: particles in isolated boundaries */
  int TestOrbitInitializeGrid(int NumberOfTestParticles,
			      FLOAT TestRadius,
			      float CentralMass,
			      float TestMass,
			      int UseBaryons);

/* Gravity Test: initialize grid. */

  int TestGravityInitializeGrid(float CentralDensity, 
				int NumberOfNewParticles, int UseBaryons);

/* Gravity Test: check results. */

  int TestGravityCheckResults(FILE *fptr, grid *TopGrid);

/* Gravity Test Motion: initialize grid. */

  int TestGravityMotionInitializeGrid(float InitialVelocity);

/* Gravity Test (Sphere): initialize grid. */

  int TestGravitySphereInitializeGrid(float InteriorDensity, 
				      float ExteriorDensity,
				      float SphereRadius, 
				      int SphereType, int UseBaryons,
				      FLOAT SphereCenter[]);

/* Gravity Test (Sphere): check results. */

  int TestGravitySphereCheckResults(FILE *fptr);

/* Spherical Infall Test: initialize grid. */

  int SphericalInfallInitializeGrid(float InitialPerturbation, int UseBaryons,
				    float SphericalInfallOmegaBaryonNow,
				    float SphericalInfallOmegaCDMNow,
				    int SubgridIsStatic);


/* Spherical Infall Test: get the profile from the center. */

  int SphericalInfallGetProfile(int level, int ReportLevel);

/* GravityEquilibriumTest: initialize grid. */

  int GravityEquilibriumTestInitializeGrid(float ScaleHeight);

/* CollapseTest: initialize grid. */

#define MAX_SPHERES 10

int CollapseTestInitializeGrid(int NumberOfSpheres,
			     FLOAT SphereRadius[MAX_SPHERES],
			     FLOAT SphereCoreRadius[MAX_SPHERES],
			     float SphereDensity[MAX_SPHERES],
			     float SphereTemperature[MAX_SPHERES],
			     float SphereMetallicity[MAX_SPHERES],
			     FLOAT SpherePosition[MAX_SPHERES][MAX_DIMENSION],
			     float SphereVelocity[MAX_SPHERES][MAX_DIMENSION],
			     float SphereFracKeplarianRot[MAX_SPHERES],
			     float SphereTurbulence[MAX_SPHERES],
			     float SphereDispersion[MAX_SPHERES],
			     float SphereCutOff[MAX_SPHERES],
			     float SphereAng1[MAX_SPHERES],
			     float SphereAng2[MAX_SPHERES],
			     int   SphereNumShells[MAX_SPHERES],
			     int   SphereType[MAX_SPHERES],
			     int   SphereUseParticles,
			     float ParticleMeanDensity,
			     float UniformVelocity[MAX_DIMENSION],
			     int   SphereUseColour,
			     int   SphereUseMetals,
			     float InitialTemperature, 
			     float InitialDensity, int level);

/* CosmologySimulation: initialize grid. */

  int CosmologySimulationInitializeGrid(
			  int   InitialGridNumber,
			  float CosmologySimulationOmegaBaryonNow,
			  float CosmologySimulationOmegaCDMNow,
			  float CosmologySimulationInitialTemperature,
			  char *CosmologySimulationDensityName,
			  char *CosmologySimulationTotalEnergyName,
			  char *CosmologySimulationGasEnergyName,
			  char *CosmologySimulationVelocityNames[],
			  char *CosmologySimulationParticlePositionName,
			  char *CosmologySimulationParticleVelocityName,
			  char *CosmologySimulationParticleMassName,
			  char *CosmologySimulationParticleTypeName,
			  char *CosmologySimulationParticleVelocityNames[],
			  int   CosmologySimulationSubgridsAreStatic,
			  int   TotalRefinement,
			  float CosmologySimulationInitialFrctionHII,
			  float CosmologySimulationInitialFrctionHeII,
			  float CosmologySimulationInitialFrctionHeIII,
			  float CosmologySimulationInitialFrctionHM,
			  float CosmologySimulationInitialFrctionH2I,
			  float CosmologySimulationInitialFrctionH2II,
			  int   CosmologySimulationUseMetallicityField,
			  int  &CurrentNumberOfParticles,
			  int CosmologySimulationManuallySetParticleMassRatio,
			  float CosmologySimulationManualParticleMassRatio,
			  int CosmologySimulationCalculatePositions);

  int CosmologyInitializeParticles(
		   char *CosmologySimulationParticleVelocityName,
		   char *CosmologySimulationParticleMassName,
		   char *CosmologySimulationParticleTypeName,
		   char *CosmologySimulationParticleVelocityNames[],
		   float CosmologySimulationOmegaBaryonNow,
		   int *Offset, int level);

/* CosmologySimulation: initialize partitioned nested grids. */

  int NestedCosmologySimulationInitializeGrid(
			  int   InitialGridNumber,
			  float CosmologySimulationOmegaBaryonNow,
			  float CosmologySimulationOmegaCDMNow,
			  float CosmologySimulationInitialTemperature,
			  char *CosmologySimulationDensityName,
			  char *CosmologySimulationTotalEnergyName,
			  char *CosmologySimulationGasEnergyName,
			  char *CosmologySimulationVelocityNames[],
			  char *CosmologySimulationParticlePositionName,
			  char *CosmologySimulationParticleVelocityName,
			  char *CosmologySimulationParticleMassName,
			  char *CosmologySimulationParticleTypeName,
			  char *CosmologySimulationParticleVelocityNames[],
			  int   CosmologySimulationSubgridsAreStatic,
			  int   TotalRefinement,
			  float CosmologySimulationInitialFrctionHII,
			  float CosmologySimulationInitialFrctionHeII,
			  float CosmologySimulationInitialFrctionHeIII,
			  float CosmologySimulationInitialFrctionHM,
			  float CosmologySimulationInitialFrctionH2I,
			  float CosmologySimulationInitialFrctionH2II,
			  int   CosmologySimulationUseMetallicityField,
			  int  &CurrentNumberOfParticles,
			  int CosmologySimulationManuallySetParticleMassRatio,
			  float CosmologySimulationManualParticleMassRatio,
			  int CosmologySimulationCalculatePositions);

/* Supernova restart initialize grid. */

  int SupernovaRestartInitialize(float EjectaDensity, float EjectaRadius,
				 float EjectaThermalEnergy, 
				 FLOAT EjectaCenter[3], int ColourField,
				 int *NumberOfCellsSet);

  /* Cooling test initialization */

  int CoolingTestInitializeGrid(float MinimumDensity,
				float MaximumDensity,
				float MinimumTemperature,
				float MaximumTemperature,
				float MinimumColour,
				float MaximumColour,
				int UseMetals,
				int UseElectronFraction);

/* Tricks for Random Forcing. */

  int ReturnNumberOfBaryonFields(){return NumberOfBaryonFields;};
  int SetNumberOfBaryonFields(int &number)
    {NumberOfBaryonFields = number; return 0;};
  int AppendForcingToBaryonFields();
  int DetachForcingFromBaryonFields();
  int RemoveForcingFromBaryonFields();
  int AddRandomForcing(float * norm, float dtTopGrid);
  int PrepareRandomForcingNormalization(float * GlobVal, int GlobNum);
  int ReadRandomForcingFields(FILE *main_file_pointer);

  inline bool isLocal () {return MyProcessorNumber == ProcessorNumber; };

 private:
//  int ReadRandomForcingFields(FILE *main_file_pointer);
 public:

/* TurbulenceSimulation: initialize grid. */

#define TURBULENCE_INIT_PARAMETERS_DECL \
     float TurbulenceSimulationInitialDensity, \
     float TurbulenceSimulationInitialTemperature, \
     char *TurbulenceSimulationDensityName, \
     char *TurbulenceSimulationTotalEnergyName, \
     char *TurbulenceSimulationGasEnergyName, \
     char *TurbulenceSimulationVelocityNames[], \
     char *TurbulenceSimulationRandomForcingNames[], \
     int   TurbulenceSimulationSubgridsAreStatic, \
     int   TotalRefinement


#define TURBULENCE_INIT_PARAMETERS \
     TurbulenceSimulationInitialDensity, \
     TurbulenceSimulationInitialTemperature, \
     TurbulenceSimulationDensityName, \
     TurbulenceSimulationTotalEnergyName, \
     TurbulenceSimulationGasEnergyName, \
     TurbulenceSimulationVelocityNames, \
     TurbulenceSimulationRandomForcingNames, \
     TurbulenceSimulationSubgridsAreStatic, \
     TotalRefinement


 int TurbulenceSimulationInitializeGrid(TURBULENCE_INIT_PARAMETERS_DECL);

  // The following are private since they should only be called by
  // TurbulenceSimulationInitializeGrid()

 private:
//  int TurbulenceSimulationInitializeGrid(TURBULENCE_INIT_PARAMETERS_DECL);
 public:


/* Comoving coordinate expansion terms. */

  int ComovingExpansionTerms();

/* Adjust the gravity source terms for comoving coordinates. */

  int ComovingGravitySourceTerm();

/* Star Particle handler routine. */

  int StarParticleHandler(int level);

/* Apply a time-action to a grid. */

  int ApplyTimeAction(int Type, float Parameter);

/* Routine to set the tracer particle velocities from the grid velocity. */

  int TracerParticleSetVelocity();

/* Creates tracer particles in this grid. */

  int TracerParticleCreateParticles(FLOAT LeftEdge[], FLOAT RightEdge[],
                                    FLOAT Spacing, int &TotalParticleCount);

// -------------------------------------------------------------------------
// Analysis functions for AnalysisBaseClass and it's derivatives.
//

  // Flag cells that overlap a subgrid (used for analysis).
  int FlagRefinedCells(grid *Subgrid);
  
#if defined(CONFIG_PFLOAT_8) || defined(CONFIG_PFLOAT_16)
  // FLOAT version of above
  inline int IsInVolume( FLOAT *LeftEdge, FLOAT *RightEdge ){
    for( int i = 0; i < GridRank; i++ ){
      if( (GridLeftEdge[i] >= RightEdge[i]) ||
	  (GridRightEdge[i] <= LeftEdge[i]) ){
	return FALSE;
      }
    }
     return TRUE;
  }
#else
  // Check to see if the grid overlaps a volume.
  // Doesn't care about periodicity. */
  inline int IsInVolume( float *LeftEdge, float *RightEdge ){
    for( int i = 0; i < GridRank; i++ ){
      if( (GridLeftEdge[i] >= RightEdge[i]) ||
	  (GridRightEdge[i] <= LeftEdge[i]) ){
	return FALSE;
      }
    }
     return TRUE;
  }

#endif

#if defined(p8) || defined(p16)
  // Check to see if a FLOAT point is in the grid.
  inline int PointInGrid( FLOAT *point ){
    for( int i = 0; i < GridRank; i++ ){
      if( ((point[i] >= GridLeftEdge[i]) &&
	  (point[i] <= GridRightEdge[i])) == FALSE )
	return FALSE;
    }
    return TRUE;
  }
#else
  // Check to see if a point is in the grid.
  inline int PointInGrid( float *point ){
    for( int i = 0; i < GridRank; i++ ){
      if( ((point[i] >= GridLeftEdge[i]) &&
	  (point[i] <= GridRightEdge[i])) == FALSE )
	return FALSE;
    }
    return TRUE;
  }
#endif

  // Flags a 3D array where the grid overlaps.
  // Very similar to the FastSib stuff. (I think.)
  void FlagGridArray( HierarchyEntry ***GridArray, int *dx,
		      float *cell_width, HierarchyEntry *my_HE );

  /* Includes for analysis tools */

#ifdef ANALYSIS_TOOLS
#   include "Grid_AnalyzeClusters.h"
#endif

#ifdef EMBEDDED_PYTHON
    void ConvertToNumpy(int GridID, PyArrayObject *container[],
                        int ParentID, int level, FLOAT WriteTime);
#endif
//------------------------------------------------------------------------
// Methods for star formation
//------------------------------------------------------------------------

  Star *ReturnStarPointer(void) { return Stars; };
  int ReturnNumberOfStars(void) { return NumberOfStars; };
  void SetNumberOfStars(int value) { NumberOfStars = value; };

  /* Calculate enclosed mass within a radius */

  int GetEnclosedMass(Star *star, float radius, float &mass,
		      float &metallicity, float &coldgas_mass, 
		      float AvgVelocity[]);

  int RemoveParticle(int ID);

  int AddFeedbackSphere(Star *cstar, int level, float radius, float VelocityUnits, 
			float TemperatureUnits, double EjectaDensity, 
			double EjectaMetalDensity, double EjectaThermalEnergy,
			int &CellsModified);

  int MoveAllStars(int NumberOfGrids, grid* FromGrid[], int TopGridDimension);

  int CommunicationSendStars(grid *ToGrid, int ToProcessor);

  int MoveSubgridStars(int NumberOfSubgrids, grid* ToGrids[], int AllLocal);
  
  int FindNewStarParticles(int level);

  int FindAllStarParticles(int level);

  int MirrorStarParticles(void);

  int UpdateStarParticles(int level);

  int AddH2Dissociation(Star *AllStars);

//------------------------------------------------------------------------
// Radiative transfer methods that don't fit in the TRANSFER define
//------------------------------------------------------------------------

  int IdentifyRadiativeTransferFields(int &kphHINum, int &gammaHINum,
				      int &kphHeINum, int &gammaHeINum,
				      int &kphHeIINum, int &gammaHeIINum,
				      int &kdissH2INum);

#ifdef TRANSFER
#include "PhotonGrid_Methods.h"
#endif /* TRANSFER */

//------------------------------------------------------------------------
// Vertex interpolation (for radiative transfer and streaming data)
//------------------------------------------------------------------------

  int ComputeVertexCenteredField(int Num);
  int ComputeCellCenteredField(int Num);
  
  float ComputeInterpolatedValue(int Num, int vci, int vcj, int vck, 
				 float mx, float my, float mz);
  
  int NeighborIndices(int index, int vi[]);
  int NeighborVertexIndices(int index, int vi[]);

  void DeleteInterpolatedFields() {
    int i;
    for (i = 0; i < NumberOfBaryonFields; i++)
      if (InterpolatedField[i] != NULL) {
	delete [] InterpolatedField[i];
	InterpolatedField[i] = NULL;
      }
  }

//-----------------------------------------------------------------------
//  Returns radiative cooling by component
//-----------------------------------------------------------------------

  int ComputeLuminosity(float *luminosity, int NumberOfLuminosityFields);

//------------------------------------------------------------------------
//  Misc.
//------------------------------------------------------------------------

  int FindMinimumParticleMass(float &min_mass, int level);

  int FindMassiveParticles(float min_mass, int level, FLOAT *pos[], int &npart,
			   int CountOnly);


};

// inline int grid::ReadRandomForcingFields (FILE *main_file_pointer);

// inline int grid::TurbulenceSimulationInitializeGrid (TURBULENCE_INIT_PARAMETERS_DECL);


#endif