Cloudy
Spectral Synthesis Code for Astrophysics
Loading...
Searching...
No Matches
h2_priv.h
Go to the documentation of this file.
1/* This file is part of Cloudy and is copyright (C)1978-2025 by Gary J. Ferland and
2 * others. For conditions of distribution and use see copyright notice in license.txt */
3
4#ifndef H2_PRIV_H_
5#define H2_PRIV_H_
6
7#include "transition.h"
8#include "container_classes.h"
9
10typedef void (*linefunc)(const TransitionProxy& t,
11 bool lgShield_this_zone,
12 realnum pestrk,
13 realnum DopplerWidth,
14 bool lgKeepLyman);
15class Parser;
16class molecule;
18
21const int N_X_COLLIDER = 5;
23const int chN_X_COLLIDER = 10;
24// colliders are (see h2.cpp diatoms_init)
25// 0 H
26// 1 He
27// 2 H2 ortho
28// 3 H2 para
29// 4 H+
30
32const int nTE_HMINUS = 7;
33
35const int N_ELEC = 7;
36
38const realnum H2_logte_hminus[nTE_HMINUS] = {1.,1.47712,2.,2.47712,3.,3.47712,4.};
39
41{
42 long n, v, j;
43};
44
46{
47public:
49 {
50 initial = a;
51 final = b;
52 energies.clear();
53 xsections.clear();
54 rate_coeff = 0.;
55 };
57 vector<double> energies;
58 vector<double> xsections;
59 double rate_coeff;
60};
61
63{
65 {
66 magic = 0;
67 filename = "";
68 };
69 long magic;
70 string filename;
71};
72
74{
75public:
76 double Abund() const
77 {
78 return *dense_total;
79 }
80 void GetIndices( long& ipHi, long& ipLo, const char* chLine, long& i ) const;
81 void CalcPhotoionizationRate(void);
82 double (*photoion_opacity_fun)( double energy );
83 long OpacityCreate( vector<double>& stack );
84 double GetHeatRate( const diss_tran& tran );
85 double GetDissociationRate( const diss_tran& tran );
86 double MolDissocOpacity( const diss_tran& tran, const double& Mol_Ene );
87 double Cont_Diss_Heat_Rate( void );
88 void Mol_Photo_Diss_Rates( void );
92 double GetExcitedElecDensity(void);
93 realnum GetXColden( long iVib, long iRot );
94
95 long int getLine( long iElecHi, long iVibHi, long iRotHi, long iElecLo, long iVibLo, long iRotLo, double *relint, double *absint );
96
97 /* compute rate coefficient for a single quenching collision */
98 realnum H2_CollidRateEvalOne( long iVibHi, long iRotHi, long iVibLo, long iRotLo, long ipHi, long ipLo, long nColl, double temp_K );
99
100 void H2_Calc_Average_Rates( void );
101
102 void H2_X_sink_and_source( void );
103 /*H2_X_coll_rate_evaluate find collisional rates within X -
104 * this is one time upon entry into H2_LevelPops */
105 void H2_X_coll_rate_evaluate( void );
106
107 /*H2_Level_low_matrix evaluate lower populations within X */
108 /* total abundance within matrix */
109 void H2_Level_low_matrix(realnum abundance );
110
111 /*H2_Read_Cosmicray_distribution read distribution function for H2 population following cosmic ray collisional excitation
112 void H2_Read_Cosmicray_distribution(void); */
113
114 // read energies for all electronic levels
115 void H2_ReadEnergies();
116 void H2_ReadEnergies( long int nelec, vector<int>& n, vector<int>& v, vector<int>&J, vector<double>& eWN );
117
121 void H2_ReadDissprob( long int nelec );
122
124 void H2_CollidRateEvalAll( void );
125
129 void H2_CollidRateRead( long int nColl );
130
134 void H2_ReadTransprob( long int nelec, TransitionList &trans );
135
138
140 void mole_H2_form( void );
141
143 void mole_H2_LTE( void );
144
147 void H2_Solomon_rate( void );
148
150 double gs_rate( void );
151
154 void H2_zero_pops_too_low( void );
155
157 void init(void);
158
160 void H2_ContPoint( void );
161
163 double H2_DR(void);
164
166 double H2_Accel(void);
167
169 void H2_RT_OTS( void );
170
172 double H2_RadPress(void);
173
175 void H2_LinesAdd(void);
176
178 void H2_Reset( void );
179
181 double H2_InterEnergy(void);
182
186 void H2_Colden( const char *chLabel );
187
189 void H2_Cooling(void);
190
193 double LTE_Cooling_per_H2();
194
200 FILE* ioPUN ,
201 bool lgDoAll );
202
208 void H2_PunchLineStuff( FILE * io , realnum xLimit , long index);
209
211 void H2_RT_diffuse(void);
212
215 void H2_RTMake( linefunc line_one );
216
218 void H2_RT_tau_inc(void);
219
221 void H2_Prt_Zone(void);
222
223 // print departure coefficients for all X levels
224 void H2_PrtDepartCoef(void);
225
227 void H2_LineZero( void );
228
230 void H2_RT_tau_reset( void );
231
233 void H2_LevelPops( bool &lgPopsConverged, double &old_value, double &new_value );
234
241 void H2_PunchDo( FILE* io , char chJOB[] , const char chTime[] , long int ipPun );
242
244 void H2_ParseSave( Parser &p, ostringstream& chHeader );
245
247 double H2_itrzn( void );
248
253 void H2_Prt_column_density( FILE *ioMEAN );
254
255 void set_numLevelsMatrix( long numLevels );
256
257 void H2_ReadDissocEnergies( void );
258
261
267
268 /* total spontaneous dissociation rate [s-1],
269 * summed over excited electronic states, weighted by pops */
271
274
278
283
288
292
296
297 double HeatDiss;
298 double HeatDexc;
302
304 double Average_A;
315
319
323
325
330
331 // single precision versions of the above
334
338
339 /* old and older ortho - para ratios, used to determine whether soln is converged */
341
342 // these remember the largest and smallest factors needed to
343 // renormalize the H2 chemistry
344 double renorm_max ,
346
347 // this will say how many times the large H2 molecule has been called in this zone -
348 // if not called (due to low H2 abundance) then not need to update its line arrays
350
351 // flag saying whether to bother with the large molecule at all,
352 // default is false, set true with atom h2 on command
354
355 // this is the number of electronic levels to include in the output - default is 1,
356 // only X. changed with PRINT LINES H2 ELECTRONIC and option on PUNCH H2 LINES commands
358
362
365
368
369 // include collision rates that come from real calculations,
370 // off with atom h2 collisions off command
372
374 bool lgLTE;
375
378
379 // which set of He - H2 collisions to use? default is ORNL, other
380 // is Le BOURlet
382
383 // flag saying whether (true) or not to use ORNL H2 - H2 collisions
386
391
393 long int nzoneEval;
394
397
400
402
405
411
412 // the number of electronic quantum states to include.
413 // To do both Lyman and Werner bands want nelec = 3
415
416 /* this is fraction of population that is within levels done with matrix */
418
419 /* used to recall the temperature used for last set of Boltzmann factors */
422
423 explicit diatomics( const string& a, const double& e_star, const double* const abund, double (*fun)(double) ) ;
424
430 vector< diss_tran > Diss_Trans;
431
432private:
433 string label;
435 string path;
436
439 /* >> chng 05 jul 15, TE, H2g = sum (v=0, J=0,1) */
440 /* >>chng 05 jul 29, to 0.5 eV, this goes up to J=8 for v=0 */
441 /* >>chng 05 aug 03, slight upward change in energy to include the J=8 level,
442 * also give energy in waveumbers for simplicity (save h2 levels give energy in ryd) */
443 /*#define ENERGY_H2_STAR (0.5/EVRYD/WAVNRYD)*/
444 /* energy of v=0, J=8 is 4051.73, J=9 is 5001.97
445 * v=1, J=0 is 4161.14 */
446public:
447 const double ENERGY_H2_STAR;
448
449private:
450 // pointer to the density of the species
451 const double* const dense_total;
452
454
455 /* these vars are private for H2 but uses same style as all other header files -
456 * the extern is extern in all except cddefines */
457
460
463
466
469
470 //int H2_nRot_add_ortho_para[N_ELEC];
473 long int nVib_hi[N_ELEC];
475 valarray<long> nRot_hi[N_ELEC];
478 long int Jlowest[N_ELEC];
485 vector<CollRateCoeffArray> RateCoefTable;
486
487 // quantities dealing with chemistry
488#if 1
489#endif
490
491 // these are quantities for each state
492#if 1
497
503#endif
504
506
507 // these are quantities for states with X
508#if 1
528
531#endif
532
533 valarray<realnum> H2_X_source;
534 valarray<realnum> H2_X_sink;
535
539
542
545
546 valarray<long> ipVib_H2_energy_sort;
547 valarray<long> ipElec_H2_energy_sort;
548 valarray<long> ipRot_H2_energy_sort;
551
561
562 long int levelAsEval;
564 long int nzone_eval;
566
569
573
575 long int nH2_pops;
576 long int nH2_zone;
577
580
585
586 /* LTE cooling per molecule */
587 vector<double> LTE_Temp, LTE_cool;
588
593
594public:
595 /* Read LTE cooling per molecule */
597
598 /* interpolate over LTE cooling array */
599 double interpolate_LTE_Cooling( double Temp );
600};
601
602/* compute H2 continuum dissociation cross sections */
603double MolDissocCrossSection( const diss_tran& tran, const double& Mol_Ene );
604
605double Yan_H2_CS( double energy_ryd /* photon energy in ryd */);
606
607/* compute H2 continuum dissoication opacities */
608//double MolDissocOpacity( const diss_tran& tran, const double& Mol_Ene );
609
610#endif /* H2_PRIV_H_ */
611
t_abund abund
Definition abund.cpp:5
float realnum
Definition cddefines.h:127
Definition atmdat.h:13
Definition parser.h:43
Definition transition.h:288
TransitionProxy::iterator iterator
Definition transition.h:293
Definition transition.h:24
long int iterationAsEval
Definition h2_priv.h:505
double Average_collH_deexcit
Definition h2_priv.h:307
void H2_PunchDo(FILE *io, char chJOB[], const char chTime[], long int ipPun)
Definition mole_h2_io.cpp:1143
long int n_elec_states
Definition h2_priv.h:414
bool lgImgMatrix
Definition h2_priv.h:592
long int nzone_eval
Definition h2_priv.h:564
double rel_pop_LTE_s
Definition h2_priv.h:291
multi_arr< realnum, 2 > H2_X_coll_rate
Definition h2_priv.h:468
char chH2ColliderLabels[N_X_COLLIDER][chN_X_COLLIDER]
Definition h2_priv.h:453
double para_density
Definition h2_priv.h:329
multi_arr< double, 2 > H2_col_rate_out
Definition h2_priv.h:511
double Solomon_elec_decay_s
Definition h2_priv.h:277
double Average_collH_dissoc_g
Definition h2_priv.h:311
double Cont_Dissoc_Rate_H2g
Definition h2_priv.h:286
multi_arr< bool, 3 > H2_lgOrtho
Definition h2_priv.h:502
double H2_renorm_chemistry
Definition h2_priv.h:465
double xSTDNoise
Definition h2_priv.h:396
double photodissoc_BigH2_H2s
Definition h2_priv.h:265
bool lgColl_gbar
Definition h2_priv.h:361
long int nzoneEval
Definition h2_priv.h:393
void Mol_Photo_Diss_Rates(void)
Definition mole_dissociate.cpp:131
bool lgREAD_DATA
Definition h2_priv.h:260
bool lgH2_He_ORNL
Definition h2_priv.h:381
long int nVib_hi[N_ELEC]
Definition h2_priv.h:473
bool lgH2_NOISECOSMIC
Definition h2_priv.h:390
double Average_collH2_deexcit
Definition h2_priv.h:306
void H2_RT_tau_inc(void)
Definition mole_h2.cpp:410
multi_arr< realnum, 3 > H2_X_hminus_formation_distribution
Definition h2_priv.h:544
void H2_Calc_Average_Rates(void)
Definition mole_h2.cpp:2445
double pops_per_elec[N_ELEC]
Definition h2_priv.h:482
molecule * sp_star
Definition h2_priv.h:426
double rate_grain_op_conserve
Definition h2_priv.h:281
void H2_ParseSave(Parser &p, ostringstream &chHeader)
Definition mole_h2_io.cpp:74
double HeatChange
Definition h2_priv.h:301
realnum mass_amu
Definition h2_priv.h:401
long int nzoneAsEval
Definition h2_priv.h:505
vector< double > stat_levn
Definition h2_priv.h:560
double Average_collH_excit
Definition h2_priv.h:309
double Abund() const
Definition h2_priv.h:76
long ip_photo_opac_offset
Definition h2_priv.h:322
long int getLine(long iElecHi, long iVibHi, long iRotHi, long iElecLo, long iVibLo, long iRotLo, double *relint, double *absint)
Definition mole_h2_io.cpp:1862
bool lgPrtMatrix
Definition h2_priv.h:590
void H2_CollidRateEvalAll(void)
Definition mole_h2_coll.cpp:13
vector< double > LTE_Temp
Definition h2_priv.h:587
bool lgEnabled
Definition h2_priv.h:353
vector< double > depart
Definition h2_priv.h:560
void GetIndices(long &ipHi, long &ipLo, const char *chLine, long &i) const
Definition mole_h2_coll.cpp:202
void H2_X_coll_rate_evaluate(void)
Definition mole_h2.cpp:199
void SolveExcitedElectronicLevels(void)
Definition mole_h2.cpp:1937
double Average_collH_dissoc_s
Definition h2_priv.h:312
double Solomon_elec_decay_g
Definition h2_priv.h:276
multi_arr< realnum, 3 > H2_dissprob
Definition h2_priv.h:494
string path
Definition h2_priv.h:435
double interpolate_LTE_Cooling(double Temp)
Definition mole_h2_etc.cpp:497
double H2_itrzn(void)
Definition mole_h2.cpp:264
multi_arr< realnum, 2 > H2_X_formation
Definition h2_priv.h:515
double ortho_para_old
Definition h2_priv.h:340
valarray< long > ipElec_H2_energy_sort
Definition h2_priv.h:547
void set_numLevelsMatrix(long numLevels)
Definition mole_h2_io.cpp:1924
double HeatDexc_old
Definition h2_priv.h:299
double TeUsedColl
Definition h2_priv.h:421
multi_arr< double, 2 > AulEscp
Definition h2_priv.h:556
void H2_Solomon_rate(void)
Definition mole_h2_etc.cpp:24
double H2_to_H_limit
Definition h2_priv.h:399
vector< double > pops
Definition h2_priv.h:560
long int nXLevelsMatrix
Definition h2_priv.h:554
double frac_matrix
Definition h2_priv.h:417
double spon_diss_tot
Definition h2_priv.h:270
void H2_ReadDissprob(long int nelec)
Definition mole_h2_io.cpp:854
t_coll_source coll_source[N_X_COLLIDER]
Definition h2_priv.h:324
double photo_heat_soft
Definition h2_priv.h:263
double ortho_colden
Definition h2_priv.h:336
double Cont_Diss_Heat_Rate(void)
Definition mole_dissociate.cpp:206
qList states
Definition h2_priv.h:427
multi_arr< realnum, 2 > H2_coll_dissoc_rate_coef
Definition h2_priv.h:527
string shortlabel
Definition h2_priv.h:434
string label
Definition h2_priv.h:433
void H2_Read_hminus_distribution(void)
Definition mole_h2_io.cpp:942
long int Jlowest[N_ELEC]
Definition h2_priv.h:478
double TeUsedBoltz
Definition h2_priv.h:420
bool lgLTE
Definition h2_priv.h:374
realnum GetXColden(long iVib, long iRot)
Definition mole_h2.cpp:2351
multi_arr< long int, 2 > ipTransitionSort
Definition h2_priv.h:550
double GetDissociationRate(const diss_tran &tran)
Definition mole_dissociate.cpp:194
int nElecLevelOutput
Definition h2_priv.h:357
double average_energy_g
Definition h2_priv.h:294
void H2_Prt_Zone(void)
Definition mole_h2_io.cpp:307
vector< double > create
Definition h2_priv.h:560
long int nH2_zone
Definition h2_priv.h:576
void H2_RT_OTS(void)
Definition mole_h2.cpp:2421
vector< double > destroy
Definition h2_priv.h:560
double(* photoion_opacity_fun)(double energy)
Definition h2_priv.h:82
double ortho_para_current
Definition h2_priv.h:340
realnum H2_CollidRateEvalOne(long iVibHi, long iRotHi, long iVibLo, long iRotLo, long ipHi, long ipLo, long nColl, double temp_K)
Definition mole_h2_coll.cpp:99
void H2_LineZero(void)
Definition mole_h2.cpp:440
bool lgH2_ORH2_ORNL
Definition h2_priv.h:384
void init(void)
Definition mole_h2_create.cpp:107
valarray< long > ipVib_H2_energy_sort
Definition h2_priv.h:546
double HeatDiss
Definition h2_priv.h:297
double ortho_para_older
Definition h2_priv.h:340
const double *const dense_total
Definition h2_priv.h:451
void SolveSomeGroundElectronicLevels(void)
Definition mole_h2.cpp:2048
void H2_ReadEnergies()
Definition mole_h2_io.cpp:634
void CalcPhotoionizationRate(void)
Definition mole_h2_etc.cpp:407
double HeatDexc_deriv
Definition h2_priv.h:300
multi_arr< realnum, 3 > CollRateErrFac
Definition h2_priv.h:484
double Average_collH2_dissoc_s
Definition h2_priv.h:314
double Solomon_dissoc_rate_g
Definition h2_priv.h:272
double H2_DissocEnergies[N_ELEC]
Definition h2_priv.h:471
void H2_Read_LTE_cooling_per_H2()
Definition mole_h2_io.cpp:1931
vector< double > excit
Definition h2_priv.h:560
void H2_Cooling(void)
Definition mole_h2.cpp:2185
int n_trace_matrix
Definition h2_priv.h:410
double rel_pop_LTE_g
Definition h2_priv.h:290
multi_arr< double, 2 > H2_col_rate_in
Definition h2_priv.h:510
double H2_Accel(void)
Definition mole_h2.cpp:295
long int levelAsEval
Definition h2_priv.h:562
multi_arr< double, 2 > CollRate_levn
Definition h2_priv.h:559
void H2_ContPoint(void)
Definition mole_h2.cpp:277
multi_arr< realnum, 6 > H2_SaveLine
Definition h2_priv.h:568
double MolDissocOpacity(const diss_tran &tran, const double &Mol_Ene)
Definition mole_dissociate.cpp:95
vector< CollRateCoeffArray > RateCoefTable
Definition h2_priv.h:485
long int loop_h2_oscil
Definition h2_priv.h:392
bool lgEvaluated
Definition h2_priv.h:318
multi_arr< double, 3 > H2_populations_LTE
Definition h2_priv.h:500
double gs_rate(void)
Definition mole_h2_etc.cpp:111
double H2_DR(void)
Definition mole_h2.cpp:2415
realnum ortho_density_f
Definition h2_priv.h:332
double H2_den_s
Definition h2_priv.h:541
long int iteration_evaluated
Definition h2_priv.h:565
double photodissoc_BigH2_H2g
Definition h2_priv.h:266
vector< diss_tran > Diss_Trans
Definition h2_priv.h:430
valarray< long > ipRot_H2_energy_sort
Definition h2_priv.h:548
multi_arr< realnum, 3 > CollRateCoeff
Definition h2_priv.h:483
void H2_X_sink_and_source(void)
Definition mole_h2.cpp:52
double LTE_Cooling_per_H2()
Definition mole_h2_etc.cpp:438
double photoionize_rate
Definition h2_priv.h:262
multi_arr< double, 2 > AulPump
Definition h2_priv.h:558
multi_arr< double, 2 > AulDest
Definition h2_priv.h:557
bool lgH2_grain_deexcitation
Definition h2_priv.h:371
bool lgColl_dissoc_coll
Definition h2_priv.h:367
bool lgH2_ortho_para_coll_on
Definition h2_priv.h:377
void H2_Prt_column_density(FILE *ioMEAN)
Definition mole_h2_io.cpp:369
double HeatChangeOld
Definition h2_priv.h:301
void H2_Colden(const char *chLabel)
Definition mole_h2.cpp:2370
void H2_LinesAdd(void)
Definition mole_h2_io.cpp:44
double GetExcitedElecDensity(void)
Definition mole_h2.cpp:2530
const double ENERGY_H2_STAR
Definition h2_priv.h:447
void Read_Mol_Diss_cross_sections(void)
Definition mole_dissociate.cpp:12
bool lgH2_NOISE
Definition h2_priv.h:388
multi_arr< double, 2 > H2_X_rate_from_elec_excited
Definition h2_priv.h:523
bool lgColl_deexec_Calc
Definition h2_priv.h:364
double xMeanNoise
Definition h2_priv.h:396
double Average_A
Definition h2_priv.h:304
double Average_collH2_dissoc_g
Definition h2_priv.h:313
double GetHeatRate(const diss_tran &tran)
Definition mole_dissociate.cpp:227
void H2_zero_pops_too_low(void)
Definition mole_h2_etc.cpp:184
double renorm_max
Definition h2_priv.h:344
vector< double > LTE_cool
Definition h2_priv.h:587
double Cont_Dissoc_Rate_H2s
Definition h2_priv.h:285
double photo_heat_hard
Definition h2_priv.h:264
int n_trace_full
Definition h2_priv.h:409
void mole_H2_LTE(void)
Definition mole_h2_etc.cpp:225
multi_arr< int, 2 > H2_ipPhoto
Definition h2_priv.h:509
diatomics(const string &a, const double &e_star, const double *const abund, double(*fun)(double))
Definition h2.cpp:13
void H2_RTMake(linefunc line_one)
Definition mole_h2.cpp:388
void H2_ReadDissocEnergies(void)
Definition mole_h2_io.cpp:793
void H2_PrtDepartCoef(void)
Definition mole_h2_io.cpp:336
long ip_photo_opac_thresh
Definition h2_priv.h:321
long int nCall_this_iteration
Definition h2_priv.h:584
void H2_Level_low_matrix(realnum abundance)
Definition mole_h2.cpp:472
void H2_ReadTransprob(long int nelec, TransitionList &trans)
Definition mole_h2_io.cpp:404
multi_arr< bool, 2 > lgH2_radiative
Definition h2_priv.h:572
multi_arr< double, 2 > pops_per_vib
Definition h2_priv.h:462
double HeatDexc
Definition h2_priv.h:298
int n_trace_iterations
Definition h2_priv.h:408
bool lgH2_PAH2_ORNL
Definition h2_priv.h:385
double H2_den_g
Definition h2_priv.h:541
double H2_RadPress(void)
Definition mole_h2.cpp:318
void H2_RT_tau_reset(void)
Definition mole_h2.cpp:456
valarray< realnum > H2_X_source
Definition h2_priv.h:533
double ortho_density
Definition h2_priv.h:327
void H2_LevelPops(bool &lgPopsConverged, double &old_value, double &new_value)
Definition mole_h2.cpp:896
multi_arr< realnum, 2 > H2_coll_dissoc_rate_coef_H2
Definition h2_priv.h:530
multi_arr< realnum, 2 > H2_X_colden_LTE
Definition h2_priv.h:521
void H2_RT_diffuse(void)
Definition mole_h2.cpp:369
void H2_PunchLineStuff(FILE *io, realnum xLimit, long index)
Definition mole_h2_io.cpp:1106
multi_arr< realnum, 2 > H2_X_colden
Definition h2_priv.h:519
double para_colden
Definition h2_priv.h:337
TransitionList trans
Definition h2_priv.h:428
multi_arr< realnum, 2 > H2_X_Hmin_back
Definition h2_priv.h:517
double average_energy_s
Definition h2_priv.h:295
multi_arr< realnum, 3 > H2_X_grain_formation_distribution
Definition h2_priv.h:538
long OpacityCreate(vector< double > &stack)
Definition mole_h2_etc.cpp:167
int n_trace_final
Definition h2_priv.h:407
void mole_H2_form(void)
Definition mole_h2_form.cpp:14
double Solomon_dissoc_rate_s
Definition h2_priv.h:273
bool lgFirst
Definition h2_priv.h:563
long int nCall_this_zone
Definition h2_priv.h:349
valarray< long > nRot_hi[N_ELEC]
Definition h2_priv.h:475
double renorm_min
Definition h2_priv.h:345
multi_arr< double, 3 > H2_old_populations
Definition h2_priv.h:499
double Average_collH2_excit
Definition h2_priv.h:308
long int ndim_allocated
Definition h2_priv.h:555
long int nH2_pops
Definition h2_priv.h:575
void H2_Reset(void)
Definition mole_h2_etc.cpp:298
realnum para_density_f
Definition h2_priv.h:333
int nTRACE
Definition h2_priv.h:404
multi_arr< double, 3 > H2_rad_rate_out
Definition h2_priv.h:496
TransitionList::iterator rad_end
Definition h2_priv.h:429
multi_arr< realnum, 3 > H2_disske
Definition h2_priv.h:495
multi_arr< double, 2 > H2_rad_rate_in
Definition h2_priv.h:512
double H2_InterEnergy(void)
multi_arr< double, 2 > H2_X_rate_to_elec_excited
Definition h2_priv.h:525
long int nLevels_per_elec[N_ELEC]
Definition h2_priv.h:480
void H2_CollidRateRead(long int nColl)
Definition mole_h2_coll.cpp:163
long int nzone_nlevel_set
Definition h2_priv.h:579
double rate_grain_J1_to_J0
Definition h2_priv.h:282
void H2_Punch_line_data(FILE *ioPUN, bool lgDoAll)
Definition mole_h2_io.cpp:1065
valarray< realnum > H2_X_sink
Definition h2_priv.h:534
molecule * sp
Definition h2_priv.h:425
multi_arr< long int, 3 > ipEnergySort
Definition h2_priv.h:549
long int nEner_H2_ground
Definition h2_priv.h:459
multi_arr< double, 3 > Cont_Dissoc_Rate
Definition h2_priv.h:287
Definition h2_priv.h:46
vector< double > xsections
Definition h2_priv.h:58
diss_tran(diss_level a, diss_level b)
Definition h2_priv.h:48
diss_level initial
Definition h2_priv.h:56
double rate_coeff
Definition h2_priv.h:59
vector< double > energies
Definition h2_priv.h:57
Definition mole.h:145
Definition container_classes.h:916
Definition quantumstate.h:36
double energy(const genericState &gs)
Definition generic_state.cpp:51
const int N_ELEC
Definition h2_priv.h:35
double Yan_H2_CS(double energy_ryd)
Definition mole_h2_etc.cpp:358
const int nTE_HMINUS
Definition h2_priv.h:32
void(* linefunc)(const TransitionProxy &t, bool lgShield_this_zone, realnum pestrk, realnum DopplerWidth, bool lgKeepLyman)
Definition h2_priv.h:10
const realnum H2_logte_hminus[nTE_HMINUS]
Definition h2_priv.h:38
double MolDissocCrossSection(const diss_tran &tran, const double &Mol_Ene)
Definition mole_dissociate.cpp:104
const int chN_X_COLLIDER
Definition h2_priv.h:23
const int N_X_COLLIDER
Definition h2_priv.h:21
STATIC long int ipPun
Definition save_do.cpp:368
Definition h2_priv.h:41
long v
Definition h2_priv.h:42
long j
Definition h2_priv.h:42
long n
Definition h2_priv.h:42
Definition h2_priv.h:63
t_coll_source()
Definition h2_priv.h:64
string filename
Definition h2_priv.h:70
long magic
Definition h2_priv.h:69