Source

enzo-3.0 / src / enzo / grid / solvers / Grid_SolveHydroEquations.C

The active_particles branch has multiple heads

Full commit
  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
/***********************************************************************
/
/  GRID CLASS (SOLVE THE HYDRO EQUATIONS, SAVING SUBGRID FLUXES)
/
/  written by: Greg Bryan
/  date:       November, 1994
/  modified1:
/
/  PURPOSE:
/
/  RETURNS:
/    SUCCESS or FAIL
/
************************************************************************/

// Solve the hydro equations with the solver, saving the subgrid fluxes
//

#include "preincludes.h"
#include "ErrorExceptions.h"
#include "EnzoTiming.h"
#include "performance.h"
#include "macros_and_parameters.h"
#include "typedefs.h"
#include "global_data.h"
#include "Fluxes.h"
#include "GridList.h"
#include "ExternalBoundary.h"
#include "Grid.h"
#include "fortran.def"

/* function prototypes */

int CosmologyComputeExpansionFactor(FLOAT time, FLOAT *a, FLOAT *dadt);
int FindField(int f, int farray[], int n);

#ifdef PPM_LR
extern "C" void FORTRAN_NAME(ppm_lr)(
			  float *d, float *E, float *u, float *v, float *w,
			    float *ge,
                          int *grav, float *gr_ax, float *gr_ay, float *gr_az,
			  float *gamma, float *dt, int *cycle_number,
                            float dx[], float dy[], float dz[],
			  int *rank, int *in, int *jn, int *kn,
                            int is[], int ie[],
			  float gridvel[], int *flatten, int *ipresfree,
			  int *diff, int *steepen, int *idual, 
                            float *eta1, float *eta2,
			  int *num_subgrids, int leftface[], int rightface[],
			  int istart[], int iend[], int jstart[], int jend[],
			  float *standard, int dindex[], int Eindex[],
			  int uindex[], int vindex[], int windex[],
			    int geindex[], float *temp,
                          int *ncolour, float *colourpt, int *coloff,
                            int colindex[], int *ifloat_size);
#endif /* PPM_LR */

