cloudy trunk
Loading...
Searching...
No Matches
iso.h File Reference
#include "two_photon.h"
#include "transition.h"
#include "freebound.h"
Include dependency graph for iso.h:

Go to the source code of this file.

Data Structures

class  t_isoCTRL
class  extra_tr
class  t_iso_sp

Macros

#define KILL_BELOW_PLASMA(E_)
#define N_(A_)
#define L_(A_)
#define S_(A_)
#define J_(A_)
#define IPRAD   0
#define IPCOLLIS   1
#define RREC_MAXN   40
#define LIKE_RREC_MAXN(A_)
#define N_ISO_TE_RECOMB   41
#define SumUpToThisN   1000
#define RECOMBMAGIC   (130216)

Enumerations

enum  {
  ipSINGLET = 1 , ipDOUBLET = 2 , ipTRIPLET = 3 , ipMULTIPLET_END ,
  ipMULTIPLET_BEGIN =ipSINGLET
}

Functions

void iso_cascade (long ipISO, long nelem)
void iso_charge_transfer_update (long nelem)
void iso_collapsed_bnl_print (long ipISO, long nelem)
void iso_collapsed_bnl_set (long ipISO, long nelem)
void iso_collapsed_Aul_update (long ipISO, long nelem)
void iso_collapsed_lifetimes_update (long ipISO, long nelem)
void iso_collide (long ipISO, long nelem)
void iso_collisional_ionization (long ipISO, long nelem)
void iso_continuum_lower (long ipISO, long nelem)
void iso_cool (long ipISO, long nelem)
void iso_create (void)
double iso_cross_section (double ERyd, double EthRyd, long n, long l, long S, long globalZ, long globalISO)
void iso_departure_coefficients (long ipISO, long nelem)
double iso_dielec_recomb_rate (long ipISO, long nelem, long ipLo)
void iso_error_generation (long ipISO, long nelem)
long iso_get_total_num_levels (long ipISO, long nmaxResolved, long numCollapsed)
void IonHydro ()
void iso_ionize_recombine (long ipISO, long nelem)
void iso_level (const long ipISO, const long nelem, double &renorm)
void iso_photo (long ipISO, long nelem)
void iso_prt_pops (long ipISO, long nelem, bool lgPrtDeparCoef)
void iso_put_error (long ipISO, long nelem, long ipHi, long ipLo, long whichData, realnum errorOpt, realnum errorPess)
void iso_radiative_recomb (long ipISO, long nelem)
void iso_radiative_recomb_effective (long ipISO, long nelem)
double iso_recomb_check (long ipISO, long nelem, long level, double temperature)
void iso_recomb_auxiliary_free (void)
void iso_recomb_malloc (void)
void iso_recomb_setup (long ipISO)
double iso_RRCoef_Te (long ipISO, long nelem, long n)
void iso_satellite_update (long nelem)
double iso_state_lifetime (long ipISO, long nelem, long n, long l)
void iso_solve (long ipISO, long nelem, double &maxerr)
void iso_suprathermal (long ipISO, long nelem)
void iso_update_num_levels (long ipISO, long nelem)
void iso_update_rates (void)
void iso_collapsed_update (void)
void iso_set_ion_rates (long ipISO, long nelem)
void iso_renorm (long nelem, long ipISO, double &renorm)

Variables

long int max_num_levels
const int ipH1s = 0
const int ipH2s = 1
const int ipH2p = 2
const int ipH3s = 3
const int ipH3p = 4
const int ipH3d = 5
const int ipH4s = 6
const int ipH4p = 7
const int ipH4d = 8
const int ipH4f = 9
const int ipHe1s1S = 0
const int ipHe2s3S = 1
const int ipHe2s1S = 2
const int ipHe2p3P0 = 3
const int ipHe2p3P1 = 4
const int ipHe2p3P2 = 5
const int ipHe2p1P = 6
const int ipHe3s3S = 7
const int ipHe3s1S = 8
const int ipHe3p3P = 9
const int ipHe3d3D = 10
const int ipHe3d1D = 11
const int ipHe3p1P = 12
const int ipH_LIKE = 0
const int ipHE_LIKE = 1
const int ipLI_LIKE = 2
const int ipBE_LIKE = 3
const int ipB_LIKE = 4
const int ipC_LIKE = 5
const int ipN_LIKE = 6
const int ipO_LIKE = 7
const int ipF_LIKE = 8
const int ipNE_LIKE = 9
const int ipNA_LIKE = 10
const int ipMG_LIKE = 11
const int ipAL_LIKE = 12
const int ipSI_LIKE = 13
const int ipP_LIKE = 14
const int ipS_LIKE = 15
const int ipCL_LIKE = 16
const int ipAR_LIKE = 17
t_isoCTRL iso_ctrl
t_iso_sp iso_sp [NISO][LIMELM]

