Carmen Code
 All Classes Files Functions Variables Macros Pages
Functions | Variables
Parameters.cpp File Reference

User parameters. More...

#include "Carmen.h"
#include "carmen.par"
Include dependency graph for Parameters.cpp:

Functions

void InitParameters ()
 Inits parameters from file carmen.par. If a parameter is not mentioned in this file, the default value is used. More...
 

Variables

real CarmenVersion =2.483
 
int StepNb =2
 
int IterationNb =0
 
int IterationNbRef =1000
 
real CFL =0.5
 
real TimeStep =0.
 
real PhysicalTime =0.
 
int Refresh =0
 
int RefreshNb =200
 
real PrintTime1 =0.
 
real PrintTime2 =0.
 
real PrintTime3 =0.
 
real PrintTime4 =0.
 
real PrintTime5 =0.
 
real PrintTime6 =0.
 
int PrintIt1 =0
 
int PrintIt2 =0
 
int PrintIt3 =0
 
int PrintIt4 =0
 
int PrintIt5 =0
 
int PrintIt6 =0
 
int PrintEvery =0
 
int ImageNb =0
 
bool UseBackup = true
 
bool Recovery =false
 
bool TimeAveraging =false
 
real StartTimeAveraging =0.
 
bool ComputeCPUTimeRef =false
 
int EquationType =7
 
int SchemeNb =1
 
int LimiterNo =5
 
int ScalarEqNb =0
 
int DivClean =2
 
int Dimension =1
 
int Coordinate =1
 
real XMin [4] ={-1.,-1.,-1.,-1.}
 
real XMax [4] ={1.,1.,1.,1.}
 
int CMin [4] ={3,3,3,3}
 
int CMax [4] ={3,3,3,3}
 
bool UseBoundaryRegions = false
 
bool Multiresolution =true
 
bool TimeAdaptivity =false
 
int TimeAdaptivityFactor = 1
 
int ScaleNb =5
 
int ScaleNbRef =2
 
int IcNb = 0
 
int PrintMoreScales =0
 
real Tolerance =0.
 
real ToleranceScale =1.
 
real PenalizeFactor =1.
 
bool ConstantTimeStep =true
 
real ElapsedTime =0.
 
real ExpectedCompression = 0.
 
bool CVS =false
 
bool LES =false
 
real SmoothCoeff =0.
 
int ThresholdNorm =0
 
real Celerity = 1.
 
real Viscosity = 1.
 
real Re =100.
 
real Pr =0.71
 
real Gamma =1.4
 
real Ma =0.3
 
real Fr =0.
 
real TRef =273.
 
real Circulation =10.
 
real ForceX =0.
 
real ForceY =0.
 
real ForceZ =0.
 
real ModelConstant =0.1
 
real ThermalConduction =0.
 
real ConstantForce =true
 
bool Resistivity =false
 
bool Diffusivity =false
 
bool ComputeTemp =false
 
real eta =0.
 
real chi =0.
 
real Alpha =0.64
 
real Ze =10.
 
real Le =1.
 
real Sigma =5.E-02
 
real DIVB =0.
 
real DIVBMax =0.
 
real Bdivergence =0.
 
real PsiGrad =0.
 
real ch =0.0
 
real auxvar =0.
 
bool debug =false
 
bool FluxCorrection =true
 
int PostProcessing =2
 
bool DataIsBinary =true
 
bool ZipData =true
 
bool WriteAsPoints =false
 
real SpaceStep =0.
 
real Eigenvalue =0.
 
real EigenvalueX =0.
 
real EigenvalueY =0.
 
real EigenvalueZ =0.
 
real EigenvalueMax =0.
 
Vector QuantityMax
 
Vector QuantityAverage
 
real pi =acos(-1.)
 
int StepNo =1
 
int IterationNo =0
 
Timer CPUTime
 
int Cluster =1
 
int QuantityNb =2
 
double FVTimeRef =0.
 
int CellNb =0
 
int LeafNb =0
 
real TotalLeafNb =0.
 
real TotalCellNb =0.
 
real ErrorMax =0.
 
real ErrorGlobalMax =0.
 
real ErrorMid =0.
 
real ErrorGlobalMid =0.
 
real ErrorL2 =0.
 
real ErrorGlobalL2 =0.
 
int ErrorNb =0
 
int ErrorGlobalNb =0
 
real RKFError =0.
 
