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