Detailed Description

  • information for isoelectronic sequences

Definition in file iso.h.

Macro Definition Documentation

◆ IPCOLLIS

#define IPCOLLIS   1

Definition at line 87 of file iso.h.

Referenced by iso_collide(), iso_collisional_ionization(), and iso_level().

◆ IPRAD

◆ J_

#define J_ ( A_)
Value:
(iso_sp[ipISO][nelem].st[A_].j())
t_iso_sp iso_sp[NISO][LIMELM]
Definition iso.cpp:8

Definition at line 23 of file iso.h.

Referenced by DGEMM(), DGER(), DGETF2(), DGETRF(), and DTRSM().

◆ KILL_BELOW_PLASMA

#define KILL_BELOW_PLASMA ( E_)
Value:
( (rfield.lgPlasNu && ((E_)<rfield.plsfrq) ) ? 0.:1. )
t_rfield rfield
Definition rfield.cpp:8

This macro is used to zero any radiative process with photon energy below the plasma frequency. The energy must be in Rydbergs!

Definition at line 17 of file iso.h.

Referenced by iso_level().

◆ L_

◆ LIKE_RREC_MAXN

#define LIKE_RREC_MAXN ( A_)
Value:
( A_ == ipHELIUM ? 40 : 20 )
const int ipHELIUM
Definition cddefines.h:306

Ions of the sequences will go up to this n, h-like He will get same as iso roots.

Definition at line 98 of file iso.h.

Referenced by iso_recomb_malloc(), iso_recomb_setup(), and ParseCompile().

◆ N_

◆ N_ISO_TE_RECOMB

#define N_ISO_TE_RECOMB   41

Definition at line 100 of file iso.h.

Referenced by iso_recomb_malloc(), iso_recomb_setup(), and iso_RRCoef_Te().

◆ RECOMBMAGIC

#define RECOMBMAGIC   (130216)

the magic number for the table of recombination coefficients, YYMMDD

Definition at line 106 of file iso.h.

Referenced by iso_recomb_setup().

◆ RREC_MAXN

#define RREC_MAXN   40

this is the number of levels used with the atom xx-like levels large command

Definition at line 95 of file iso.h.

Referenced by iso_recomb_malloc(), iso_recomb_setup(), ParseAtomISO(), and ParseCompile().

◆ S_

◆ SumUpToThisN

#define SumUpToThisN   1000

This is the n to go up to when calculating total recombination. Any change here will not be reflected in total recomb until "compile xxlike" is run

Definition at line 104 of file iso.h.

Referenced by iso_radiative_recomb(), and iso_recomb_setup().

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
Enumerator
ipSINGLET 
ipDOUBLET 
ipTRIPLET 
ipMULTIPLET_END 
ipMULTIPLET_BEGIN 

Definition at line 81 of file iso.h.

Function Documentation

◆ IonHydro()

void IonHydro ( )

IonHydro this controls hydrogen atomic and molecular crosstalk

Definition at line 149 of file iso_solve.cpp.

References atmdat, colliders, DEBUG_ENTRY, dense, findspecieslocal(), fixit(), fnzone, hmi, hydro, ion_solver(), ioQQQ, ipH1s, ipH2p, ipH_LIKE, ipHYDROGEN, iso_renorm(), iso_sp, MAX2, PrintE82(), secondaries, SMALLDOUBLE, SMALLFLOAT, and trace.

Referenced by ion_wrapper().

Here is the call graph for this function:

◆ iso_cascade()

void iso_cascade ( long ipISO,
long nelem )

iso_cascade - calculate cascade probabilities, branching ratios, and associated errors

Parameters
ipISO
nelem

Cascade probabilities are as defined in Robbins 68, generalized here for cascade probability for any iso sequence.

>refer He triplets Robbins, R.R. 1968, ApJ 151, 497R
>refer He triplets Robbins, R.R. 1968a, ApJ 151, 511R

Definition at line 1101 of file iso_create.cpp.

References ASSERT, DEBUG_ENTRY, ex, ioQQQ, ipH_LIKE, ipHe2p3P0, ipHe2p3P1, ipHe2p3P2, ipHE_LIKE, ipHELIUM, IPRAD, iso_ctrl, iso_sp, L_, MALLOC, N_, opac, S, and S_.

Referenced by iso_collapsed_update(), and iso_create().

◆ iso_charge_transfer_update()

void iso_charge_transfer_update ( long nelem)

iso_charge_transfer_update - update rate coefficients for CT of H and He with everything else

Definition at line 20 of file iso_ionize_recombine.cpp.

References atmdat, DEBUG_ENTRY, dense, ipHYDROGEN, iso_sp, LIMELM, and t_atmdat::NCX.

Referenced by ion_solver().