real RKFAccuracyFactor =1.E-03
 
real RKFSafetyFactor =0.01
 
real cr =0.
 
real Helicity =0.
 
real GlobalMomentum =0.
 
real GlobalEnergy =0.
 
real GlobalEnstrophy =0.
 
real ExactMomentum =0.
 
real ExactEnergy =0.
 
real FlameVelocity =0.
 
real GlobalMomentumOld =0.
 
real GlobalVolume =0.
 
real GlobalReactionRate =0.
 
real AverageRadius =0.
 
real PreviousAverageRadius =0.
 
real PreviousAverageRadius2 =0.
 
real XCenter [4] ={0.,0.,0.,0.}
 
real ReactionRateMax =0.
 
FILE * GlobalFile
 
int ChildNb =0
 
real IntVorticity =0.
 
real IntDensity =0.
 
Vector IntMomentum
 
real IntEnergy =0.
 
real BaroclinicEffect =0.
 
int rank
 
int size
 
int CPUScales
 
real AllXMin [4]
 
real AllXMax [4]
 
int AllTaskScaleNb
 
int NeighbourNb
 
int coords [3]
 
int CartDims [3]
 
int CellElementsNb
 
int MaxCellElementsNb
 
int one_D
 
int two_D
 
int MPISendType
 
int MPIRecvType
 
Timer CommTimer
 
int rank_il
 
int rank_iu
 
int rank_jl
 
int rank_ju
 
int rank_kl
 
int rank_ku
 
int WhatSend
 
int SendD = 1 << 0
 
int SendGrad = 1 << 1
 
int SendQ = 1 << 2
 
int SendQs = 1 << 3
 
int SendX = 1 << 4
 
int SenddX = 1 << 5
 
char BackupName [255]
 

Detailed Description

User parameters.

Function Documentation

void InitParameters ( )

Inits parameters from file carmen.par. If a parameter is not mentioned in this file, the default value is used.

Returns
void

— Compute ch --------------------------------------------------------—

