cloudy  trunk
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
parse_commands.cpp
Go to the documentation of this file.
1 /* This file is part of Cloudy and is copyright (C)1978-2022 by Gary J. Ferland and
2  * others. For conditions of distribution and use see copyright notice in license.txt */
3 /*ParseCommands main command line parser, decode command, then call other routines to read */
4 #include "cddefines.h"
5 #include "parse.h"
6 #include "stopcalc.h"
7 #include "abund.h"
8 #include "geometry.h"
9 #include "dense.h"
10 #include "plot.h"
11 #include "grid.h"
12 #include "grainvar.h"
13 #include "dynamics.h"
14 #include "magnetic.h"
15 #include "trace.h"
16 #include "atmdat.h"
17 #include "h2.h"
18 #include "rt.h"
19 #include "thermal.h"
20 #include "opacity.h"
21 #include "called.h"
22 #include "wind.h"
23 #include "hextra.h"
24 #include "iterations.h"
25 #include "radius.h"
26 #include "input.h"
27 #include "monitor_results.h"
28 #include "phycon.h"
29 #include "fudgec.h"
30 #include "version.h"
31 #include "conv.h"
32 #include "cosmology.h"
33 #include "pressure.h"
34 #include "parser.h"
35 #include "dark_matter.h"
36 #include "iso.h"
37 #include "mole.h"
38 #include "parse_species.h"
39 #include "doppvel.h"
40 #include "rfield.h"
41 #include "prt.h"
42 
43 void ParseAperture(Parser &p);
45 void ParseCExtra(Parser &p);
46 void ParseCMBOuter(Parser &p);
47 void ParseCosm(Parser &p);
48 void ParseCovering(Parser &p);
49 void ParseCylinder(Parser &p);
50 void ParseDarkMatter(Parser &p);
51 void ParseDatabase(Parser &p);
53 void ParseDiffuse(Parser &p);
54 void ParseDistance(Parser &p);
55 void ParseDoubleTau(Parser &);
56 void ParseEden(Parser &p);
57 void ParseEnergy(Parser &p);
58 void ParseFail(Parser &p);
59 void ParseFill(Parser &p);
60 void ParseF_nuSpecific(Parser &p);
62 void ParseFudge(Parser &p);
63 void ParsePGrains(Parser &);
64 void ParseGravity(Parser &p);
65 void ParseHeLike(Parser &);
66 void ParseHelp(Parser &);
67 void ParseHExtra(Parser &p);
68 void ParseConvHighT(Parser &);
69 void ParseHydrogen(Parser &);
70 void ParseInitCount(Parser &p);
71 void ParseIntensity(Parser &p);
72 void ParseIterations(Parser &p);
73 void ParseL_nu(Parser &p);
74 void ParseLaser(Parser &p);
75 void ParseLuminosity(Parser &p);
76 void ParseNeutrons(Parser &p);
77 void ParseNuF_nu(Parser &p);
78 void ParseNuL_nu(Parser &p);
79 void ParsePhi(Parser &p);
80 void ParseQH(Parser &p);
81 void ParseRoberto(Parser &);
82 void ParseSpecial(Parser &);
83 void ParseTauMin(Parser &p);
84 void ParseTitle(Parser &);
85 void ParseTolerance(Parser &);
86 void ParseVLaw(Parser &p);
87 void ParseTurbulence(Parser &p);
88 
89 void ParseCommands(void)
90 {
91  bool lgStop ,
92  lgStop_not_enough_info;
93 
94  DEBUG_ENTRY( "ParseCommands()" );
95 
96  /* following says abundances are solar */
97  abund.lgAbnSolar = true;
98 
99  /* there are no plots desired yet */
100  plotCom.lgPlotON = false;
101  plotCom.nplot = 0;
102 
103  /* this flag remembers whether grains have ever been turned on. It is passed
104  * to routine ParseAbundances - there grains will be turned on with commands
105  * such as abundances ism, unless grains were previously set
106  * with a grains command. this way abundances will not clobber explicitly set
107  * grains. */
108 
109  rfield.nShape = 0;
110 
111  /* initialize some code variables in case assert command used in input stream */
113 
114  for( long int i=0; i < LIMSPC; i++ )
115  {
116  strcpy( rfield.chRSpec[i], "UNKN" );
117  }
118  optimize.nparm = 0;
119 
120  /* this is an option to turn on trace printout on the nth
121  * call from the optimizer */
122  if( optimize.lgTrOpt )
123  {
124  /* nTrOpt was set with "optimize trace" command,
125  * is iteration to turn on trace */
127  {
128  trace.lgTrace = true;
129  called.lgTalk = cpu.i().lgMPI_talk();
130  trace.lgTrOvrd = true;
131  fprintf( ioQQQ, " READR turns on trace from optimize option.\n" );
132  }
133  }
134 
135  if ( called.lgTalk && prt.lgPrintHTML )
136  {
137  fprintf( ioQQQ, "<!DOCTYPE html PUBLIC "
138  "\"-//W3C//DTD XHTML 1.0 Transitional//EN\" "
139  "\"http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd\">\n");
140  fprintf( ioQQQ,"<html xmlns=\"http://www.w3.org/1999/xhtml\">\n");
141  fprintf( ioQQQ,"<head>\n");
142  fprintf( ioQQQ,"<meta http-equiv=\"Content-Type\" "
143  "content=\"text/html; charset=UTF-8\" />\n");
144  fprintf( ioQQQ,"<title>Cloudy output</title>\n");
145  fprintf( ioQQQ,"</head>\n");
146  fprintf( ioQQQ,"<body>\n");
147  fprintf( ioQQQ,"<pre>\n");
148  }
149 
150  /* say this is a beta version if we are talking and it is the truth */
151  if( t_version::Inst().nBetaVer > 0 && called.lgTalk )
152  {
153  fprintf( ioQQQ,
154  "\n This is a beta release of Cloudy, and is intended for testing only.\n" );
155  fprintf( ioQQQ,
156  "Please help make Cloudy better by posing problems or suggestions on cloudyastrophysics.groups.io.\n\n" );
157  }
158 
159  if( called.lgTalk )
160  {
161  /* this code prints pretty lines at top of output box */
162  int indent = (int)((122 - strlen(t_version::Inst().chVersion))/2);
163  fprintf( ioQQQ, "%*cCloudy %s\n", indent, ' ', t_version::Inst().chVersion );
164 
165  fprintf( ioQQQ, "%57cwww.nublado.org\n\n", ' ' );
166 
167  /* now print box and date of version, before printing commands */
168  fprintf( ioQQQ, "%23c", ' ' );
169  fprintf( ioQQQ, "**************************************");
170  fprintf( ioQQQ, "%7.7s", t_version::Inst().chDate);
171  fprintf( ioQQQ, "**************************************\n");
172 
173  fprintf( ioQQQ, "%23c*%81c*\n", ' ', ' ' );
174  }
175 
176  /* read in commands and print them */
177 
178  /* initialize array reader, this sub does nothing but set
179  * initial value of a variable */
180  input.init();
181 
182  static const CloudyCommand commands[] = {
183  {"ABSO",ParseAbsMag},
184  /* enter luminosity in absolute magnitudes, in reads2 */
185  {"AGE",ParseAge},
186  /* enter age of cloud so we can check for time-steady reads2 */
187  {"AGN ",ParseAgn},
188  /* enter generic style AGN continuum, in reads2 */
189  {"ABUN",ParseAbundances},
190  {"APER",ParseAperture},
191  {"BACK", ParseBackgrd},
192  /* cosmic background, in parse_backgrd */
193  {"BLAC", ParseBlackbody},
194  /* black body, in reads2 */
195  {"BREM", ParseBremsstrahlung},
196  {"CASE",ParseCaseB},
197  /* do Case A or Case B (usually Case B since most realistic) */
198  {"CEXT",ParseCExtra},
199  {"CHEM",ParseChemistry},
200  /* CMB command */
201  {"CLUM",ParseFill},
202  {"CMB", ParseCMBOuter},
203  /* cosmic thermal background radiation, argument is redshift */
204  /* if no number on line then (realnum)FFmtRead returns z=0; i.e., now */
205  {"COMP", ParseCompile},
206  /* compile ascii version of stellar atmosphere continua in volk */
207  {"CONS",ParseConstant},
208  /* constant temperature, pressure, density, or gas pressure
209  * in readsun */
210  {"CORO",ParseCoronal},
211  /* coronal equilibrium; set constant temperature to number on line
212  * in readsun */
213  {"COSMOLOGY",ParseCosmology},
214  {"COSMI", ParseCosmicRays},
215  {"COSM",ParseCosm}, // Backstop for ambiguity
216  {"COVE",ParseCovering},
217  {"CRAS",ParseCrashDo},
218  /* any of several tests to cause the code to crash */
219  {"CYLI",ParseCylinder},
220  {"DARK",ParseDarkMatter},
221  {"DATA",ParseDatabase},
222  {"ATOM",ParseDatabase},// accept the old ATOM command as a pseudonym for DATABASE
223  {"DIEL",ParseDielectronic},
224  {"DIFF",ParseDiffuse},
225  {"DIST",ParseDistance},
226  {"DLAW",ParseDLaw},
227  /* either use dense_fabden routine, or read in table of points */
228  {"DOUB",ParseDoubleTau},
229  /* double optical depth scale after each iteration */
230  {"DRIV",ParseDrive},
231  /* call one of several drivers, readsun */
232  {"EDEN",ParseEden},
233  /* option to add "extra" electrons, to test Compton limit
234  * for very low T(star) - option is log of eden */
235  {"ELEM",ParseElement},
236  /* element command;
237  * scale or abundance options, to change abundance of specific element
238  * read option to change order of elements
239  * in reads2.f */
240  {"ENER",ParseEnergy},
241  {"EXTI",ParseExtinguish},
242  /* extinguish ionizing continuum by absorbing column AFTER
243  * setting luminosity or Q(H). First number is the column
244  * density (log), second number is leakage (def=0%)
245  * last number is lowest energy (ryd), last two may be omitted
246  * from right to left
247  *
248  * extinction is actually done in extin, which is called by ContSetIntensity */
249  {"FAIL",ParseFail},
250  /* reset number of temp failures allowed, default=20 */
251  {"FILL",ParseFill},
252  {"FLUC",ParseFluc},
253  /* rapid density fluctuations, in readsun */
254  {"F(NU",ParseF_nuSpecific},
255  /* this is the specific flux at nu
256  * following says F(nu) not nuF(nu) */
257  {"FORC",ParseForceTemperature},
258  /* force temperature of first zone, but don't keep constant
259  * allow to then go to nearest equilibrium
260  * log if < 10 */
261  {"FUDG",ParseFudge},
262  {"GLOB",ParseGlobule},
263  /* globule with density increasing inward
264  * parameters are outer density, radius of globule, and density power */
265  {"GRAI",ParseGrain},
266  /* read parameters dealing with grains, in reads2 */
267  {"PGRA",ParsePGrains},
268  {"GRAV",ParseGravity},
269  /* (self-)Gravity forces: Yago Ascasibar (UAM, Spring 2009) */
270  {"GRID",ParseGrid},
271  /* option to run grid of models by varying certain parameters
272  * in reads2 */
273  {"HDEN",ParseHDEN},
274  /* parse the hden command to set the hydrogen density, in reads2 */
275  {"HELI",ParseHeLike},
276  {"HELP",ParseHelp},
277  {"HEXT",ParseHExtra},
278  /* "extra" heating rate, so that first= log(erg(cm-3, s-1),
279  * second optional number is scale radius, so that HXTOT = TurbHeat*SEXP(DEPTH/SCALE)
280  * if missing then constant heating.
281  * third option is depth from shielded face, to mimic irradiation from both sides*/
282  {"HIGH",ParseConvHighT},
283  /* approach equilibrium from high te */
284  {"HYDROGEN",ParseHydrogen},
285  {"ILLU",ParseIlluminate},
286  // illuminate command
287  {"INIT",ParseInitCount},
288  {"INTEN",ParseIntensity},
289  {"INTER",ParseInterp},
290  /* interpolate on input tables for continuum, set of power laws used
291  * input ordered pairs nu( ryd or log(Hz),, log(fnu)
292  * additional lines begin CONTINUE
293  * first check that this is the one and only INTERP command
294  * in readsun */
295  {"IONI",ParseIonParI},
296  /* inter ionization parameter U=Q/12 R*R N C;
297  * defined per hydrogen, not per electron (as before)
298  * radius must also be entered if spherical, not needed if plane */
299  {"ITER",ParseIterations},
300  {"L(NU",ParseL_nu},
301  {"LASE",ParseLaser},
302  {"LUMI",ParseLuminosity},
303  {"MAGN", ParseMagnet},
304  /* parse the magnetic field command, routine in magnetic.c */
305  {"MAP ",ParseMap},
306  /* do cooling space map for specified zones, in reads2 */
307  {"META",ParseMetal},
308  /* read depletion for metals, all elements heavier than He
309  * in reads2 */
310  {"MONI",ParseMonitorResults},
311  /* monitor that code predicts certain results, in monitor_results.h */
312  {"NEUT",ParseNeutrons},
313  {"NO ", ParseDont},
314  /* don't do something, in readsun */
315  {"NORM",ParseNorm},
316  /* normalize lines to this rather than h-b, sec number is scale factor */
317  {"NUF(",ParseNuF_nu},
318  {"NUL(",ParseNuL_nu},
319  {"OPTI",ParseOptimize},
320  /* option to optimize the model by varying certain parameters
321  * in reads2 */
322  {"PHI(",ParsePhi},
323  {"PLOT", ParsePlot},
324  /* make plot of log(nu.f(nu), vs log(nu) or opacity
325  * in file plot */
326  {"POWE", ParsePowerlawContinuum},
327  /* power law with cutoff, in reads2 */
328  {"PRIN",ParsePrint},
329  /* adjust print schedule, in readsun */
330  {"PUNC",ParseSave},
331  /* save something, in save */
332  {"Q(H)",ParseQH},
333  {"RATI",ParseRatio},
334  /* enter a continuum luminosity as a ratio of
335  * nuFnu for this continuum relative to a previous continuum
336  * format; first number is ratio of second to first continuum
337  * second number is energy for this ratio
338  * if third number on line, then 2nd number is energy of
339  * first continuum, while 3rd number is energy of second continuum
340  * in reads2 */
341  {"RADI", ParseRadius},
342  /* log of inner and outer radii, default second=infinity,
343  * if R2<R1 then R2=R1+R2
344  * there is an optional keyword, "PARSEC" on the line
345  * to use PC as units, reads2 */
346  {"ROBE",ParseRoberto},
347  {"SAVE",ParseSave},
348  {"SET ",ParseSet},
349  {"SPECIAL", ParseSpecial},
350  {"SPECIES",ParseSpecies},
351  {"SPHE", ParseSphere},
352  /* compute a spherical model, diffuse field from other side in
353  * in reads2 */
354  {"STAT",ParseState},
355  /* either get or put the code's state as a file on disk */
356  {"STOP",ParseStop},
357  /* stop model at desired zone, temperature, column density or tau-912
358  * in readsun */
359  {"TABL",ParseTable},
360  /* interpolate on input tables for continuum, set of power laws used
361  * input stored in big BLOCK data
362  * first check that this is the one and only INTERP command
363  * in readsun */
364  {"TAUM",ParseTauMin},
365  {"TEST",ParseTest},
366  /* parse the test command and its options */
367  {"TIME",ParseDynaTime},
368  /* parse the time dependent command, in dynamics.c */
369  {"TITL",ParseTitle},
370  {"TLAW", ParseTLaw},
371  /* some temperature vs depth law */
372  {"TOLE", ParseTolerance},
373  {"TRAC", ParseTrace},
374  /* turn on trace output, in reads2 */
375  {"VLAW",ParseVLaw},
376  {"TURB",ParseTurbulence},
377  {"WIND",ParseDynaWind},
378  /* NB - advection and wind commands are now a single command */
379  /* parse the wind command, in dynamics.c */
380  {"XI",ParseIonParX},
381  {NULL,NULL}}; // {NULL,NULL} sentinel must come last
382 
383  Parser p(commands);
384 
385  p.m_nqh = 0;
386  p.m_lgDSet = false;
387  p.m_lgEOF = false;
388 
389  // set default solar abundances
390  char chDUMMY[INPUT_LINE_LENGTH];
391  sprintf(chDUMMY, "abundances \"default.abn\"" );
392  p.setline(chDUMMY);
393  ParseAbundances( p );
394 
395  // set default isotopic abundances
396  sprintf( chDUMMY, "abundances isotopes \"default-iso.abn\"" );
397  p.setline( chDUMMY );
398  ParseAbundances( p );
399 
400  input.lgVisibilityStatus = true;
401 
402  /* read until eof or blank line, then return control back to main program */
403  while (p.getline())
404  {
405  /* this checks for in-file EOF markers */
406  if ( p.last( ) )
407  break;
408 
409  /* echo the line but only if it does not contain the keyword HIDE */
410  p.echo( );
411 
412  /* check whether VARY is on line */
413  if( p.nMatch("VARY") )
414  {
415  optimize.lgVarOn = true;
416  if( optimize.nparm + 1 > LIMPAR )
417  {
418  fprintf( ioQQQ, " Too many VARY lines entered; the limit is%4ld\n",
419  LIMPAR );
421  }
422  }
423  else
424  {
425  optimize.lgVarOn = false;
426  }
427 
428  if( p.isComment() )
429  {
430  ((void)0); // do nothing for comments
431  }
432  else if( p.isVar() )
433  {
434  p.doSetVar();
435  }
436  else
437  {
438  long int i;
439  for (i=0; commands[i].name != NULL; ++i)
440  {
441  if (p.Command(commands[i].name,commands[i].action))
442  break;
443  }
444  if (commands[i].name == NULL)
445  {
446  p.CommandError(); // No command was recognized
447  }
448  }
449  }
450 
451  /*------------------------------------------------------------------- */
452  /* fall through - hit lgEOF or blank line */
453 
454  if( called.lgTalk )
455  {
456  fprintf( ioQQQ, "%23c*%81c*\n", ' ', ' ' );
457  fprintf( ioQQQ, "%23c***********************************************************************************\n\n\n\n", ' ' );
458  }
459 
460  if( input.lgDeprecatedComment ) {
461  fprintf( ioQQQ, "\n NOTE: one of these comment characters was found in the input: //, %%, *, or ;.\n" );
462  fprintf( ioQQQ, " Using these is deprecated and they will be removed in the next major release.\n" );
463  fprintf( ioQQQ, " Please convert your input deck using scripts/ccc.pl.\n\n" );
464  }
465 
466  /* set R to large value if U specified but R is not */
467  /* set radius to very large value if not already set */
468  /* >>chng 01 jul 24, from Radius == 0 to this, as per PvH comments */
469  if( !radius.lgRadiusKnown )
472 
473  /* this is an option to turn on trace printout on the nth
474  * call from the optimizer - only allow trace if
475  * this is the case and nOptimiz 1 below nTrOpt */
476  if( optimize.lgTrOpt )
477  {
478  /* nTrOpt was set with "optimize trace" command,
479  * is iteration to turn on trace */
481  {
482  trace.lgTrace = false;
483  /* following overrides turning on trace elsewhere */
484  trace.lgTrOvrd = false;
485  }
486  else
487  {
488  trace.lgTrace = true;
489  called.lgTalk = cpu.i().lgMPI_talk();
490  trace.lgTrOvrd = true;
491  fprintf( ioQQQ, " READR turns on trace from optimize option.\n" );
492  }
493  }
494 
495  /* set density from DLAW command, must be done here since it may depend on later input commands */
496  if( strcmp(dense.chDenseLaw,"DLW1") == 0 )
497  {
499 
500  if( dense.gas_phase[ipHYDROGEN] <= 0. )
501  {
502  fprintf( ioQQQ, " PROBLEM DISASTER Hydrogen density set by DLAW must be > 0.\n" );
504  }
505  }
506  else if( strcmp(dense.chDenseLaw,"DLW2") == 0 )
507  {
509 
510  if( dense.gas_phase[ipHYDROGEN] <= 0. )
511  {
512  fprintf( ioQQQ, " PROBLEM DISASTER Hydrogen density set by DLAW must be > 0.\n" );
514  }
515  }
516  else if( strcmp(dense.chDenseLaw,"DLW3") == 0 )
517  {
519 
520  if( dense.gas_phase[ipHYDROGEN] <= 0. )
521  {
522  fprintf( ioQQQ, " PROBLEM DISASTER Hydrogen density set by DLAW must be > 0.\n" );
524  }
525  }
526 
527  /* start checks on parameters set properly - this begins with same line saying start .. */
528 
529  /* lgStop_not_enough_info says that not enough info for model, so stop
530  * set true in following tests if anything missing */
531  lgStop_not_enough_info = false;
532  lgStop = false;
533 
534  /* check whether hydrogen density has been set - this value was set to 0 in zero */
535  if( dense.gas_phase[ipHYDROGEN] <= 0. )
536  {
537  fprintf( ioQQQ, " PROBLEM DISASTER Hydrogen density MUST be specified.\n" );
538  lgStop_not_enough_info = true;
539  lgStop = true;
540  /* need to set something since used below - will abort
541  * since lgStop is set */
543  }
544 
545  /* the SAVE XSPEC command cannot be combined with negative increments on the GRID command */
547  {
548  if( called.lgTalk )
549  {
550  fprintf( ioQQQ, " PROBLEM DISASTER The SAVE XSPEC command cannot be combined with negative grid increments.\n" );
551  fprintf( ioQQQ, " PROBLEM DISASTER Please check your GRID commands.\n\n\n" );
553  }
554  }
555 
556  if( geometry.covaper < 0.f || geometry.iEmissPower == 2 )
558 
559  /* mass flux for wind model - used for mass conservation */
561 
562  /* set converge criteria - limit number of iterations and zones */
564  {
566  for( long int j=0; j < iterations.iter_malloc; j++ )
567  {
569  }
570  }
571 
572  /* NSAVE is number of lines saved, =0 if no commands entered */
573  if( input.nSave < 0 )
574  {
575  fprintf( ioQQQ, " PROBLEM DISASTER No commands were entered - whats up?\n" );
577  }
578 
579  /* iterate to convergence and wind models are mutually exclusive */
580  if( wind.lgBallistic() && conv.lgAutoIt )
581  {
582  if( called.lgTalk )
583  {
584  fprintf( ioQQQ, " NOTE PROBLEM Due to the nature of the Sobolev approximation, it makes no sense to converge a windy model.\n" );
585  fprintf( ioQQQ, " NOTE Iterate to convergence is turned off\n\n\n" );
586  }
587  conv.lgAutoIt = false;
588  iterations.itermx = 0;
589  }
590 
591  /* iterate to convergence and case b do not make sense together */
592  /* WJH 22 May 2004: unless we are doing i-front dynamics (negative wind.windv0) */
593  if( opac.lgCaseB && conv.lgAutoIt && (wind.lgBallistic() || wind.lgStatic()) )
594  {
595  if( called.lgTalk )
596  {
597  fprintf( ioQQQ, " NOTE Case B is an artificial test, it makes no sense to converge this model.\n" );
598  fprintf( ioQQQ, " NOTE Iterate to convergence is turned off.\n\n\n" );
599  }
600  conv.lgAutoIt = false;
601  iterations.itermx = 0;
602  }
603 
604  /* specifying a density power and constant pressure makes no sense */
605  if( dense.DensityPower!=0. && strcmp( dense.chDenseLaw, "CPRE" )==0 )
606  {
607  if( called.lgTalk )
608  {
609  fprintf( ioQQQ, " NOTE Specifying both a density power law and constant pressure is impossible.\n" );
610  }
611  lgStop = true;
612  }
613 
614  if( !rfield.lgIonizReevaluate && strcmp( dense.chDenseLaw, "CDEN" )!=0 )
615  {
616  if( called.lgTalk )
617  {
618  fprintf( ioQQQ, " NOTE NO REEVALUATE IONIZATION can only be used with constant density.\n" );
619  fprintf( ioQQQ, " NOTE Resetting to reevaluate ionization.\n\n" );
620  }
621  rfield.lgIonizReevaluate = true;
622  }
623 
624  if( !rfield.lgOpacityReevaluate && strcmp( dense.chDenseLaw, "CDEN" )!=0 )
625  {
626  if( called.lgTalk )
627  {
628  fprintf( ioQQQ, " NOTE NO REEVALUATE OPACITY can only be used with constant density.\n" );
629  fprintf( ioQQQ, " NOTE Resetting to reevaluate opacity.\n\n" );
630  }
632  }
633 
634  /* check that a symmetry is specified if gravity from an external mass has been added */
635  if( pressure.external_mass[0].size()>0 && pressure.gravity_symmetry==-1 )
636  {
637  if( called.lgTalk )
638  {
639  fprintf( ioQQQ, " NOTE Gravity from an external mass has been added, but no symmetry (spherical/mid-plane) was specified.\n" );
640  fprintf( ioQQQ, " NOTE It will be ignored.\n\n\n" );
641  }
642  }
643 
644  /* check if the combination of stopping column density and H density are physically plausible */
645  double thickness = min( MIN3( StopCalc.tauend, StopCalc.colpls, StopCalc.colnut ),
647  if( thickness < COLUMN_INIT )
648  {
649  /* a stop column density was specified - check on physical thickness this corresponds to */
650  thickness /= (dense.gas_phase[ipHYDROGEN]*geometry.FillFac);
651  /* don't complain if outer radius set small with `stop thickness' command. */
652  if( thickness > 1e25 && iterations.StopThickness[0] > 1e25 )
653  {
654  fprintf( ioQQQ,
655  "NOTE The specified column density and hydrogen density correspond to a thickness of %.2e cm.\n",
656  thickness);
657  fprintf( ioQQQ,
658  "NOTE This seems large to me.\n");
659  fprintf(ioQQQ,"NOTE a very large radius may cause overflow.\n\n");
660  }
661  }
662 
664  {
665  /* warn if constant grain temperature but gas-grain thermal effects
666  * are still included */
667  fprintf( ioQQQ,
668  "NOTE The grain temperatures are set to a constant value with the "
669  "CONSTANT GRAIN TEMPERATURE command, but "
670  "energy exchange \n");
671  fprintf( ioQQQ,
672  "NOTE is still included. The grain-gas heating-cooling will be incorrect. "
673  "Consider turning off gas-grain collisional energy\n");
674  fprintf( ioQQQ,
675  "NOTE exchange with the NO GRAIN GAS COLLISIONAL ENERGY EXCHANGE command.\n\n\n");
676  }
677 
679  {
680  if( called.lgTalk )
681  {
682  fprintf( ioQQQ, " NOTE NO LINE TRANSER set but fine opacities still computed.\n" );
683  fprintf( ioQQQ, " NOTE Turning off fine opacities.\n\n" );
684  }
685  rfield.lgOpacityFine = false;
686  }
687 
689  {
690  if( called.lgTalk )
691  {
692  fprintf( ioQQQ, " NOTE Large H2 molecule turned on but line transfer and fine opacities are not.\n" );
693  fprintf( ioQQQ, " NOTE Turning on line transfer and fine opacities.\n\n" );
694  }
695  rfield.lgOpacityFine = true;
696  rfield.lgDoLineTrans = true;
697  }
698 
700  {
701  /* one of the input continua needs to have H-ionizing radiation
702  * blocked with extinguish command, but it was not done */
703  fprintf( ioQQQ, "\n NOTE\n"
704  " NOTE One of the incident continuum is a form used when no H-ionizing radiation is produced.\n" );
705  fprintf( ioQQQ, " NOTE You must also include the EXTINGUISH command to make sure this is done.\n" );
706  fprintf( ioQQQ, " NOTE The EXTINGUISH command was not included.\n" );
707  fprintf( ioQQQ, " NOTE YOU MAY BE MAKING A BIG MISTAKE!!\n NOTE\n\n\n\n" );
708  }
709 
710  /* if stop temp set below default then we are going into cold and possibly molecular
711  * gas - check some parameters in this case */
713  /* thermal.ConstTemp def is zero, set pos when constant temperature is set */
715  {
716  /* print warning if temperature set below default but cosmic rays not turned on
717  * do not print if molecules are off */
718  if( (hextra.cryden == 0.) && !mole_global.lgNoMole )
719  {
720  fprintf( ioQQQ, "\n NOTE\n"
721  " NOTE The simulation is going into neutral gas but cosmic rays are not included.\n" );
722  fprintf( ioQQQ, " NOTE Ion-molecule chemistry will not occur without a source of ionization.\n" );
723  fprintf( ioQQQ, " NOTE The chemistry network may collapse deep in molecular regions.\n" );
724  fprintf( ioQQQ, " NOTE Consider adding galactic background cosmic rays with the COSMIC RAYS BACKGROUND command.\n" );
725  fprintf( ioQQQ, " NOTE You may be making a BIG mistake.\n NOTE\n\n\n\n" );
726  }
727  }
728 
729  /* dense.gas_phase[ipHYDROGEN] is linear hydrogen density (cm-3) */
730  /* test for hydrogen density properly set has already been done above */
731  if( called.lgTalk && dense.gas_phase[ipHYDROGEN] < 1e-4 )
732  {
733  fprintf( ioQQQ, " NOTE Is the entered value of the hydrogen density (%.2e) reasonable?\n",
735  fprintf( ioQQQ, " NOTE It seems pretty low to me.\n\n\n" );
736  }
737  else if( called.lgTalk && dense.gas_phase[ipHYDROGEN] > 1e15 )
738  {
739  fprintf( ioQQQ, " NOTE Is this value of the hydrogen density reasonable?\n" );
740  fprintf( ioQQQ, " NOTE It seems pretty high to me.\n\n\n" );
741  }
742 
743  /* is the model going to crash because of extreme density? */
744  if( called.lgTalk && !lgStop && !lgStop_not_enough_info )
745  {
746  if( dense.gas_phase[ipHYDROGEN] < 1e-6 || dense.gas_phase[ipHYDROGEN] > 1e19 )
747  {
748  fprintf( ioQQQ, " NOTE Simulation may crash because of extreme "
749  "density. The value was %.2e\n\n" ,
751  }
752  }
753 
754  if( rfield.nShape == 0 && (p.m_nqh) == 0 )
755  {
756  // no incident radiation field, at all
757  // This is ok if temperature or heating is specified
758  if( thermal.ConstTemp <=0 && hextra.TurbHeat<=0. )
759  {
760  fprintf( ioQQQ, " PROBLEM DISASTER No incident radiation field was specified - "
761  "at least put in the CMB.\n" );
762  lgStop = true;
763  lgStop_not_enough_info = true;
764  }
765 
766  }
767  else if( rfield.nShape == 0 )
768  {
769  fprintf( ioQQQ, " PROBLEM DISASTER No incident radiation field was specified - "
770  "at least put in the CMB.\n" );
771  lgStop = true;
772  lgStop_not_enough_info = true;
773  }
774  else if( (p.m_nqh) == 0 )
775  {
776  fprintf( ioQQQ, " PROBLEM DISASTER Luminosity of continuum MUST be specified.\n" );
777  lgStop = true;
778  lgStop_not_enough_info = true;
779  }
780 
781  radius.lgPredLumin = false;
782  for( long i=0; i < rfield.nShape; i++ )
783  {
784  if( strcmp(rfield.chRSpec[i],"4 PI") == 0 )
785  {
786  radius.lgPredLumin = true;
787  break;
788  }
789  }
790 
792  {
793  fprintf(ioQQQ," PROBLEM DISASTER A continuum source was specified as a luminosity,"
794  " but the inner radius of the cloud was not set.\n");
795  lgStop = true;
796  lgStop_not_enough_info = true;
797  }
798 
799  for( long i=0; i < rfield.nShape; i++ )
800  {
801  /* do spherical dilution for TABLE READ SCALE command */
803  {
804  if( rfield.TableRadius[i] > radius.rinner )
805  fprintf(ioQQQ, " WARNING: outer radius of TABLE READ sim is larger than inner radius of this sim.\n\n");
806  rfield.totpow[i] += 2.*log10(rfield.TableRadius[i]/radius.rinner);
807  }
809  fprintf(ioQQQ, " CAUTION: no outer radius set in TABLE READ sim, but this sim uses inner radius.\n\n");
811  fprintf(ioQQQ, " CAUTION: outer radius set in TABLE READ sim, but inner radius not set in this sim.\n\n");
812  }
813 
814  if( rfield.nShape != p.m_nqh )
815  {
816  fprintf( ioQQQ, " PROBLEM DISASTER There were not the same number of continuum shapes and luminosities entered.\n" );
817  lgStop = true;
818  }
819 
820  /* we only want to do this test on the first call to the command
821  * parser - it will be called many more times but with no grid command
822  * during the actual grid calculation */
823  static bool lgFirstPass = true;
824 
825  /* the NO VARY option sets this flag, and can be used to turn off
826  * the grid command, as well as the optimizer */
827  if( optimize.lgNoVary && grid.lgGrid )
828  {
829  /* ignore grids */
830  grid.lgGrid = false;
831  optimize.nparm = 0;
832  grid.nGridCommands = 0;
833  }
834 
835  if( lgFirstPass && grid.lgGrid && (optimize.nparm!=grid.nGridCommands) )
836  {
837  /* number of grid vary options do match */
838  fprintf( ioQQQ, " PROBLEM DISASTER The GRID command was entered "
839  "but there were %li GRID commands and %li commands with a VARY option.\n" ,
841  fprintf( ioQQQ, " There must be the same number of GRIDs and VARY.\n" );
842  lgStop = true;
843  }
844  lgFirstPass = false;
845 
846  if( lgStop_not_enough_info )
847  {
848  fprintf( ioQQQ, " PROBLEM DISASTER I do not have enough information to do the simulation, I cannot go on.\n" );
849  fprintf( ioQQQ, "\n\n Sorry.\n\n\n" );
851  }
852 
853  {
854  bool lgParserTest = false;
855  if ( lgParserTest )
856  {
857  // Quit after parse phase, to allow quick verification that
858  // there are no immediate syntax handling failures.
859  fprintf(ioQQQ,"Parser phase PASSED\n");
860  exit(0);
861  }
862  }
863 
864  if( lgStop )
865  {
867  }
868 
869  /* end checks on parameters set properly - this begins with same line saying start .. */
870  return;
871 }
872 
874 {
875  DEBUG_ENTRY( "ParseAperture()" );
876  /* aperture command to simulate pencil beam or long slit */
877 
878  /* if the "BEAM" or "SLIT" option is specified then only part
879  * of the geometry is observed, and intensities
880  * should not be weighted by r^2. There are two limiting cases, SLIT in which
881  * the slit is longer than the diameter of the nebula and the contribution to the
882  * detected luminosity goes as r^1, and BEAM when the contribution is r^0,
883  * the same as plane parallel
884  *
885  * default value of geometry.iEmissPower is 2 (set in zero.c) for full geometry
886  */
887  if( p.nMatch("SLIT") )
888  {
889  /* long slit is case where slit is longer than diameter, so emissivity contributes
890  * r^1 to the observed luminosity */
891  geometry.iEmissPower = 1;
892  }
893  else if( p.nMatch("BEAM") )
894  {
895  /* pencil beam is case where we view the nebula through a narrow square
896  * centered on the central source, this gives r^0 dependence */
897  geometry.iEmissPower = 0;
898  }
899  else if( p.nMatch("SIZE") )
900  {
901  /* set the aperture size: slit width or suface area of the pencil beam */
902  /* units are arcsec for slit width, or arcsec^2 for pencil beam */
903  geometry.size = realnum(p.FFmtRead());
904  if( p.lgEOL() )
905  {
906  p.NoNumb("aperture size");
907  }
908  if( geometry.size <= 0.f )
909  {
910  fprintf( ioQQQ, " The aperture size must be positive. Sorry.\n" );
912  }
913  geometry.lgSizeSet = true;
914  }
915  else if( p.nMatch("COVE") )
916  {
917  /* set the aperture covering factor, see Hazy 1 for a discussion
918  * this is a dimensionless fraction between 0 and 1 */
920  if( p.lgEOL() )
921  {
922  p.NoNumb("aperture covering factor");
923  }
924  if( geometry.covaper <= 0.f || geometry.covaper > 1.f )
925  {
926  fprintf( ioQQQ, " The aperture covering factor must be > 0 and <= 1. Sorry.\n" );
928  }
929  }
930  else
931  {
932  fprintf( ioQQQ, " One of the keywords SLIT, BEAM, SIZE or COVEring factor must appear.\n" );
933  fprintf( ioQQQ, " Sorry.\n" );
935  }
936 }
937 
939 {
940  DEBUG_ENTRY( "ParseDatabase()" );
941 
942  string chString_quotes_original;
943  bool lgQuotesFound = true;
944  if (p.GetQuote(chString_quotes_original))
945  lgQuotesFound = false;
946 
947  /* accept both forms of feii */
948  if( p.nMatch("FEII") || p.nMatch("FE II") )
949  {
950  fprintf( ioQQQ, " Warning: The 'atom feii' command is obsolete. "
951  " Instead, please use 'species \"Fe+\" levels=all'.\n Sorry.\n\n" );
952  cdEXIT( EXIT_FAILURE );
953  }
954 
955  else if( p.nMatch("H-LI") )
956  {
957  /* parse the atom h-like command */
959  }
960 
961  else if( p.nMatch("HE-L") )
962  {
963  /* parse the atom he-like command */
965  }
966 
967  else if( p.nMatch("ROTO") || p.nMatch(" CO ") )
968  {
969  fprintf(ioQQQ," The old CO models no longer exist, and this command is no longer supported.\n" );
970  fprintf( ioQQQ, " Sorry.\n" );
972  }
973 
974  else if( p.nMatch(" H2 ") )
975  {
976  ParseDatabaseH2(p);
977  }
978 
979  /* enable models from CHIANTI database */
980  else if (p.nMatch("CHIANTI"))
981  {
982  // option to specify different CloudyChianti.ini file, was initialized
983  // with string CloudyChianti.ini
984  if (lgQuotesFound == true)
985  strcpy(atmdat.chCloudyChiantiFile, chString_quotes_original.c_str());
986 
987  atmdat.lgChiantiOn = true;
988  if (p.nMatch(" OFF"))
989  {
990  atmdat.lgChiantiOn = false;
991  atmdat.lgChiantiHybrid = false;
992  }
993 
994  // hybrid, chianti with OP for higher energy lines
995  // Turn off hybrid, use Chianti only
996  if (p.nMatch("NO HYBR"))
997  atmdat.lgChiantiHybrid = false;
998 
999  // Print which species are being used in output and # of levels
1000  if (p.nMatch("PRINT"))
1001  atmdat.lgChiantiPrint = true;
1002 
1003  // Use Experimental energies exclusively. Default use experimental.
1004  if (p.nMatch("THEOR"))
1005  atmdat.lgChiantiExp = false;
1006 
1007  // Input the maximum number of Chianti levels to use
1008  if (p.nMatch("LEVEL"))
1009  {
1010  if (p.nMatch(" MAX"))
1011  {
1012  atmdat.nChiantiMaxLevels = LONG_MAX;
1013  atmdat.nChiantiMaxLevelsFe = LONG_MAX;
1014  }
1015  else
1016  {
1017  atmdat.nChiantiMaxLevelsFe = (long)p.FFmtRead();
1018  atmdat.nChiantiMaxLevels = (long)p.FFmtRead();
1019  if (p.lgEOL())
1020  {
1021  p.NoNumb("two numbers, the maximum number of levels in Fe, and in other elements, or the keyword MAX,");
1022  }
1024  {
1025  fprintf(ioQQQ,
1026  " \nPROBLEM The maximum number of chianti levels should be two or greater.\n");
1027  fprintf(ioQQQ, " To turn off the Chianti data use \"atom Chianti off\" instead.\n");
1028  fprintf(ioQQQ, " See Hazy 1 for details.\n");
1029  cdEXIT( EXIT_FAILURE );
1030  }
1031  }
1032  atmdat.lgChiantiLevelsSet = true;
1033  }
1034  }
1035 
1036  /* enable models from STOUT database */
1037  else if (p.nMatch("STOUT"))
1038  {
1039 
1040  // option to specify different Stout.ini file, was initialized
1041  // with string Stout.ini
1042  if (lgQuotesFound == true)
1043  strcpy(atmdat.chStoutFile, chString_quotes_original.c_str());
1044 
1045  atmdat.lgStoutOn = true;
1046  // Print which species are being used in output and # of levels
1047  if (p.nMatch("PRINT"))
1048  atmdat.lgStoutPrint = true;
1049 
1050  if (p.nMatch(" OFF"))
1051  {
1052  atmdat.lgStoutOn = false;
1053  atmdat.lgStoutHybrid = false;
1054  }
1055 
1056  // hybrid, Stout with OP for higher energy lines
1057  // Turn off hybrid, use Stout only
1058  if (p.nMatch("NO HYBR"))
1059  atmdat.lgStoutHybrid = false;
1060 
1061  // Input the maximum number of Stout levels to use
1062  if (p.nMatch("LEVEL"))
1063  {
1064  if (p.nMatch(" MAX"))
1065  {
1066  atmdat.nStoutMaxLevels = LONG_MAX;
1067  atmdat.nStoutMaxLevelsFe = LONG_MAX;
1068  }
1069  else
1070  {
1071  atmdat.nStoutMaxLevelsFe = (long)p.FFmtRead();
1072  atmdat.nStoutMaxLevels = (long)p.FFmtRead();
1073  if (p.lgEOL())
1074  {
1075  p.NoNumb("two numbers, the maximum number of levels in Fe, and in other elements, or the keyword MAX,");
1076  }
1078  {
1079  fprintf(ioQQQ,
1080  " \nPROBLEM The maximum number of stout levels should be two or greater.\n");
1081  fprintf(ioQQQ, " To turn off the Stout data use \"atom Stout off\" instead.\n");
1082  fprintf(ioQQQ, " See Hazy 1 for details.\n");
1083  cdEXIT( EXIT_FAILURE );
1084  }
1085  }
1086  atmdat.lgStoutLevelsSet = true;
1087  }
1088  }
1089 
1090  /* enable models from LAMDA (Leiden Atomic and Molecular Database) */
1091  else if (p.nMatch("LAMDA"))
1092  {
1093  // option to specify different Lamda.ini file, was initialized
1094  // with string Lamda.ini
1095  if (lgQuotesFound == true)
1096  strcpy(atmdat.chLamdaFile, chString_quotes_original.c_str());
1097 
1098  atmdat.lgLamdaOn = true;
1099  // Print which species are being used in output and # of levels
1100  if (p.nMatch("PRINT"))
1101  atmdat.lgLamdaPrint = true;
1102 
1103  if(p.nMatch(" OFF"))
1104  {
1105  atmdat.lgLamdaOn = false;
1106  }
1107 
1108  // Input the maximum number of Lamda levels to use
1109  if (p.nMatch("LEVEL"))
1110  {
1111  if(p.nMatch(" MAX"))
1112  {
1113  atmdat.nLamdaMaxLevels = LONG_MAX;
1114  }
1115  else
1116  {
1117  atmdat.nLamdaMaxLevels = (long)p.FFmtRead();
1118  if( p.lgEOL())
1119  {
1120  p.NoNumb("the maximum number of levels,");
1121  }
1122  if( atmdat.nLamdaMaxLevels < 2 )
1123  {
1124  fprintf(ioQQQ,
1125  " \nPROBLEM The maximum number of Lamda levels should be two or greater.\n");
1126  fprintf(ioQQQ, " To turn off the Lamda data use \"atom lamda off\" instead.\n");
1127  fprintf(ioQQQ, " See Hazy 1 for details.\n");
1128  cdEXIT( EXIT_FAILURE );
1129  }
1130  }
1131  atmdat.lgLamdaLevelsSet = true;
1132  }
1133  }
1134 
1135  /* no single species, but print summaries of all */
1136  else if (p.nMatch("PRINT"))
1137  {
1138  atmdat.lgChiantiPrint = true;
1139  atmdat.lgLamdaPrint = true;
1140  atmdat.lgStoutPrint = true;
1142  }
1143 
1144  else
1145  {
1146  fprintf( ioQQQ, " I could not recognize a keyword on this species command.\n");
1147  fprintf( ioQQQ, " The available keys are FeII, H-Like, He-like, H2, Chianti, Lamda, and Stout.\n");
1148  fprintf( ioQQQ, " Sorry.\n" );
1150  }
1151 }
1153 {
1154  DEBUG_ENTRY( "ParseBremsstrahlung()" );
1155  /* bremsstrahlung continuum from central object */
1156  strcpy( rfield.chSpType[rfield.nShape], "BREMS" );
1158  (realnum)p.FFmtRead();
1159  if( p.lgEOL() )
1160  {
1161  p.NoNumb("temperature");
1162  }
1163 
1164  /* temperature is interpreted as log if <= 10, or if keyword is present */
1165  if( rfield.slope[rfield.nShape] <= 10. || p.nMatch(" LOG") )
1166  {
1169  }
1170  rfield.cutoff[rfield.nShape][0] = 0.;
1171 
1172  /* option for vary keyword */
1173  if( optimize.lgVarOn )
1174  {
1175  /* only one parameter */
1177  strcpy( optimize.chVarFmt[optimize.nparm], "BREMS, T=%f LOG" );
1178  /* pointer to where to write */
1180  /* log of temp will be pointer */
1182  optimize.vincr[optimize.nparm] = 0.5;
1183  ++optimize.nparm;
1184  }
1185  ++rfield.nShape;
1186  if( rfield.nShape >= LIMSPC )
1187  {
1188  /* too many continua were entered */
1189  fprintf( ioQQQ, " Too many continua entered; increase LIMSPC\n" );
1191  }
1192 }
1194 {
1195  /* add "extra" cooling, power law temp dependence */
1196  thermal.lgCExtraOn = true;
1198  if( p.lgEOL() )
1199  {
1200  p.NoNumb("extra cooling");
1201  }
1202  thermal.cextpw = (realnum)p.FFmtRead();
1203 }
1205 {
1208 
1209  /* >>chng 06 mar 22, add time option to vary only some continua with time */
1210  if( p.nMatch( "TIME" ) )
1211  rfield.lgTimeVary[(p.m_nqh)] = true;
1212 
1214 
1215  /* option to also set the hydrogen density at given redshift. */
1216  if( p.nMatch("DENS") && !p.lgEOL() )
1217  {
1218  char chStuff[INPUT_LINE_LENGTH];
1219 
1220  /* hydrogen density */
1221  sprintf( chStuff , "HDEN %.2e LINEAR", GetDensity( cosmology.redshift_start ) );
1222  p.setline(chStuff);
1223  p.set_point(4);
1224  ParseHDEN(p);
1225  }
1226 }
1228 {
1229  DEBUG_ENTRY( "ParseCosm()" );
1230  fprintf(ioQQQ,
1231  "This command is now ambiguous -- please specify either COSMIC RAYS or COSMOLOGY.\nSorry.\n");
1233 }
1235 {
1236  DEBUG_ENTRY( "ParseCovering()" );
1237  /* covering factor for gas */
1238  /* The geometric covering factor accounts for how much of 4\pi is
1239  * covered by gas, and so linearly multiplies the predicted intensities */
1241 
1242  if( p.lgEOL() )
1243  {
1244  p.NoNumb("covering factor");
1245  }
1246 
1247  /* if negative then log, convert to linear */
1248  if( geometry.covgeo <= 0. )
1249  {
1251  }
1252 
1253  if( geometry.covgeo > 1. )
1254  {
1255  fprintf( ioQQQ, " A covering factor greater than 1 makes no physical sense. Sorry.\n" );
1257  }
1258 
1259  /* radiative transfer covering factor will be equal to the geometric one */
1261 }
1263 {
1264  /* height of cylinder in cm */
1265  radius.lgCylnOn = true;
1267  if( p.lgEOL() )
1268  {
1269  p.NoNumb("height");
1270  }
1271 }
1273 {
1274  DEBUG_ENTRY( "ParseDarkMatter()" );
1275 
1276  if( p.nMatch(" NFW") )
1277  {
1278  /* specify an NFW profile */
1279  /* >>refer dark matter Navarro, Frenk, & White, 1996, ApJ, 462, 563 */
1280 
1281  dark.r_200 = (realnum)p.getNumberCheckAlwaysLog("NFW r_200");
1282  /* Let r_s have a default corresponding to c_200 = 10. */
1283  dark.r_s = (realnum)p.getNumberDefaultAlwaysLog("NFW r_s", log10(dark.r_200)-1. );
1284  dark.lgNFW_Set = true;
1285 
1286  /* option for vary keyword */
1287  if( optimize.lgVarOn )
1288  {
1289  /* only one parameter */
1291  strcpy( optimize.chVarFmt[optimize.nparm], "DARK NFW %f" );
1292  /* pointer to where to write */
1294  /* log of temp will be pointer */
1295  optimize.vparm[0][optimize.nparm] = (realnum)log10(dark.r_200);
1296  optimize.vincr[optimize.nparm] = 0.5;
1297  ++optimize.nparm;
1298  }
1299  }
1300  else
1301  {
1302  fprintf( ioQQQ, " Did not recognize a valid option for this DARK command.\nSorry.\n\n" );
1304  }
1305 
1306  return;
1307 }
1309 {
1310  DEBUG_ENTRY( "ParseDielectronic()" );
1311  /* >>chng 05 dec 21, change dielectronic command to set dielectronic recombination command */
1312  fprintf( ioQQQ, " The DIELectronic command has been replaced with the SET DIELectronic recombination command.\n" );
1313  fprintf( ioQQQ, " Please have a look at Hazy.\n Sorry.\n\n" );
1315 }
1317 {
1318  DEBUG_ENTRY( "ParseDiffuse()" );
1319  /* set method of transferring diffuse fields,
1320  * default is outward only, cdDffTrns label "OU2", set in zero.c */
1321  if( p.nMatch(" OTS") )
1322  {
1323  if( p.nMatch("SIMP") )
1324  {
1325  /* this is simple ots, which never allows photons to escape */
1326  strcpy( rfield.chDffTrns, "OSS" );
1327  }
1328  else
1329  {
1330  /* this is regular ots, which allows photons to escape */
1331  strcpy( rfield.chDffTrns, "OTS" );
1332  }
1333  rfield.lgOutOnly = false;
1334  }
1335  else if( p.nMatch(" OUT") )
1336  {
1337  rfield.lgOutOnly = true;
1338  long int j = (long int)p.FFmtRead();
1339  if( p.lgEOL() )
1340  {
1341  /* this is the default set in zero */
1342  strcpy( rfield.chDffTrns, "OU2" );
1343  }
1344  else
1345  {
1346  if( j > 0 && j < 10 )
1347  {
1348  sprintf( rfield.chDffTrns, "OU%1ld", j );
1349  }
1350  else
1351  {
1352  fprintf( ioQQQ, " must be between 1 and 9 \n" );
1354  }
1355  }
1356  }
1357 
1358  else
1359  {
1360  fprintf( ioQQQ, " There should have been OUTward or OTS on this line. Sorry.\n" );
1362  }
1363 }
1365 {
1366  /* distance to the object */
1367  radius.distance = p.FFmtRead();
1368  if( p.lgEOL() )
1369  {
1370  p.NoNumb("distance");
1371  }
1372  /* default is for quantity to be log of distance, linear keyword
1373  * overrides this - is LINEAR is not on line then exp */
1374  if( !p.nMatch("LINE" ) )
1375  {
1377  }
1378  /* default is radius in cm - if parsecs appears then must
1379  * convert to cm */
1380  if( p.nMatch("PARS") )
1381  {
1382  radius.distance *= PARSEC;
1383  }
1384 
1385  /* vary option */
1386  if( optimize.lgVarOn )
1387  {
1388  strcpy( optimize.chVarFmt[optimize.nparm], "DISTANCE %f LOG" );
1391  optimize.vincr[optimize.nparm] = 0.3f;
1393  ++optimize.nparm;
1394  }
1395 }
1397 {
1398  rt.DoubleTau = 2.;
1399 }
1401 {
1403  if( p.lgEOL() )
1404  {
1405  p.NoNumb("electron density");
1406  }
1407  /* warn that this model is meaningless */
1408  phycon.lgPhysOK = false;
1409 
1410  /* vary option */
1411  if( optimize.lgVarOn )
1412  {
1413  strcpy( optimize.chVarFmt[optimize.nparm], "EDEN %f LOG " );
1414  /* pointer to where to write */
1417  optimize.vincr[optimize.nparm] = 0.1f;
1419  optimize.varang[optimize.nparm][0] = -FLT_MAX;
1420  optimize.varang[optimize.nparm][1] = 20.f;
1421  ++optimize.nparm;
1422  }
1423 
1424 }
1426 {
1427  DEBUG_ENTRY( "ParseEnergy()" );
1428  /* energy density (degrees K) of this continuum source */
1429  if( (p.m_nqh) >= LIMSPC )
1430  {
1431  /* too many continua were entered */
1432  fprintf( ioQQQ, " Too many continua entered; increase LIMSPC\n" );
1434  }
1435  /* silly, but calms down the lint */
1436  ASSERT( (p.m_nqh) < LIMSPC );
1437  strcpy( rfield.chRSpec[(p.m_nqh)], "SQCM" );
1438  realnum teset = (realnum)p.FFmtRead();
1439  if( p.lgEOL() )
1440  {
1441  p.NoNumb("energy density");
1442  }
1443 
1444  /* convert temp to log, recognize linear option */
1445  if( !p.nMatch(" LOG") && (p.nMatch("LINE") || teset > 10.) )
1446  {
1447  /* option for linear temperature, must store log */
1448  teset = (realnum)log10(teset);
1449  }
1450 
1451  if( teset > 5. )
1452  {
1453  fprintf( ioQQQ, " This intensity may be too large. The code may crash due to overflow. Was log intended?\n" );
1454  }
1455 
1456  /* teset is now log of temp, now get log of total luminosity */
1457  strcpy( rfield.chSpNorm[(p.m_nqh)], "LUMI" );
1458 
1459  /* full range of continuum will be used */
1460  rfield.range[(p.m_nqh)][0] = rfield.emm();
1461  rfield.range[(p.m_nqh)][1] = rfield.egamry();
1462  rfield.totpow[(p.m_nqh)] = (4.*(teset) - 4.2464476 + 0.60206);
1463 
1464  /* >>chng 06 mar 22, add time option to vary only some continua with time */
1465  if( p.nMatch( "TIME" ) )
1466  rfield.lgTimeVary[(p.m_nqh)] = true;
1467 
1468  /* vary option */
1469  if( optimize.lgVarOn )
1470  {
1471  strcpy( optimize.chVarFmt[optimize.nparm], "ENERGY DENSITY %f LOG" );
1472  if( rfield.lgTimeVary[(p.m_nqh)] )
1473  strcat( optimize.chVarFmt[optimize.nparm], " TIME" );
1474  /* pointer to where to write */
1476  optimize.vparm[0][optimize.nparm] = teset;
1477  optimize.vincr[optimize.nparm] = 0.1f;
1479  ++optimize.nparm;
1480  }
1481 
1482  /* finally increment number of continua */
1483  ++(p.m_nqh);
1484 }
1486 {
1487  /* save previous value */
1488  long int j = conv.LimFail;
1489  conv.LimFail = (long int)p.FFmtRead();
1490  if( p.lgEOL() )
1491  {
1492  p.NoNumb("limit");
1493  }
1494 
1495  /* >>chng 01 mar 14, switch logic on maps */
1496  /* ' map' flag, make map when failure, default is no map,
1497  * second check is so that no map does not trigger a map */
1498  if( p.nMatch(" MAP") && !p.nMatch(" NO ") )
1499  {
1500  conv.lgMap = true;
1501  }
1502 
1503  /* complain if failures was increased above default */
1504  if( conv.LimFail > j )
1505  {
1506  fprintf( ioQQQ, " This command should not be necessary.\n" );
1507  fprintf( ioQQQ, " Please show this input stream to Gary Ferland if this command is really needed for this simulation.\n" );
1508  }
1509 }
1511 {
1512  /* filling factor, power law exponent (default=1., 0.) */
1513  realnum a = (realnum)p.FFmtRead();
1514  if( p.lgEOL() )
1515  {
1516  p.NoNumb("filling factor");
1517  }
1518 
1519  if( a <= 0. || p.nMatch(" LOG") )
1520  {
1521  /* number less than or equal to 0, must have been entered as log */
1522  geometry.FillFac = exp10(a);
1523  }
1524  else
1525  {
1526  /* number greater than zero, must have been the real thing */
1527  geometry.FillFac = a;
1528  }
1529  if( geometry.FillFac > 1.0 )
1530  {
1531  if( called.lgTalk )
1532  fprintf( ioQQQ, " Filling factor > 1, reset to 1\n" );
1533  geometry.FillFac = 1.;
1534  }
1536 
1537  /* option to have power law dependence, filpow set to 0 in zerologic */
1539 
1540  /* vary option */
1541  if( optimize.lgVarOn )
1542  {
1543  strcpy( optimize.chVarFmt[optimize.nparm], "FILLING FACTOR= %f LOG power= %f" );
1544 
1545  /* pointer to where to write */
1548  optimize.vincr[optimize.nparm] = 0.5f;
1549 
1550  /* power law dependence here, but cannot be varied */
1553 
1554  /* do not allow filling factor to go positive */
1555  optimize.varang[optimize.nparm][0] = -FLT_MAX;
1556  optimize.varang[optimize.nparm][1] = 0.f;
1557  ++optimize.nparm;
1558  }
1559 }
1561 {
1562  bool lgNu2 = false;
1563  /* in reads2 */
1564  ParseF_nu(p,"SQCM",lgNu2);
1565 }
1567 {
1569  if( p.lgEOL() )
1570  {
1571  p.NoNumb("temperature");
1572  }
1573 
1574  if( p.nMatch(" LOG") || (thermal.ConstTemp <= 10. &&
1575  !p.nMatch("LINE")) )
1576  {
1578  }
1579 
1580  /* low energy bounds of continuum array do not permit too-low a temp */
1581  if( thermal.ConstTemp < 3. )
1582  {
1583  fprintf( ioQQQ, " TE reset to 3K: entered number too small.\n" );
1584  thermal.ConstTemp = 3.;
1585  }
1586 }
1588 {
1589  /* enter a fudge factor */
1590  fudgec.nfudge = 0;
1591  for( long int j=0; j < NFUDGC; j++ )
1592  {
1593  fudgec.fudgea[j] = (realnum)p.FFmtRead();
1594  /* fudgec.nfudge is number of factors on the line */
1595  if( !p.lgEOL() )
1596  fudgec.nfudge = j+1;
1597  }
1598  if( fudgec.nfudge == 0 )
1599  p.NoNumb("fudge factor");
1600 
1601  /* vary option */
1602  if( optimize.lgVarOn )
1603  {
1604  /* no luminosity options on vary */
1606  strcpy( optimize.chVarFmt[optimize.nparm], "FUDGE= %f" );
1607  /* enter remaining parameters as constants */
1608  char chHold[1000];
1609  for( long int j=1; j<fudgec.nfudge; ++j )
1610  {
1611  sprintf( chHold , " %f" , fudgec.fudgea[j] );
1612  strcat( optimize.chVarFmt[optimize.nparm] , chHold );
1613  }
1615 
1616  /* pointer to where to write */
1618  /* fudge factor stored here */
1620  /* the increment in the first steps away from the original value */
1621  optimize.vincr[optimize.nparm] = abs(0.2f*fudgec.fudgea[0]);
1622  /* we have no clue what to use when initial estimate is 0... */
1623  if( optimize.vincr[optimize.nparm] == 0.f )
1624  optimize.vincr[optimize.nparm] = 1.f;
1625  ++optimize.nparm;
1626  }
1627 }
1629 {
1630  DEBUG_ENTRY( "ParsePGrains()" );
1631  fprintf(ioQQQ," Sorry, this command is obsolete, you can now use the normal GRAINS command.\n");
1633 }
1635 {
1636  DEBUG_ENTRY( "ParseGravity()" );
1637 
1638  if( p.nMatch("EXTE") )
1639  {
1640  /* external mass: M/Msun or Sigma/(Msun/pc^2) depending on symmetry */
1641  /* if no number is read, FFmtRead returns 0 */
1642  /* default is linear, unless "log" keyword is present */
1643  double M_i = p.FFmtRead();
1644 
1645  if( !p.lgEOL() && p.nMatch("LOG") )
1646  M_i = exp10( M_i );
1647  pressure.external_mass[0].push_back( M_i );
1648 
1649  /* extent of the mass distribution, in pc */
1650  /* default is linear, unless "log" keyword is present */
1651  double x_i = p.FFmtRead();
1652 
1653  if( !p.lgEOL() && p.nMatch("LOG") )
1654  x_i = exp10( x_i );
1655  pressure.external_mass[1].push_back( x_i * PARSEC );
1656 
1657  /* exponential index of the mass distribution */
1658  pressure.external_mass[2].push_back( p.FFmtRead() );
1659  }
1660  else
1661  {
1662  /* choose spherical or mid-plane symmetry for the gas mass distribution */
1663  if( p.nMatch("SPHE") )
1664  {
1666  }
1667  else if( p.nMatch("PLAN") )
1668  {
1670  }
1671  else
1672  {
1673  fprintf( ioQQQ, " The symmetry of the gravitational mass must be specified explicitly. Sorry.\n" );
1675  }
1676 
1677  /* external mass, proportional to the gas density */
1678  /* default is linear, unless "log" keyword is present */
1680 
1681  if( p.lgEOL() )
1683  else if( p.nMatch("LOG") )
1684  {
1686  }
1687  }
1688 }
1690 {
1691  DEBUG_ENTRY( "ParseHeLike()" );
1692  fprintf(ioQQQ,"Sorry, this command is replaced with SPECIES HE-LIKE\n");
1694 }
1696 {
1697  DEBUG_ENTRY( "ParseHelp()" );
1698  p.help(ioQQQ);
1699 }
1701 {
1702  DEBUG_ENTRY( "ParseHExtra()" );
1704  if( p.lgEOL() )
1705  p.NoNumb("extra heating first parameter" );
1706 
1707  /* save initial value in case reset with time dependent command */
1709 
1710  /* remember which type of scale dependence we will use */
1711  const char *chHextraScale;
1712  /* these are optional numbers on depth or density option */
1713  realnum par1=0. , par2=0.;
1714  long int nPar;
1715  if( p.nMatch("DEPT") )
1716  {
1717  /* depend on depth */
1718  hextra.lgHextraDepth = true;
1719  chHextraScale = "DEPTH";
1720  /* optional scale radius */
1721  realnum a = (realnum)p.FFmtRead();
1722  if( p.lgEOL() )
1723  p.NoNumb("depth");
1724  else
1725  hextra.turrad = exp10(a);
1726 
1727  /* depth from shielded face, to mimic illumination from both sides
1728  * may not be specified */
1729  realnum b = (realnum)p.FFmtRead();
1730  if( p.lgEOL() )
1731  {
1732  hextra.turback = 0.;
1733  nPar = 2;
1734  }
1735  else
1736  {
1737  hextra.turback = exp10(b);
1738  nPar = 3;
1739  }
1740  par1 = a;
1741  par2 = b;
1742  }
1743  else if( p.nMatch("DENS") )
1744  {
1745  /* depends on density */
1746  chHextraScale = "DENSITY";
1747  hextra.lgHextraDensity = true;
1748 
1749  /* the log of the density */
1750  realnum a = (realnum)p.FFmtRead();
1751  /* if number not entered then unity is used, heating
1752  * proportional to density */
1754  par1 = a;
1755  par2 = 0;
1756  nPar = 2;
1757  }
1758  else if( p.nMatch("SS") )
1759  {
1760  /* alpha disk model, specify alpha, mass of black hole, and distance */
1761  chHextraScale = "SS";
1762  hextra.lgHextraSS = true;
1763 
1764  /* the parameter alpha of alpha-disk model */
1766 
1767  /* mass in solar masses of center black hole */
1768  realnum a = (realnum)p.FFmtRead();
1769  if( p.lgEOL() )
1770  p.NoNumb("hextraSS Mass");
1771  hextra.HextraSS_M = a*SOLAR_MASS;
1772 
1773  /* radius (cm) from center */
1774  realnum b = (realnum)p.FFmtRead();
1775  if( p.lgEOL() )
1776  p.NoNumb("hextraSS radius");
1777  hextra.HextraSSradius = b;
1778  par1 = a;
1779  par2 = 0;
1780  nPar = 2;
1781  }
1782  else
1783  {
1784  /* constant heating */
1785  chHextraScale = "";
1786  nPar = 1;
1787  }
1788 
1789  /* option to have heating vary with time */
1790  if( p.nMatch( "TIME" ) )
1791  hextra.lgTurbHeatVaryTime = true;
1792 
1793  /* vary option */
1794  if( optimize.lgVarOn )
1795  {
1796  if( hextra.lgHextraSS )
1797  {
1798  fprintf(ioQQQ,"Sorry, HEXTRA SS command does not now support vary option.\n");
1800  }
1801  /* 1 to 3 options on hextra vary */
1802  optimize.nvarxt[optimize.nparm] = nPar;
1803  optimize.vparm[0][optimize.nparm] = log10(hextra.TurbHeat);
1804  optimize.vparm[1][optimize.nparm] = par1;
1805  optimize.vparm[2][optimize.nparm] = par2;
1806 
1807  /* the keyword LOG is not used above, but is checked elsewhere */
1808  strcpy( optimize.chVarFmt[optimize.nparm], "HEXTra %f LOG " );
1809  strcat( optimize.chVarFmt[optimize.nparm], chHextraScale );
1810  while( nPar > 1 )
1811  {
1812  /* add extra spots to write parameters */
1813  --nPar;
1814  strcat( optimize.chVarFmt[optimize.nparm], " %f" );
1815  }
1817  strcat( optimize.chVarFmt[optimize.nparm], " TIME" );
1818  /* pointer to where to write */
1820  /* the increment in the first steps away from the original value */
1821  optimize.vincr[optimize.nparm] = 0.1f;
1822  ++optimize.nparm;
1823  }
1824 }
1825 
1827 {
1828  thermal.lgTeHigh = true;
1829 }
1831 {
1832  DEBUG_ENTRY( "ParseHydrogen()" );
1833  fprintf(ioQQQ," Sorry, this command has been replaced with the SPECIES H-LIKE command.\n");
1835 }
1837 {
1838  // parsing the init file was already done in cdRead()
1839  // so there is nothing left to do here...
1840  (void)0;
1841 }
1843 {
1844  DEBUG_ENTRY( "ParseIntensity()" );
1845  /* intensity of this continuum source */
1846  if( (p.m_nqh) >= LIMSPC )
1847  {
1848  /* too many continua were entered */
1849  fprintf( ioQQQ, " Too many continua entered; increase LIMSPC\n" );
1851  }
1852 
1853  /* silly, but calms down the lint */
1854  ASSERT( (p.m_nqh) < LIMSPC );
1855  strcpy( rfield.chRSpec[(p.m_nqh)], "SQCM" );
1856  rfield.totpow[(p.m_nqh)] = p.FFmtRead();
1857  if( p.lgEOL() )
1858  p.NoNumb("intensity");
1859 
1860  if( p.nMatch("LINE") )
1861  {
1862  /* silly, but calms down the lint */
1863  ASSERT( (p.m_nqh) < LIMSPC );
1864  /* option for linear input parameter */
1865  rfield.totpow[(p.m_nqh)] = log10(rfield.totpow[(p.m_nqh)]);
1866  }
1867  strcpy( rfield.chSpNorm[(p.m_nqh)], "LUMI" );
1868  /* ParseRangeOption in readsun */
1869  ParseRangeOption(p);
1870 
1871  /* >>chng 06 mar 22, add time option to vary only some continua with time */
1872  if( p.nMatch( "TIME" ) )
1873  rfield.lgTimeVary[(p.m_nqh)] = true;
1874 
1875  /* vary option */
1876  if( optimize.lgVarOn )
1877  {
1878  /* create the format we will use to write the command */
1879  strcpy( optimize.chVarFmt[optimize.nparm], "INTENSITY %f LOG range %f %f" );
1880  if( rfield.lgTimeVary[(p.m_nqh)] )
1881  strcat( optimize.chVarFmt[optimize.nparm], " TIME" );
1882  /* array index for this command */
1885  optimize.vincr[optimize.nparm] = 0.5;
1886  /* range option, but cannot be varied */
1887  optimize.vparm[1][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][0]);
1888  optimize.vparm[2][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][1]);
1890  ++optimize.nparm;
1891  }
1892  ++(p.m_nqh);
1893 }
1895 {
1896  /* iterate to get optical depths and diffuse field properly */
1897  iterations.itermx = (long int)p.FFmtRead() - 1;
1899  /* >>chng 06 mar 19, malloc itrdim arrays
1900  * iterations.iter_malloc is -1 before space allocated and set to
1901  * itermx if not previously set */
1903  {
1904  long int iter_malloc_save = iterations.iter_malloc;
1905  /* add one for mismatch between array dim and iterations defn,
1906  * another few for safety */
1908  iterations.IterPrnt.resize((size_t)iterations.iter_malloc);
1909  iterations.nend.resize((size_t)iterations.iter_malloc);
1911  iterations.StopRadius.resize((size_t)iterations.iter_malloc);
1912  /* >>chng 06 jun 07, need to set IterPrnt, thickness, and nend */
1913  for(long int j=iter_malloc_save; j<iterations.iter_malloc; ++j )
1914  {
1915  iterations.IterPrnt[j] = iterations.IterPrnt[iter_malloc_save-1];
1916  iterations.nend[j] = iterations.nend[iter_malloc_save-1];
1917  iterations.StopThickness[j] = iterations.StopThickness[iter_malloc_save-1];
1918  iterations.StopRadius[j] = iterations.StopRadius[iter_malloc_save-1];
1919  }
1920  }
1921 
1922  if( p.nMatch("CONV") )
1923  {
1924  /* option to keep iterating until it converges, checks on convergence
1925  * are done in update, and checked again in prtcomment */
1926  conv.lgAutoIt = true;
1927  /* above would have been legitimate setting of ITERMX, set to default 10 */
1928  if( p.lgEOL() )
1929  {
1930  iterations.itermx = 10 - 1;
1931  }
1932  realnum a = (realnum)p.FFmtRead();
1933  /* change convergence criteria, preset in zero */
1934  if( !p.lgEOL() )
1935  {
1936  conv.autocv = a;
1937  }
1938  if( p.nMatch("ALL") )
1939  {
1940  conv.lgAllTransitions = true;
1941  }
1942  }
1943 }
1945 {
1946  /* this is the specific luminosity at nu
1947  * following says n(nu) not nuF(nu) */
1948  bool lgNu2 = false;
1949  /* in reads2 */
1950  ParseF_nu(p,"4 PI",lgNu2);
1951 }
1953 {
1954  DEBUG_ENTRY( "ParseLaser()" );
1955  /* mostly one frequency (a laser) to check gamma's,*/
1956  /* say the continuum type */
1957  strcpy( rfield.chSpType[rfield.nShape], "LASER" );
1958 
1959  /* scan off the laser's energy ion Rydbergs */
1961 
1962  /* negative energies are logs */
1963  if( rfield.slope[rfield.nShape] <= 0. )
1964  {
1967  }
1968  if( p.lgEOL() )
1969  {
1970  p.NoNumb("energy");
1971  }
1972 
1973  /* next number is optional relative width of laser */
1974  rfield.cutoff[rfield.nShape][0] = p.FFmtRead();
1975  if( p.lgEOL() )
1976  {
1977  /* default width is 0.05 */
1978  rfield.cutoff[rfield.nShape][0] = 0.05;
1979  }
1980 
1981  /* increment counter making sure we still have space enough */
1982  ++rfield.nShape;
1983  if( rfield.nShape >= LIMSPC )
1984  {
1985  /* too many continua were entered */
1986  fprintf( ioQQQ, " Too many continua entered; increase LIMSPC\n" );
1988  }
1989 }
1991 {
1992  DEBUG_ENTRY( "ParseLuminosity()" );
1993  /* luminosity of this continuum source */
1994  if( (p.m_nqh) >= LIMSPC )
1995  {
1996  /* too many continua were entered */
1997  fprintf( ioQQQ, " Too many continua entered; increase LIMSPC\n" );
1999  }
2000 
2001  strcpy( rfield.chRSpec[(p.m_nqh)], "4 PI" );
2002  rfield.totpow[(p.m_nqh)] = p.FFmtRead();
2003  if( p.lgEOL() )
2004  {
2005  p.NoNumb("luminosity");
2006  }
2007  if( p.nMatch("LINE") )
2008  {
2009  /* option for linear input parameter */
2010  rfield.totpow[(p.m_nqh)] = log10(rfield.totpow[(p.m_nqh)]);
2011  }
2012  if( rfield.totpow[(p.m_nqh)] > 200. )
2013  {
2014  fprintf( ioQQQ, " The log of the luminosity is very large: %g\n", rfield.totpow[(p.m_nqh)] );
2015  fprintf( ioQQQ, " Did you omit the keyword LINEAR?\n" );
2017  }
2018 
2019  strcpy( rfield.chSpNorm[(p.m_nqh)], "LUMI" );
2020 
2021  /* the solar option - number is log total solar luminosity */
2022  if( p.nMatch("SOLA") )
2023  {
2024  /* option to use log of total luminosity in solar units */
2025  rfield.range[(p.m_nqh)][0] = rfield.emm();
2026  rfield.range[(p.m_nqh)][1] = rfield.egamry();
2027  /* log of solar luminosity in erg s-1 */
2028  rfield.totpow[(p.m_nqh)] += 33.5827f;
2029  }
2030  else
2031  {
2032  /* check if range is present and parse it if it is - ParseRangeOption in readsun
2033  * this includes TOTAL keyword for total luminosity */
2034  ParseRangeOption(p);
2035  }
2036 
2037  /* >>chng 06 mar 22, add time option to vary only some continua with time */
2038  if( p.nMatch( "TIME" ) )
2039  rfield.lgTimeVary[(p.m_nqh)] = true;
2040 
2041  /* vary option */
2042  if( optimize.lgVarOn )
2043  {
2044  strcpy( optimize.chVarFmt[optimize.nparm], "LUMINOSITY %f LOG range %f %f" );
2045  if( rfield.lgTimeVary[(p.m_nqh)] )
2046  strcat( optimize.chVarFmt[optimize.nparm], " TIME" );
2047  /* pointer to where to write */
2050  optimize.vincr[optimize.nparm] = 0.5;
2051  /* range option in, but cannot be varied */
2052  optimize.vparm[1][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][0]);
2053  optimize.vparm[2][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][1]);
2055  ++optimize.nparm;
2056  }
2057  ++(p.m_nqh);
2058 }
2060 {
2061  /* heating and ionization due to fast neutrons */
2062  hextra.lgNeutrnHeatOn = true;
2063 
2064  /* first number is neutron luminosity
2065  * relative to bolometric luminosity */
2066  hextra.frcneu = (realnum)p.FFmtRead();
2067  if( p.lgEOL() )
2068  p.NoNumb("neutron luminosity");
2069 
2070  /* save as log of fraction */
2071  if( hextra.frcneu > 0. )
2072  hextra.frcneu = (realnum)log10(hextra.frcneu);
2073 
2074  /* second number is efficiency */
2075  hextra.effneu = (realnum)p.FFmtRead();
2076  if( p.lgEOL() )
2077  {
2078  hextra.effneu = 1.0;
2079  }
2080  else
2081  {
2082  if( hextra.effneu <= 0. )
2084  }
2085 }
2087 {
2088  /* flux density of this continuum source, at optional frequency
2089  * expressed as product nu*f_nu */
2090  bool lgNu2 = true;
2091  /* in reads2 */
2092  ParseF_nu(p,"SQCM",lgNu2);
2093 }
2095 {
2096  /* specific luminosity density of this continuum source, at opt freq
2097  * expressed as product nu*f_nu */
2098  bool lgNu2 = true;
2099  /* in reads2 */
2100  ParseF_nu(p,"4 PI",lgNu2);
2101 }
2103 {
2104  DEBUG_ENTRY( "ParsePhi()" );
2105  /* enter phi(h), the number of h-ionizing photons/cm2 */
2106  if( (p.m_nqh) >= LIMSPC )
2107  {
2108  /* too many continua were entered */
2109  fprintf( ioQQQ, " Too many continua entered; increase LIMSPC\n" );
2111  }
2112  /* silly, but calms down the lint */
2113  ASSERT( (p.m_nqh) < LIMSPC );
2114  strcpy( rfield.chRSpec[(p.m_nqh)], "SQCM" );
2115  strcpy( rfield.chSpNorm[(p.m_nqh)], "PHI " );
2116  rfield.totpow[(p.m_nqh)] = p.FFmtRead();
2117  if( p.lgEOL() )
2118  {
2119  p.NoNumb("number of h-ionizing photons");
2120  }
2121  /* check if continuum so intense that crash is likely */
2122  if( rfield.totpow[(p.m_nqh)] > 29. )
2123  {
2124  fprintf( ioQQQ, " Is the flux for this continuum correct?\n" );
2125  fprintf( ioQQQ, " It appears too bright to me.\n" );
2126  }
2127  /* ParseRangeOption in readsun xx parse_rangeoption*/
2128  ParseRangeOption(p);
2129 
2130  /* >>chng 06 mar 22, add time option to vary only some continua with time */
2131  if( p.nMatch( "TIME" ) )
2132  rfield.lgTimeVary[(p.m_nqh)] = true;
2133 
2134  /* vary option */
2135  if( optimize.lgVarOn )
2136  {
2137  strcpy( optimize.chVarFmt[optimize.nparm], "phi(h) %f LOG range %f %f" );
2138  if( rfield.lgTimeVary[(p.m_nqh)] )
2139  strcat( optimize.chVarFmt[optimize.nparm], " TIME" );
2140  /* pointer to where to write */
2143  optimize.vincr[optimize.nparm] = 0.5;
2144  /* range option in, but cannot be varied */
2145  optimize.vparm[1][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][0]);
2146  optimize.vparm[2][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][1]);
2148  ++optimize.nparm;
2149  }
2150  ++(p.m_nqh);
2151 }
2152 void ParseQH(Parser &p)
2153 {
2154  DEBUG_ENTRY( "ParseQH()" );
2155  /* log of number of ionizing photons */
2156  if( (p.m_nqh) >= LIMSPC )
2157  {
2158  /* too many continua were entered */
2159  fprintf( ioQQQ, " Too many continua entered; increase LIMSPC\n" );
2161  }
2162 
2163  /* silly, but calms down the lint */
2164  ASSERT( (p.m_nqh) < LIMSPC );
2165  strcpy( rfield.chRSpec[(p.m_nqh)], "4 PI" );
2166  strcpy( rfield.chSpNorm[(p.m_nqh)], "Q(H)" );
2167  rfield.totpow[(p.m_nqh)] = p.FFmtRead();
2168  if( rfield.totpow[(p.m_nqh)] > 100. && called.lgTalk )
2169  {
2170  fprintf( ioQQQ, " Is this reasonable?\n" );
2171  }
2172  if( p.lgEOL() )
2173  {
2174  p.NoNumb("number of ionizing photons");
2175  }
2176  /* ParseRangeOption in readsun */
2177  ParseRangeOption(p);
2178 
2179  /* >>chng 06 mar 22, add time option to vary only some continua with time */
2180  if( p.nMatch( "TIME" ) )
2181  rfield.lgTimeVary[(p.m_nqh)] = true;
2182 
2183  /* vary option */
2184  if( optimize.lgVarOn )
2185  {
2186  strcpy( optimize.chVarFmt[optimize.nparm], "Q(H) %f LOG range %f %f" );
2187  if( rfield.lgTimeVary[(p.m_nqh)] )
2188  strcat( optimize.chVarFmt[optimize.nparm], " TIME" );
2189  /* pointer to where to write */
2192  optimize.vincr[optimize.nparm] = 0.5;
2193  /* range option in, but cannot be varied */
2194  optimize.vparm[1][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][0]);
2195  optimize.vparm[2][optimize.nparm] = (realnum)log10(rfield.range[(p.m_nqh)][1]);
2197  ++optimize.nparm;
2198  }
2199  /* increment number of luminosity sources specified */
2200  ++(p.m_nqh);
2201 }
2203 {
2204  /* this is the Roberto Terlivich command, no telling if it still works */
2205  radius.dRadSign = -1.;
2206 }
2208 {
2209  DEBUG_ENTRY( "ParseSpecial()" );
2210  /* special key, can do anything */
2212 }
2214 {
2215  /* taumin command minimum optical depths for lines dafault 1e-20 */
2216  opac.taumin = (realnum)exp10(p.FFmtRead());
2217  if( p.lgEOL() )
2218  p.NoNumb("minimum optical depth");
2219 }
2221 {
2222  /* read in title of model starting in col 5 -- prefer to get string
2223  in quotes, but for the moment if it's not present take what you
2224  can get */
2225  string title;
2226  if (p.GetQuote(title) != 0)
2227  strcpy( input.chTitle , p.getRawTail().c_str() );
2228  else
2229  strcpy( input.chTitle , title.c_str() );
2230 }
2232 {
2233  DEBUG_ENTRY( "ParseTolerance()" );
2234  fprintf(ioQQQ,
2235  "Sorry, this command has been replaced with the SET TEMPERATURE TOLERANCE command.\n");
2237 }
2239 {
2240  /* velocity power law as a function of radius. */
2242 
2243  DoppVel.lgTurbLawOn = true;
2245  /* for now, insist upon non-positive power,
2246  * so that velocity decreases with increasing radius. */
2247  ASSERT( DoppVel.TurbVelLaw <= 0.f );
2248 }
2250 {
2251  DEBUG_ENTRY( "ParseTurbulence()" );
2252  string ExtraPars;
2253 
2254  if( p.nMatch("EQUIPART") )
2255  {
2256  /* turbulence equipartition option */
2257  DoppVel.lgTurbEquiMag = true;
2258 
2259  /* this is optional F parameter from equation 34 of
2260  *>>refer pressure turb Heiles, C. & Troland, T.H. 2005, ApJ, 624, 773
2261  * turbulent energy density will be rho F v^2 / 2 */
2263  if( p.lgEOL() )
2264  {
2265  /* this is the default value of 3 for isotropic turbulence */
2266  DoppVel.Heiles_Troland_F = 3.f;
2267  }
2268  }
2269  else
2270  {
2271  /* line has turbulent velocity in km/s */
2272  DoppVel.lgTurbEquiMag = false;
2274  if( p.lgEOL() )
2275  p.NoNumb("microturbulent velocity");
2276 
2277  if( p.nMatch(" LOG") )
2278  {
2279  if( DoppVel.TurbVel > 32. )
2280  {
2281  fprintf( ioQQQ, "PROBLEM the log of the turbulence is "
2282  "%.2e - I cannot handle a number this big.\n",
2283  DoppVel.TurbVel );
2284  fprintf( ioQQQ, " The line image was\n" );
2285  p.PrintLine(ioQQQ);
2286  fprintf( ioQQQ, " Sorry.\n" );
2288  }
2289 
2291  }
2292  /* now convert from km/s to cm/s */
2293  DoppVel.TurbVel *= 1e5;
2294 
2295  if( DoppVel.TurbVel <= 0. )
2296  {
2297  fprintf( ioQQQ, " PROBLEM: the turbulent velocity needs to be > 0, but this was entered: %e\n",
2298  DoppVel.TurbVel );
2299  fprintf( ioQQQ, " Bailing out. Sorry.\n" );
2301  }
2302  else if( DoppVel.TurbVel >= SPEEDLIGHT )
2303  {
2304  fprintf( ioQQQ, " PROBLEM: A turbulent velocity greater than speed of light is not allowed, this was entered: %e\n",
2305  DoppVel.TurbVel );
2306  fprintf( ioQQQ, " Bailing out. Sorry.\n" );
2308  }
2309 
2310  /* this is optional F parameter from equation 34 of
2311  *>>refer pressure turb Heiles, C. & Troland, T.H. 2005, ApJ, 624, 773
2312  * turbulent energy density will be rho F v^2 / 2 */
2314  if( p.lgEOL() )
2315  {
2316  /* this is the default value of 3 for isotropic turbulence */
2317  DoppVel.Heiles_Troland_F = 3.f;
2318  }
2319 
2320  /* option to have turbulence be dissipative, keyword is dissipate,
2321  * argument is log of scale length in cm */
2322  if( p.nMatch("DISS") )
2323  {
2325  if( p.lgEOL() )
2326  p.NoNumb("turbulence dissipation scale");
2327  ExtraPars += " DISSIPATE %f";
2328  }
2329  }
2330 
2331  /* include turbulent pressure in equation of state?
2332  * >>chng 06 mar 24, include turbulent pressure in gas equation of state by default,
2333  * but not if NO PRESSURE occurs */
2334  if( p.nMatch(" NO ") && p.nMatch("PRES") )
2335  {
2336  DoppVel.lgTurb_pressure = false;
2337  ExtraPars += " NO PRESSURE";
2338  }
2339  else
2340  {
2341  /* default is to include turbulent pressure in equation of state */
2342  DoppVel.lgTurb_pressure = true;
2343  }
2344 
2345  /* vary option */
2346  if( optimize.lgVarOn && !p.nMatch("EQUIPART") )
2347  {
2348  /* only one parameter to vary */
2350  // the line image
2351  strcpy( optimize.chVarFmt[optimize.nparm], "TURBULENCE= %f LOG %f" );
2352  strcat( optimize.chVarFmt[optimize.nparm], ExtraPars.c_str() );
2353 
2354  /* pointer to where to write */
2356  /* turbulent velocity */
2357  optimize.vparm[0][optimize.nparm] = log10(DoppVel.TurbVel/1e5f);
2359  if( p.nMatch("DISS") )
2360  {
2363  }
2364  optimize.vincr[optimize.nparm] = 0.1f;
2365  ++optimize.nparm;
2366  }
2367 
2368  /* set initial turbulence */
2370 }
void ParseState(Parser &p)
Definition: parse_state.cpp:9
long int nSave
Definition: input.h:102
realnum col_h2
Definition: stopcalc.h:74
void ParseF_nuSpecific(Parser &p)
void ParseL_nu(Parser &p)
long int lim_iter
Definition: iterations.h:62
bool lgStoutLevelsSet
Definition: atmdat.h:415
realnum cextpw
Definition: thermal.h:152
double emm() const
Definition: mesh.h:93
bool nMatch(const char *chKey) const
Definition: parser.h:150
void ParseHDEN(Parser &p)
Definition: parse_hden.cpp:10
t_fudgec fudgec
Definition: fudgec.cpp:5
void ParseAperture(Parser &p)
t_mole_global mole_global
Definition: mole.cpp:7
long int iter_malloc
Definition: iterations.h:40
char chLamdaFile[FILENAME_PATH_LENGTH]
Definition: atmdat.h:395
double Radius
Definition: radius.h:31
void ParseAgn(Parser &p)
Definition: parse_agn.cpp:10
realnum colnut
Definition: stopcalc.h:69
double depth
Definition: radius.h:31
void ParseHydrogen(Parser &)
realnum GetDensity(realnum z)
Definition: cosmology.cpp:28
t_atmdat atmdat
Definition: atmdat.cpp:6
void InitMonitorResults(void)
t_thermal thermal
Definition: thermal.cpp:6
void SetGasPhaseDensity(const long nelem, const realnum density)
Definition: dense.cpp:106
void ParseDont(Parser &p)
Definition: parse_dont.cpp:28
bool lgTurbLawOn
Definition: doppvel.h:33
bool lgStoutHybrid
Definition: atmdat.h:404
void ParseInitCount(Parser &p)
realnum DensityPower
Definition: dense.h:250
bool lgChiantiLevelsSet
Definition: atmdat.h:388
realnum size
Definition: geometry.h:77
void echo(void) const
Definition: parser.cpp:189
double FFmtRead(void)
Definition: parser.cpp:472
double exp10(double x)
Definition: cddefines.h:1368
void ParseDynaWind(Parser &p)
Definition: dynamics.cpp:1818
const int ipHE_LIKE
Definition: iso.h:65
void ParseAge(Parser &p)
Definition: parse_age.cpp:38
char chStoutFile[FILENAME_PATH_LENGTH]
Definition: atmdat.h:408
void ParsePlot(Parser &p)
Definition: parse_plot.cpp:16
void ParseDatabaseISO(long ipISO, Parser &p)
t_input input
Definition: input.cpp:12
void ParseCylinder(Parser &p)
bool lgGrid
Definition: grid.h:41
void ParseBlackbody(Parser &p)
t_opac opac
Definition: opacity.cpp:5
realnum Heiles_Troland_F
Definition: doppvel.h:37
void setline(const char *const card)
Definition: parser.h:82
realnum turrad
Definition: hextra.h:51
long int nvfpnt[LIMPAR]
Definition: optimize.h:198
bool lgSphericalDilution[LIMSPC]
Definition: rfield.h:321
void ParseBackgrd(Parser &p)
realnum turback
Definition: hextra.h:51
bool Command(const char *name, OptionParser doOpts)
Definition: parser.h:193
bool lgStoutOn
Definition: atmdat.h:402
long int m_nqh
Definition: parser.h:54
t_isoCTRL iso_ctrl
Definition: iso.cpp:9
t_cpu_i & i()
Definition: cpu.h:419
void ParseFail(Parser &p)
double totpow[LIMSPC]
Definition: rfield.h:284
void ParseEnergy(Parser &p)
bool lgCylnOn
Definition: radius.h:127
char chDffTrns[4]
Definition: rfield.h:217
vector< long int > nend
Definition: iterations.h:71
void ParseDatabaseH2(Parser &p)
realnum windv0
Definition: wind.h:11
char chRSpec[LIMSPC][5]
Definition: rfield.h:335
void ParseMonitorResults(Parser &p)
void ParseNeutrons(Parser &p)
double CylindHigh
Definition: radius.h:126
void ParseExtinguish(Parser &p)
realnum redshift_start
Definition: cosmology.h:26
long int iEmissPower
Definition: geometry.h:71
void ParseLuminosity(Parser &p)
void ParseTable(Parser &p)
Definition: parse_table.cpp:94
void ParseSet(Parser &p)
Definition: parse_set.cpp:38
realnum HextraSSradius
Definition: hextra.h:73
bool lgTrOpt
Definition: optimize.h:256
long nStoutMaxLevelsFe
Definition: atmdat.h:410
t_StopCalc StopCalc
Definition: stopcalc.cpp:7
t_conv conv
Definition: conv.cpp:5
long int nOptimiz
Definition: optimize.h:250
t_hextra hextra
Definition: hextra.cpp:5
bool lgNeutrnHeatOn
Definition: hextra.h:81
long nLamdaMaxLevels
Definition: atmdat.h:397
bool lgOpacityFine
Definition: rfield.h:404
double distance
Definition: radius.h:71
int GetQuote(string &chLabel)
Definition: parser.cpp:213
realnum varang[LIMPAR][2]
Definition: optimize.h:201
realnum fiscal
Definition: geometry.h:29
long int nRead
Definition: input.h:105
t_phycon phycon
Definition: phycon.cpp:6
double TableRadius[LIMSPC]
Definition: rfield.h:319
realnum effneu
Definition: hextra.h:85
realnum TurbVel
Definition: doppvel.h:21
void ParseFill(Parser &p)
void ParseRatio(Parser &p)
Definition: parse_ratio.cpp:10
void ParseHExtra(Parser &p)
void ParseConstant(Parser &p)
realnum TurbHeatSave
Definition: hextra.h:42
void ParseTitle(Parser &)
double getNumberDefaultAlwaysLog(const char *chDesc, double fdef)
Definition: parser.cpp:446
bool lgIonizReevaluate
Definition: rfield.h:110
bool lgTeHigh
Definition: thermal.h:72
void ParseHelp(Parser &)
bool lgChiantiPrint
Definition: atmdat.h:378
#define NFUDGC
Definition: fudgec.h:9
void ParseFudge(Parser &p)
void ParseCrashDo(Parser &p)
char chVarFmt[LIMPAR][FILENAME_PATH_LENGTH_2]
Definition: optimize.h:267
vector< double > StopThickness
Definition: iterations.h:77
realnum covgeo
Definition: geometry.h:45
FILE * ioQQQ
Definition: cddefines.cpp:7
void ParseStop(Parser &p)
Definition: parse_stop.cpp:17
realnum vparm[LIMEXT][LIMPAR]
Definition: optimize.h:192
realnum FillFac
Definition: geometry.h:29
bool lgAllTransitions
Definition: conv.h:252
char chTitle[INPUT_LINE_LENGTH]
Definition: input.h:80
bool lgCExtraOn
Definition: thermal.h:151
bool lgTalk
Definition: called.h:12
void ParseCosm(Parser &p)
void ParseSave(Parser &p)
Definition: parse_save.cpp:85
void ParseGrid(Parser &p)
Definition: parse_grid.cpp:10
t_DoppVel DoppVel
Definition: doppvel.cpp:5
void ParseDoubleTau(Parser &)
bool lgConverge_set
Definition: iterations.h:60
const char * name
Definition: parser.h:28
void ParseDLaw(Parser &p)
Definition: parse_dlaw.cpp:10
#define MIN2(a, b)
Definition: cddefines.h:803
Definition: parser.h:43
bool lgTurbEquiMag
Definition: doppvel.h:46
realnum colpls
Definition: stopcalc.h:69
realnum TurbVelZero
Definition: doppvel.h:25
void ParseRangeOption(Parser &p)
double cutoff[LIMSPC][3]
Definition: rfield.h:284
bool lgDeprecatedComment
Definition: input.h:119
realnum frcneu
Definition: hextra.h:83
void ParseVLaw(Parser &p)
void ParsePrint(Parser &p)
bool lgVarOn
Definition: optimize.h:207
bool lgPhysOK
Definition: phycon.h:111
t_dense dense
Definition: global.cpp:15
void ParseSphere(Parser &p)
Definition: parse_sphere.cpp:9
void ParseTurbulence(Parser &p)
bool lgTimeVary[LIMSPC]
Definition: rfield.h:290
realnum autocv
Definition: conv.h:258
void ParseMap(Parser &p)
Definition: parse_map.cpp:9
static t_version & Inst()
Definition: cddefines.h:209
double range[LIMSPC][2]
Definition: rfield.h:331
void ParseIlluminate(Parser &p)
void ParseGravity(Parser &p)
const int LIMSPC
Definition: rfield.h:21
void ParseQH(Parser &p)
void ParseCommands(void)
Wind wind
Definition: wind.cpp:5
t_trace trace
Definition: trace.cpp:5
bool lgAbnSolar
Definition: abund.h:202
vector< long int > IterPrnt
Definition: iterations.h:43
void ParseCosmology(Parser &p)
bool lgBallistic(void) const
Definition: wind.h:31
OptionParser action
Definition: parser.h:29
void ParseCompile(Parser &p)
double rinner
Definition: radius.h:31
t_abund abund
Definition: abund.cpp:5
t_geometry geometry
Definition: geometry.cpp:5
void ParseTest(Parser &p)
Definition: parse_test.cpp:9
bool lgOutOnly
Definition: rfield.h:222
const double TEMP_STOP_DEFAULT
Definition: phycon.h:119
bool lgChiantiHybrid
Definition: atmdat.h:376
t_dark_matter dark
Definition: dark_matter.cpp:5
bool lgMap
Definition: conv.h:231
realnum CoolExtra
Definition: thermal.h:152
void ParseSpecies(Parser &p)
realnum HColStop
Definition: stopcalc.h:69
realnum redshift_current
Definition: cosmology.h:26
double slope[LIMSPC]
Definition: rfield.h:284
bool lgLamdaPrint
Definition: atmdat.h:393
bool lgSaveXspec
Definition: grid.h:38
void ParseNuL_nu(Parser &p)
void ParseHeLike(Parser &)
void ParseGlobule(Parser &p)
void ParseCaseB(Parser &p)
Definition: parse_caseb.cpp:9
void ParseTolerance(Parser &)
string getRawTail()
Definition: parser.h:217
void ParseMagnet(Parser &p)
Definition: magnetic.cpp:156
bool lgChiantiOn
Definition: atmdat.h:374
bool lgTrace
Definition: trace.h:12
void ParsePGrains(Parser &)
bool m_lgDSet
Definition: parser.h:55
void set_point(long int ipnt)
Definition: parser.h:95
void ParseDiffuse(Parser &p)
bool lgEnabled
Definition: h2_priv.h:352
void ParseF_nu(Parser &p, const char *chType, bool lgNU2)
Definition: parse_f_nu.cpp:9
void ParseIonParX(Parser &p)
char chCloudyChiantiFile[FILENAME_PATH_LENGTH]
Definition: atmdat.h:382
realnum TurbHeat
Definition: hextra.h:42
long int nparm
Definition: optimize.h:204
t_plotCom plotCom
Definition: plot.cpp:20
t_pressure pressure
Definition: pressure.cpp:9
t_rfield rfield
Definition: rfield.cpp:9
long nStoutMaxLevels
Definition: atmdat.h:413
realnum TurbVelLaw
Definition: doppvel.h:29
bool lgCaseB
Definition: opacity.h:174
void ParseAbsMag(Parser &p)
Definition: parse_absmag.cpp:9
realnum DispScale
Definition: doppvel.h:51
long int nTrOpt
Definition: optimize.h:255
float realnum
Definition: cddefines.h:124
#define EXIT_FAILURE
Definition: cddefines.h:168
bool lgPrintNumberOfLevels
Definition: iso.h:346
bool last(void) const
Definition: parser.cpp:196
const realnum COLUMN_INIT
Definition: stopcalc.h:14
const int INPUT_LINE_LENGTH
Definition: cddefines.h:301
realnum HextraScaleDensity
Definition: hextra.h:60
void ParseFluc(Parser &p)
Definition: parse_fluc.cpp:9
void ParseNorm(Parser &p)
Definition: parse_norm.cpp:10
bool lgTurbHeatVaryTime
Definition: hextra.h:76
#define cdEXIT(FAIL)
Definition: cddefines.h:482
DepthTable DLW
Definition: dense.h:198
bool lgMPI_talk() const
Definition: cpu.h:397
void ParseCMB(double z, long int *nqh)
Definition: parse_CMB.cpp:10
double tabval(double r0, double depth) const
Definition: depth_table.cpp:8
NORETURN void NoNumb(const char *chDesc) const
Definition: parser.cpp:345
bool lgRadiusKnown
Definition: radius.h:122
realnum covrt
Definition: geometry.h:61
long min(int a, long b)
Definition: cddefines.h:762
void ParsePhi(Parser &p)
bool lgMustBlockHIon
Definition: rfield.h:92
diatomics h2("h2", 4100.,&hmi.H2_total, Yan_H2_CS)
const long LIMPAR
Definition: optimize.h:61
long nChiantiMaxLevels
Definition: atmdat.h:386
int gravity_symmetry
Definition: pressure.h:84
t_iterations iterations
Definition: iterations.cpp:6
void ParseLaser(Parser &p)
void ParseDynaTime(Parser &p)
Definition: dynamics.cpp:1671
void ParseConvHighT(Parser &)
void ParseIntensity(Parser &p)
t_optimize optimize
Definition: optimize.cpp:6
void ParseDrive(Parser &p)
Definition: parse_drive.cpp:28
char chSpNorm[LIMSPC][5]
Definition: rfield.h:335
t_grid grid
Definition: grid.cpp:5
realnum cryden
Definition: hextra.h:24
bool isVar(void) const
Definition: parser.cpp:144
t_radius radius
Definition: radius.cpp:5
double dRadSign
Definition: radius.h:74
realnum vincr[LIMPAR]
Definition: optimize.h:195
realnum TempLoStopZone
Definition: stopcalc.h:42
long int lim_zone
Definition: iterations.h:61
void ParseTrace(Parser &p)
Definition: parse_trace.cpp:11
t_prt prt
Definition: prt.cpp:14
void ParseTauMin(Parser &p)
bool lgTrOvrd
Definition: trace.h:121
bool lgDoLineTrans
Definition: rfield.h:99
bool lgPredLumin
Definition: radius.h:145
void ParseGrain(Parser &p)
Definition: parse_grain.cpp:12
long int nfudge
Definition: fudgec.h:17
bool lgHextraSS
Definition: hextra.h:64
realnum gas_phase[LIMELM]
Definition: dense.h:76
long int itermx
Definition: iterations.h:37
void ParsePowerlawContinuum(Parser &p)
double HextraSS_M
Definition: hextra.h:70
bool lgOptimizeAsLinear[LIMPAR]
Definition: optimize.h:184
void help(FILE *fp) const
Definition: parser.cpp:294
long int LimFail
Definition: conv.h:228
bool lgNegativeIncrements
Definition: grid.h:37
void ParseIonParI(Parser &p)
#define ASSERT(exp)
Definition: cddefines.h:613
realnum covaper
Definition: geometry.h:54
bool lgVisibilityStatus
Definition: input.h:99
char chDenseLaw[5]
Definition: dense.h:176
bool lgOpacityReevaluate
Definition: rfield.h:103
bool lgPrintHTML
Definition: prt.h:277
realnum ConstTemp
Definition: thermal.h:56
bool getline()
Definition: parser.cpp:273
const int ipH_LIKE
Definition: iso.h:64
realnum emdot
Definition: wind.h:39
bool lgHextraDensity
Definition: hextra.h:57
#define DEBUG_ENTRY(funcname)
Definition: cddefines.h:723
void ParseCovering(Parser &p)
t_cosmology cosmology
Definition: cosmology.cpp:8
void ParseSpecial(Parser &)
void ParseNuF_nu(Parser &p)
bool lgAutoIt
Definition: conv.h:249
double getNumberCheckAlwaysLog(const char *chDesc)
Definition: parser.cpp:427
#define MIN3(a, b, c)
Definition: cddefines.h:808
void ParseRadius(Parser &p)
realnum ConstGrainTemp
Definition: thermal.h:59
double egamry() const
Definition: mesh.h:97
void ParseTLaw(Parser &p)
Definition: parse_tlaw.cpp:13
void ParseRoberto(Parser &)
realnum tauend
Definition: stopcalc.h:23
bool lgEOL(void) const
Definition: parser.h:113
#define MAX2(a, b)
Definition: cddefines.h:824
void ParseIterations(Parser &p)
void ParseCosmicRays(Parser &p)
bool lgStoutPrint
Definition: atmdat.h:406
int fprintf(const Output &stream, const char *format,...)
Definition: service.cpp:1121
NORETURN void CommandError(void) const
Definition: parser.cpp:260
double dense_fabden(double radius, double depth)
bool lgTurb_pressure
Definition: doppvel.h:42
void ParseBremsstrahlung(Parser &p)
vector< double > external_mass[3]
Definition: pressure.h:86
realnum DoubleTau
Definition: rt.h:178
long int nShape
Definition: rfield.h:306
void ParseCMBOuter(Parser &p)
double dense_parametric_wind(double rad)
realnum HextraSSalpha
Definition: hextra.h:67
void ParseOptimize(Parser &p)
int PrintLine(FILE *fp) const
Definition: parser.h:206
bool lgStatic(void) const
Definition: wind.h:24
vector< double > StopRadius
Definition: iterations.h:80
double self_mass_factor
Definition: pressure.h:85
bool lgSizeSet
Definition: geometry.h:80
bool lgPlotON
Definition: plot.h:21
double rdfalt
Definition: radius.h:130
GrainVar gv
Definition: grainvar.cpp:5
bool lgLamdaOn
Definition: atmdat.h:391
static t_cpu cpu
Definition: cpu.h:427
bool m_lgEOF
Definition: parser.h:55
bool lgNoMole
Definition: mole.h:325
bool lgNoVary
Definition: optimize.h:175
void ParseMetal(Parser &p)
Definition: parse_metal.cpp:12
long nChiantiMaxLevelsFe
Definition: atmdat.h:384
const int ipHYDROGEN
Definition: cddefines.h:349
long int nGridCommands
Definition: grid.h:56
bool lgChiantiExp
Definition: atmdat.h:380
void doSetVar(void)
Definition: parser.cpp:161
void ParseCExtra(Parser &p)
bool lgBlockHIon
Definition: rfield.h:96
void ParseDielectronic(Parser &)
void ParseDatabase(Parser &p)
realnum EdenExtra
Definition: dense.h:217
realnum col_h2_nut
Definition: stopcalc.h:80
double r_200
Definition: dark_matter.h:19
void ParseCoronal(Parser &p)
long int nvarxt[LIMPAR]
Definition: optimize.h:198
char chSpType[LIMSPC][6]
Definition: rfield.h:335
void ParseElement(Parser &p)
realnum taumin
Definition: opacity.h:167
t_called called
Definition: called.cpp:4
void ParseForceTemperature(Parser &p)
bool lgLamdaLevelsSet
Definition: atmdat.h:399
bool isComment(void) const
Definition: parser.cpp:140
long int nplot
Definition: plot.h:27
realnum filpow
Definition: geometry.h:29
void ParseInterp(Parser &p)
void ParseEden(Parser &p)
void ParseDistance(Parser &p)
void ParseDarkMatter(Parser &p)
realnum fudgea[NFUDGC]
Definition: fudgec.h:15
void init(void)
Definition: input.cpp:186
void ParseAbundances(Parser &p)
void ParseChemistry(Parser &p)
Definition: mole.cpp:144
bool lgDColOn
Definition: grainvar.h:494
bool lgHextraDepth
Definition: hextra.h:49
t_rt rt
Definition: rt.cpp:5