◆ iso_collapsed_Aul_update()

void iso_collapsed_Aul_update ( long ipISO,
long nelem )

iso_collapsed_Aul_update - update decays from collapsed levels

Parameters
ipISO
nelem

Definition at line 1760 of file iso_create.cpp.

References ASSERT, DEBUG_ENTRY, ipH_LIKE, ipHE_LIKE, iso_sp, L_, and TotalInsanity().

Referenced by iso_collapsed_update(), and iso_zero().

Here is the call graph for this function:

◆ iso_collapsed_bnl_print()

void iso_collapsed_bnl_print ( long ipISO,
long nelem )

iso_collapsed_bnl_print - print departure coefficients for collapsed levels

Parameters
ipISO
nelem

Definition at line 1714 of file iso_create.cpp.

References DEBUG_ENTRY, elementnames, ioQQQ, ipH_LIKE, ipHE_LIKE, iso_ctrl, and iso_sp.

◆ iso_collapsed_bnl_set()

void iso_collapsed_bnl_set ( long ipISO,
long nelem )

iso_collapsed_bnl_set - set departure coefficients for collapsed levels

Parameters
ipISO
nelem

these are just sanity checks, the interpolated value should be between values at interpolation points

Definition at line 1519 of file iso_create.cpp.

References ASSERT, DEBUG_ENTRY, dense, hunt_bisect(), ipH_LIKE, ipHE_LIKE, ipHELIUM, ipHYDROGEN, iso_sp, MAX2, MAX4, MIN2, MIN4, phycon, and TotalInsanity().

Referenced by iso_collapsed_update(), and iso_zero().

Here is the call graph for this function:

◆ iso_collapsed_lifetimes_update()

void iso_collapsed_lifetimes_update ( long ipISO,
long nelem )

iso_collapsed_lifetimes_update - update lifetimes of collapsed levels

Parameters
ipISO
nelem

Definition at line 1807 of file iso_create.cpp.

References ASSERT, DEBUG_ENTRY, iso_ctrl, iso_sp, PI4, and SMALLFLOAT.

Referenced by iso_collapsed_update(), and iso_zero().

◆ iso_collapsed_update()

void iso_collapsed_update ( void )

Definition at line 27 of file iso_solve.cpp.

References conv, dense, ipH_LIKE, ipHELIUM, ipHYDROGEN, iso_cascade(), iso_collapsed_Aul_update(), iso_collapsed_bnl_set(), iso_collapsed_lifetimes_update(), MIN2, and NISO.

Referenced by ConvBase().

Here is the call graph for this function:

◆ iso_collide()

◆ iso_collisional_ionization()

void iso_collisional_ionization ( long ipISO,
long nelem )

iso_collisional_ionization - calculate collisional ionization rate for ipISO, nelem

Parameters
ipISO
nelem

Definition at line 25 of file iso_collide.cpp.

References ASSERT, t_ADfA::coll_ion_wrapper(), DEBUG_ENTRY, t_iso_sp::fb, Hion_coll_ioniz_ratecoef(), hydro_vs_ioniz(), Singleton< t_ADfA >::Inst(), IPCOLLIS, iso_ctrl, iso_put_error(), iso_sp, t_iso_sp::lgLevelsLowered, N_, NISO, t_iso_sp::numLevels_max, and phycon.

Referenced by iso_collide().

Here is the call graph for this function:

◆ iso_continuum_lower()

void iso_continuum_lower ( long ipISO,
long nelem )

◆ iso_cool()

void iso_cool ( long ipISO,
long nelem )

iso_cool compute net heating/cooling due to hydrogenc atom species

Parameters
ipISOthe isoelectronic sequence, 0 for H
nelemis element, so 0 for H itself

References globalISO, globalZ, and S.

Referenced by CoolEvaluate().

◆ iso_create()

◆ iso_cross_section()

double iso_cross_section ( double ERyd,
double EthRyd,
long n,
long l,
long S,
long globalZ,
long globalISO )

iso_cross_section get cross section for a particular level of an iso sequence ion

Parameters
ERyd
EthRyd
n
l
S
Z
ipISO

Definition at line 132 of file iso_radiative_recomb.cpp.

References cross_section(), DEBUG_ENTRY, globalISO, globalZ, H_cross_section(), He_cross_section(), ipH_LIKE, ipHE_LIKE, S, and TotalInsanity().

Referenced by iso_recomb_integrand().

Here is the call graph for this function:

◆ iso_departure_coefficients()

void iso_departure_coefficients ( long ipISO,
long nelem )

iso_departure_coefficients - calculate departure coefficients

Parameters
ipISO
nelem

Definition at line 363 of file iso_solve.cpp.

References DEBUG_ENTRY, dense, iso_sp, safe_div(), and SMALLFLOAT.

Referenced by ion_solver().