284 {
285  // --- Local variables ---------------------------------------------------------------------
286 
287  int i; // Counter
288 
289  // --- Set global variables from file "carmen.par" -----------------------------------------
290 
291  #include "carmen.par"
292 
293  // --- Adapt IterationNbRef to the dimension -----------------------------------------------
294 
295  IterationNbRef=(int)(exp((4.-Dimension)*log(10.)));
296 
297  // --- Compute the number of children of a given parent cell ---
298 
299  ChildNb = (1<<Dimension);
300 
301 #if defined PARMPI
302 
304  for (i=0;i<4;i++)
305  {
306  AllXMax[i]=XMax[i];
307  AllXMin[i]=XMin[i];
308  }
309 
310  //some combinations give deadlock...
311  MPISendType = 10; //0 - Ibsend; 10 - Isend; 20 - Issend;
312  MPIRecvType = 1; //0 - Recv; 1 - Irecv;
313 
314  CPUScales=0;
315  int tmp=size;
316  while ((tmp=(tmp>>1))>0) CPUScales++;
318 
319  one_D=1; two_D=1;
320  if (Dimension >= 2) one_D=1<<ScaleNb;
321  if (Dimension == 3) two_D=1<<ScaleNb;
322 
323 //#if defined PARMPI
324 
325  NeighbourNb=2;
327 
328  // -- Create memory arrays thats are needs for the MPI Type creation ---
329  disp = new MPI_Aint[NeighbourNb*MaxCellElementsNb*one_D*two_D];
330  blocklen = new int [NeighbourNb*MaxCellElementsNb*one_D*two_D];
331 
332  // --- Allocate additional memory for MPI buffer send---
333  Cell tc;
334  int CellElNb,bufsize;
335  CellElNb=tc.size().dimension()+tc.center().dimension()+tc.average().dimension()+tc.tempAverage().dimension()+tc.divergence().dimension();
336 
337  if (EquationType==6)
338  CellElNb += tc.gradient().lines()*tc.gradient().columns();
339 
340  bufsize=(CellElNb*one_D*two_D*NeighbourNb+MPI_BSEND_OVERHEAD)*2*Dimension+1024;
341  MPIbuffer=new real[bufsize];
342  MPI_Buffer_attach(MPIbuffer,bufsize*sizeof(real));
343 
344 #else
345  NeighbourNb=0;
346 #endif
347 
348 #if defined PARMPI
350 
351  // --- Compute domain coordinates for the processors ---
352  XMin[1] = AllXMin[1] + coords[0]*(AllXMax[1]-AllXMin[1])/CartDims[0];
353  XMax[1] = AllXMin[1] + (coords[0]+1)*(AllXMax[1]-AllXMin[1])/CartDims[0];
354 
355  if (Dimension >= 2)
356  {
357  XMin[2] = AllXMin[2] + coords[1]*(AllXMax[2]-AllXMin[2])/CartDims[1];
358  XMax[2] = AllXMin[2] + (coords[1]+1)*(AllXMax[2]-AllXMin[2])/CartDims[1];
359  }
360 
361  if (Dimension == 3)
362  {
363  XMin[3] = AllXMin[3] + coords[2]*(AllXMax[3]-AllXMin[3])/CartDims[2];
364  XMax[3] = AllXMin[3] + (coords[2]+1)*(AllXMax[3]-AllXMin[3])/CartDims[2];
365  }
366 
367  // --- Set the backup file name for the current processor
368  sprintf(BackupName,"%d_%d_%d_%s",coords[0],coords[1],coords[2],"carmen.bak");
369 #else
370  sprintf(BackupName,"%s","carmen.bak");
371 #endif
372 
373 
374  // --- Use CVS only if Dimension > 1 -------------------------------------------------------
375 
376  if (Dimension == 1)
377  CVS = false;
378 
379  // --- TimeAveraging always false if not Navier-Stokes -------------------------------------
380 
381  if (EquationType != 6)
382  TimeAveraging = false;
383 
384  // --- If there is no file "carmen.bak", set Recovery=false --------------------------------
385 
386  if (!fopen(BackupName,"r"))
387  Recovery = false;
388 
389  // --- If PrintMoreScales != 0 or 1 with Multiresolution = false, print error and stop ---
390 
391  if (!Multiresolution && (!(PrintMoreScales == 0 || PrintMoreScales == -1)) )
392  {
393  cout << "Parameters.cpp: In method `void InitParameters()':\n";
394  cout << "Parameters.cpp: value of PrintMoreScales incompatible with FV computations\n";
395  cout << "Parameters.cpp: must be 0 or -1\n";
396  cout << "carmen: *** [Parameters.o] Execution error\n";
397  cout << "carmen: abort execution.\n";
398  exit(1);
399  }
400 
401  // --- Compute global volume ---------------------------------------------------------------
402 
403  GlobalVolume = fabs(XMax[1]-XMin[1]);
404 
405  if (Dimension > 1)
406  GlobalVolume *= fabs(XMax[2]-XMin[2]);
407 
408  if (Dimension > 2)
409  GlobalVolume *= fabs(XMax[3]-XMin[3]);
410 
411  // --- Compute PostProcessing and DataIsBinary ---------------------------------------------
412 
413  // In 1D, use Gnuplot instead of Data Explorer
414 
415  if (Dimension == 1 && PostProcessing == 2)
416  PostProcessing = 1;
417 
418  // In 2D-3D, use Data Explorer instead of Gnuplot
419 
420  if (Dimension != 1 && PostProcessing == 1)
421  PostProcessing = 2;
422 
423  // --- Compute number of conservative quantities -------------------------------------------
424 
425  QuantityNb = 9;
426 
427  // --- Set the dimension of QuantityMax to QuantityNb ---------------------------------------
428 
430 
431  // --- Set the dimension of QuantityAverage to 4 (pressure, vorticity, entropy, volume)
432 
434 
435  // --- Set the dimension of IntMomentum to dimension ----------------------------------------
436 
438 
439  // --- Compute minimal space step -----------------------------------------------------------
440 
441  SpaceStep = fabs(XMax[1]-XMin[1]);
442 
443  for (i = 2; i <= Dimension; i++)
444  SpaceStep = Min(SpaceStep, fabs(XMax[i]-XMin[i]));
445 
446  SpaceStep /= (1<<ScaleNb);
447 
448  // --- Compute time step from CFL if TimeStep = 0 -------------------------------------------
449 
450  if (TimeStep == 0.)
451  {
452  if (fabs(Eigenvalue)>0.0e-20)
454  else
455  TimeStep = 0.0001;
456  }
457  else
458  ConstantTimeStep = true;
459 
462 
463 }
int one_D
Definition: Parameters.cpp:236
An object Cell contains all the informations of a cell for both multiresolution and finite volume com...
Definition: Cell.h:41
real ch
Definition: Parameters.cpp:140
int ScaleNb
Definition: Parameters.cpp:87
int CPUScales
Definition: Parameters.cpp:225
int IterationNbRef
Definition: Parameters.cpp:38
real center(const int AxisNo) const
Returns the component no. AxisNo of the cell-center position.
Definition: Cell.h:1112
bool TimeAveraging
Definition: Parameters.cpp:60
int QuantityNb
Definition: Parameters.cpp:171
real XMax[4]
Definition: Parameters.cpp:77
int MPIRecvType
Definition: Parameters.cpp:239
int PrintMoreScales
Definition: Parameters.cpp:90
real Eigenvalue
Definition: Parameters.cpp:157
int ChildNb
Definition: Parameters.cpp:213
int coords[3]
Definition: Parameters.cpp:230
real SpaceStep
Definition: Parameters.cpp:156
int Dimension
Definition: Parameters.cpp:74
real GlobalVolume
Definition: Parameters.cpp:203
real AllXMax[4]
Definition: Parameters.cpp:226
int AllTaskScaleNb
Definition: Parameters.cpp:227
#define Min(x, y)
Definition: Carmen.h:62
Vector IntMomentum
Definition: Parameters.cpp:217
int MPISendType
Definition: Parameters.cpp:238
real divergence(const int QuantityNo) const
Returns the component no. QuantityNo of the divergence vector.
Definition: Cell.h:1192
int PostProcessing
Definition: Parameters.cpp:144
real size(const int AxisNo) const
Returns the cell size in the direction AxisNo.
Definition: Cell.h:1095
bool CVS
Definition: Parameters.cpp:97
real average(const int QuantityNo) const
Returns the component no. QuantityNo of the cell-average value.
Definition: Cell.h:1128
void CreateMPITopology()
Parallel function DOES NOT WORK!
Definition: Parallel.cpp:22
real AllXMin[4]
Definition: Parameters.cpp:226
real XMin[4]
Definition: Parameters.cpp:76
real tempAverage(const int QuantityNo) const
Returns the component no. QuantityNo of the temporary cell-average value.
Definition: Cell.h:1144
char BackupName[255]
Definition: Parameters.cpp:257
real TimeStep
Definition: Parameters.cpp:40
int two_D
Definition: Parameters.cpp:236
Vector QuantityAverage
Definition: Parameters.cpp:163
int CartDims[3]
Definition: Parameters.cpp:231
bool ConstantTimeStep
Definition: Parameters.cpp:94
int NeighbourNb
Definition: Parameters.cpp:228
int size
Definition: Parameters.cpp:224
bool Recovery
Definition: Parameters.cpp:59
Vector QuantityMax
Definition: Parameters.cpp:162
int EquationType
Definition: Parameters.cpp:66
bool Multiresolution
Definition: Parameters.cpp:84
real CFL
Definition: Parameters.cpp:39
real gradient(const int i, const int j) const
Returns the component no. i, j of the velocity gradient. Does not work for MHD!
Definition: Cell.h:1225
#define real
Definition: PreProcessor.h:31
int MaxCellElementsNb
Definition: Parameters.cpp:234
void setDimension(const int n)
Sets the dimension of the vector to n and reset values to zero.
Definition: Vector.cpp:298

