cloudy  trunk
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
conv_base.cpp
Go to the documentation of this file.
1 /* This file is part of Cloudy and is copyright (C)1978-2017 by Gary J. Ferland and
2  * others. For conditions of distribution and use see copyright notice in license.txt */
3 /*ConvBase main routine to drive ionization solution for all species, find total opacity
4  * called by ConvIoniz */
5 /*lgConverg check whether ionization of element nelem has converged */
6 #include "cddefines.h"
7 #include "dynamics.h"
8 #include "trace.h"
9 #include "elementnames.h"
10 #include "save.h"
11 #include "phycon.h"
12 #include "secondaries.h"
13 #include "stopcalc.h"
14 #include "grainvar.h"
15 #include "highen.h"
16 #include "hmi.h"
17 #include "pressure.h"
18 #include "taulines.h"
19 #include "rt.h"
20 #include "grains.h"
21 #include "atmdat.h"
22 #include "ionbal.h"
23 #include "ion_trim.h"
24 #include "opacity.h"
25 #include "cooling.h"
26 #include "thermal.h"
27 #include "iso.h"
28 #include "conv.h"
29 #include "h2.h"
30 #include "deuterium.h"
31 #include "mole.h"
32 #include "rfield.h"
33 #include "dense.h"
34 
35 STATIC bool lgNetEdenSrcSmall( void );
36 
37 void UpdateUTAs( void );
38 
39 /*lgIonizConverg check whether ionization of element nelem has converged, called by ionize,
40  * returns true if element is converged, false if not */
41 namespace
42 {
43  static double MIN_CONV_FRAC=1e-4;
44 
45  class IonizConverg
46  {
47  /* this is fractions [ion stage][nelem], ion stage = 0 for atom, nelem=0 for H*/
48  double OldFracs[LIMELM+1][LIMELM+1];
49  public:
50  IonizConverg()
51  {
52  for (long nelem = ipHYDROGEN; nelem < LIMELM; ++nelem)
53  {
54  if (!dense.lgElmtOn[nelem])
55  continue;
56  for (long ion = 0; ion <= nelem+1; ++ion)
57  {
58  OldFracs[nelem][ion] = dense.xIonDense[nelem][ion];
59  }
60  }
61  }
62  void operator()(
63  /* atomic number on the C scale, 0 for H, 25 for Fe */
64  long loop_ion ,
65  /* this is allowed error as a fractional change. Most are 0.15 */
66  double delta ,
67  /* option to print abundances */
68  bool lgPrint )
69  {
70  bool lgConverg_v = false;
71  double Abund,
72  bigchange ,
73  change ,
74  one;
75 
76  DEBUG_ENTRY( "IonizConverg::operator()" );
77 
78  for (long nelem =ipHYDROGEN; nelem<LIMELM; ++nelem )
79  {
80  if( !dense.lgElmtOn[nelem] )
81  continue;
82 
83  double abundold=0. , abundnew=0.;
84  long ionchg=-1;
85 
86  /* arguments are atomic number, ionization stage
87  * OldFracs[nelem][0] is abundance of nelem (cm^-3) */
88 
89  /* this function returns true if ionization of element
90  * with atomic number nelem has not changed by more than delta,*/
91 
92  /* check whether abundances exist yet, only do this for after first zone */
93  /*if( nzone > 0 )*/
94  /* >>chng 03 sep 02, check on changed ionization after first call through,
95  * to insure converged constant eden / temperature models */
96  if( conv.nPres2Ioniz )
97  {
98  /* >>chng 04 aug 31, this had been static, caused last hits on large changes
99  * in ionization */
100  bigchange = 0.;
101  change = 0.;
102  Abund = dense.gas_phase[nelem];
103 
104  /* loop over all ionization stages, loop over all ions, not just active ones,
105  * since this also sets old values, and so will zero out non-existant ones */
106  for( long ion=0; ion <= (nelem+1); ++ion )
107  {
108  /*lint -e727 OlsdFracs not initialized */
109  if( OldFracs[nelem][ion]/Abund > MIN_CONV_FRAC &&
110  dense.xIonDense[nelem][ion]/Abund > MIN_CONV_FRAC )
111  /*lint +e727 OlsdFracs not initialized */
112  {
113  /* check change in old to new value */
114  one = fabs(dense.xIonDense[nelem][ion]-OldFracs[nelem][ion])/
115  OldFracs[nelem][ion];
116  change = MAX2(change, one );
117  /* remember abundances for largest change */
118  if( change>bigchange )
119  {
120  bigchange = change;
121  abundold = OldFracs[nelem][ion]/Abund;
122  abundnew = dense.xIonDense[nelem][ion]/Abund;
123  ionchg = ion;
124  }
125  }
126  /* now set new value */
127  OldFracs[nelem][ion] = dense.xIonDense[nelem][ion];
128  }
129 
130  if( change >= delta )
131  {
132  char chConvIoniz[INPUT_LINE_LENGTH];
133  sprintf( chConvIoniz , "%2s ion" , elementnames.chElementSym[nelem] );
134  conv.setConvIonizFail(chConvIoniz,abundold,abundnew);
135  ASSERT( abundold>0. && abundnew>0. );
136  }
137  else
138  lgConverg_v = true;
139  }
140  else
141  {
142  bigchange = 0.;
143  for( long ion=0; ion <= (nelem+1); ++ion )
144  {
145  OldFracs[nelem][ion] = dense.xIonDense[nelem][ion];
146  }
147 
148  }
149 
150  /* option to print abundances */
151  if( lgPrint )
152  {
153  fprintf(ioQQQ," nz %ld loop %ld element %li converged? %c worst %ld change %g\n",
154  nzone, loop_ion, nelem, TorF(lgConverg_v),ionchg,bigchange);
155  for( long ion=0; ion<(nelem+1); ++ion )
156  {
157  fprintf(ioQQQ,"\t%.5e", dense.xIonDense[nelem][ion]/dense.gas_phase[nelem]);
158  }
159  fprintf(ioQQQ,"\n");
160  }
161  }
162  }
163  };
164 
165  double dground(long nelem, long ion)
166  {
167  DEBUG_ENTRY( "IonizConverg::dground()" );
168  if( ! ( nelem >= ipHYDROGEN && nelem < LIMELM &&
169  ion >= 0 && ion <= nelem+1 ) )
170  {
171  fprintf( ioQQQ,
172  "ERROR: Invalid ionization stage in dground()."
173  " nelem= %ld\t ion= %ld\n",
174  nelem+1, ion );
175  cdEXIT( EXIT_FAILURE );
176  }
177  long ipISO = nelem - ion;
178  if (ipISO >= 0 && ipISO < NISO)
179  return iso_sp[ipISO][nelem].st[0].Pop();
180  else
181  return dense.xIonDense[nelem][ion];
182  }
183 }
184 
185 /*ConvBase main routine to drive ionization solution for all species, find total opacity
186  * called by ConvIoniz
187  * return 0 if ok, 1 if abort */
188 int ConvBase(
189  /* this is zero the first time ConvBase is called by convIoniz,
190  * counts number of call thereafter */
191  long loopi )
192 {
193  double HeatOld,
194  EdenTrue_old,
195  EdenFromMolecOld,
196  EdenFromGrainsOld,
197  HeatingOld ,
198  CoolingOld;
199  static double SecondOld;
200  static long int nzoneOTS=-1;
201  const int LOOP_ION_LIMIT = 10;
202  long int loop_ion;
203  static double SumOTS=0. , OldSumOTS[2]={0.,0.};
204  double save_iso_grnd[NISO][LIMELM];
205  valarray<realnum> mole_save(mole_global.num_calc);
206  IonizConverg lgIonizConverg;
207 
208  DEBUG_ENTRY( "ConvBase()" );
209 
210  bool lgNewTrim = false;
211  bool lgIonizWidened = false;
212  // Adaptive ionization trim widening approach -- currently disables
213  if (lgNewTrim)
214  {
215  static double wide_te_used = -1.;
216  static long int nZoneIonizWidenCalled = 0;
217 
218  bool lgWiden = ( nZoneIonizWidenCalled != nzone ||
219  fabs(phycon.te-wide_te_used) > 0.1*phycon.te );
220 
221  if( lgWiden )
222  {
223  nZoneIonizWidenCalled = nzone;
224  wide_te_used = phycon.te;
225  lgIonizWidened = true;
226 
227  for( long nelem=ipHELIUM; nelem<LIMELM; ++nelem )
228  {
229  if( dense.lgElmtOn[nelem] )
230  {
231  /* ion_trim will set conv.lgIonStageTrimed true is any ion has its
232  * lowest stage of ionization dropped or trimmed */
233  ion_widen(nelem);
234  }
235  }
236  }
237 
238  /* following block only set of asserts */
239 # if !defined(NDEBUG)
240  /* check that proper abundances are either positive or zero */
241  for( long nelem=ipHELIUM; nelem<LIMELM; ++nelem)
242  {
243  if( dense.lgElmtOn[nelem] )
244  {
245  ion_trim_validate(nelem, false);
246  }
247  }
248 # endif
249  }
250 
251 
252  /* this is set to phycon.te in tfidle, is used to insure that all temp
253  * vars are properly updated when conv_ionizeopacitydo is called
254  * NB must be same type as phycon.te */
256 
257  /* this allows zone number to be printed with slight increment as zone converged
258  * conv.nPres2Ioniz is incremented at the bottom of this routine */
259  fnzone = (double)nzone + (double)conv.nPres2Ioniz/100.;
260 
261  /* reevaluate pressure */
262  /* this sets values of pressure.PresTotlCurr, also calls tfidle,
263  * and sets the total energy content of gas, which may be important for acvection */
264  PresTotCurrent();
265 
266  /* >>chng 04 sep 15, move EdenTrue_old up here, and will redo at bottom
267  * to find change
268  * find and save current true electron density - if this changes by more than the
269  * tolerance then ionization solution is not converged */
270  /* >>chng 04 jul 27, move eden_sum here from after this loop, so that change in EdenTrue
271  * can be monitored */
272  /* update EdenTrue, eden itself is actually changed in ConvEdenIoniz */
273  /* must not call eden_sum on very first time since for classic PDR total
274  * ionization may still be zero on first call */
275  if( conv.nTotalIoniz )
276  {
277  if( eden_sum() )
278  {
279  /* non-zero return indicates abort condition */
280  ++conv.nTotalIoniz;
281  return 1;
282  }
283  }
284 
285  /* the following two were evaluated in eden_sum
286  * will confirm that these are converged */
287  EdenTrue_old = dense.EdenTrue;
288  EdenFromMolecOld = mole.elec;
289  EdenFromGrainsOld = gv.TotalEden;
290  HeatingOld = thermal.htot;
291  CoolingOld = thermal.ctot;
292 
293  /* remember current ground state population - will check if converged */
294  for( long ipISO=ipH_LIKE; ipISO<NISO; ++ipISO )
295  {
296  for( long nelem=ipISO; nelem<LIMELM;++nelem )
297  {
298  if( dense.lgElmtOn[nelem] )
299  {
300  /* save the ground state population */
301  save_iso_grnd[ipISO][nelem] = iso_sp[ipISO][nelem].st[0].Pop();
302  }
303  }
304  }
305 
306  if( loopi==0 )
307  {
308  /* these will be used to look for oscillating ots rates */
309  OldSumOTS[0] = 0.;
310  OldSumOTS[1] = 0.;
311  conv.lgOscilOTS = false;
312  }
313 
314  if( trace.lgTrace )
315  {
316  fprintf( ioQQQ,
317  " ConvBase called. %.2f Te:%.3e HI:%.3e HII:%.3e H2:%.3e Ne:%.3e htot:%.3e CSUP:%.2e Conv?%c\n",
318  fnzone,
319  phycon.te,
322  hmi.H2_total,
323  dense.eden,
324  thermal.htot,
326  TorF(conv.lgConvIoniz()) );
327  }
328  /* want this flag to be true when we exit, various problems will set false */
330 
331  /* this routine is in heatsum.c, and zeros out the heating array */
332  HeatZero();
333 
334  /* if this is very first call, say not converged, since no meaningful way to
335  * check on changes in quantities. this counter is false only on very first
336  * call, when equal to zero */
337  if( !conv.nTotalIoniz )
338  {
339  conv.setConvIonizFail( "first call", 0.0, 0.0 );
340  }
341 
342  /* must redo photoionization rates for all species on second and later tries */
343  /* always reevaluate the rates when . . . */
344  /* this flag says not to redo opac and photo rates, and following test will set
345  * true if one of several tests not done*/
346  opac.lgRedoStatic = false;
347  if(
348  /* opac.lgOpacStatic (usually true), is set false with no static opacity command */
349  !opac.lgOpacStatic ||
350  /* we are in search mode */
351  conv.lgSearch ||
352  /* this is the first call to this zone */
353  conv.nPres2Ioniz == 0 )
354  {
355  /* we need to redo ALL photoionization rates */
356  opac.lgRedoStatic = true;
357  }
358 
359  /* calculate radiative and dielectronic recombination rate coefficients */
361 
362  /* this adjusts the lowest and highest stages of ionization we will consider,
363  * only safe to call when lgRedoStatic is true since this could lower the
364  * lowest stage of ionization, which needs all its photo rates */
365 
366  /* this will be flag to check whether any ionization stages
367  * were trimmed down */
368  conv.lgIonStageTrimed = false;
369  static long int nZoneIonizTrimCalled = 0;
370  if( ! conv.nTotalIoniz )
371  {
372  nZoneIonizTrimCalled = 0;
373  }
374 
375  /* conv.nTotalIoniz is only 0 (false) on the very first call to here,
376  * when the ionization distribution is not yet done */
377  if( conv.nTotalIoniz )
378  {
379 #ifndef NDEBUG
380  bool lgIonizTrimCalled = false;
381 #endif
382 
383  /* ionization trimming only used for He and heavier, not H */
384  /* only do this one time per zone since up and down cycle can occur */
385  /* >>chng 05 jan 15, increasing temperature above default first conditions, also
386  * no trim during search - this fixed major logical error when sim is
387  * totally mechanically heated to coronal temperatures -
388  * problem discovered by Ronnie Hoogerwerf */
389  /* do not keep changing the trim after the first call within
390  * this zone once we are deep into layer - doing so introduced very
391  * small level of noise as some stages
392  * went up and down - O+2 - O+3 was good example, when small H+3 after He+ i-front
393  * limit to one increase per element per zone */
394  if( !lgNewTrim && conv.nTotalIoniz>2 &&
395  /* only call one time per zone except during search phase,
396  * when only call after 20 times only if temperature has changed */
397  ( conv.lgSearch || nZoneIonizTrimCalled!=nzone ) )
398  {
399 #ifndef NDEBUG
400  lgIonizTrimCalled = true;
401 #endif
402  for( long nelem=ipHELIUM; nelem<LIMELM; ++nelem )
403  {
404  if( dense.lgElmtOn[nelem] )
405  {
406  /* ion_trim will set conv.lgIonStageTrimed true is any ion has its
407  * lowest stage of ionization dropped or trimmed */
408  ion_trim(nelem);
409  }
410  }
411  nZoneIonizTrimCalled = nzone;
412  }
413 
414  /* following block only set of asserts */
415 # if !defined(NDEBUG)
416  /* check that proper abundances are either positive or zero */
417  for( long nelem=ipHELIUM; nelem<LIMELM; ++nelem)
418  {
419  if( dense.lgElmtOn[nelem] )
420  {
421  ion_trim_validate(nelem, lgIonizTrimCalled);
422  }
423  }
424 # endif
425  }
426 
427  /* now check if anything trimmed down */
428  if( conv.lgIonStageTrimed )
429  {
430  /* something was trimmed down, so say that ionization not yet stable */
431  /* say that ionization has not converged, secondaries changed by too much */
432  conv.setConvIonizFail( "IonTrimmed", 0.0, 0.0 );
433  }
434 
436 
437 #ifndef NDEBUG
438  long ionLowCache[LIMELM], ionHighCache[LIMELM];
439 #endif
440  for( long nelem=ipHYDROGEN; nelem<LIMELM; nelem++ )
441  {
442  if (!dense.lgElmtOn[nelem])
443  continue;
444 #ifndef NDEBUG
445  ionLowCache[nelem] = dense.IonLow[nelem];
446  ionHighCache[nelem] = dense.IonHigh[nelem];
447 #endif
448  }
449 
450  /* reevaluate advective terms if turned on */
452  DynaIonize();
453 
454  /* evaluate Compton heating, bound E Compton, cosmic rays */
455  highen();
456 
457  if (!lgConvBaseHeatTest)
458  {
459  SecIoniz();
460  }
461 
462  // Depends on dense.eden, phycon.te and trimming level of H, He
464 
465  iso_update_rates( );
466 
467  for( diatom_iter diatom = diatoms.begin(); diatom != diatoms.end(); ++diatom )
468  {
469  (*diatom)->CalcPhotoionizationRate();
470  }
471 
472  /* >>chng 04 feb 15, add loop over ionization until converged. Non-convergence
473  * in this case means ionization ladder pop changed, probably because of way
474  * that Auger ejection is treated - loop exits when conditions tested at end are met */
475 
476  static ConvergenceCounter cctr=conv.register_("CONV_BASE_CALLS");
477  ++cctr;
478 
479  const int nconv = 3;
480  double xIonDense0[nconv][LIMELM][LIMELM+1];
481  double xIonDense_prev[LIMELM][LIMELM+1];
482  bool lgOscill[LIMELM];
483  for (long nelem=0; nelem < LIMELM; ++nelem)
484  {
485  if (!dense.lgElmtOn[nelem])
486  continue;
487  lgOscill[nelem] = false;
488  }
489 
490  loop_ion = 0;
491  do
492  {
493  static ConvergenceCounter cctrl=conv.register_("CONV_BASE_LOOPS");
494  ++cctrl;
496 
497  /* set the convergence flag to true,
498  * if anything changes in ConvBase, it will be set false */
499  if( loop_ion )
501 
502  /* charge transfer evaluation needs to be here so that same rate
503  * coefficient used for H ion and other recombination */
504  /* fill in master charge transfer array, and zero out some arrays that track effects */
505 
506  ChargTranEval();
507 
508  /* find grain temperature, charge, and gas heating rate */
509  /* >>chng 04 apr 15, moved here from after call to HeLike(), PvH */
510  GrainDrive();
511 
512  /* evaluate Compton heating, bound E Compton, cosmic rays */
513  highen();
514 
515  /* find corrections for three-body rec - collisional ionization */
516  atmdat_3body();
517 
518  /* update UTA inner-shell ionization rates */
519  UpdateUTAs();
520 
522 
523  for( long i=0; i < mole_global.num_calc; i++ )
524  {
525  mole_save[i] = (realnum) mole.species[i].den;
526  }
527 
528  if (loop_ion != 0)
529  {
530  for (long nelem = ipHYDROGEN; nelem < LIMELM; ++nelem)
531  {
532  if (!dense.lgElmtOn[nelem] || lgOscill[nelem])
533  continue;
534  for (long ion = 0; ion <= nelem+1; ++ion)
535  {
536  xIonDense_prev[nelem][ion] = xIonDense0[0][nelem][ion];
537  }
538  }
539  }
540 
541  long ion_loop = 0;
542  bool lgShortCircuit = false;
543  for( ion_loop=0; ion_loop<nconv && !lgShortCircuit; ++ion_loop)
544  {
545  static ConvergenceCounter cctra=conv.register_("CONV_BASE_ACCELS");
546  ++cctra;
547  for (long nelem = ipHYDROGEN; nelem < LIMELM; ++nelem)
548  {
549  if (!dense.lgElmtOn[nelem])
550  continue;
551  for (long ion = 0; ion <= nelem+1; ++ion)
552  {
553  xIonDense0[ion_loop][nelem][ion] = dense.xIonDense[nelem][ion];
554  }
555  }
556 
557  // netion has net ionization for pairs of ions.
558  double netion[t_atmdat::NCX][LIMELM];
559  for ( long nlo=0; nlo<t_atmdat::NCX; ++nlo)
560  for ( long nhi=0; nhi<LIMELM; ++nhi)
561  netion[nlo][nhi] = 0.0;
562 
563  for( long nelem=ipHYDROGEN; nelem<LIMELM; ++nelem )
564  {
565  if (!dense.lgElmtOn[nelem])
566  continue;
567 
568  ASSERT(dense.IonLow[nelem] >= ionLowCache[nelem]);
569  ASSERT(dense.IonHigh[nelem] <= ionHighCache[nelem]);
570  ion_wrapper( nelem );
571  for (long n2=ipHYDROGEN; n2<LIMELM; ++n2)
572  {
573  if (!dense.lgElmtOn[n2] || nelem == n2)
574  continue;
575  long nlo = MIN2(nelem, n2);
576  long nhi = MAX2(nelem, n2);
577  if (nlo >= t_atmdat::NCX)
578  continue;
579 
580  double dl1 = dground(nlo,0);
581  double ul1 = dground(nlo,1);
582  double dl2 = dground(nhi,0);
583  double ul2 = dground(nhi,1);
584  double hion =
585  atmdat.CharExcRecTo[nlo][nhi][0]*dl1*ul2-
586  atmdat.CharExcIonOf[nlo][nhi][0]*ul1*dl2;
587  if (nelem < n2)
588  {
589  netion[nlo][nhi] += hion;
590  }
591  else
592  {
593  netion[nlo][nhi] -= hion;
594  }
595  }
596 
599  }
600 
601  if (ion_loop == 1)
602  {
603  for (long nlo = ipHYDROGEN; nlo < t_atmdat::NCX; ++nlo)
604  {
605  if ( !dense.lgElmtOn[nlo] )
606  continue;
607  for (long nhi = nlo+1; nhi < LIMELM; ++nhi)
608  {
609  if ( !dense.lgElmtOn[nhi] )
610  continue;
611  double dl1 = dground(nlo,0);
612  double ul1 = dground(nlo,1);
613  double dl2 = dground(nhi,0);
614  double ul2 = dground(nhi,1);
615  // Comparison rate is whether error is either a small fraction of the minimum ionization rate,
616  // or a tighter check on whether the fluxes are balanced as well as is feasible numerically
617  double ion_cmp = MAX2(0.01*MIN2(ionbal.RateIonizTot(nlo,0)*dl1, ionbal.RateIonizTot(nhi,0)*dl2),
618  1e-4*atmdat.CharExcRecTo[nlo][nhi][0]*dl1*ul2);
619 
620  bool lgCheckAll = false;
621  if ( (lgCheckAll || (nlo == ipHYDROGEN && nhi == ipOXYGEN)) &&
622  fabs(netion[nlo][nhi]) > ion_cmp && ul1*ul2 > 1e-8*dl1*dl2 )
623  {
624  ostringstream oss;
625  oss << elementnames.chElementSym[nlo] << "/" << elementnames.chElementSym[nhi]
626  << " CX inconsistency";
627  conv.setConvIonizFail( oss.str().c_str(), ion_cmp, netion[nlo][nhi]);
628  }
629  }
630  }
631  }
632 
633  // If not going to end anyhow, check whether changes were small
634  bool lgCanShortCircuit = (ion_loop+1 < nconv);
635  for( long nelem=ipHYDROGEN; nelem<LIMELM && lgCanShortCircuit; ++nelem )
636  {
637  if (!dense.lgElmtOn[nelem])
638  continue;
639  for (long ion = dense.IonLow[nelem];
640  ion <= dense.IonHigh[nelem]; ++ion)
641  {
642  double x0 = xIonDense0[ion_loop][nelem][ion];
643  double x1 = dense.xIonDense[nelem][ion];
644  if (fabs(x0-x1) > 1e-6*(x0+x1) && x1 > 1e-12*dense.gas_phase[nelem])
645  {
646  lgCanShortCircuit = false;
647  break;
648  }
649  }
650  }
651  lgShortCircuit = lgCanShortCircuit;
652  }
653 
654  if (!lgShortCircuit)
655  {
656  // Apply convergence acceleration
657  for (long nelem = ipHYDROGEN; nelem < LIMELM; ++nelem)
658  {
659  if (!dense.lgElmtOn[nelem])
660  continue;
661  double tot0 = 0., tot1 = 0.;
662  double xIonNew[LIMELM+1];
663  for (long ion = 0; ion <= nelem+1; ++ion)
664  {
665  double x0 = xIonDense0[nconv-2][nelem][ion];
666  double x1 = xIonDense0[nconv-1][nelem][ion];
667  double x2 = dense.xIonDense[nelem][ion];
668  xIonNew[ion] = x2;
669  tot0 += x2;
670  // Richardson extrapolation formula to accelerate convergence
671  // Assumes convergence to x^ is geometric, i.e.
672  // x_i = x^ + a d^i
673 
674  double step0 = x1-x0, step1 = x2-x1, abs1 = fabs(step1);
675  // Protect against roundoff/noise in inner solver
676  if ( abs1 > 1000.0*((double)DBL_EPSILON)*x2 )
677  {
678  double denom = fabs(step1-step0);
679  double sgn = (step1*step0 > 0)? 1.0 : -1.0;
680  // Greatest acceleration allowed is MAXACC*latest step length
681  // Can we do better than static tuning this parameter?
682  const double MAXACC=100.0;
683  double extfac = 1.0/(denom/abs1 + 1.0/MAXACC);
684  double extstep = sgn*extfac*step1;
685  //extstep = sgn*MIN2(extfac*step1,0.01*x2);
686  double predict = x2+extstep;
687  if (predict > 0.0)
688  xIonNew[ion] = predict;
689  if ( 0 )
690  if ( nelem == ipHYDROGEN || (nelem == ipOXYGEN && ion <=2 ) )
691  fprintf(ioQQQ,"Extrap %3ld %3ld %13.6g %13.6g %13.6g %13.6g %13.6g %13.6g\n",
692  nelem,ion,
693  x0,x0-xIonDense0[nconv-3][nelem][ion],x1-x0,x2-x1,extstep,predict);
694  }
695  tot1 += xIonNew[ion];
696  }
697  if ( tot1 > SMALLFLOAT )
698  {
699  double scal = tot0/tot1;
700  for (long ion = 0; ion <= nelem+1; ++ion)
701  {
702  dense.xIonDense[nelem][ion] = scal*xIonNew[ion];
703  }
704  for ( long ipISO=ipH_LIKE; ipISO<NISO; ++ipISO )
705  {
706  long ion = nelem-ipISO;
707  if (ion < 0 || ion < dense.IonLow[nelem] || ion > dense.IonHigh[nelem])
708  continue;
709  double thiserror;
710  iso_solve( ipISO, nelem, thiserror );
711  double renorm;
712  iso_renorm(nelem, ipISO, renorm);
713 
714  if ( dense.xIonDense[nelem][ion] > SMALLFLOAT &&
715  fabs(renorm - 1.0) > conv.IonizErrorAllowed)
716  {
717  // fprintf(ioQQQ,"Renorm failed nelem %ld iso %ld renorm %g\n",nelem,ipISO,renorm);
718  conv.setConvIonizFail( "iso/ion accel", 0.0, renorm);
719  }
720  }
721  }
722  }
723 
724  bool lgPostExtrapSolve = true;
725  if (lgPostExtrapSolve)
726  {
727  for (long nelem = ipHYDROGEN; nelem < LIMELM; ++nelem)
728  {
729  if (!dense.lgElmtOn[nelem])
730  continue;
731  ion_wrapper( nelem );
732  }
733  }
734  }
735 
736  // If solver appears to be struggling, take smaller steps where
737  // there is evidence of oscillations
738  if (loop_ion > 3)
739  {
740  for (long nelem = ipHYDROGEN; nelem < LIMELM; ++nelem)
741  {
742  if (!dense.lgElmtOn[nelem] || !lgOscill[nelem])
743  continue;
744  for (long ion = 0; ion <= nelem+1; ++ion)
745  {
746  dense.xIonDense[nelem][ion] = 0.5*(xIonDense0[0][nelem][ion]+dense.xIonDense[nelem][ion]);
747  }
748  for ( long ipISO=ipH_LIKE; ipISO<NISO; ++ipISO )
749  {
750  long ion = nelem-ipISO;
751  if (ion < 0 || ion < dense.IonLow[nelem] || ion > dense.IonHigh[nelem])
752  continue;
753  double thiserror;
754  iso_solve( ipISO, nelem, thiserror );
755  double renorm;
756  iso_renorm(nelem, ipISO, renorm);
757  }
758  }
759  }
760 
763 
764  /*>>chng 04 may 09, add option to abort here, inspired by H2 pop failures
765  * which cascaded downstream if we did not abort */
766  /* return if one of above solvers has declared abort condition */
767  if( lgAbort )
768  {
769  ++conv.nTotalIoniz;
770  return 1;
771  }
772 
774 
775  /* drive chemistry network*/
776  mole_drive();
777 
778  if (0)
779  {
780  fprintf(ioQQQ,"DEBUG H2 %g\n",findspecieslocal("H2")->den);
781  vector<const molecule*>debug_list;
782  debug_list.push_back(findspecies("H2"));
783  mole_dominant_rates( debug_list, ioQQQ,true,3,0.);
784  }
785 
786  if (loop_ion != 0)
787  {
788  for (long nelem = ipHYDROGEN; nelem < LIMELM; ++nelem)
789  {
790  if (!dense.lgElmtOn[nelem] || lgOscill[nelem])
791  continue;
792  for (long ion = 0; ion <= nelem+1; ++ion)
793  {
794  if ((dense.xIonDense[nelem][ion]-xIonDense0[0][nelem][ion]) *
795  (xIonDense0[0][nelem][ion]-xIonDense_prev[nelem][ion]) < 0)
796  {
797  lgOscill[nelem] = true;
798  break;
799  }
800  }
801  }
802  }
803 
805 
806  /* all elements have now had ionization reevaluated - in some cases we may have upset
807  * the ionization that was forced with an "element ionization" command - here we will
808  * re-impose that set ionization */
809  /* >>chng 04 oct 13, add this logic */
811 
812  /* redo ct ion rate for reasons now unclear */
813  ChargTranEval();
814 
815  /* evaluate molecular hydrogen level populations */
816  for( diatom_iter diatom = diatoms.begin(); diatom != diatoms.end(); ++diatom )
817  {
818  bool lgPopsConverged = true;
819  double old_val, new_val;
820  (*diatom)->H2_LevelPops( lgPopsConverged, old_val, new_val );
821  if( !lgPopsConverged )
822  {
823  conv.setConvIonizFail( "H2 pops", old_val, new_val);
824  }
825  }
826 
828 
829  /* lgIonizConverg is a function to check whether ionization has converged
830  * check whether ionization changed by more than relative error
831  * given by second number */
832  /* >>chng 04 feb 14, loop over all elements rather than just a few */
833  lgIonizConverg(loop_ion, conv.IonizErrorAllowed , false );
834 
835  if( deut.lgElmtOn )
836  {
837  static double OldDeut[2] = {0., 0.};
838  for( long ion=0; ion<2; ++ion )
839  {
840  if( !conv.nTotalIoniz && !loop_ion )
841  OldDeut[ion] = deut.xIonDense[ion];
842  if( deut.xIonDense[ion] > MIN_CONV_FRAC*deut.gas_phase &&
843  fabs(deut.xIonDense[ion] - OldDeut[ion] ) > 0.2*conv.IonizErrorAllowed*fabs(OldDeut[ion]) )
844  {
845  conv.setConvIonizFail( "D ion" , OldDeut[ion], deut.xIonDense[ion]);
846  }
847  OldDeut[ion] = deut.xIonDense[ion];
848  }
849  }
850 
851  if( fabs(EdenTrue_old - dense.EdenTrue) > conv.EdenErrorAllowed/2.*fabs(dense.EdenTrue) )
852  {
853  conv.setConvIonizFail( "EdTr cng" , EdenTrue_old, dense.EdenTrue);
854  }
855 
856  for( long ipISO=ipH_LIKE; ipISO<NISO; ++ipISO )
857  {
858  for( long nelem=ipISO; nelem<LIMELM; ++nelem )
859  {
860  if (!dense.lgElmtOn[nelem])
861  continue;
862  lgStatesConserved( nelem, nelem-ipISO, iso_sp[ipISO][nelem].st, iso_sp[ipISO][nelem].numLevels_local, 1e-3f, loop_ion );
863  }
864  }
865 
866  {
867  long iworst = -1;
868  double xworst = 0.0;
869  /* check whether molecular abundances are stable */
870  for( long i=0; i < mole_global.num_calc; ++i )
871  {
872  // Should probably make floor more specific to the species
873  double xval = fabs(mole.species[i].den-mole_save[i])/MAX2(mole_save[i],1e-10*dense.xNucleiTotal);
874  if( xval > xworst )
875  {
876  xworst = xval;
877  iworst = i;
878  }
879  }
880  if (iworst != -1 && xworst > conv.EdenErrorAllowed/2.)
881  {
882  string chConvIoniz = "change " + mole_global.list[iworst]->label;
883  conv.setConvIonizFail( chConvIoniz.c_str(),
884  mole_save[iworst],
885  mole.species[iworst].den);
886  }
887  }
888 
889  if ( loop_ion == 0 && lgIonizWidened )
890  {
891  for( long nelem=ipHELIUM; nelem<LIMELM; ++nelem )
892  {
893  if( dense.lgElmtOn[nelem] )
894  {
895  ion_trim2(nelem);
896  }
897  }
898  }
899 
900 
901  if( trace.nTrConvg>=4 )
902  {
903  /* trace ionization */
904  fprintf( ioQQQ,
905  " ConvBase4 ionization driver loop_ion %li converged? %c reason not converged %s\n" ,
906  loop_ion ,
907  TorF( conv.lgConvIoniz()) ,
908  conv.chConvIoniz() );
909  }
910 
911  ++loop_ion;
912  }
913  /* loop is not converged, less than loop limit, and we are reevaluating */
914  while( !conv.lgConvIoniz() && loop_ion < LOOP_ION_LIMIT && rfield.lgIonizReevaluate);
915 
916  if( conv.lgConvIoniz() )
918 
919  /* >>chng 05 oct 29, move CT heating here from heat_sum since this sometimes contributes
920  * cooling rather than heat and so needs to be sorted out before either heating or cooling
921  * are derived first find net heating - */
923  /* net is cooling if negative */
926 
927  /* evaluate current opacities
928  * rfield.lgOpacityReevaluate normally true,
929  * set false with no opacity reevaluate command, an option to only
930  * evaluate opacity one time per zone */
932  OpacityAddTotal();
933 
934  /* >>chng 02 jun 11, call even first time that this routine is called -
935  * this seemed to help convergence */
936 
937  /* do OTS rates for all lines and all continua since
938  * we now have ionization balance of all species. Note that this is not
939  * entirely self-consistent, since destruction probabilities here are not the same as
940  * the ones used in the model atoms. Problems?? if near convergence
941  * then should be nearly identical */
943  conv.lgIonStageTrimed || conv.lgSearch || nzone!=nzoneOTS )
944  {
945  RT_OTS();
946  nzoneOTS = nzone;
947 
948  /* remember old ots rates */
949  OldSumOTS[0] = OldSumOTS[1];
950  OldSumOTS[1] = SumOTS;
951  /*fprintf(ioQQQ," calling RT_OTS zone %.2f SumOTS is %.2e\n",fnzone,SumOTS);*/
952 
953  /* now update several components of the continuum, this only happens after
954  * we have gone through entire solution for this zone at least one time.
955  * there can be wild ots oscillation on first call */
956  /* the rel change of 0.2 was chosen by running hizqso - smaller increased
957  * itrzn but larger did not further decrease it. */
958  RT_OTS_Update(&SumOTS);
959  /*fprintf(ioQQQ,"RT_OTS_Updateee\t%.3f\t%.2e\t%.2e\n", fnzone,SumOTS , OldSumOTS[1] );*/
960  }
961  else
962  SumOTS = 0.;
963 
964  /* now check whether the ots rates changed */
965  if( SumOTS> 0. )
966  {
967  /* the ots rate must be converged to the error in the electron density */
968  /* how fine should this be converged?? originally had above, 10%, but take
969  * smaller ratio?? */
970  if( fabs(1.-OldSumOTS[1]/SumOTS) > conv.EdenErrorAllowed )
971  {
972  /* this branch, ionization not converged due to large change in ots rates.
973  * check whether ots rates are oscillating, if loopi > 1 so we have enough info*/
974  if( loopi > 1 )
975  {
976  /* here we have three values, are they changing sign? */
977  if( (OldSumOTS[0]-OldSumOTS[1]) * ( OldSumOTS[1] - SumOTS ) < 0. )
978  {
979  /* ots rates are oscillating */
980  conv.lgOscilOTS = true;
981  }
982  }
983 
984  conv.setConvIonizFail( "OTSRatChng" , OldSumOTS[1], SumOTS);
985  }
986 
987  /* produce info on the ots fields if either "trace ots" or
988  * "trace convergence xxx ots " was entered */
989  if( ( trace.lgTrace && trace.lgOTSBug ) ||
990  ( trace.nTrConvg && trace.lgOTSBug ) )
991  {
992  RT_OTS_PrtRate(SumOTS*0.05 , 'b' );
993  }
994  /*fprintf(ioQQQ,"DEBUG opac\t%.2f\t%.3e\t%.3e\n",fnzone,
995  dense.xIonDense[ipNICKEL][0] ,
996  dense.xIonDense[ipZINC][0] );*/
997  {
998  /* DEBUG OTS rates - turn this on to debug line, continuum or both rates */
999  enum {DEBUG_LOC=false};
1000  if( DEBUG_LOC && (nzone>110) )
1001  {
1002 # if 0
1003 # include "lines_service.h"
1004  DumpLine( &iso_sp[ipH_LIKE][ipHYDROGEN].trans(2,0) );
1005 # endif
1006  /* last par 'l' for line, 'c' for continua, 'b' for both,
1007  * the numbers printed are:
1008  * cell i, [i], so 1 less than ipoint
1009  * anu[i],
1010  * otslin[i],
1011  * opacity_abs[i],
1012  * otslin[i]*opacity_abs[i],
1013  * rfield.chLineLabel[i] ,
1014  * rfield.line_count[i] */
1015  }
1016  }
1017  }
1018 
1019  /* option to print OTS continuum with TRACE OTS */
1020  if( trace.lgTrace && trace.lgOTSBug )
1021  {
1022  /* find ots rates here, so we only print fraction of it,
1023  * SumOTS is both line and continuum contributing to ots, and is multiplied by opacity */
1024  /* number is weakest rate to print */
1025  RT_OTS_PrtRate(SumOTS*0.05 , 'b' );
1026  }
1027 
1028  // all RT routines called
1029  realnum rt_err = 0.0;
1030  bool lgCheckEsc = false;
1031  RT_line_all_escape( lgCheckEsc ? &rt_err : NULL );
1032 
1033  if ( rt_err > 0.1 )
1034  {
1035  conv.setConvIonizFail( "escp change", rt_err, 0.);
1036  }
1037 
1038  /* get total heating rate - first save old quantities to check how much it changes */
1039  HeatOld = thermal.htot;
1040 
1041  /* HeatSum will update ElecFrac,
1042  * secondary electron ionization and excitation efficiencies,
1043  * and sum up current secondary rates - remember old values before we enter */
1044  SecondOld = secondaries.csupra[ipHYDROGEN][0];
1045 
1046  if (lgConvBaseHeatTest)
1047  {
1048  /* get total cooling, thermal.ctot = does not occur since passes as pointer. This can add heat.
1049  * it calls coolSum at end to sum up the total cooling */
1051 
1052 
1053  /* update the total heating - it was all set to zero in HeatZero at top of this routine,
1054  * occurs before secondaries bit below, since this updates electron fracs */
1055  HeatSum();
1056  }
1057 
1058  /* test whether we can just set the rate to the new one, or whether we should
1059  * take average and do it again. secondaries.sec2total was set in hydrogenic, and
1060  * is ratio of secondary to total hydrogen destruction rates */
1061  /* >>chng 02 nov 20, add test on size of csupra - primal had very close to underflow */
1062  if( (secondaries.csupra[ipHYDROGEN][0]>SMALLFLOAT) && (secondaries.sec2total>0.001) &&
1063  fabs( 1. - SecondOld/SDIV(secondaries.csupra[ipHYDROGEN][0]) ) > 0.1 &&
1064  SecondOld > 0. && secondaries.csupra[ipHYDROGEN][0] > 0.)
1065  {
1066  /* say that ionization has not converged, secondaries changed by too much */
1067  conv.setConvIonizFail( "SecIonRate", SecondOld,
1069  }
1070 
1071 # if 0
1072  static realnum hminus_old=0.;
1073  /* >>chng 04 apr 15, add this convergence test */
1074  if( conv.nTotalIoniz )
1075  {
1076  realnum hminus_den = findspecieslocal("H-")->den;
1077  if( fabs( hminus_old-hminus_den ) > fabs( hminus_den * conv.EdenErrorAllowed ) )
1078  {
1079  conv.setConvIonizFail( "Big H- chn", hminus_old, hminus_den );
1080  }
1081  hminus_old = hminus_den;
1082  }
1083 # endif
1084 
1085  if( lgConvBaseHeatTest && HeatOld > 0. && !thermal.lgTemperatureConstant )
1086  {
1087  /* check if heating has converged - tolerance is final match */
1088  if( fabs(1.-thermal.htot/HeatOld) > conv.HeatCoolRelErrorAllowed*.5 )
1089  {
1090  conv.setConvIonizFail( "Big d Heat", HeatOld, thermal.htot);
1091  }
1092  }
1093 
1094  /* abort flag may have already been set - if so bail */
1095  if( lgAbort )
1096  {
1097 
1098  return 1;
1099  }
1100 
1101  /* >>chng 01 mar 16, evaluate pressure here since changing and other values needed */
1102  /* reevaluate pressure */
1103  /* this sets values of pressure.PresTotlCurr, also calls tfidle */
1104  PresTotCurrent();
1105 
1107 
1108  /* update some counters that keep track of how many times this routine
1109  * has been called */
1110 
1111  if( trace.lgTrace )
1112  {
1113  fprintf( ioQQQ,
1114  " ConvBase return. fnzone %.2f nPres2Ioniz %li Te:%.3e HI:%.3e HII:%.3e H2:%.3e Ne:%.3e htot:%.3e CSUP:%.2e Conv?%c reason:%s\n",
1115  fnzone,
1116  conv.nPres2Ioniz,
1117  phycon.te,
1118  dense.xIonDense[ipHYDROGEN][0],
1119  dense.xIonDense[ipHYDROGEN][1],
1120  hmi.H2_total,
1121  dense.eden,
1122  thermal.htot,
1124  TorF(conv.lgConvIoniz()) ,
1125  conv.chConvIoniz());
1126  }
1127 
1128  /* this counts number of times we are called by ConvPresTempEdenIoniz,
1129  * number of calls in this zone so first call is zero
1130  * reset to zero each time ConvPresTempEdenIoniz is called */
1131  ++conv.nPres2Ioniz;
1132 
1133  /* this is abort option set with SET PRES IONIZ command,
1134  * test on nzone since init can take many iterations
1135  * this is seldom used except in special cases */
1136  if( conv.nPres2Ioniz > conv.limPres2Ioniz && nzone > 0)
1137  {
1138  fprintf(ioQQQ," PROBLEM ConvBase sets lgAbort since nPres2Ioniz exceeds limPres2Ioniz. ");
1139  fprintf(ioQQQ,"Their values are %li and %li. ",conv.nPres2Ioniz , conv.limPres2Ioniz);
1140  fprintf(ioQQQ,"Search phase? %c T:%.3e\n",TorF(conv.lgSearch) , phycon.te );
1141  fprintf(ioQQQ," Reset this limit with the SET PRES IONIZ command.\n");
1142  lgAbort = true;
1143  return 1;
1144  }
1145 
1146  /* various checks on the convergence of the current solution */
1147  if( eden_sum() )
1148  {
1149  /* non-zero return indicates abort condition */
1150  return 1;
1151  }
1152 
1153  /* is electron density converged? */
1155  if( fabs(EdenTrue_old - dense.EdenTrue) > fabs(dense.EdenTrue * conv.EdenErrorAllowed/2.) )
1156  {
1157  conv.setConvIonizFail( "eden chng", EdenTrue_old, dense.EdenTrue);
1158  }
1159 
1160  /* check on molecular electron den */
1161  if( fabs(EdenFromMolecOld - mole.elec) > fabs(dense.EdenTrue * conv.EdenErrorAllowed/2.) )
1162  {
1163  conv.setConvIonizFail( "edn chnCO", EdenFromMolecOld, dense.EdenTrue);
1164  }
1165 
1166  if( gv.lgGrainElectrons )
1167  {
1168  /* check on grain electron den */
1169  if( fabs(EdenFromGrainsOld - gv.TotalEden) > fabs(dense.EdenTrue * conv.EdenErrorAllowed/2.) )
1170  {
1171  conv.setConvIonizFail( "edn grn e", EdenFromGrainsOld, gv.TotalEden);
1172  }
1173 
1174  /* check on sum of grain and molecular electron den - often two large numbers that cancel */
1175  if( fabs( (EdenFromMolecOld-EdenFromGrainsOld) - (mole.elec-gv.TotalEden) ) >
1176  fabs(dense.EdenTrue * conv.EdenErrorAllowed/4.) )
1177  {
1178  conv.setConvIonizFail( "edn mole-grn",
1179  (EdenFromMolecOld-EdenFromGrainsOld),
1180  (mole.elec-gv.TotalEden));
1181  }
1182  }
1183 
1184  /* check on heating and cooling if vary temp model
1185  * >>chng 08 jul 01, over the code's entire history it had tested whether
1186  * this is a constant temperature simulation and did not do this test if
1187  * the thermal solution was not done. There are some cases where we do
1188  * want to specify the temperature and then find the heating or cooling -
1189  * this is done in calculations of cooling curves for instance. With this
1190  * change the heating/cooling are converged even in a constant temperature
1191  * sim. this does make CT sims run more slowly but with greater accuracy
1192  * if heating or cooling is reported */
1194  {
1195  if( fabs(HeatingOld - thermal.htot)/thermal.htot > conv.HeatCoolRelErrorAllowed/2. )
1196  {
1197  conv.setConvIonizFail( "heat chg", HeatingOld, thermal.htot);
1198  }
1199 
1200  if( fabs(CoolingOld - thermal.ctot)/thermal.ctot > conv.HeatCoolRelErrorAllowed/2. )
1201  {
1202  conv.setConvIonizFail( "cool chg", CoolingOld, thermal.ctot);
1203  }
1204  if (0)
1205  {
1206  // Print heating and coolants at current iteration if required
1207  SaveHeat(ioQQQ);
1208  CoolSave(ioQQQ,"COOL");
1209  }
1210  }
1211 
1212  /* >>chng 05 mar 26, add this convergence test - important for molecular or advective
1213  * sims since iso ion solver must sync up with chemistry */
1214  /* remember current ground state population - will check if converged */
1215  for( long ipISO=ipH_LIKE; ipISO<NISO; ++ipISO )
1216  {
1217  for( long nelem=ipISO; nelem<LIMELM;++nelem )
1218  {
1219  if( dense.lgElmtOn[nelem] )
1220  {
1221  /* only do this check for "significant" levels of ionization */
1222  /*lint -e644 var possibly not init */
1223  if( dense.xIonDense[nelem][nelem-ipISO]/dense.gas_phase[nelem] > 1e-5 )
1224  {
1225  if( fabs(iso_sp[ipISO][nelem].st[0].Pop()-save_iso_grnd[ipISO][nelem])/SDIV(iso_sp[ipISO][nelem].st[0].Pop())-1. >
1227  {
1228  char chConvIoniz[INPUT_LINE_LENGTH];
1229  sprintf( chConvIoniz,"iso %2li %2li",ipISO, nelem );
1230  conv.setConvIonizFail(chConvIoniz,
1231  save_iso_grnd[ipISO][nelem],
1232  iso_sp[ipISO][nelem].st[0].Pop());
1233  }
1234  }
1235  /*lint +e644 var possibly not init */
1236  }
1237  }
1238  }
1239 
1240  /* this counts how many times ConvBase has been called in this iteration,
1241  * located here at bottom of routine so that number is false on first
1242  * call, set to 0 in when iteration starts - used to create itr/zn
1243  * number in printout often used to tell whether this is the very
1244  * first attempt at solution in this iteration */
1245  ++conv.nTotalIoniz;
1246 
1247  /* first sweep is over if this is not search phase */
1248  if( !conv.lgSearch )
1249  conv.lgFirstSweepThisZone = false;
1250 
1251  /* this was set with STOP NTOTALIONIZ command - only a debugging aid
1252  * by default is zero so false */
1254  {
1255  /* non-zero return indicates abort condition */
1256  fprintf(ioQQQ , "ABORT flag set since STOP nTotalIoniz was set and reached.\n");
1257  return 1;
1258  }
1259 
1261  {
1263  {
1264  static int iter_punch=-1;
1265  if( iteration !=iter_punch )
1267  iter_punch = iteration;
1268  }
1269 
1271  "%li\t%.4e\t%.4e\t%.4e\n",
1273  }
1274 
1275  return 0;
1276 }
1277 
1278 void UpdateUTAs( void )
1279 {
1280  DEBUG_ENTRY( "UpdateUTAs()" );
1281 
1282  /* only reevaluate this on first pass through on this zone */
1284  {
1285  /* inner shell ionization */
1286  for( long nelem=0; nelem< LIMELM; ++nelem )
1287  {
1288  for( long ion=0; ion<nelem+1; ++ion )
1289  {
1290  ionbal.UTA_ionize_rate[nelem][ion] = 0.;
1291  ionbal.UTA_heat_rate[nelem][ion] = 0.;
1292  }
1293  }
1294  /* inner shell ionization by UTA lines */
1295  /* this flag turned off with no UTA command */
1296  for( size_t i=0; i < UTALines.size(); ++i )
1297  {
1298  /* rateone is inverse lifetime of level against autoionization */
1299  double rateone = UTALines[i].Emis().pump() * UTALines[i].Emis().AutoIonizFrac();
1300  ionbal.UTA_ionize_rate[(*UTALines[i].Hi()).nelem()-1][(*UTALines[i].Hi()).IonStg()-1] += rateone;
1301  /* heating rate in erg atom-1 s-1 */
1302  ionbal.UTA_heat_rate[(*UTALines[i].Hi()).nelem()-1][(*UTALines[i].Hi()).IonStg()-1] += rateone*UTALines[i].Coll().heat();
1303  {
1304  /* DEBUG OTS rates - turn this on to debug line, continuum or both rates */
1305  /*@-redef@*/
1306  enum {DEBUG_LOC=false};
1307  /*@+redef@*/
1308  if( DEBUG_LOC /*&& UTALines[i].nelem==ipIRON+1 && (UTALines[i].IonStg==15||UTALines[i].IonStg==14)*/ )
1309  {
1310  fprintf(ioQQQ,"DEBUG UTA %3i %3i %.3f %.2e %.2e %.2e\n",
1311  (*UTALines[i].Hi()).nelem() , (*UTALines[i].Hi()).IonStg() , UTALines[i].WLAng() ,
1312  rateone, UTALines[i].Coll().heat(),
1313  UTALines[i].Coll().heat()*dense.xIonDense[(*UTALines[i].Hi()).nelem()-1][(*UTALines[i].Hi()).IonStg()-1] );
1314  }
1315  }
1316  }
1317  }
1318 
1319  return;
1320 }
1321 
1323 {
1324  DEBUG_ENTRY( "lgNetEdenSrcSmall()" );
1325 
1326  fixit("lgNetEdenSrcSmall needs to be enabled");
1327  return true;
1328 
1329  if( conv.lgSearch )
1330  return true;
1331  fixit("grain rates not well tested below");
1332  if( gv.lgDustOn() )
1333  return true;
1334 
1335  // Check for consistency of explicit source and sink rates for
1336  // electrons with population derived from neutrality
1337  double ionsrc = 0., ionsnk = 0.;
1338  for( long nelem=0; nelem < LIMELM; ++nelem )
1339  {
1340  if( !dense.lgElmtOn[nelem] )
1341  continue;
1342  ionsrc += ionbal.elecsrc[nelem];
1343  ionsnk += ionbal.elecsnk[nelem];
1344  for( long ion_from = 0; ion_from <= nelem + 1; ++ion_from )
1345  {
1346  for( long ion_to = 0; ion_to <= nelem + 1; ++ion_to )
1347  {
1348  if( ion_to-ion_from > 0 )
1349  {
1350  ionsrc += gv.GrainChTrRate[nelem][ion_from][ion_to] *
1351  dense.xIonDense[nelem][ion_from] * (ion_to-ion_from);
1352  }
1353  else if( ion_to-ion_from < 0 )
1354  {
1355  ionsnk += gv.GrainChTrRate[nelem][ion_from][ion_to] *
1356  dense.xIonDense[nelem][ion_from] * (ion_from-ion_to);
1357  }
1358  }
1359  }
1360  }
1361  long ipMElec = findspecies("e-")->index;
1362  const double totsrc = ionsrc + mole.species[ipMElec].src;
1363  const double totsnk = ionsnk + mole.species[ipMElec].snk*dense.EdenTrue;
1364  const double diff = (totsrc - totsnk);
1365  const double ave = ( fabs(totsrc) + fabs(totsnk) )/2.;
1366  fixit("Need to tighten up e- population convergence criterion");
1367  const double error_allowed = 0.05 * ave; //conv.EdenErrorAllowed * ave;
1368  if( fabs(diff) > error_allowed )
1369  {
1370  enum {DEBUG_LOC=false};
1371  if( DEBUG_LOC )
1372  {
1373  fprintf(ioQQQ,"PROBLEM large NetEdenSrc nzone %li\t%e\t%e\t%e\t%e\n",
1374  nzone,
1375  totsrc/SDIV(totsnk),
1376  dense.EdenTrue,
1377  ionsrc + mole.species[ipMElec].src,
1378  ionsnk + mole.species[ipMElec].snk*dense.EdenTrue);
1379  }
1380  conv.setConvIonizFail( "NetEdenSrc", diff, error_allowed);
1381  return false;
1382  }
1383  else
1384  return true;
1385 }
void ion_trim(long int nelem)
Definition: ion_trim.cpp:176
t_mole_global mole_global
Definition: mole.cpp:7
void RT_OTS_Update(double *SumOTS)
Definition: rt_ots.cpp:476
double htot
Definition: thermal.h:169
void DumpLine(const TransitionProxy &t)
Definition: transition.cpp:138
t_atmdat atmdat
Definition: atmdat.cpp:6
t_thermal thermal
Definition: thermal.cpp:6
void SaveHeat(FILE *io)
Definition: heat_save.cpp:22
static double x2[63]
void CoolSave(FILE *io, const char chJob[])
Definition: cool_save.cpp:20
size_t size(void) const
Definition: transition.h:331
bool lgElmtOn
Definition: deuterium.h:21
qList st
Definition: iso.h:482
double EdenErrorAllowed
Definition: conv.h:265
TransitionList UTALines("UTALines",&AnonStates)
double ** UTA_heat_rate
Definition: ionbal.h:173
bool lgGrainElectrons
Definition: grainvar.h:500
double te_update
Definition: thermal.h:148
t_opac opac
Definition: opacity.cpp:5
int num_calc
Definition: mole.h:362
static double x1[83]
void mole_drive(void)
Definition: mole_drive.cpp:29
int ConvBase(long loopi)
Definition: conv_base.cpp:188
const realnum SMALLFLOAT
Definition: cpu.h:246
realnum xNucleiTotal
Definition: dense.h:114
bool lgFirstSweepThisZone
Definition: conv.h:148
const int NISO
Definition: cddefines.h:310
ConvergenceCounter register_(const string name)
Definition: conv.cpp:88
void ion_widen(long nelem)
Definition: ion_trim.cpp:864
long int IonHigh[LIMELM+1]
Definition: dense.h:130
char TorF(bool l)
Definition: cddefines.h:753
void RT_OTS_PrtRate(double weak, int chFlag)
Definition: rt_ots.cpp:687
const int ipOXYGEN
Definition: cddefines.h:355
int eden_sum(void)
Definition: eden_sum.cpp:17
bool lgAdvection
Definition: dynamics.h:66
bool lgTimeDependentStatic
Definition: dynamics.h:102
t_StopCalc StopCalc
Definition: stopcalc.cpp:7
t_conv conv
Definition: conv.cpp:5
STATIC bool lgNetEdenSrcSmall(void)
Definition: conv_base.cpp:1322
double elecsrc[LIMELM]
Definition: ionbal.h:245
double ChargTranSumHeat(void)
long int limPres2Ioniz
Definition: conv.h:154
t_phycon phycon
Definition: phycon.cpp:6
void CoolEvaluate(double *tot)
Definition: cool_eval.cpp:58
void resetConvIoniz()
Definition: conv.h:93
char chHashString[INPUT_LINE_LENGTH]
Definition: save.h:398
bool lgIonizReevaluate
Definition: rfield.h:112
double char_tran_heat
Definition: thermal.h:166
bool lgTemperatureConstant
Definition: thermal.h:44
double char_tran_cool
Definition: thermal.h:166
FILE * ioQQQ
Definition: cddefines.cpp:7
molezone * findspecieslocal(const char buf[])
long int nzone
Definition: cddefines.cpp:14
bool lgConvIoniz() const
Definition: conv.h:108
t_dynamics dynamics
Definition: dynamics.cpp:42
#define MIN2(a, b)
Definition: cddefines.h:807
void RT_line_all_escape(realnum *error)
Definition: rt_line_all.cpp:21
void lgStatesConserved(long nelem, long ionStage, qList states, long numStates, realnum err_tol, long loop_ion)
Definition: dense.cpp:195
bool lgOscilOTS
Definition: conv.h:188
void mole_save(FILE *punit, const char speciesname[], const char args[], bool lgHeader, bool lgData, bool lgCoef, double depth)
t_dense dense
Definition: global.cpp:15
t_elementnames elementnames
Definition: elementnames.cpp:5
t_iso_sp iso_sp[NISO][LIMELM]
Definition: iso.cpp:11
void PresTotCurrent(void)
bool lgIonStageTrimed
Definition: conv.h:184
STATIC void ion_trim_from_set(long nelem)
Definition: ion_trim.cpp:103
double xIonDense[LIMELM][LIMELM+1]
Definition: dense.h:135
void UpdateUTAs(void)
Definition: conv_base.cpp:1278
long int iteration
Definition: cddefines.cpp:16
bool lgSearch
Definition: conv.h:168
t_trace trace
Definition: trace.cpp:5
bool fp_equal(sys_float x, sys_float y, int n=3)
Definition: cddefines.h:858
void iso_multiplet_opacities(void)
Definition: iso_level.cpp:756
void iso_collapsed_update(void)
Definition: iso_solve.cpp:26
const char * chConvIoniz() const
Definition: conv.h:112
t_ionbal ionbal
Definition: ionbal.cpp:8
double TotalEden
Definition: grainvar.h:531
FILE * ipTraceConvergeBase
Definition: save.h:454
long int IonLow[LIMELM+1]
Definition: dense.h:129
void ion_trim2(long int nelem)
Definition: ion_trim.cpp:612
long int nPres2Ioniz
Definition: conv.h:145
realnum sec2total
Definition: secondaries.h:39
#define STATIC
Definition: cddefines.h:118
bool lgTrace
Definition: trace.h:12
t_mole_local mole
Definition: mole.cpp:8
molecule * findspecies(const char buf[])
double EdenTrue
Definition: dense.h:232
t_rfield rfield
Definition: rfield.cpp:9
void DynaIonize(void)
Definition: dynamics.cpp:160
double ** UTA_ionize_rate
Definition: ionbal.h:171
void OpacityAddTotal(void)
static double x0[83]
float realnum
Definition: cddefines.h:124
valarray< class molezone > species
Definition: mole.h:468
#define EXIT_FAILURE
Definition: cddefines.h:168
realnum IonizErrorAllowed
Definition: conv.h:278
const int INPUT_LINE_LENGTH
Definition: cddefines.h:301
vector< diatomics * > diatoms
Definition: h2.cpp:8
void ion_recom_calculate(void)
#define cdEXIT(FAIL)
Definition: cddefines.h:484
int index
Definition: mole.h:194
void iso_solve(long ipISO, long nelem, double &maxerr)
Definition: iso_solve.cpp:99
int nTrConvg
Definition: trace.h:27
realnum HeatCoolRelErrorAllowed
Definition: conv.h:276
realnum gas_phase
Definition: deuterium.h:22
long int nTotalIonizStop
Definition: stopcalc.h:127
void iso_update_rates(void)
Definition: iso_solve.cpp:48
double xIonDense[2]
Definition: deuterium.h:23
long int nTotalIoniz
Definition: conv.h:159
bool lgElmtOn[LIMELM]
Definition: dense.h:160
char chElementSym[LIMELM][CHARS_ELEMENT_SYM]
Definition: elementnames.h:25
realnum gas_phase[LIMELM]
Definition: dense.h:76
bool lgOpacStatic
Definition: opacity.h:152
void mole_update_sources(void)
Definition: mole_drive.cpp:55
void SetDeuteriumIonization(const double &xNeutral, const double &xIonized)
Definition: deuterium.cpp:37
#define ASSERT(exp)
Definition: cddefines.h:617
bool lgOpacityReevaluate
Definition: rfield.h:105
void mole_dominant_rates(const vector< const molecule * > &debug_list, FILE *ioOut, bool lgPrintReagents, size_t NPRINT, double fprint)
const int ipH_LIKE
Definition: iso.h:64
const int LIMELM
Definition: cddefines.h:307
double elecsnk[LIMELM]
Definition: ionbal.h:245
double den
Definition: mole.h:421
void SecIoniz(void)
Definition: heat_sum.cpp:100
#define DEBUG_ENTRY(funcname)
Definition: cddefines.h:729
const int ipHELIUM
Definition: cddefines.h:349
void GrainDrive()
Definition: grains.cpp:1597
bool lgTraceConvergeBase
Definition: save.h:452
bool lgRedoStatic
Definition: opacity.h:159
double H2_total
Definition: hmi.h:25
bool lgUpdateCouplings
Definition: conv.h:257
t_deuterium deut
Definition: deuterium.cpp:7
double eden
Definition: dense.h:201
void iso_renorm(long nelem, long ipISO, double &renorm)
Definition: iso_solve.cpp:310
double elec
Definition: mole.h:460
#define MAX2(a, b)
Definition: cddefines.h:828
double RateIonizTot(long nelem, long ion) const
Definition: ionbal.cpp:223
int fprintf(const Output &stream, const char *format,...)
Definition: service.cpp:1217
MoleculeList list
Definition: mole.h:365
void ion_wrapper(long nelem)
realnum ** csupra
Definition: secondaries.h:33
bool lgOTSBug
Definition: trace.h:100
sys_float SDIV(sys_float x)
Definition: cddefines.h:1006
void setConvIonizFail(const char *reason, double oldval, double newval)
Definition: conv.h:100
void HeatZero(void)
Definition: heat_sum.cpp:1071
#define fixit(a)
Definition: cddefines.h:416
bool lgElemsConserved(void)
Definition: dense.cpp:119
GrainVar gv
Definition: grainvar.cpp:5
t_secondaries secondaries
Definition: secondaries.cpp:5
t_hmi hmi
Definition: hmi.cpp:5
bool lgTraceConvergeBaseHash
Definition: save.h:452
double fnzone
Definition: cddefines.cpp:15
void ChargTranEval(void)
double CharExcIonOf[NCX][LIMELM][LIMELM+1]
Definition: atmdat.h:300
t_save save
Definition: save.cpp:5
double te
Definition: phycon.h:21
void highen(void)
Definition: highen.cpp:20
void ion_trim_validate(long nelem, bool lgIonizTrimCalled)
Definition: ion_trim.cpp:909
void atmdat_3body(void)
bool lgDustOn() const
Definition: grainvar.h:477
const int ipHYDROGEN
Definition: cddefines.h:348
realnum GrainChTrRate[LIMELM][LIMELM+1][LIMELM+1]
Definition: grainvar.h:543
double CharExcRecTo[NCX][LIMELM][LIMELM+1]
Definition: atmdat.h:300
void RT_OTS(void)
Definition: rt_ots.cpp:41
void HeatSum(void)
Definition: heat_sum.cpp:498
EmissionList & Emis()
Definition: transition.h:363
vector< diatomics * >::iterator diatom_iter
Definition: h2.h:13
bool lgAbort
Definition: cddefines.cpp:10
double ctot
Definition: thermal.h:130
const bool lgConvBaseHeatTest
Definition: cooling.h:7