Here is the call graph for this function:

◆ iso_dielec_recomb_rate()

double iso_dielec_recomb_rate ( long ipISO,
long nelem,
long ipLo )

iso_dielec_recomb_rate - get state-specific dielectronic recombination rate

Parameters
ipISO
nelem
ipLo

Definition at line 1149 of file iso_radiative_recomb.cpp.

References ASSERT, DEBUG_ENTRY, freeBound::DielecRecombVsTemp, hunt_bisect(), ipHe1s1S, ipHE_LIKE, iso_ctrl, iso_sp, NUM_DR_TEMPS, and phycon.

Referenced by iso_radiative_recomb().

Here is the call graph for this function:

◆ iso_error_generation()

void iso_error_generation ( long ipISO,
long nelem )

iso_error_generation generate gaussian errors

Parameters
ipISO
nelem

Definition at line 39 of file iso_error.cpp.

References ASSERT, DEBUG_ENTRY, ex, iso_sp, and MyGaussRand().

Referenced by iso_update_rates().

Here is the call graph for this function:

◆ iso_get_total_num_levels()

long iso_get_total_num_levels ( long ipISO,
long nmaxResolved,
long numCollapsed )

iso_get_total_num_levels - get total number of levels with the given number of resolved and collapsed

Parameters
ipISO
nmaxResolved
numCollapsed

Definition at line 1465 of file iso_create.cpp.

References DEBUG_ENTRY, ipH_LIKE, ipHE_LIKE, and TotalInsanity().

Referenced by iso_continuum_lower(), iso_recomb_malloc(), and iso_update_num_levels().

Here is the call graph for this function:

◆ iso_ionize_recombine()

void iso_ionize_recombine ( long ipISO,
long nelem )

iso_ionize_recombine evaluate state specific creation and destruction processes

Parameters
ipISO
nelem

Referenced by iso_update_rates().

◆ iso_level()

void iso_level ( const long ipISO,
const long nelem,
double & renorm )

iso_level solve for iso-sequence ionization balance

Parameters
ipISO
nelem

Referenced by iso_solve().

◆ iso_photo()

void iso_photo ( long ipISO,
long nelem )

iso_photo do photoionization rates for element nelem on the ipISO isoelectronic sequence

Parameters
ipISO
nelem

Referenced by iso_update_rates().

◆ iso_prt_pops()

void iso_prt_pops ( long ipISO,
long nelem,
bool lgPrtDeparCoef )

iso_prt_pops routine to print level pops or departure coefficients for iso sequences

Parameters
ipISO
nelem
lgPrtDeparCoef

Definition at line 386 of file iso_solve.cpp.

References ASSERT, DEBUG_ENTRY, elementnames, ioQQQ, ipH_LIKE, ipHe2p3P0, ipHe2p3P1, ipHe2p3P2, ipHE_LIKE, iso_ctrl, iso_sp, ITEM_TO_PRINT, and NISO.

Referenced by PrtZone().

◆ iso_put_error()

void iso_put_error ( long ipISO,
long nelem,
long ipHi,
long ipLo,
long whichData,
realnum errorOpt,
realnum errorPess )

iso_put_error put an error bar on a piece of data, to be used with Gaussian random noise gen

Parameters
ipISO
nelem
ipHi
ipLo
whichData
errorOpt
errorPess

Referenced by ForbiddenAuls(), he_1trans(), helike_transprob(), hydro_transprob(), iso_collide(), iso_collisional_ionization(), iso_put_recomb_error(), and iso_radiative_recomb().

◆ iso_radiative_recomb()

void iso_radiative_recomb ( long ipISO,
long nelem )

◆ iso_radiative_recomb_effective()

void iso_radiative_recomb_effective ( long ipISO,
long nelem )

iso_radiative_recomb_effective - get effective recomb rate coefficients into each level (including indirect)

Parameters
ipISO
nelem

Definition at line 602 of file iso_radiative_recomb.cpp.

References ASSERT, DEBUG_ENTRY, dprintf(), EN1RYD, ex, ioQQQ, IPRAD, ipRecRad, iso_ctrl, iso_sp, L_, MAX2, N_, phycon, prt_wl(), RYDLAM, S_, and wavelength.

Referenced by iso_update_rates().

Here is the call graph for this function:

◆ iso_recomb_auxiliary_free()

void iso_recomb_auxiliary_free ( void )

iso_recomb_auxiliary_free - free up some auxiliary space associated with iso recombination tables.

Definition at line 825 of file iso_radiative_recomb.cpp.

References DEBUG_ENTRY, NISO, and NumLevRecomb.

Referenced by iso_create().

◆ iso_recomb_check()

double iso_recomb_check ( long ipISO,
long nelem,
long level,
double temperature )