Here is the caller graph for this function:

Variable Documentation

int AllTaskScaleNb

Global Scale number

real AllXMax[4]
real AllXMin[4]

Global domain parameters

real Alpha =0.64

Temperature ratio

real auxvar =0.

Auxiliar variable

real AverageRadius =0.

Average radius of the flame ball

char BackupName[255]
real BaroclinicEffect =0.

Intensity of the baroclinic effects

real Bdivergence =0.

Auxilian divergence variable

real CarmenVersion =2.483

Version release

int CartDims[3]
real Celerity = 1.

Advection-diffusion celerity (0, 1, -1)

int CellElementsNb
int CellNb =0

Number of cells

real CFL =0.5

Courant-Friedrich-Levy number

real ch =0.0

Divergence cleaning ch parameter

real chi =0.

Artificial diffusion constant

int ChildNb =0

Number of children for a given parent (equal to 2**Dimension)

real Circulation =10.

Circulation parameter

int Cluster =1

0 for local execution, 1 for cluster

int CMax[4] ={3,3,3,3}

Max. boundary condition (1 = Boundary, 2 = Symetric, 3 = Periodic)

int CMin[4] ={3,3,3,3}

Min. boundary condition (1 = Boundary, 2 = Symetric, 3 = Periodic)

Timer CommTimer