int grid::SolveHydroEquations(int CycleNumber, int NumberOfSubgrids, 
			      fluxes *SubgridFluxes[], int level)
{

  /* Return if this doesn't concern us. */
 
  if (ProcessorNumber != MyProcessorNumber || !UseHydro)
    return SUCCESS;

  LCAPERF_START("grid_SolveHydroEquations");
  TIMER_START("SolveHydroEquations");

  this->DebugCheck("SolveHydroEquations");

  if (NumberOfBaryonFields > 0) {

    /* initialize */

    // MAX_COLOR is defined in fortran.def
    int dim, i, j, field, size, subgrid, n, colnum[MAX_COLOR];
    Elong_int GridGlobalStart[MAX_DIMENSION];
    FLOAT a = 1, dadt;

    /* Compute size (in floats) of the current grid. */

    size = 1;
    for (dim = 0; dim < GridRank; dim++)
      size *= GridDimension[dim];

    /* If multi-species being used, then treat them as colour variables
       (note: the solver has been modified to treat these as density vars). */

    int NumberOfColours = 0, ColourNum;

    // use different color fields for RadiativeTransferFLD problems
    //   first, the standard Enzo color field advection
    if (MultiSpecies > 0 && RadiativeTransferFLD != 2) {
      NumberOfColours = 6 + 3*(MultiSpecies-1);

      if ((ColourNum =
           FindField(ElectronDensity, FieldType, NumberOfBaryonFields)) < 0) {
        ENZO_FAIL("Could not find ElectronDensity.");
      }

      /* Generate an array of field numbers corresponding to the colour fields
	 (here assumed to start with ElectronDensity and continue in order). */

      for (i = 0; i < NumberOfColours; i++)
        colnum[i] = ColourNum+i;

    }
    // second, the color field advection if using RadiativeTransferFLD for 
    // a radiation propagation problem (i.e. not using ray-tracing)
    if (RadiativeTransferFLD == 2) {
      if (ImplicitProblem < 4)  {  // grey radiation problem
	
	// set the grey radiation field (required)
	if ((ColourNum =
	     FindField(RadiationFreq0, FieldType, NumberOfBaryonFields)) < 0) 
	  ENZO_FAIL("Could not find RadiationFreq0.");
	colnum[0] = ColourNum;

	// check for other chemistry fields; add if they're present
	//   ElectronDensity
	if ((ColourNum =
	     FindField(ElectronDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   HIDensity
	if ((ColourNum =
	     FindField(HIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   HIIDensity
	if ((ColourNum =
	     FindField(HIIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   HeIDensity
	if ((ColourNum =
	     FindField(HeIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   HeIIDensity
	if ((ColourNum =
	     FindField(HeIIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   HeIIIDensity
	if ((ColourNum =
	     FindField(HeIIIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   HMDensity
	if ((ColourNum =
	     FindField(HMDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   H2IDensity
	if ((ColourNum =
	     FindField(H2IDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   H2IIDensity
	if ((ColourNum =
	     FindField(H2IIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   DIDensity
	if ((ColourNum =
	     FindField(DIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   DIIDensity
	if ((ColourNum =
	     FindField(DIIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
	//   HDIDensity
	if ((ColourNum =
	     FindField(HDIDensity, FieldType, NumberOfBaryonFields)) >= 0) 
	  colnum[++NumberOfColours] = ColourNum;
      }
    }

    /* Add "real" colour fields (metallicity, etc.) as colour variables. */

    int SNColourNum, MetalNum, MBHColourNum, Galaxy1ColourNum, Galaxy2ColourNum,
      MetalIaNum; 

    if (this->IdentifyColourFields(SNColourNum, MetalNum, MetalIaNum, MBHColourNum, 
				   Galaxy1ColourNum, Galaxy2ColourNum) == FAIL)
      ENZO_FAIL("Error in grid->IdentifyColourFields.\n");

    if (MetalNum != -1) {
      colnum[NumberOfColours++] = MetalNum;
      if (MultiMetals || TestProblemData.MultiMetals) {
	colnum[NumberOfColours++] = MetalNum+1; //ExtraType0
	colnum[NumberOfColours++] = MetalNum+2; //ExtraType1
      }
    }

    if (MetalIaNum       != -1) colnum[NumberOfColours++] = MetalIaNum;
    if (SNColourNum      != -1) colnum[NumberOfColours++] = SNColourNum;
    if (MBHColourNum     != -1) colnum[NumberOfColours++] = MBHColourNum;
    if (Galaxy1ColourNum != -1) colnum[NumberOfColours++] = Galaxy1ColourNum;
    if (Galaxy2ColourNum != -1) colnum[NumberOfColours++] = Galaxy2ColourNum;


    /* Add Simon Glover's chemistry species as color fields */

    if(TestProblemData.GloverChemistryModel){

      // Declarations for Simon Glover's cooling.
      int DeNum, HINum, HIINum, HeINum, HeIINum, HeIIINum, HMNum, H2INum, H2IINum,
	DINum, DIINum, HDINum;

      int CINum,CIINum,OINum,OIINum,SiINum,SiIINum,SiIIINum,CHINum,CH2INum,
	CH3IINum,C2INum,COINum,HCOIINum,OHINum,H2OINum,O2INum;

      int GCM = TestProblemData.GloverChemistryModel;  // purely for convenience

      if (IdentifyGloverSpeciesFields(HIINum,HINum,H2INum,DINum,DIINum,HDINum,
				      HeINum,HeIINum,HeIIINum,CINum,CIINum,OINum,
				      OIINum,SiINum,SiIINum,SiIIINum,CHINum,CH2INum,
				      CH3IINum,C2INum,COINum,HCOIINum,OHINum,H2OINum,
				      O2INum) == FAIL) {
	ENZO_FAIL("Error in IdentifyGloverSpeciesFields.");
      }

      colnum[NumberOfColours++] = HIINum;
      colnum[NumberOfColours++] = HINum;
      colnum[NumberOfColours++] = H2INum;

      if( (GCM==1) || (GCM==2) || (GCM==3) || (GCM==7) ){
	colnum[NumberOfColours++] = DINum;
	colnum[NumberOfColours++] = DIINum;
	colnum[NumberOfColours++] = HDINum;
	colnum[NumberOfColours++] = HeINum;
	colnum[NumberOfColours++] = HeIINum;
	colnum[NumberOfColours++] = HeIIINum;
      }

      if( (GCM==3) || (GCM==5) || (GCM==7) ){
	colnum[NumberOfColours++] = COINum;
      }

      if( (GCM==2) || (GCM==3) || (GCM==7) ){
	colnum[NumberOfColours++] = CINum;
	colnum[NumberOfColours++] = CIINum;
	colnum[NumberOfColours++] = OINum;
	colnum[NumberOfColours++] = OIINum;
      }

      if( (GCM==2) || (GCM==3) ){
	colnum[NumberOfColours++] = SiINum;
	colnum[NumberOfColours++] = SiIINum;
	colnum[NumberOfColours++] = SiIIINum;
      }

      if( (GCM==3) || (GCM==7) ){
	colnum[NumberOfColours++] = CHINum;
	colnum[NumberOfColours++] = CH2INum;
	colnum[NumberOfColours++] = CH3IINum;
	colnum[NumberOfColours++] = C2INum;
	colnum[NumberOfColours++] = HCOIINum;
	colnum[NumberOfColours++] = OHINum;
	colnum[NumberOfColours++] = H2OINum;
	colnum[NumberOfColours++] = O2INum;
      }
      
    } // if(TestProblemData.GloverChemistryModel)


    /* Add shock/cosmic ray variables as a colour variable. */

    if(ShockMethod){
      int MachNum, PSTempNum,PSDenNum;
      
      if (IdentifyShockSpeciesFields(MachNum,PSTempNum,PSDenNum) == FAIL) {
	ENZO_FAIL("Error in IdentifyShockSpeciesFields.")
      }
      
      colnum[NumberOfColours++] = MachNum;
      if(StorePreShockFields){
	colnum[NumberOfColours++] = PSTempNum;
	colnum[NumberOfColours++] = PSDenNum;
      }
    }
    /* Determine if Gamma should be a scalar or a field. */
    
    int UseGammaField = FALSE;
    float *GammaField = NULL;
    if (HydroMethod == Zeus_Hydro && MultiSpecies > 1) {
      UseGammaField = TRUE;
      GammaField = new float[size];
      if (this->ComputeGammaField(GammaField) == FAIL) {
	ENZO_FAIL("Error in grid->ComputeGammaField.");
      }
    } else {
      GammaField = new float[1];
      GammaField[0] = Gamma;

    }
    
    /* Set lowest level flag (used on Zeus hydro). */

    int LowestLevel = (level > MaximumRefinementLevel-1) ? TRUE : FALSE;

    /* Set minimum support (used natively in zeus hydro). */

    float MinimumSupportEnergyCoefficient = 0;
    if (UseMinimumPressureSupport == TRUE && level > MaximumRefinementLevel-1)
      if (this->SetMinimumSupport(MinimumSupportEnergyCoefficient) == FAIL) {
	ENZO_FAIL("Error in grid->SetMinimumSupport,");
      }

    /* allocate space for fluxes */

    /* Set up our restart dump fluxes container */
    this->SubgridFluxStorage = SubgridFluxes;
    this->NumberOfSubgrids = NumberOfSubgrids;

    for (i = 0; i < NumberOfSubgrids; i++) {
      for (dim = 0; dim < GridRank; dim++)  {

	/* compute size (in floats) of flux storage */

        size = 1;
        for (j = 0; j < GridRank; j++)
          size *= SubgridFluxes[i]->LeftFluxEndGlobalIndex[dim][j] -
                  SubgridFluxes[i]->LeftFluxStartGlobalIndex[dim][j] + 1;

	/* set unused dims (for the solver, which is hardwired for 3d). */

        for (j = GridRank; j < 3; j++) {
          SubgridFluxes[i]->LeftFluxStartGlobalIndex[dim][j] = 0;
          SubgridFluxes[i]->LeftFluxEndGlobalIndex[dim][j] = 0;
          SubgridFluxes[i]->RightFluxStartGlobalIndex[dim][j] = 0;
          SubgridFluxes[i]->RightFluxEndGlobalIndex[dim][j] = 0;
        }

	/* Allocate space (if necessary). */

        for (field = 0; field < NumberOfBaryonFields; field++) {
	  //
	  if (SubgridFluxes[i]->LeftFluxes[field][dim] == NULL)
	    SubgridFluxes[i]->LeftFluxes[field][dim]  = new float[size];
	  if (SubgridFluxes[i]->RightFluxes[field][dim] == NULL)
	    SubgridFluxes[i]->RightFluxes[field][dim] = new float[size];
	  for (n = 0; n < size; n++) {
	    SubgridFluxes[i]->LeftFluxes[field][dim][n] = 0;
	    SubgridFluxes[i]->RightFluxes[field][dim][n] = 0;
	  }
        }

	for (field = NumberOfBaryonFields; field < MAX_NUMBER_OF_BARYON_FIELDS;
	     field++) {
          SubgridFluxes[i]->LeftFluxes[field][dim] = NULL;
          SubgridFluxes[i]->RightFluxes[field][dim] = NULL;
	}

      }  // next dimension

      /* make things pretty */

      for (dim = GridRank; dim < 3; dim++)
        for (field = 0; field < MAX_NUMBER_OF_BARYON_FIELDS; field++) {
          SubgridFluxes[i]->LeftFluxes[field][dim] = NULL;
          SubgridFluxes[i]->RightFluxes[field][dim] = NULL;
	}

    } // end of loop over subgrids

    /* compute global start index for left edge of entire grid 
       (including boundary zones) */

    for (dim = 0; dim < GridRank; dim++)
      GridGlobalStart[dim] = nlongint((GridLeftEdge[dim]-DomainLeftEdge[dim])/(*(CellWidth[dim]))) -
	GridStartIndex[dim];

    /* fix grid quantities so they are defined to at least 3 dims */

    for (i = GridRank; i < 3; i++) {
      GridDimension[i]   = 1;
      GridStartIndex[i]  = 0;
      GridEndIndex[i]    = 0;
      GridVelocity[i]    = 0.0;
      GridGlobalStart[i] = 0;
    }

    /* If using comoving coordinates, multiply dx by a(n+1/2).
       In one fell swoop, this recasts the equations solved by solver
       in comoving form (except for the expansion terms which are taken
       care of elsewhere). */

    if (ComovingCoordinates)
      if (CosmologyComputeExpansionFactor(Time+0.5*dtFixed, &a, &dadt) 
	  == FAIL) {
	ENZO_FAIL("Error in CosmologyComputeExpansionFactors.");
      }

    /* Create a cell width array to pass (and convert to absolute coords). */

    float *CellWidthTemp[MAX_DIMENSION];
    for (dim = 0; dim < MAX_DIMENSION; dim++) {
      CellWidthTemp[dim] = new float[GridDimension[dim]];
      for (i = 0; i < GridDimension[dim]; i++)
	if (dim < GridRank)
	  CellWidthTemp[dim][i] = float(a*CellWidth[dim][i]);
	else
	  CellWidthTemp[dim][i] = 1.0;
    }

    /* Prepare Gravity. */

    int GravityOn = 0, FloatSize = sizeof(float);
    if (SelfGravity || UniformGravity || PointSourceGravity)
      GravityOn = 1;
#ifdef TRANSFER
    if (RadiationPressure)
      GravityOn = 1;
#endif    

    /* Call Solver on this grid.
       Notice that it is hard-wired for three dimensions, but it does
       the right thing for < 3 dimensions. */
    /* note: Start/EndIndex are zero based */

    printf("Before PPM: %"GOUTSYM" \n", this->SumGasMass());
    if (HydroMethod == PPM_DirectEuler)
      this->SolvePPM_DE(CycleNumber, NumberOfSubgrids, SubgridFluxes, 
			CellWidthTemp, GridGlobalStart, GravityOn, 
			NumberOfColours, colnum);
    printf("After PPM:  %"GOUTSYM" \n", this->SumGasMass());
    /* PPM LR has been withdrawn. */

    if (HydroMethod == PPM_LagrangeRemap) {
#ifdef PPM_LR
      FORTRAN_NAME(ppm_lr)(
			density, totalenergy, velocity1, velocity2, velocity3,
                          gasenergy,
			&GravityOn, AccelerationField[0],
                           AccelerationField[1],
                           AccelerationField[2],
			&Gamma, &dtFixed, &CycleNumber,
                          CellWidthTemp[0], CellWidthTemp[1], CellWidthTemp[2],
			&GridRank, &GridDimension[0], &GridDimension[1],
                           &GridDimension[2], GridStartIndex, GridEndIndex,
			GridVelocity, &PPMFlatteningParameter,
                           &PressureFree,
			&PPMDiffusionParameter, &PPMSteepeningParameter,
                           &DualEnergyFormalism, &DualEnergyFormalismEta1,
                           &DualEnergyFormalismEta2,
			&NumberOfSubgrids, leftface, rightface,
			istart, iend, jstart, jend,
			standard, dindex, Eindex, uindex, vindex, windex,
			  geindex, temp,
                        &NumberOfColours, colourpt, coloff, colindex);
#else /* PPM_LR */
      ENZO_FAIL("PPM LR is not supported.");
#endif /* PPM_LR */
    }

    if (HydroMethod == Zeus_Hydro)
      if (this->ZeusSolver(GammaField, UseGammaField, CycleNumber, 
			   CellWidthTemp[0], CellWidthTemp[1], CellWidthTemp[2],
			   GravityOn, NumberOfSubgrids, GridGlobalStart,
			   SubgridFluxes,
			   NumberOfColours, colnum, LowestLevel,
			   MinimumSupportEnergyCoefficient) == FAIL)
	ENZO_FAIL("ZeusSolver() failed!\n");
	


    /* Clean up allocated fields. */

    delete [] GammaField;   

    for (dim = 0; dim < MAX_DIMENSION; dim++)
      delete [] CellWidthTemp[dim];

  /* If we're supposed to be outputting on Density, we need to update
     the current maximum value of that Density. */

    if(OutputOnDensity == 1){
      int DensNum = FindField(Density, FieldType, NumberOfBaryonFields);
      for(i = 0; i < size; i++)
        CurrentMaximumDensity =
            max(BaryonField[DensNum][size], CurrentMaximumDensity);
    }

  }  // end: if (NumberOfBaryonFields > 0)


  this->DebugCheck("SolveHydroEquations (after)");

  TIMER_STOP("SolveHydroEquations");
  LCAPERF_STOP("grid_SolveHydroEquations");
  return SUCCESS;

}