iso_recomb_check - called by SanityCheck to confirm that recombination coef are ok, return value is relative error between new calculation of recom, and interp value

Parameters
ipISO
nelemthe chemical element, 1 for He
levelthe level, 0 for ground
temperaturethe temperature to be used

Definition at line 734 of file iso_radiative_recomb.cpp.

References DEBUG_ENTRY, iso_radrecomb_from_cross_section(), iso_RRCoef_Te(), MAX2, phycon, and TempChange().

Referenced by SanityCheckBegin().

Here is the call graph for this function:

◆ iso_recomb_malloc()

void iso_recomb_malloc ( void )

iso_recomb_malloc - malloc space needed for iso recombination tables.

Definition at line 765 of file iso_radiative_recomb.cpp.

References DEBUG_ENTRY, dense, iso_get_total_num_levels(), iso_sp, LIKE_RREC_MAXN, LIMELM, MALLOC, MAX2, N_ISO_TE_RECOMB, NISO, NumLevRecomb, RRCoef, RREC_MAXN, TeRRCoef, and TotalRecomb.

Referenced by iso_create().

Here is the call graph for this function:

◆ iso_recomb_setup()

void iso_recomb_setup ( long ipISO)

iso_recomb_setup - read in or compile iso recombination tables.

Parameters
ipISO

Establish radiative recombination rate coefficients - RRC

Definition at line 838 of file iso_radiative_recomb.cpp.

References AS_LOCAL_ONLY, ASSERT, cdEXIT, chLine_LENGTH, DEBUG_ENTRY, dense, elementnames, EXIT_FAILURE, EXIT_SUCCESS, FFmtRead(), t_ADfA::H_rad_rec(), Singleton< t_ADfA >::Inst(), ioQQQ, ipHe1s1S, iso_ctrl, iso_radrecomb_from_cross_section(), iso_sp, LIKE_RREC_MAXN, LIMELM, N_ISO_TE_RECOMB, NHYDRO_MAX_LEVEL, NISO, NumLevRecomb, open_data(), read_whole_line(), Recomb_Seaton59(), RECOMBMAGIC, RRCoef, RREC_MAXN, SumUpToThisN, TeRRCoef, TotalRecomb, and trace.

Referenced by iso_create().

Here is the call graph for this function:

◆ iso_renorm()

void iso_renorm ( long nelem,
long ipISO,
double & renorm )

iso_renorm - renormalize H-like so that it agrees with the ionization balance

Definition at line 272 of file iso_solve.cpp.

References ASSERT, conv, DEBUG_ENTRY, dense, fp_equal(), ioQQQ, iso_ctrl, iso_sp, and SMALLFLOAT.

Referenced by ConvBase(), IonHydro(), iso_level(), lgTrivialSolution(), ScaleIonDensities(), and store_new_densities().

Here is the call graph for this function:

◆ iso_RRCoef_Te()

double iso_RRCoef_Te ( long ipISO,
long nelem,
long n )

iso_RRCoef_Te - interpolate iso recomb coeff as function of temperature

Parameters
ipISO
nelem
n

Definition at line 708 of file iso_radiative_recomb.cpp.

References ASSERT, DEBUG_ENTRY, iso_ctrl, iso_sp, N_ISO_TE_RECOMB, RRCoef, TempInterp(), TeRRCoef, and TotalRecomb.

Referenced by iso_radiative_recomb(), and iso_recomb_check().

Here is the call graph for this function:

◆ iso_satellite_update()

void iso_satellite_update ( long nelem)

iso_satellite_update - update iso satellite line information

Definition at line 1381 of file iso_create.cpp.

References abscf(), ATOMIC_MASS_UNIT, DEBUG_ENTRY, dense, dsexp(), ELECTRON_MASS, ERG1CM, GetGF(), HION_LTE_POP, ipHE_LIKE, ipSatelliteLines, iso_ctrl, iso_sp, max(), MIN2, NISO, phycon, PI4, SatelliteLines, SMALLDOUBLE, and SMALLFLOAT.

Referenced by ion_solver(), iso_create(), and lines().

Here is the call graph for this function:

◆ iso_set_ion_rates()

void iso_set_ion_rates ( long ipISO,
long nelem )

◆ iso_solve()

void iso_solve ( long ipISO,
long nelem,
double & maxerr )

iso_solve - main routine to call iso_level and determine iso level balances

Parameters
ipISO

Definition at line 102 of file iso_solve.cpp.

References ASSERT, DEBUG_ENTRY, dense, HydroLevel(), ipH_LIKE, iso_ctrl, iso_level(), and iso_sp.

Referenced by ion_solver().

Here is the call graph for this function:

◆ iso_state_lifetime()

double iso_state_lifetime ( long ipISO,
long nelem,
long n,
long l )

◆ iso_suprathermal()