Communication timer for perfomance analyse

bool ComputeCPUTimeRef =false

True = the reference CPU time is being computed

bool ComputeTemp =false
real ConstantForce =true

False = adapt force to maintain constant energy

bool ConstantTimeStep =true

true = constant TimeStep

int Coordinate =1

1 = Cartesian, 2 = Spherical in x

int coords[3]

Current CPU coodinates in the virtual CPU processors cart

int CPUScales

CPUScales=log2(Number of processors)

Timer CPUTime

Timer for CPU time

real cr =0.

Alpha parameter divergence cleaning

bool CVS =false

True = use Donoho thresholding to perform CVS.

bool DataIsBinary =true

true = write data in binary format, false = write data in ASCII format

bool debug =false

true = check if tree is graded

bool Diffusivity =false

True = artificial diffusion. False = no artificial diffusion

int Dimension =1

Dimension (1,2,3)

real DIVB =0.

Divergence of B

real DIVBMax =0.

Maximum divergence of B

int DivClean =2

Divergence cleaning: 1-EGLM 2-GLM

real Eigenvalue =0.

Maximal eigenvalue

real EigenvalueMax =0.

Maximal eigenvalue

real EigenvalueX =0.

Eigenvalue at x direction

real EigenvalueY =0.

Eigenvalue at y direction

real EigenvalueZ =0.

Eigenvalue at z direction

real ElapsedTime =0.

ElapsedTime

int EquationType =7

Type of equation

real ErrorGlobalL2 =0.

Global L2 error on space and time

real ErrorGlobalMax =0.

Global error max on space and time

real ErrorGlobalMid =0.

Global mean error on space and time

int ErrorGlobalNb =0

Number of points for the computation of the global mean error

real ErrorL2 =0.

L2 error on the grid

real ErrorMax =0.

Error Max on the grid

real ErrorMid =0.

Mean error on the grid

int ErrorNb =0

Number of points for the computation of the mean error

real eta =0.

Resistivity function

real ExactEnergy =0.

Global energy for the exact solution (only for EquationType = 1 or 2)

real ExactMomentum =0.

Global momentum for the exact solution (only for EquationType = 1 or 2)

real ExpectedCompression = 0.

Expected memory compression

real FlameVelocity =0.

Flame velocity

bool FluxCorrection =true

true = conservative flux correction

real ForceX =0.
real ForceY =0.
real ForceZ =0.
real Fr =0.

Froude number

double FVTimeRef =0.

FV reference CPU time for 1 iteration

real Gamma =1.4

Adiabatic function

real GlobalEnergy =0.

Global energy (only for EquationType = 1 or 2)

real GlobalEnstrophy =0.

Global enstrophy (only for EquationType = 6)

FILE* GlobalFile

Global file

real GlobalMomentum =0.

Global momentum (only for EquationType = 1 or 2)

real GlobalMomentumOld =0.

Old global momentum (only for EquationType = 6)

real GlobalReactionRate =0.

Global reaction rate

real GlobalVolume =0.

Global volume

real Helicity =0.

Time derivative of helicity (must be zero)

int IcNb = 0

Initial condition suavization

int ImageNb =0

Print ImageNb images

real IntDensity =0.

Integral of the density

real IntEnergy =0.

Integral of the energy

Vector IntMomentum

Integral of the modulus of the momentum

real IntVorticity =0.

Integral of the modulus of the vorticity

int IterationNb =0

Number of iterations

int IterationNbRef =1000

Number of iterations for the FV reference computation

int IterationNo =0

Current iteration number

real Le =1.

Lewis number

int LeafNb =0

Number of leaves

bool LES =false

True = use eddy-viscosity.

int LimiterNo =5

Limiter number

real Ma =0.3

Mach number

int MaxCellElementsNb
real ModelConstant =0.1

Constant used in the turbulence model

int MPIRecvType
int MPISendType

Type of calling MPI communication functions. See Parameters.cpp for the more information

bool Multiresolution =true

true = Multiresolution, false = FV on fine mesh

int NeighbourNb

Important parameter: The deep of the inter-CPU domain overlapping