void iso_suprathermal ( long ipISO,
long nelem )

iso_suprathermal - calculate secondary excitation by suprathermal electrons for iso sequences

Parameters
ipISO
nelem

Definition at line 79 of file iso_collide.cpp.

References ASSERT, TransitionProxy::Coll(), DEBUG_ENTRY, TransitionProxy::EnergyWN(), TransitionProxy::ipCont(), ipH2p, ipH_LIKE, ipHYDROGEN, iso_ctrl, iso_sp, LIMELM, NISO, t_iso_sp::numLevels_max, CollisionProxy::rate_lu_nontherm_set(), secondaries, and t_iso_sp::trans().

Referenced by iso_collide().

Here is the call graph for this function:

◆ iso_update_num_levels()

void iso_update_num_levels ( long ipISO,
long nelem )

iso_update_num_levels - update level informations for iso sequences

Parameters
ipISO
nelem

Definition at line 1488 of file iso_create.cpp.

References ASSERT, cdEXIT, DEBUG_ENTRY, EXIT_FAILURE, ioQQQ, iso_get_total_num_levels(), iso_sp, MAX2, and max_num_levels.

Referenced by InitCoreload(), InitCoreloadPostparse(), InitSimPostparse(), ParseAtomISO(), and ParseCompile().

Here is the call graph for this function:

◆ iso_update_rates()

void iso_update_rates ( void )

iso_update_rates routine to set up iso rates, level balance is done elsewhere

Todo
2 the indices for the two-photon rates must be changed for further iso sequences.

Definition at line 51 of file iso_solve.cpp.

References ASSERT, CalcTwoPhotonRates(), conv, dense, ionbal, ipH_LIKE, ipHE_LIKE, ipHYDROGEN, iso_collide(), iso_continuum_lower(), iso_ctrl, iso_error_generation(), iso_ionize_recombine(), iso_photo(), iso_radiative_recomb(), iso_radiative_recomb_effective(), iso_sp, LIMELM, MIN2, NISO, nzone, rfield, and t_iso_sp::TwoNu.

Referenced by ConvBase().

Here is the call graph for this function:

Variable Documentation

◆ ipAL_LIKE

const int ipAL_LIKE = 12

Definition at line 74 of file iso.h.

Referenced by atmdat_readin().

◆ ipAR_LIKE

const int ipAR_LIKE = 17

Definition at line 79 of file iso.h.

◆ ipB_LIKE

const int ipB_LIKE = 4

Definition at line 66 of file iso.h.

◆ ipBE_LIKE

const int ipBE_LIKE = 3

Definition at line 65 of file iso.h.

◆ ipC_LIKE

const int ipC_LIKE = 5

Definition at line 67 of file iso.h.

◆ ipCL_LIKE

const int ipCL_LIKE = 16

Definition at line 78 of file iso.h.

◆ ipF_LIKE

const int ipF_LIKE = 8

Definition at line 70 of file iso.h.

◆ ipH1s

◆ ipH2p

◆ ipH2s

◆ ipH3d

const int ipH3d = 5

Definition at line 32 of file iso.h.

Referenced by iso_create(), lines_continuum(), OpacityAddTotal(), and RT_tau_init().

◆ ipH3p

◆ ipH3s

const int ipH3s = 3

Definition at line 30 of file iso.h.

Referenced by iso_create(), lines_continuum(), lines_hydro(), OpacityAddTotal(), PrtAllTau(), and RT_tau_init().

◆ ipH4d

const int ipH4d = 8

Definition at line 35 of file iso.h.

Referenced by lines_general(), and OpacityAddTotal().

◆ ipH4f

const int ipH4f = 9

Definition at line 36 of file iso.h.

Referenced by OpacityAddTotal().

◆ ipH4p

const int ipH4p = 7

Definition at line 34 of file iso.h.

Referenced by lines_general(), lines_hydro(), OpacityAddTotal(), and PrtAllTau().

◆ ipH4s

const int ipH4s = 6

Definition at line 33 of file iso.h.

Referenced by lines_general(), lines_hydro(), OpacityAddTotal(), and PrtAllTau().

◆ ipH_LIKE

const int ipH_LIKE = 0

these are array indices for isoelectronic sequences, same as element but used for array addressing to make context totally clear

Definition at line 62 of file iso.h.

Referenced by atmdat_2phot_shapefunction(), atom_oi_calc(), t_fe2ovr_la::atoms_fe2ovr(), Badnell_rec_init(), cdTemp(), ChargTranSumHeat(), chkCaHeps(), collision_strength_VF01(), ContCreatePointers(), ContRate(), ContSetIntensity(), ConvBase(), ConvIterCheck(), CoolCalc(), CoolEvaluate(), CS_l_mixing_VF01(), CS_PercivalRichards78(), dBase_solve(), DynaCreateArrays(), DynaIonize(), DynaNewStep(), DynaSaveLast(), DynaStartZone(), DynaZero(), eden_sum(), Fe2_cooling(), FeIILevelPops(), FeIILyaPump(), FillExtraLymanLine(), GrainRateDr(), GrainTemperature(), H21_cm_pops(), t_ADfA::h_coll_str(), HeatSum(), hydro_transprob(), HydroCSInterp(), HydroLevel(), HydroRecCool(), InitCoreload(), InitCoreloadPostparse(), InitDefaultsPreparse(), InitSimPostparse(), ion_photo(), ion_recomb(), ion_solver(), IonHelium(), IonHydro(), iso_allocate(), iso_assign_quantum_numbers(), iso_cascade(), iso_collapsed_Aul_update(), iso_collapsed_bnl_print(), iso_collapsed_bnl_set(), iso_collapsed_update(), iso_collide(), iso_create(), iso_cross_section(), iso_get_total_num_levels(), iso_ionize_recombine(), iso_level(), iso_prt_pops(), iso_rad_rec_cooling_approx(), iso_rad_rec_cooling_extra(), iso_radrecomb_from_cross_section(), iso_solve(), iso_suprathermal(), iso_update_rates(), iso_zero(), IterRestart(), IterStart(), lgCheckMonitors(), lgTrivialSolution(), lines(), lines_continuum(), lines_general(), lines_helium(), lines_hydro(), lines_lv1_k_zn(), mole_h_rate_diagnostics(), Opacity_iso_photo_cs(), OpacityAddTotal(), OpacityCreateAll(), ParseAtom(), ParseAtomISO(), ParseCompile(), ParseDont(), ParseElement(), ParseMonitorResults(), ParsePrint(), ParseTest(), ParseTrace(), PresTotCurrent(), PrtAllTau(), PrtComment(), PrtFinal(), PrtHeader(), PrtHydroTrace1(), PrtHydroTrace1a(), prtmet(), PrtZone(), radius_first(), radius_increment(), radius_next(), RT_continuum(), RT_DestProb(), RT_diffuse(), RT_line_all(), RT_line_driving(), RT_line_pumping(), RT_OTS(), RT_OTS_Update(), RT_stark(), RT_tau_inc(), RT_tau_init(), RT_tau_reset(), SanityCheckBegin(), save_opacity(), SaveDo(), SaveLineData(), SaveLineStuff(), SaveSpecial(), state_get_put(), store_new_densities(), tfidle(), and zero().

◆ ipHe1s1S

◆ ipHe2p1P

const int ipHe2p1P = 6

Definition at line 49 of file iso.h.

Referenced by ForbiddenAuls(), he_1trans(), PrtAllTau(), PrtZone(), RT_tau_init(), RT_tau_reset(), SaveDo(), and zero().

◆ ipHe2p3P0

◆ ipHe2p3P1

◆ ipHe2p3P2

◆ ipHe2s1S

const int ipHe2s1S = 2

Definition at line 45 of file iso.h.

Referenced by ForbiddenAuls(), he_1trans(), IonCSInterp(), PrtAllTau(), and PrtZone().

◆ ipHe2s3S

◆ ipHe3d1D

const int ipHe3d1D = 11

Definition at line 56 of file iso.h.

Referenced by AGN_He1_CS(), and PrtZone().

◆ ipHe3d3D

const int ipHe3d3D = 10

Definition at line 55 of file iso.h.

Referenced by AGN_He1_CS(), lines_helium(), PrtAllTau(), and PrtZone().

◆ ipHe3p1P

const int ipHe3p1P = 12

Definition at line 57 of file iso.h.

Referenced by ForbiddenAuls(), PrtAllTau(), and PrtZone().

◆ ipHe3p3P

const int ipHe3p3P = 9

Definition at line 54 of file iso.h.

Referenced by AGN_He1_CS(), ForbiddenAuls(), lines_helium(), PrtAllTau(), and PrtZone().

◆ ipHe3s1S

const int ipHe3s1S = 8

Definition at line 53 of file iso.h.

Referenced by ForbiddenAuls(), and PrtZone().

◆ ipHe3s3S

const int ipHe3s3S = 7

Definition at line 52 of file iso.h.

Referenced by AGN_He1_CS(), ForbiddenAuls(), he_1trans(), lines_helium(), and PrtZone().

◆ ipHE_LIKE

◆ ipLI_LIKE

const int ipLI_LIKE = 2

Definition at line 64 of file iso.h.

Referenced by atmdat_readin().

◆ ipMG_LIKE

const int ipMG_LIKE = 11

Definition at line 73 of file iso.h.

Referenced by atmdat_readin().

◆ ipN_LIKE

const int ipN_LIKE = 6

Definition at line 68 of file iso.h.