int one_D
real PenalizeFactor =1.

Factor of penalization (obsolete)

real PhysicalTime =0.

Physical elapsed time

real pi =acos(-1.)

Pi constant

int PostProcessing =2

1 = Gnuplot, 2 = Data Explorer, 3 = Tecplot

real Pr =0.71

Prandtl number

real PreviousAverageRadius =0.

Previous average radius

real PreviousAverageRadius2 =0.

Previous average radius

int PrintEvery =0

Print every PrintEvery iteration

int PrintIt1 =0

Iteration for print 1

int PrintIt2 =0

Iteration for print 2

int PrintIt3 =0

Iteration for print 3

int PrintIt4 =0

Iteration for print 4

int PrintIt5 =0

Iteration for print 5

int PrintIt6 =0

Iteration for print 6

int PrintMoreScales =0

More scales to print

real PrintTime1 =0.

Time for print 1

real PrintTime2 =0.

Time for print 2

real PrintTime3 =0.

Time for print 3

real PrintTime4 =0.

Time for print 4

real PrintTime5 =0.

Time for print 5

real PrintTime6 =0.

Time for print 6

real PsiGrad =0.

Psi gradient

Vector QuantityAverage

Vector containing the average quantities

Vector QuantityMax

Vector containing the maximal quantities

int QuantityNb =2

Number of conservative quantites

int rank

Current CPU

int rank_il

axis X, direction - low

int rank_iu

axis X, direction - high

int rank_jl

axis Y, direction - low

int rank_ju
int rank_kl
int rank_ku
real Re =100.

Reynolds number

real ReactionRateMax =0.

Maximum of the reaction rate

bool Recovery =false

true = restore data from previous computation

int Refresh =0

/ Refresh rate

int RefreshNb =200

Number of refreshments

bool Resistivity =false

True = resistive model. False = ideal model

real RKFAccuracyFactor =1.E-03

Desired value of RKFError (only if TimeAdaptivity=true)

real RKFError =0.

Maximum of the relative errors between RK2 and RK3 (only if TimeAdaptivity=true)

real RKFSafetyFactor =0.01

Safety factor for the computation of the time step (only if TimeAdaptivity=true)

int ScalarEqNb =0

Number of additional scalar equations

int ScaleNb =5

Maximal number of scales allowed

int ScaleNbRef =2

Number of scales for the FV reference computation

int SchemeNb =1

Scheme number

int SendD = 1 << 0
int SenddX = 1 << 5
int SendGrad = 1 << 1
int SendQ = 1 << 2
int SendQs = 1 << 3
int SendX = 1 << 4
real Sigma =5.E-02

Radiation coefficient

int size

Number of processors

real SmoothCoeff =0.

Smoothing coefficient

real SpaceStep =0.

Space step

real StartTimeAveraging =0.

Time where the time-averaging procedure must start

int StepNb =2

Number of steps for the time integration

int StepNo =1

Current step number for the time integration

real ThermalConduction =0.

Dimensionless thermal conduction

int ThresholdNorm =0

Normalization of the wavelet basis for the threshold

bool TimeAdaptivity =false

true = use time adaptivity (only when Multiresolution = true, dummy else)

int TimeAdaptivityFactor = 1

The factor of time step between two scales is 2^TimeAdaptivityFactor

bool TimeAveraging =false

true = use a time-averaging grid (only for turbulence)

real TimeStep =0.

Time step

real Tolerance =0.

Prescribed tolerance

real ToleranceScale =1.

Scale factor for tolerance (when Tolerance = 0).

real TotalCellNb =0.

Total number of cells for all iterations

real TotalLeafNb =0.

Total number of leaves for all iterations

real TRef =273.

Reference temperature for Sutherland's law

int two_D
bool UseBackup = true

true = use Backup procedure.

bool UseBoundaryRegions = false

true = use file carmen.bc

real Viscosity = 1.

0 or 1. 0 means no viscosity.

int WhatSend
bool WriteAsPoints =false

true = write data as point-values, false = write data as cell-averages

real XCenter[4] ={0.,0.,0.,0.}

Coordinate of the center of the flame ball

real XMax[4] ={1.,1.,1.,1.}

Maximal values of coordinates;

real XMin[4] ={-1.,-1.,-1.,-1.}

Minimal values of coordinates;

real Ze =10.

Equivalent to Beta

bool ZipData =true

true = zip data files