◆ ipNA_LIKE

const int ipNA_LIKE = 10

Definition at line 72 of file iso.h.

Referenced by atmdat_readin().

◆ ipNE_LIKE

const int ipNE_LIKE = 9

Definition at line 71 of file iso.h.

◆ ipO_LIKE

const int ipO_LIKE = 7

Definition at line 69 of file iso.h.

◆ ipP_LIKE

const int ipP_LIKE = 14

Definition at line 76 of file iso.h.

◆ ipS_LIKE

const int ipS_LIKE = 15

Definition at line 77 of file iso.h.

◆ ipSI_LIKE

const int ipSI_LIKE = 13

Definition at line 75 of file iso.h.

◆ iso_ctrl

◆ iso_sp

t_iso_sp iso_sp[NISO][LIMELM]
extern

Definition at line 8 of file iso.cpp.

Referenced by atmdat_readin(), atom_oi_calc(), AtomCSInterp(), t_fe2ovr_la::atoms_fe2ovr(), Badnell_rec_init(), cdTemp(), ChargTranSumHeat(), chkCaHeps(), collision_strength_VF01(), ContCreatePointers(), ContRate(), ContSetIntensity(), ConvBase(), ConvIterCheck(), CoolCalc(), CoolEvaluate(), cross_section(), CS_l_mixing_S62(), dBase_solve(), DoFSMixing(), DoSatelliteLines(), DynaCreateArrays(), DynaIonize(), DynaNewStep(), DynaSaveLast(), DynaStartZone(), eden_sum(), Fe2_cooling(), FeIILevelPops(), FeIILyaPump(), fill_array(), FillExtraLymanLine(), ForbiddenAuls(), GetStandardHeLines(), GrainRateDr(), GrainTemperature(), H21_cm_pops(), he_1trans(), He_cross_section(), HeatSum(), HeCollidSetup(), HeCSInterp(), helike_energy(), helike_quantum_defect(), helike_transprob(), hydro_transprob(), hydro_vs_coll_str(), hydro_vs_deexcit(), HydroCSInterp(), HydroLevel(), HydroRecCool(), InitCoreload(), InitCoreloadPostparse(), InitSimPostparse(), ion_CX(), ion_photo(), ion_trim(), IonCSInterp(), IonHelium(), IonHydro(), iso_allocate(), iso_assign_quantum_numbers(), iso_cascade(), iso_charge_transfer_update(), iso_collapsed_Aul_update(), iso_collapsed_bnl_print(), iso_collapsed_bnl_set(), iso_collapsed_lifetimes_update(), iso_collide(), iso_collisional_ionization(), iso_continuum_lower(), iso_cool(), iso_create(), iso_departure_coefficients(), iso_dielec_recomb_rate(), iso_error_generation(), iso_ionize_recombine(), iso_level(), iso_photo(), iso_prt_pops(), iso_put_error(), iso_put_recomb_error(), iso_rad_rec_cooling_approx(), iso_rad_rec_cooling_extra(), iso_radiative_recomb(), iso_radiative_recomb_effective(), iso_radrecomb_from_cross_section(), iso_recomb_malloc(), iso_recomb_setup(), iso_renorm(), iso_RRCoef_Te(), iso_satellite(), iso_satellite_update(), iso_set_ion_rates(), iso_solve(), iso_suprathermal(), iso_update_num_levels(), iso_update_rates(), iso_zero(), IterRestart(), IterStart(), lgCheckMonitors(), lines(), lines_continuum(), lines_general(), lines_helium(), lines_hydro(), lines_lv1_k_zn(), mole_h_rate_diagnostics(), mole_h_reactions(), Opacity_iso_photo_cs(), OpacityAdd1Element(), OpacityAddTotal(), OpacityCreateAll(), ParseAtomISO(), ParseCompile(), ParseElement(), ParsePrint(), PresTotCurrent(), PrintRates(), PrtAllTau(), PrtComment(), PrtFinal(), PrtHeader(), PrtHydroTrace1(), PrtHydroTrace1a(), PrtLinePres(), prtmet(), PrtZone(), radius_first(), radius_increment(), radius_next(), RT_continuum(), RT_DestProb(), RT_diffuse(), RT_iso_integrate_RRC(), RT_line_all(), RT_line_driving(), RT_line_pumping(), RT_OTS(), RT_OTS_Update(), RT_stark(), RT_tau_inc(), RT_tau_init(), RT_tau_reset(), SanityCheckBegin(), Save_Line_RT(), save_opacity(), SaveDo(), SaveLineData(), SaveLineStuff(), SaveSpecial(), state_get_put(), store_new_densities(), and tfidle().

◆ max_num_levels

long int max_num_levels
extern

Definition at line 10 of file iso.cpp.

Referenced by InitCoreload(), and iso_update_num_levels().