API
The following sections describe the call, and the input and output parameters of each XCLASS function in detail.
myXCLASS
# import myXCLASS function
from xclass.task_myXCLASS import myXCLASSCore
myXCLASSCore(MolfitFile_Class = None,
MolfitFileName = "",
MolfitsFileName = None,
ObsXMLFile_Class = None,
ObsXMLFileName = None,
JobDir = None,
ParameterMapDir = None,
regionFileName = "", \
FITSImageMaskFileName = "", \
Threshold = None, \
FreqMin = 1.0,
FreqMax = 1.e8,
FreqStep = 1.0,
TelescopeSize = None,
BMIN = None,
BMAJ = None,
BPA = None,
Redshift = None,
Inter_Flag = False,
t_back_flag = True,
tBack = 0.0,
tSlope = 0.0,
nH_flag = False,
N_H = 0.0,
beta_dust = 0.0,
kappa_1300 = 0.0,
DustFileName = "",
BackgroundFileName = "",
Te_ff = None,
EM_ff = None,
kappa_sync = None,
B_sync = None,
p_sync = None,
l_sync = None,
ContPhenFuncID = None,
ContPhenFuncParam1 = None,
ContPhenFuncParam2 = None,
ContPhenFuncParam3 = None,
ContPhenFuncParam4 = None,
ContPhenFuncParam5 = None,
iso_flag = False,
IsoTableFileName = "",
IsoRatioFile_Class = None,
CollisionFileName = "",
CollPartnerFile_Class = None,
NumModelPixelXX = 100,
NumModelPixelYY = 100,
LocalOverlapFlag = False,
NoSubBeamFlag = False,
DB_Class = None,
DBFileName = "",
dbFileName = "",
CheckMoleculeInDB = True,
RestFreq = 0.0,
vLSR = 0.0,
EmAbsFlag = False,
EmAbsPATH = None,
ChannelIntegrationFlag = True,
OpacityFlag = False,
printFlag = None,
NoWarningFlag = False,
SpectrumOnlyFlag = False,
UsemyXCLASSProgramFlag = True,
BibTeXFlag = False,
NumberCores = 1,
ParallelizationMethod = "process")
The function models a spectrum by solving the radiative transfer
equation for an isothermal object in one dimension, called detection
equation. Furthermore, Non-LTE, recombination-lines, extended sources
and dust, free-free, and synchrotron continuum contribution can be
taken into account as well. The myXCLASS function can be use to
generate single spectra or entire cubes. For cubes, the input parameter
for each pixel can be described by parameter maps, which are located in
a directory described by parameter ParameterMapDir.
The files produced by a myXCLASS function call are copied to a
so-called “job-directory”. If input parameter JobDir is not
set, a new job directory located in the XCLASS interface working
directory path-of-XCLASS/run/myXCLASS/ is created! The name of a
job directory for a myXCLASS run is made up of four components: The
first part consists of the phrase “job_” whereas the second and
third part describe the date and the time stamp (hours, minutes,
seconds) of the function execution, respectively. The last part
indicates a so-called job ID which is composed of the so-called PID
followed by a four digit random integer number to create a really
unambiguous job number, e.g.
path-of-XCLASS/run/myXCLASS/job__25-07-2013__12-02-03__189644932/
Input parameters:
MolfitsFileName: path and name of the molfit file, see Sect. “The molfit file”.
MolfitFileName: (optional) alias ofMolfitsFileName, (default:None).
MolfitFile_Class: (optional) molfit file class object, (default:None).
ObsXMLFileName: (optional) path and name of an obs. xml file, (default:None), see Sect. “Observational xml file” for a detailed description. The obs. xml file can be used to create a synthetic spectrum for each frequency range.Note, the name of the xml file has to end with
".xml".
ObsXMLFile_Class: (optional) obs. xml file class object, (default:None).
JobDir: (optional) job directory used to store output files, if None, a new directory is created (default:None).
ParameterMapDir: (optional) path of subdirectory containing FITS images for different molfit parameters, (default:None).
regionFileName: (optional) path and name of a so-called ds9 [2] region file, (default:"").
FITSImageMaskFileName: (optional) path and name of a FITS image used for masking. Here, each pixel with a real numerical value (except 0), i.e. which is not NaN or 0, is fitted, (default:""). The value of the selecting pixel is used to weight the contribution of the corresponding pixel to the total \(\chi^2\) value.
Threshold: defines a threshold for a pixel. If the spectrum of a pixel has an max. intensity lower than the value defined by this parameter the pixel is not fitted (ignored), (default:None).Please note, the value for the
Thresholdparameter has to be given in the same unit than the spectrum in the FITS file.NOTE, if an observational xml-file is used, the user can define different threshold values for each frequency range by using the tag
<Threshold>. Please note, that the tag<Threshold>has to be given for each frequency range. If tag<Threshold>is defined in the xml file, the value of the input parameterThresholdis ignored.
FreqMin: start frequency of simulated spectrum (in MHz), (default: 1).
FreqMax: end frequency of simulated spectrum (in MHz), (default: \(10^8\)).
FreqStep: step frequency of simulated spectrum (in MHz), (default: 1).
vLSR: velocity (local standard of rest) in km s\(^{-1}\) (default: 0) used in the calculation of the synthetic spectra, i.e. all velocity offsets described in the molfit file have to be defined relative to the v\(_{\rm LSR}\) parameter. For example, the user defines v\(_{\rm LSR}\) = 20 km/s, and the molfit file describes CH\(_3\)CN with one component and v\(_{\rm off}^{\rm molfit}\) = 2 km/s, XCLASS uses v\(_{\rm off}\) = v\(_{\rm LSR}\) + v\(_{\rm off}^{\rm molfit}\) = 20 km/s + 2 km/s = 22 km/s for the calculation of the spectra.Additionally, XCLASS considers the v\(_{\rm LSR}\) parameter for the determination of the transition frequencies contained in the database: For each frequency range, XCLASS shifts the lowest and highest frequency by v\(_{\rm LSR}\) (e.g. to describe lines in high-red-shifted galaxies) and expands each range by taking the lowest and highest velocity offset defined in the molfit file into account. This offers the possibility to describe lines which are located at the edges of the frequency ranges.
TelescopeSize: for single dish observations (Inter_Flag = False):TelescopeSizedescribes the size of telescope (in m), (default: 1); for interferometric observations (Inter_Flag = True):TelescopeSizedescribes the interferometric beam FWHM size (in arcsec), (default: 1).
BMIN: (optional) Beam minor axis length (in arsec) (default:None), see Sect. “Sub-beam description”.
BMAJ: (optional) Beam major axis length (in arsec) (default:None), see Sect. “Sub-beam description”.
BPA: (optional) Beam position angle (in degree) (default:None), see Sect. “Sub-beam description”.
Redshift: (optional) red shift (dimensionless), (default:None) NOT YET AVAILABLE!
Inter_Flag(True/False): defines, if single dish (False) or interferometric observations (True) are described, (default:False).
t_back_flag(True/False): defines, if the user defined background temperature \(T_{\rm bg}\) and temperature slope \(T_{\rm slope}\) given by the input parameterstBackandtSlopedescribe the continuum contribution completely (t_back_flag = True) or not (t_back_flag = False) (default:True).
tBackbackground temperature (in K), (default: 0).
tSlopetemperature slope (dimensionless), (default: 0).
BackgroundFileName(optional) path and name of file describing the background as function of frequency (in MHz), (default:None).
nH_flag: is deprecatedN_H(optional) Hydrogen column density (in cm\(^{-2}\)), (default: 0.0), see Sect. “Dust contribution”.
beta_dust(optional) spectral index for dust (dimensionless), (default: 0.0), see Sect. “Dust contribution”.
kappa_1300(optional) kappa at 1.3 mm (cm\(^2 \,\)g\(^{-1}\)), (default: 0.0), see Sect. “Dust contribution”.
DustFileName(optional) path and name of file describing optical depth of dust as function of frequency (in MHz), (default:None), see Sect. “Dust contribution”.
Te_ff: (optional) electronic temperature (in K) used for free-free continuum (default:None), see Sect. “Free-free contribution”.
EM_ff: (optional) emission measure (in pc cm\(^{-6}\)) used for free-free continuum (default:None), see Sect. “Free-free contribution”.
kappa_sync: (optional) kappa of energy spectrum of electrons for synchrotron continuum (electrons m\(^{-3}\) GeV\(^{-1}\)), (default:None), see Sect. “Synchrotron contribution”.
B_sync: (optional) magnetic field (in Gauss) used for synchrotron continuum (default:None), see Sect. “Synchrotron contribution”.
p_sync: (optional) energy spectral index (dimensionless) used for synchrotron continuum (default:None), see Sect. “Synchrotron contribution”.
l_sync: (optional) thickness of slab (in AU) used for synchrotron continuum (default:None), see Sect. “Synchrotron contribution”.
ContPhenFuncID: (optional) describes which phenomenological function is used to describe the continuum (default:None), see Sect. “Phenomenological continuum description”.
ContPhenFuncParam1: (optional) first parameter for phenomenological function describing the continuum (default:None), see Sect. “Phenomenological continuum description”.
ContPhenFuncParam2: (optional) second parameter for phenomenological function describing the continuum (default:None), see Sect. “Phenomenological continuum description”.
ContPhenFuncParam3: (optional) third parameter for phenomenological function describing the continuum (default:None), see Sect. “Phenomenological continuum description”.
ContPhenFuncParam4: (optional) fourth parameter for phenomenological function describing the continuum (default:None), see Sect. “Phenomenological continuum description”.
ContPhenFuncParam5: (optional) fifth parameter for phenomenological function describing the continuum (default:None), see Sect. “Phenomenological continuum description”.
iso_flag: use isotopologues (True/False). Ifiso_flagis set toTrueisotopologues defined in the iso ratio file are used (default:False).
IsoTableFileName: (has to be given only ifiso_flagis set toTrue): path and file name of an ASCII file including the iso ratios between certain molecules. The so-called “iso ratio file” defines the iso ratios between molecules. The ASCII file consists of three columns, where the first two columns indicates the molecules, respectively. The third column defines the ratio for both molecules. The columns are separated by blanks or tabs. So, the names of the molecules must not contain blanks. Comments have to be marked with the “%” sign, i.e. all characters on the right side of a “%” sign are ignored.The myXCLASSFit function offers the possibility to optimize the ratios between isotopologues as well. For that purpose, the user has to add two additional columns on the right indicating the lower and the upper limit of a certain ratio, respectively. A detailed description of the iso ratio file is given in Sect. “The iso ratio file”.
IsoRatioFile_Class: (optional) iso ratio file class object (default:None).
CollisionFileName: (has to be defined only for Non-LTE calculation of molecules): path and name of file describing collision partners (default:None), see Sect. “Non-LTE description of molecules”.
CollPartnerFile_Class: (optional) collision partner file class object, (default:None).
NumModelPixelXX: (optional) used for sub-beam modeling, describes the number of model pixels used in x-direction (default: 100), see Sect. “Sub-beam description”.
NumModelPixelYY: (optional) used for sub-beam modeling, describes the number of model pixels used in y-direction (default: 100), see Sect. “Sub-beam description”.
LocalOverlapFlag(True/False): (optional) take local overlap into account (default:False), see Sect. “Local-overlap of neighboring lines”.
NoSubBeamFlag(True/False): (optional) do not use sub-beam description into (default:False), see Sect. “Sub-beam description”, i.e. XCLASS uses the old source description, see Sect. “Simple 1-d description”.
DBFileName: (optional) path and name of a database file, which is used instead of the default filecdms_sqlite.db, (default:""").
dbFilename: (optional) alias ofDBFileName, (default:""").
DB_Class: (optional) database class object, (default:None).
RestFreq: rest frequency in MHz (default: 0). (If this parameter is unequal zero, the intensity is also plotted against velocity (in km s\(^{-1}\)) using parametervLSR, see above.
EmAbsFlag(True/False): (optional) flag indicating that emission and absorption as function of frequency are stored (default:False).
EmAbsPATH: (optional) path containing files describing functions emission / absorption (default:None).
CheckMoleculeInDB(True/False): (optional) check, if molecules included in molfit file have transitions within given frequency range (default:True).
ChannelIntegrationFlag(True/False): (optional, for full fitting method only) flag indicating channel integration (default:True).
OpacityFlag(True/False): (optional, for full fitting method only) flag indicating storage of opacities and intensities for each distance (ifLocalOverlapFlagis set toTrue) or for each component (ifLocalOverlapFlagis set toFalse) (default:False).
printFlag(True/False): (optional) flag for printing messages to the screen (default:True).
NoWarningFlag(True/False): flag for deactivating warning messages, (default:False).
UsemyXCLASSProgramFlag(True/False): flag for using old myXCLASS program, (default:False). Is set toFalseif cubes are computed. The flag is deprecated and will be removed in one of the next releases.
SpectrumOnlyFlag(True/False): (optional, only used with old myXCLASS program) flag indicating that only the spectra (and no optical depth etc.) are calculated, (default:False). The flag is deprecated and will be removed in one of the next releases.
BibTeXFlag(True/False): flag indicating that a BibTeX file calledreferences_db.bibis generated within the job directory describing all papers used to generate the molecular parameters of transitions within the selected frequency ranges.
NumberCores: (optional, only used if new old myXCLASS program is not used and if cubes are calculated) number of cores used for parallelization, (default:1).
ParallelizationMethod: (optional, only used if new old myXCLASS program is not used and if cubes are calculated) method used for parallelization, (default:process). Methods other thanprocesswill be available in one of the next versions.
Output parameters:
myXCLASSspectrum: contains the calculated myXCLASS spectrumIf
RestFreqis unequal zero, the myXCLASS function adds a column to the output parametermyXCLASSspectrumwhich contains the velocities. So, for a rest frequency unequal zero, the parametermyXCLASSspectrumrepresents a python array with three columns, where the first column describes the frequencies, the second column describes the velocities and the third column the corresponding intensities.
myXCLASSlog: (only generated by the old myXCLASS program) contains the corresponding log file.
myXCLASSTrans: (only generated by the old myXCLASS program, python) list containing the transition frequencies (in MHz) from the last myXCLASS run, the Doppler-shifted transition frequencies (in MHz), the corresponding absolute and integrated intensities (in K), the energy of the lower level (in K and K MHz), the upper state degeneracy, the Einstein A coefficient (in s\(^{-1}\)), and the molecule names within the defined range.
myXCLASSIntOptical: contains intensities and optical depths for each molecule and component. Depending on the local-overlap flag, i.e. if local-overlap is taken into account or not, the intensities and optical depths are given for each component or for each distance, see Sect. “Local-overlap of neighboring lines”.For the new calculation option (
SpectrumOnlyFlagis set toFalse), the flagOpacityFlaghas to be set toTrueas well.The intensities and optical depths are stored as follows:
myXCLASSIntOptical[0], contains all information about the intensitiesmyXCLASSIntOptical[0][i][0]contains the name of the \(i\)th moleculemyXCLASSIntOptical[0][i][1]contains the index of the component of the \(i\)th moleculemyXCLASSIntOptical[0][i][2]contains the intensities of the \(i\)th molecule as a 2D numpy array. Please note, the intensities for each core component \(\left[T_{\rm mb}^{i=1}\right]^{m,c}(\nu)\) and \(\left[T_{\rm mb}^{i>1}\right]^{m,c}(\nu)\) are given by different expression, depending if local-overlap is taken into account or not.If local-overlap is not taken into account, the intensities \(\left[T_{\rm mb}^{i=1}\right]^{m,c}(\nu)\) for components contained in the layer which is closest to the background is
\[\begin{split}\left[T_{\rm mb}^{i=1}\right]^{m,c}(\nu) &= \Bigg[\eta \left(\theta^{m,c}\right) \left[S^{m,c}(\nu) \left(1 - e^{-\tau_{\rm total}^{m,c}(\nu)}\right) + I_{\rm bg} (\nu) \left(e^{-\tau_{\rm total}^{m,c}(\nu)} - 1\right)\right]\Bigg] \\ & \qquad + \frac{1}{N_{\rm all}^{i=1}} \, \left[I_{\rm bg} (\nu) - J_\mathrm{CMB} \right],\end{split}\]where \(N_{\rm all}^{i=1}\) indicates the number of all components of all molecules located in this layer. The intensities of components \(\left[T_{\rm mb}^{\rm i > 1}\right]^{m,c}\) with smaller distances are given by
\[\left[T_{\rm mb}^{i > 1}\right]^{m,c}(\nu) = \Big[S^{m,c}(\nu) \left(1 - e^{-\tau_{\rm total}^{m,c}(\nu)}\right) + \left[T_{\rm mb}^{(i - 1)}\right]^{m,c} e^{-\tau_{\rm total}^{m,c}(\nu)}\Big].\]If local-overlap is taken into account, the expression for \(\left[T_{\rm mb}^{i=1}\right]^{m,c}(\nu)\) is slightly different
(87)\[\begin{split}\left[T_{\rm mb}^{i=1}\right]^{m,c}(\nu) &= \sum_{m,c \in i} \Bigg[\eta \left(\theta^{m,c}\right) \left[S^{m,c}(\nu) \left(1 - e^{-\tau_{\rm total}^{m,c}(\nu)}\right) + I_{\rm bg} (\nu) \left(e^{-\tau_{\rm total}^{m,c}(\nu)} - 1\right)\right]\Bigg] \\ & \qquad + \frac{1}{N_{\rm all}^{i=1}} \, \left[I_{\rm bg} (\nu) - J_\mathrm{CMB} \right].\end{split}\]Additionally, the intensities for components at smaller distances \(\left[T_{\rm mb}^{\rm i > 1}\right]^{m,c}\) are given by
(88)\[\left[T_{\rm mb}^{i > 1}\right]^{m,c}(\nu) = \sum_{m,c \in i} \Big[S^{m,c}(\nu) \left(1 - e^{-\tau_{\rm total}^{m,c}(\nu)}\right) + \left[T_{\rm mb}^{(i - 1)}\right]^{m,c} e^{-\tau_{\rm total}^{m,c}(\nu)}\Big].\]myXCLASSIntOptical[0][i][3]contains the integrated intensities of the \(i\)th molecule.
myXCLASSIntOptical[1], contains all information about the optical depths.myXCLASSIntOptical[1][i][0]contains the name of the \(i\)th moleculemyXCLASSIntOptical[1][i][1]contains the index of the component of the \(i\)th moleculemyXCLASSIntOptical[1][i][2]contains the optical depth of the \(i\)th molecule as a 2D numpy array. In general, the optical depths are given by\[\tau_l^{m,c}(\nu) = \sum_t \tau_t^{m,c}(\nu),\]where the sum runs over all transition \(t\) of a certain component (if local-overlap is not taken into account) or over all transitions of components, which belongs to the same distances (if local-overlap is taken into account).
myXCLASSJobDir: absolute path of the job directory created for the current run
Note, the user is free to define different names for the output
parameters.
myXCLASSFit
# import myXCLASSFit function
from xclass.task_myXCLASSFit import myXCLASSFitCore
myXCLASSFitCore(MolfitFile_Class = None,
MolfitFileName = "",
MolfitsFileName = "",
ObsXMLFile_Class = None,
ObsXMLFileName = "",
ObsDataFileName = "",
experimentalData = "",
observationalData = None,
IsoRatioFile_Class = None,
CollPartnerFile_Class = None,
AlgXMLFile_Class = None,
AlgorithmXMLFileName = "",
AlgorithmXMLFile = "",
Optimizer = "magix",
NumberIteration = 10,
NumberProcessors = 4,
FastFlag = False,
Parallelization_Class = None,
ClusterFileName = "",
clusterdef = "",
daskFlag = False,
ParallelizationMethod = "process",
JobDir = None,
ChannelIntegrationFlag = True,
OpacityFlag = False,
CheckMoleculeInDB = True,
TelescopeSize = None,
BMIN = None,
BMAJ = None,
BPA = None,
Inter_Flag = False,
Redshift = 0.0,
vLSR = 0.0,
RestFreq = 0.0,
FreqMin = None,
FreqMax = None,
t_back_flag = True,
tBack = 0.0,
tSlope = 0.0,
BackgroundFileName = "",
nH_flag = None, N_H = 0.0,
beta_dust = 0.0,
kappa_1300 = 0.0,
DustFileName = "",
Te_ff = None,
EM_ff = None,
kappa_sync = None,
B_sync = None,
p_sync = None,
l_sync = None,
ContPhenFuncID = None,
ContPhenFuncParam1 = None,
ContPhenFuncParam2 = None,
ContPhenFuncParam3 = None,
ContPhenFuncParam4 = None,
ContPhenFuncParam5 = None,
iso_flag = False,
IsoTableFileName = "",
CollisionFileName = None,
CopyRatesFileFlag = True,
EmAbsFlag = False,
EmAbsPATH = None,
NoSubBeamFlag = None,
NumModelPixelXX = 100,
NumModelPixelYY = 100,
LocalOverlapFlag = None,
DBFileName = "",
dbFileName = "",
PrintFlag = True,
ExtendedPlotFlag = False,
vel_low = -50.0,
vel_up = 50.0,
BibTeXFlag = False,
DictOutFlag = False,
DicOutFlag = None)
If input parameter JobDir is not defined, the myXCLASSFit
function creates for each run a so-called job directory located in
the run directory (path-of-XCLASS-run/myXCLASSFit/) where all
files created by the myXCLASSFit function are stored in. The
name of this job directory is made up of four components: The first
part consists of the phrase “job_” whereas the second and third
part describe the date (day, month, year) and the time stamp (hours,
minutes, seconds) of the function execution, respectively. The last
part indicates a so-called job ID which is composed of the so-called
PID followed by a four digit random integer number to create a really
unambiguous job number, e.g.
/path-of-XCLASS-run/myXCLASSFit/job__25-07-2013__12-02-03__189644/
Before the fit procedure starts, the function copies the molfit, the observational data, and xml-files to the myXCLASSFit job directory. The path(s) defined in the observational xml-file are adjusted so that these paths point to the current myXCLASSFit job directory.
Please note, that the original observational xml- and data files are not modified. Only the copy of the observational xml file located in the myXCLASSFit job directory is modified!
If no observational xml-file is defined, the myXCLASSFit function
creates an observational xml-file containing the settings for the
different parameters in the myXCLASSFit job directory. Additionally, the
function creates an ASCII file located in the myXCLASSFit job directory
containing the observational data given by the parameter
"observationalData". If the parameter "observationalData"
defines the path and the name of an observational data file, then the
data file is copied to the myXCLASSFit job directory as well. The path
and the name of this ASCII file is pasted in the created observational
xml file. Furthermore, the function creates automatically the algorithm
and i/o control xml files, need by MAGIX, so that the user does not need
to edit any xml-file.
Additionally, the myXCLASSFit function converts the column density \(N_{\rm tot}\) (and if the hydrogen column density \(N_H\) is given for each component) the hydrogen column density \(N_H\) as well to a log scale, i.e. these two densities are converted automatically to their log10 values to get a better fit. At the end of the fitting process, the log10 values are converted back to the original linear values. So, the MAGIX log files contain the log10 values of these parameters, whereas the input and output molfit file contain the linear values, respectively.
Input parameters:
MolfitFile_Class: (optional) molfit file class object, (default:None).
MolfitsFileName: path and name of the (extended) molfit file, see Sect. “The molfit file”
MolfitFileName: (optional) alias ofMolfitsFileName, (default:None).
AlgorithmXMLFileName: (optional) path and name of an algorithm xml-file describing settings for an algorithm or algorithm chain, see Sect. “Algorithm xml file”.NOTE, if the user specify a xml file, the number of iterations given by the parameter
NumberIterationis ignored. The number of iteration is then given by the xml file.
AlgorithmXMLFile: (optional) alias ofMolfitsFileName, (default:"").
NumberIteration: max. number of iterations (default: 10). (This parameter is ignored if a algorithm xml file is defined). If parameterOptimizer, see below, is set to'magix', the Levenberg-Marquardt algorithm, otherwise the trf algorithm is used.
AlgXMLFile_Class: (optional) alg. xml file class object, (default:None).
ObsXMLFile_Class: (optional) obs. xml file class object, (default:None).
ObsXMLFileName: (optional) path and name of an obs. xml file (recommended procedure), (default:""), see Sect. “Observational xml file” for a detailed description.Note, the name of the xml file has to end with
".xml".
ObsDataFileName: (optional) describes the observational data, if no obs. xml file is given, (default:None). There are two options for transferring the data to XCLASS:the parameter describes path and name of an ASCII file called observational data file, where the first column describe the frequency (in MHz) and the second column the brightness temperature (in Kelvin)
Please note, without using an obs. xml file, it is not possible to fit more than one data cube, simultaneously.
the parameter describes a 2d numpy array, where the first column describes the frequency (in MHz) and the second column the brightness temperature (in Kelvin).
experimentalData: (optional) alias ofObsDataFileName, (default:None)
observationalData: is deprecated.
The following parameters are needed, if parameter ObsDataFileName is
used instead of parameter ObsXMLFileName:
vLSR: source velocity (local standard of rest) in km s\(^{-1}\) (default: 0).
TelescopeSize: for single dish observations (Inter_Flag = False):TelescopeSizedescribes the size of telescope (in m), (default: 1); for interferometric observations (Inter_Flag = True):TelescopeSizedescribes the interferometric beam FWHM size (in arcsec), (default: 1).
BMIN: (optional) Beam minor axis length (in arcsec) (default:None)
BMAJ: (optional) Beam major axis length (in arcsec) (default:None)
BPA: (optional) Beam position angle (in degree) (default:None)
Redshift: (optional) red shift (dimensionless), (default:None)
Inter_Flag(True/False): defines, if single dish (False) or interferometric observations (True) are described, (default:False).
t_back_flag(True/False): defines, if the user defined background temperature \(T_{\rm bg}\) and temperature slope \(T_{\rm slope}\) given by the input parameterstBackandtSlopedescribe the continuum contribution completely (t_back_flag = True) or not (t_back_flag = False) (default:True).
tBackbackground temperature (in K), (default: 0).
tSlopetemperature slope (dimensionless), (default: 0).
BackgroundFileName(optional) path and name of file describing the background as function of frequency (in MHz), (default:None).
nH_flag: is deprecatedN_H(optional) Hydrogen column density (in cm\(^{-2}\)), (default:None).
beta_dust(optional) spectral index for dust (dimensionless), (default:None).
kappa_1300(optional) kappa at 1.3 mm (cm\(^2 \,\)g\(^{-1}\)), (default:None).
DustFileName(optional) path and name of file describing dust contribution (default:"").
Te_ff: (optional) electronic temperature (in K) used for free-free continuum (default:None).
EM_ff: (optional) emission measure (in pc cm\(^{-6}\)) used for free-free continuum (default:None).
kappa_sync: (optional) kappa of energy spectrum of electrons for synchrotron continuum (electrons m\(^{-3}\) GeV\(^{-1}\)), (default:None)
B_sync: (optional) magnetic field (in Gauss) used for synchrotron continuum (default:None).
p_sync: (optional) energy spectral index (dimensionless) used for synchrotron continuum (default:None).
l_sync: (optional) thickness of slab (in AU) used for synchrotron continuum (default:None).
ContPhenFuncID: (optional) describes which phenomenological function is used to describe the continuum (default:None).
ContPhenFuncParam1: (optional) first parameter for phenomenological function describing the continuum (default:None).
ContPhenFuncParam2: (optional) second parameter for phenomenological function describing the continuum (default:None).
ContPhenFuncParam3: (optional) third parameter for phenomenological function describing the continuum (default:None).
ContPhenFuncParam4: (optional) fourth parameter for phenomenological function describing the continuum (default:None).
ContPhenFuncParam5: (optional) fifth parameter for phenomenological function describing the continuum (default:None).
iso_flag: use isotopologues (True/False). If the flag is set toTrueisotopologues defined in the iso ratio file are taken into account, (default:False).
IsoTableFileName(has to be given only ifiso_flagis set toTrue): path and name of an ASCII file including the iso ratios between certain molecules. A detailed description of the iso ratio file is given in Sect. “The iso ratio file”. If no file name is given (default), the so-called iso-flag is set toFalse.
NOTE, if the parameter
observationalDatadefines a relative path, the path has to be defined relative to the current working directory!
IsoRatioFile_Class: (optional) iso ratio file class object (default:None).
CollisionFileName: (has to be defined only for Non-LTE calculation): path and name of file describing collision partners (default:"").
CollPartnerFile_Class: (optional) collision partner file class object, (default:None).
DBFileName: (optional) path and name of a database file, which is used instead of the default filecdms_sqlite.db, (default:""").
dbFilename: (optional) alias ofDBFileName, (default:""").
NumModelPixelXX: (optional) used for sub-beam modeling, describes the number of pixels used in x-direction (default: 100).
NumModelPixelYY: (optional) used for sub-beam modeling, describes the number of pixels used in y-direction (default: 100).
LocalOverlapFlag(True/False): (optional) take local overlap into account (default:False), see Sect. “Local-overlap of neighboring lines”.
NoSubBeamFlag(True/False): (optional) do not use sub-beam description (default:False), see Sect. “Sub-beam description”.
EmAbsFlag(True/False): (optional) flag indicating that emission and absorption as function of frequency are stored (default:False).
EmAbsPATH: (optional) path containing files describing functions emission / absorption (default:None).
The following parameter is needed to add a velocity axis to the output
array model_values, see below:
RestFreq: rest frequency in MHz (default: 0).
The myXCLASSFit function offers the possibility to use two different
optimization packages: 'magix' or 'scipy' [Virtanen_et_al_2020],
providing different optimization algorithms, see Sect. .
Optimizer: package used for optimization ('magix'and'scipy') (default:'magix').Note, the optimization package used different fit algorithms. A list of the fit algorithms for each optimization package and the corresponding settings in the algorithm xml file is given in Sect. “Algorithm xml file”.
FastFlag(True/False): (optional) fast flag for usage of lite XCLASS routines (default:False).
NumberProcessors: (optional) max. number of cores, only used if no algorithm xml file is applied (default: 4).
ParallelizationMethod: (optional) method used for parallelization, (default:process). Methods other thanprocesswill be available in one of the next versions.
daskFlag: will be used in later releases.
ClusterFileName: (optional, only used for parallelization methods other thanprocess) path and name of a file describing a cluster, (default:"").In order to start the myXCLASSFit function on a multi-computer cluster, a so-called “cluster configuration file” is required. The following requirements are necessary for all nodes which are included in the cluster:
Password-less ssh access from the controller (user) machine to all other nodes included in the cluster.
NOTE: This is not necessary when using only “localhost”, i.e. if the cluster is deployed only on the machine where XCLASS is running.
XCLASS interface and gfortran (with OpenMP/OpenMPI) is installed on all nodes of the cluster.
The XCLASS job directory must be located in a shared file-system, accessible from all nodes comprising the cluster, and mounted in the same path of the file-system.
NOTE: If the cluster contains computers with different operating systems, e.g. Linux and MAC OS, the user has to define the different paths of the XCLASS interface directories.
Example of a “cluster configuration file” used for the myXCLASS function:
# cluster configuration file # node number of cores path of XCLASS interface anu, 2 lugal, 3
In the example described above, the nodes “anu” and “lugal” will be used with two (three) cores, respectively.
clusterdef: (optional) alias ofClusterFileName, (default:"")
Parallelization_Class: (optional) parallel class object, (default:None).
ChannelIntegrationFlag(True/False): (only for optimizer'scipy') flag indicating usage of channel integration (default:False).
CheckMoleculeInDB(True/False): (optional) check, if molecules included in molfit file have transitions within given frequency range (default:True).
OpacityFlag(True/False): (optional, for full fitting method only) flag indicating storage of opacities and intensities for each distance (ifLocalOverlapFlagis set toTrue) or for each component (ifLocalOverlapFlagis set toFalse) (default:False). If the OpacityFlag is set toTrue, the job directory contains an additional subdirectory calledintensities_and_opacities, which includes files describing the intensities and opacities per component. In addition, if theDictOutFlagis set toTrue, the returned dictionary contains two different items"IntensityPerComponent"and"OpacityPerComponent", which describe the intensities and opacities as python lists.
printFlag(True/False): (optional) flag for printing messages to the screen (default:True).
ExtendedPlotFlag(True/False): (optional) flag for creating transitions and extended plots (default:False). The extended plots, which describe for each obs. data file, the observational data together with the final fit, are stored within a directory calledplot_myXCLASSwithin the job directory. The so-called transition diagrams, describe for each transition of each molecule used in the molfit file and its isotopologues, the observed data together with the overall model, i.e. the fit that considers all molecules and their isotopologues, and the model for the current molecule. These transition plots are contained in the directoryplot_molecule-contributionlocated in the current job directory. In addition to the corresponding pdf file, each subdirectory also contains the molfit and, if applicable, the collision partner file, which was used to create the molecule model in the respective transition plots .
vel_low: (optional, only used, ifExtendedPlotFlagis setto
True) lower velocity limit used for the transition plots, (default: -50.0).
vel_up: (optional, only used, ifExtendedPlotFlagis setto
True) upper velocity limit used for the transition plots, (default: 50.0).
BibTeXFlag(True/False): flag indicating that a BibTeX file calledreferences_db.bibis generated within the job directory describing all papers used to generate the molecular parameters of transitions within the selected frequency ranges.
Finally, the results of the myXCLASSFit function can be exported in two different ways
DictOutFlag(True/False): flag indicating, that all return parameters are collected in a dictionary (default:True), see below.
DicOutFlag: (optional) alias ofDictOutFlag, (default:"")
Output parameters:
Please note, the user is free to define different names for the output parameters.
If the input parameter DictOutFlag is set to False, the
myXCLASSFit function returns the results to the following three
output parameters:
input_file: the contents of the molfit file containing the parameters for the best fit.
model_values: the model function values for each data point, which correspond to the best fit.If
RestFreqis unequal zero, the myXCLASSFit function adds a column to the output parametermodel_valueswhich contains the velocities. So, for a rest frequency unequal zero, the parametermodel_valuesrepresents a python array with three columns, where the first column describes the frequencies, the second column describes the velocities and the third column the corresponding intensities.
JobDir: absolute path of the job directory created for the current run.
Please note, if more than one algorithm is used, the output parameters
input_file and model_values contain all input files and model
function values of all results, i.e. input_file[0] contains the
molfit file and model_values[0] the model function values for the
first applied fit algorithm.
The order of the results, i.e. the name of the molfit file, which
correspond to input_file[0], is printed to the screen.
If the input parameter DictOutFlag is set to True, the
myXCLASSFit function exports the results to a python dictionary,
e.g. OutDict, containing the following items:
"JobDir": describes path of the current job directory, i.e.OutDict["JobDir"]is the path and name of the current job directory.
"MolfitFileNames": is a python list indicating all used molfit file names, describing the optimized fit parameters after the application of a certain optimization algorithm.
"MolfitFilesContent": list of the corresponding molfit file contents,
"IsoRatioFileNames": list of iso ratio file names indicating the optimized fit parameters after the application of a certain optimization algorithm,
"IsoRatioFilesContent": list of the corresponding iso ratio file contents,
"CollisionPartnerFileNames": list of collision partner file names after the application of a certain optimization algorithm,
"CollisionPartnerFilesContent": list of the corresponding collision partner file contents,
"ListOfSpectraFileNames": list of spectra file names after the application of a certain optimization algorithm,
"ListOfSpectraFilesContent": list of numpy arrays describing the optimized spectra after the application of a certain optimization algorithm,
"BestMolfitFileName": name of the molfit file of the best fit
"BestMolfitFileContent": content of the molfit file corresponding to the best fit,
"BestIsoRatioFileName": name of the iso ratio file of the best fit,
"BestIsoRatioFileContent": content of the iso ratio file which belongs to the best fit,
"BestCollisionPartnerFileName": name of the collision partner file, which corresponds to the best fit,
"BestCollisionPartnerFileContent": content of the collision partner file for the best fit,
"BestSpectraFileNames": path and name of file describing the spectra of the best result,
"BestSpectraFilesContent": spectra of best fit.
Note, if only one algorithm is applied, the content of the items
"MolfitFileNames" and "BestMolfitFileName" etc. are identical.
myXCLASSMapFit
# import myXCLASSMapFit function
from xclass.task_myXCLASSMapFit import myXCLASSMapFitCore
myXCLASSMapFitCore(MolfitFile_Class = None,
MolfitFileName = "",
MolfitsFileName = "",
ObsXMLFile_Class = None,
ObsXMLFileName = "",
ObsDataFileName = "",
experimentalData = "",
observationalData = None,
IsoRatioFile_Class = None,
CollPartnerFile_Class = None,
AlgXMLFile_Class = None,
AlgorithmXMLFileName = "",
AlgorithmXMLFile = "",
NumberIteration = 10,
NumberProcessors = 4,
FastFitFlag = False,
FullFitFlag = None,
Parallelization_Class = None,
ClusterFileName = "",
clusterdef = "",
daskFlag = False,
ParallelizationMethod = "process",
ChannelIntegrationFlag = True,
OpacityFlag = False,
CheckMoleculeInDB = True,
RemoveContinuumFlag = False,
ZipDirFlag = False,
PrintFitStatusFlag = True,
ParameterMapDir = None,
ParamMapIterations = 1,
ParamSmoothMethodParam = 1.0,
ParamSmoothMethod = "uniform",
DoParameterEstimationFlag = False,
InternalParameterList = {},
JobDir = None,
NameOfFunction = None,
regionFileName = "",
FITSImageMaskFileName = "",
Threshold = None,
EstimateContinuumFlag = False,
EstimateConNumParts = 20,
UsePreviousResults = None,
DataSmoothValue = 0.0,
TelescopeSize = None,
BMIN = None,
BMAJ = None,
BPA = None,
Inter_Flag = False,
Redshift = None,
vLSR = 0.0,
RestFreq = 0.0,
FreqMin = None,
FreqMax = None,
t_back_flag = True,
tBack = 0.0,
tSlope = 0.0,
nH_flag = None,
N_H = 0.0,
beta_dust = 0.0,
kappa_1300 = 0.0,
DustFileName = "",
BackgroundFileName = "",
Te_ff = None,
EM_ff = None,
kappa_sync = None,
B_sync = None,
p_sync = None,
l_sync = None,
ContPhenFuncID = None,
ContPhenFuncParam1 = None,
ContPhenFuncParam2 = None,
ContPhenFuncParam3 = None,
ContPhenFuncParam4 = None,
ContPhenFuncParam5 = None,
iso_flag = False,
IsoTableFileName = "",
CollisionFileName = "",
NumModelPixelXX = 100,
NumModelPixelYY = 100,
LocalOverlapFlag = False,
NoSubBeamFlag = False,
EmAbsFlag = False,
EmAbsPATH = None,
DictOutFlag = False,
BibTeXFlag = False,
DBFileName = "",
dbFileName = "",
PrintFlag = True)
In contrast to the myXCLASSFit function, the myXCLASSMapFit
function fits a complete (FITS) data cube [1] Here, the
myXCLASSMapFit function reads in the data cube(s), extracts the
spectra for each pixel and fits these spectra separately using an
user-defined optimization algorithm. The myXCLASSMapFit function
provides two different fitting procedures: "fast" and "full".
The used optimization algorithms have to be described in an algorithm
xml file, see Sect. “Algorithm xml file”. It is also
possible to combine algorithms, i.e. to start with a global optimization
algorithm and use the results of this algorithm as input for a local
optimization algorithm such as "trf".
If an algorithm xml file contains a definition for the number of
processors (cores) for an optimization algorithm, this number now
defines the number of pixels that are fitted simultaneously. If
more processor cores are defined than are available on the computer,
XCLASS reduces the number to the number of available processors
(cores) minus one. If no algorithm xml file is used, the
myXCLASSMapFit function uses the trf algorithm with four
cores and performs no more than 50 iterations.
With the "fast" method, the myXCLASSMapFit function uses
simplified versions of the XCLASS core functions to compute the
synthetic spectra used for the optimization algorithms. These
simplified functions are able to compute the contribution of
molecules and recombination lines in LTE in conjunction with
dust, free-free, and synchrotron continuum contributions.
Additionally, it is possible to arrange the components
according to the distance parameter along the line-of-sight, see
Sect. “Stacking”. In contrast to the
"full" method it is not possible to describe molecules and
recombination lines in Non-LTE Sect. “Non-LTE description of molecules”
as well as sub-beam description Sect. “Sub-beam description.
Compared to the "fast" method mentioned above, the "full"
fitting procedure is computationally more complex, so the fitting
takes a bit more time.
The files produced by a myXCLASSMapFit function call are copied
to a so-called “job-directory”. If input parameter JobDir is
not set, a new job directory located in the XCLASS interface working
directory path-of-XCLASS/run/myXCLASS/ is created, where the name
of a job directory is made up of four components: The first part of
the name consists of the phrase "job_" whereas the second part
describes the date (day, month, year), the third part the time stamp
(hours, minutes, seconds) of the function execution. The last part
describes a so-called job ID which is composed of the so-called PID
followed by a four digit random integer number to create a really
unambiguous job number, e.g.
/path-of-XCLASS-run/myXCLASSMapFit/job__25-07-2013__12-02-03__189644932/.
At the end of the whole fit procedure, the myXCLASSMapFit function
creates FITS images for each optimized parameter, where each pixel
corresponds to the value of the optimized parameter taken from the best
fit for this pixel. The name of each parameter FITS image consists of
the phrase "BestResult" followed by the name of the parameter,
e.g. T_rot for the rotation temperature etc, the name of the
corresponding molecule, and finally by the index of the related
component. For example, the FITS file
BestResult___parameter__T_rot___molecule__CH3OH_v=0____component__1.fits,
describes the optimized excitation (rotation) temperature of the first component of molecule “CH\(_3\)OH,v=0” for each pixel of the selected fit region. In addition to the parameter maps, the myXCLASSMapFit function creates one FITS image, where each pixel corresponds to the \(\chi^2\) value of the best fit, i.e. the quality of the fit, for the corresponding pixel.
Furthermore, the myXCLASSMapFit function creates FITS cubes for
each used algorithm and fitted data cube, where each pixel contains
the modeled spectrum. The name of these FITS cubes ends with the
expression "__model.out.fits". (If more than one fit algorithm was
used, this expression is extended by the respective name of the
optimization algorithm). If more than one frequency range was used for
a given FITS cube, the name of the corresponding FITS cube describing
the synthetic spectra also contains the frequency limits (in MHz),
i.e.
name-of-FITS-file__MinFreq_-_MaxFreq__MHz__model.out.fits.
The FITS files, which end with ".chi2.out.fits", describe the
\(\chi^2\) function for each pixel.
In addition, the myXCLASSMapFit function converts the following
parameters N_tot, EM_ff, EM_RRL, N_e,
nHcolumn_cont_dust, nHcolumn, kappa_sync, B_sync, and
l_sync to a log scale, i.e. these parameter values are converted
automatically to their log10 values to get a better fit. At the end of
the fitting process, the log10 values are converted back to the
original linear values. So, the MAGIX log files contain the log10
values of these parameters, whereas the input and output molfit files
contain the linear values, respectively.
The myXCLASSMapFit function offers the possibility to use a
so-called cluster consisting of at least two computers. The nodes
(computers) of the cluster are defined by the input parameter
clusterdef, see below. Unfortunately, due to problems with the
used parallelization packages, cluster computation with more than one
computer is currently not possible.
The myXCLASSMapFit function can handle FITS cubes (and images) for background and dust files as well, i.e. the user can specify the path and name of a FITS cube describing the background (or dust) spectrum for each pixel and frequency range. If a FITS image is provided instead of a cube, the corresponding spectrum is automatically generated for the required frequency points. If a background spectrum is given in Jy/beam, XCLASS automatically converts the intensities to Kelvin. The FITS cubes (or images) need not to cover the entire region. Where no background (or dust) file is given, no background (or dust) spectrum is used.
Input parameters:
MolfitFile_Class: (optional) molfit file class object, (default:None).
MolfitsFileName: path and name of the (extended) molfit file, see Sect. “The molfit file”
MolfitFileName: (optional) alias ofMolfitsFileName, (default:None).
FastFitFlag(True/False): (optional) flag indicating fast fitting, (default:True).
FullFitFlag: is deprecated.
AlgorithmXMLFileName: only necessary, if the user wants to use another fit algorithm (than"trf") for fitting. Therefore, the path and name of an algorithm xml-file describing settings for an algorithm or algorithm tree has to be given, see Sect. “Algorithm xml file”NOTE, if the user specify an xml file, the number of iterations given by the parameter
"NumberIteration"and the number of cores (NumberProcessors) are ignored. The number of iteration is then given by the algorithm xml file.
AlgorithmXMLFile: (optional) alias ofAlgorithmXMLFileName.
NumberIteration: (only necessary if no algorithm xml and no cluster file is given) max. number of iterations (default: 50).
AlgXMLFile_Class: (optional) alg. xml file class object, (default:None).
NumberProcessors: (optional) max. number of cores, only used if no algorithm xml file is applied (default: 4).
ParallelizationMethod: (optional) method used for parallelization, (default:process). Methods other thanprocesswill be available in one of the next versions.
daskFlag: will be used in later releases.
ClusterFileName: (optional, only used for parallelization methods other thanprocess) path and name of a file describing a cluster, (default:"").
clusterdef: (optional) alias ofClusterFileName, (default:"")
Parallelization_Class: (optional) parallel class object, (default:None).
"JobDir": describes path of the current job directory, i.e.OutDict["JobDir"]is the path and name of the current job directory.
NameOfFunction: (optional) name of calling function (default:None).
regionFileName: (optional) path and name of a so-called ds9 [2] region file, (default:"").
FITSImageMaskFileName: (optional) path and name of a FITS image used for masking. Here, each pixel with a real numerical value (except 0), i.e. which is not NaN or 0, is fitted, (default:""). The value of the selecting pixel is used to weight the contribution of the corresponding pixel to the total \(\chi^2\) value.
EmAbsFlag(True/False): (optional) flag indicating that emission and absorption as function of frequency are stored (default:False).
EmAbsPATH: (optional) path containing files describing functions emission / absorption (default:None).
ChannelIntegrationFlag(True/False): (optional, for full fitting method only) flag indicating channel integration (default:False).
CheckMoleculeInDB(True/False): (optional) check, if molecules included in molfit file have transitions within given frequency range (default:True).
OpacityFlag(True/False): (optional, for full fitting method only) flag indicating storage of opacities and intensities for each distance (ifLocalOverlapFlagis set toTrue) or for each component (ifLocalOverlapFlagis set toFalse) (default:False). If the OpacityFlag is set toTrue, the job directory contains an additional subdirectory calledintensities_and_opacities, which includes files describing the intensities and opacities per component. In addition, if theDictOutFlagis set toTrue, the returned dictionary contains two different items"IntensityPerComponent"and"OpacityPerComponent", which describe the intensities and opacities as python lists.
RemoveContinuumFlag(True/False): (optional) flag forremoving continuum, (default:
False).
ZipDirFlag(True/False): (optional) flag indicatingif all files in the job directory are compresses, (default:
False).
PrintFitStatusFlag(True/False): (optional) flagindicating output of fit status message, (default:
True).
printFlag(True/False): (optional) flag for printing messages to the screen (default:True).
Threshold: defines a threshold for a pixel. If the spectrum of a pixel has an max. intensity lower than the value defined by this parameter the pixel is not fitted (ignored), (default:None).Please note, the value for the
Thresholdparameter has to be given in the same unit than the spectrum in the FITS file.NOTE, if an observational xml-file is used, the user can define different threshold values for each frequency range by using the tag
<Threshold>. Please note, that the tag<Threshold>has to be given for each frequency range. If tag<Threshold>is defined in the xml file, the value of the input parameterThresholdis ignored.
EstimateContinuumFlag(True/False): (optional) flag indicating continuum estimation using thestatcont[SanchezMonge_et_al_2018] package (default:False).
EstimateConNumParts: (optional) number of sub-parts used for continuum estimation to describe non-constant continuum levels as well, (default: 20). In order to estimate the continuum level, the obs. data (of each pixel) are divided into \(n\) (=EstimateConNumParts) equidistant parts to whichstatcontis applied. Now, the obs. data within each part is replaced by the corresponding continuum level located at the center of each part. Afterwards, the well-known expression\[I_{\rm bg} (\nu) = T_{\rm bg} \cdot \left(\frac{\nu}{\nu_{\rm min}} \right)^{T_{\rm slope}}\]is fitted to these data points to achieve the parameters \(T_{\rm bg}\) and \(T_{\rm slope}\).
EstimateConNumEntireCubeFlag(True/False): (optional) flag indicating that the continuum estimation is done for the entire spectral range of the corresponding FITS cube and not only for the local frequency range, (default:True).
DataSmoothValue: (optional) smooth factor (=0, no smoothing is performed) used for spectral smoothing, (default: 0.0). Here each pixel is convolved with a Gaussian whose width is given by this parameter.
UsePreviousResults: is deprecated. used.
ParamMapIterations: (optional) number of iterations to smooth parameter maps, (default: 1).A value greater 1 indicates a parameter map smoothing procedure, i.e. after finishing the first myXCLASSMapFit run, the created parameters maps are smoothed using the method and parameter described by
ParamSmoothMethodParamandParamSmoothMethod. The values of the smoothed parameter maps are used in a second myXCLASSMapFit run, where the parameters in the molfit files at each fitted pixel position are updated by the values of the smoothed maps. The whole procedure is repeated as many times as indicated byParamMapIterations. The names of the FITS files describing the parameter maps are extended by the phrase"_____Iter__n.fits", where \(n\) indicates the current iteration, The FITS files ending with \(n=\)ParamMapIterationsdescribe the final parameter maps.
ParamSmoothMethodParam: (optional) parameter used for smoothing the parameter map(s), see input parameterParamSmoothMethod, (default: 1.0)
ParamSmoothMethod: (optional) SciPy method ("gaussian","uniform","median") used for parameter map smoothing, (default:"uniform").For
"gaussian", SciPy’s multidimensional Gaussian filter [3] is used, whereParamSmoothMethodParamindicates the standard deviation for Gaussian kernel (in pixel) which is assumed to be equal for both axes.For
"uniform", SciPy’s multidimensional uniform filter [4] is applied, whereParamSmoothMethodParamdescribes the sizes (in pixel) of the uniform filter which is equal for both axes.For
"median", SciPy’s multidimensional median filter [5] is utilized, whereParamSmoothMethodParamdefines the sizes (in pixel) of the median filter which is equal for both axes.
DoParameterEstimationFlag: is deprecated
InternalParameterList: is deprecated
ParameterMapDir: (optional) path of a directory containing FITS images describing parameter maps to update the parameter defined in the molfit file (or iso ratio file or file describing collision partners) for each pixel (default:None). By using this parameter the user can modify the content of the molfit file, defined by input parameterMolfitsFileName, for each pixel. Therefore, the name of the FITS image file describing a parameter map has to consists of the name of the molfit parameter (e.g.N_tot), which has to be updated, the name of the molecule, and the corresponding component number. The different settings has to be separated by"___". The definition of the parameter name has to start with the phrase"parameter__"followed by the name of the parameter, i.e. the name has to be"source_size","source_radius","s_radius","source_center_x","s_off_x","source_center_y","s_off_y","T_rot","T_e","Te_ff","T_cont_dust","T_dOff","T_dSlope","N_tot","EM_RRL","N_e","EM_ff","nHcolumn_cont_dust","nHcolumn","V_width,"V_width_Gauss","V_width_Lorentz","V_width_Voigt_G","V_width_Voigt_L","V_width_Horn_W","V_width_Horn_H","V_off","beta_cont_dust","kappa_cont_dust","kappa","beta","kappa_sync","B_sync","p_sync","l_sync","ContFuncID_param_1","ContPhenFuncParam1","ContFuncID_param_2","ContPhenFuncParam2","ContFuncID_param_3","ContPhenFuncParam3","ContFuncID_param_4","ContPhenFuncParam4","ContFuncID_param_5","ContPhenFuncParam5". Additionally, the name of the molecule has to start with phrase"molecule__"followed by the name of the molecule, where the characters (";",",","'","(",")","#") has to be replaced by"_". Finally, the definition of the component number has to start with"component__"followed by an integer number. For example, we want to use a FITS image describing the velocity offsets for the first component of CH\(_3\)CN,v=0 for each pixel of the selected region. The name of the corresponding FITS file is"parameter__V_off___molecule__CH3CN_v=0____component__1.fits".In order to define iso ratios for a specific isotopologue, the name of the parameter has to be set to
"isoratio"and the name of the molecule has to indicate the corresponding isotopologue. The part, which describes the component has to be removed, e.g."parameter__isoratio___molecule__C-13H3CN_v=0_.fits".Parameter maps describing the density of a certain collision partner have a more complex naming scheme: The name of the parameter has to be set to
"collpartdensity". Additionally, the name of the molecule has to indicate the name of the molecule in the molfit file, which is described in Non-LTE. In contrast to other parameters, the phrase"coll-partner__"followed by the name of the collision partner ("H2","pH2","oH2","e-","H","He", or"H+") has to be added after the definition of the molecule. If the density is defined for a specific distance, this distance has to be defined by the phrase"dist__"followed by distance, e.g.parameter__collpartdensity___molecule__CH3CN_v=0____coll-partner__H2___dist__1.0.fits
ObsXMLFileName: (optional) path and name of an obs. xml file (recommended procedure), (default:""), see Sect. “Observational xml file” for a detailed description.Note, the name of the xml file has to end with
".xml".
ObsDataFileName: (optional) path and name of a file describing a 3d or 4d FITS data cube.Note, the name of the data file has to end with
".fits".If the intensities of a data cube are given in Jy/beam, the myXCLASSMapFit function converts the intensities to Kelvin using the following expression [6] [Rohlfs_Wilson_1996]
\[T_{\rm mb} (\nu) = \frac{c^2 \cdot S_\nu}{2 \, k_B \, \nu^2 \, \Omega} \cdot 10^{-26},\]where \(S_\nu\) describes the intensity at a certain frequency \(\nu\) in Jy/beam and \(c\) and \(k_B\) represents the speed of light and the Boltzmann constant, respectively. Additionally, \(\Omega\) indicates the solid angle of the beam and can be expressed under the assumption that the beam can be approximated by a two-dimensional, elliptical Gaussian function \(G(\alpha, \beta)\)
\[\Omega = \iint G(\alpha, \beta) \, d\alpha \, d \beta = {\rm BMAJ} \cdot {\rm BMIN} \cdot \frac{\pi}{4 \, \ln(2)},\]where BMAJ and BMIN indicate the FWHM of the beam major and minor axis length in arcsec, respectively [7].
Please note, without using an obs. xml file, it is not possible to fit more than one data cube, simultaneously.
experimentalData: (optional) alias ofObsDataFileName, (default:None)
observationalData: is deprecated.
ObsXMLFile_Class: (optional) obs. xml file class object, (default:None).
IsoRatioFile_Class: (optional) iso ratio file class object (default:None).
CollPartnerFile_Class: (optional) collision partner file class object, (default:None).
DictOutFlag(True/False): flag indicating, that all return parameters are collected in a dictionary (default:True), see below.
DicOutFlag: (optional) alias ofDictOutFlag, (default:"")
BibTeXFlag(True/False): flag indicating that a BibTeX file calledreferences_db.bibis generated within the job directory describing all papers used to generate the molecular parameters of transitions within the selected frequency ranges.
The following parameters are needed, if parameter ObsDataFileName is
used instead of parameter ObsXMLFileName:
FreqMin: start frequency of simulated spectrum (default: 0).Please note, in contrast to previous XCLASS releases, the min. and max. frequency has to be defined! Additionally, the step size, i.e. the difference between two frequencies is taken from the FITS header.
FreqMax: end frequency of simulated spectrum (default: 0).
vLSR: velocity (local standard of rest) in km s\(^{-1}\) (default: 0) used in the calculation of the synthetic spectra
TelescopeSize: for single dish observations (Inter_Flag = False):TelescopeSizedescribes the size of telescope (in m), (default: 1); for interferometric observations (Inter_Flag = True):TelescopeSizedescribes the interferometric beam FWHM size (in arcsec), (default: 1).
BMIN: (optional) beam minor axis length (in arcsec), (default:None).Note, in the header of FITS files, the beam minor and major axis lengths are often given in degrees and not in arcsec, i.e. the values have to be divided by 3600.
BMAJ: (optional) beam major axis length (in arcsec), (default:None).
BPA: (optional) beam position angle (in degree), (default:None).
Redshift: (optional) red shift (dimensionless), (default:None).
Inter_Flag(True/False): defines, if single dish (False) or interferometric observations (True) are described, (default:False).
t_back_flag(True/False): defines, if the user defined background temperature \(T_{\rm bg}\) and temperature slope \(T_{\rm slope}\) given by the input parameterstBackandtSlopedescribe the continuum contribution completely (t_back_flag = True) or not (t_back_flag = False) (default:True).
tBackbackground temperature (in K), (default: 0.0).
tSlopetemperature slope (dimensionless), (default: 0.0).
BackgroundFileName(optional) path and name of file describing the background as function of frequency (in MHz), (default:None).
nH_flag: is deprecatedN_H(optional) Hydrogen column density (in cm\(^{-2}\)), (default: 0.0).
beta_dust(optional) spectral index for dust (dimensionless), (default: 0.0).
kappa_1300(optional) kappa at 1.3 mm (cm\(^2 \,\)g\(^{-1}\)), (default: 0.0).
DustFileName(optional) path and name of file describing dust contribution, (default:None).
Te_ff: (optional) electronic temperature (in K) used for free-free continuum, (default:None).
EM_ff: (optional) emission measure (in pc cm\(^{-6}\)) used for free-free continuum, (default:None).
kappa_sync: (optional) kappa of energy spectrum of electrons for synchrotron continuum (electrons m\(^{-3}\) GeV\(^{-1}\)), (default:None).
B_sync: (optional) magnetic field (in Gauss) used for synchrotron continuum, (default:None).
p_sync: (optional) energy spectral index (dimensionless) used for synchrotron continuum, (default:None).
l_sync: (optional) thickness of slab (in AU) used for synchrotron continuum, (default:None).
ContPhenFuncID: (optional) describes which phenomenological function is used to describe the continuum, (default:None).
ContPhenFuncParam1: (optional) first parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam2: (optional) second parameter for phenomenological function describing the continuum,(default:None).
ContPhenFuncParam3: (optional) third parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam4: (optional) fourth parameter for phenomenological function describing the continuum,(default:None).
ContPhenFuncParam5: (optional) fifth parameter for phenomenological function describing the continuum, (default:None).
iso_flag: use isotopologues (True/False). Ifiso_flagis set toTrueisotopologues defined in the iso ratio file are used, (default:False).
IsoTableFileName(has to be given only ifiso_flagis set toTrue): path and name of an ASCII file including the iso ratios between certain molecules. A detailed description of the iso ratio file is given in Sect. “The iso ratio file”. If no file name is given (default), the so-called iso-flag is set toFalse.NOTE, if the given path is a relative path, i.e. the path does not start with
/, this path has to be defined relative to the current working directory!
CollisionFileName: (has to be defined only for Non-LTE calculation): path and name of file describing collision partners, (default:None).
NumModelPixelXX: (optional) used for sub-beam modeling, Sect. “Sub-beam description”, describes the number of pixels used in x-direction, (default: 100).
NumModelPixelYY: (optional) used for sub-beam modeling, Sect. “Sub-beam description”, describes the number of pixels used in y-direction, (default: 100).
LocalOverlapFlag(True/False): (optional) flag indicates if local-overlap description is used or not, see Sect. “Local-overlap of neighboring lines”, (default:False).
NoSubBeamFlag(True/False): (optional) do not use sub-beam description (default:False), see Sect. “Sub-beam description”, i.e. XCLASS uses the old source description, see Sect. “Simple 1-d description”.
Output parameters:
JobDir: absolute path of the job directory created for the current run.
myXCLASSMapRedoFit
# import myXCLASSMapRedoFit function
from xclass.task_myXCLASSMapRedoFit import myXCLASSMapRedoFitCore
myXCLASSMapRedoFitCore(JobNumber = None,
PixelList = [],
JobPath = "",
MaskFileName = "",
Threshold = None,
MolfitsFileName = "",
MolfitFileName = None,
EstimateContinuumFlag = False,
EstimateConNumParts = 20,
EstimateConNumEntireCubeFlag = True,
DataSmoothValue = 0.0,
NumberIteration = 50,
NumberProcessors = 2,
AlgorithmXMLFile = "",
AlgorithmXMLFileName = "",
ClusterFileName = "",
clusterdef = "",
daskFlag = False,
ParallelizationMethod = "process",
FastFitFlag = True,
FullFitFlag = False,
ChannelIntegrationFlag = False,
FreqMin = None,
FreqMax = None,
t_back_flag = True,
tBack = 0.0,
tSlope = 0.0,
nH_flag = False,
N_H = 0.0,
beta_dust = 0.0,
kappa_1300 = 0.0,
DustFileName = "",
BackgroundFileName = "",
Te_ff = None,
EM_ff = None,
kappa_sync = None,
B_sync = None,
p_sync = None,
l_sync = None,
ContPhenFuncID = None,
ContPhenFuncParam1 = None,
ContPhenFuncParam2 = None,
ContPhenFuncParam3 = None,
ContPhenFuncParam4 = None,
ContPhenFuncParam5 = None,
iso_flag = False,
IsoTableFileName = "",
CollisionFileName = "",
DBFileName = "",
dbFileName = "",
NumModelPixelXX = 100,
NumModelPixelYY = 100,
BibTeXFlag = False,
LocalOverlapFlag = False,
NoSubBeamFlag = True,
EmAbsPATH = None)
This function offers the possibility to redo one or more so called
pixel fits of a previous myXCLASSMapFit run. The function performs
fits for the selected pixels and recreates the different parameter
maps using the new optimized parameter values together with the
corresponding FITS cube(s) describing the synthetic spectra for all
selected pixels. The names of these maps are created in the same way
as described in Sect. “myXCLASSMapFit”. In
addition, the myXCLASSMapRedoFit function renames the files
describing the parameter maps from the previous call of the
myXCLASSMapFit function by replacing the end of the each filename
".fits" by the phrase "__REDO-OLD.fits". The user has to
define the pixel coordinates, which should be re-fitted, the job
number of the previous myXCLASSMapFit run, and the new molfit
file.
Input parameters:
JobNumber: job number of a previous myXCLASSMapFit run which should be improved.
PixelList: list of all pixel coordinates which should be re-fitted. (These pixel coordinates are used to identify the pixel directories created by the selected myXCLASSMapFit call, see Sect. “myXCLASSMapFit”). Each coordinate consists of two numbers separated by a comma, where the first number corresponds to the right ascension and second to the declination. Different pixel coordinates are enclosed by squared brackets and separated by commas. For example, we want to refit the pixels 83.2013422325, -15.1345324232 and 83.1111111111, -15.2211111111. ThePixelListparameter has to be defined as followPixelList = [[83.2013422325, -15.1345324232], [83.1111111111, -15.2211111111]]
Please note, that these coordinates are used to identify the pixel directories in the selected myXCLASSMapFit job directory, i.e. the values defined in the
PixelListparameter have to be included in the names of the corresponding pixel directories. For example, we want to refit the pixel with pixel directory83.2013422325__-15.1345324232__1_-_3. In order to identify this pixel thePixelListparameter can be defined asPixelList = [[83.2013422325, -15.1345324232]]
Additionally, the
myXCLASSMapRedoFitfunction offers the possibility to define the indices of the pixels which should be refitted. In the example described above, thePixelListparameter can be defined asPixelList = [[1, 3]]
to identify the pixel 83.2013422325, -15.1345324232.
JobPath: (optional) path of job directory for application of the myXCLASSMapFit function, (default:"").
MaskFileName: (optional) path and name of a FITS image used for masking. Here, each pixel with a real numerical value, i.e. which is not NaN, is fitted, (default:"").
Threshold: defines a threshold for a pixel (default:None). If the spectrum of a pixel has a max. intensity lower than the value defined by this parameter the pixel is not fitted (ignored).
MolfitsFileName: path and name of the extended molfit file.
MolfitFileName: (optional) alias ofMolfitsFileName, (default:None).
EstimateContinuumFlag(True/False): (optional) flag indicating continuum estimation, (default:False).
EstimateConNumParts: (optional) number of sub parts used for continuum estimation, (default: 20).
EstimateConNumEntireCubeFlag(True/False): (optional) flag indicating that continuum estimation is done for the entire spectral range, (default:True).
DataSmoothValue: (optional) smooth factor (=0, no smoothing is performed) used for smoothing obs. data, (default: 0.0).
NumberIteration: max. number of iterations, (default: 50).
NumberProcessors: (optional, only used if no alg. xml file is defined) number of cores, (default: 2).
AlgorithmXMLFileName: path and name of an algorithm xml-file, (default:"").. (A relative path has to be defined relative to the current working directory!)
AlgorithmXMLFile: (optional) alias ofAlgorithmXMLFileName.
FastFitFlag(True/False): (optional) flag indicating usage of fast-fitting method, (default:True).
FullFitFlag(True/False): (optional) flag indicating usage of full-fitting method, (default:False).
ParallelizationMethod: (optional) method used for parallelization, (default:process). Methods other thanprocesswill be available in one of the next versions.
daskFlag: will be used in later releases.
ClusterFileName: (optional, only used for parallelization methods other thanprocess) path and name of a file describing a cluster, (default:"").
clusterdef: (optional) alias ofClusterFileName, (default:"")
ChannelIntegrationFlag(True/False): (optional, forfull-fitting method only) flag indicating channel integration, (default:False).
FreqMin: new start frequency of simulated spectrum, (default: 0).
FreqMax: new end frequency of simulated spectrum, (default: 0).
t_back_flag(True/False): defines, if the user defined background temperature \(T_{\rm bg}\) and temperature slope \(T_{\rm slope}\) given by the input parameterstBackandtSlopedescribe the continuum contribution completely (t_back_flag = True) or not (t_back_flag = False), (default:True).
tBackbackground temperature (in K), (default: 0).
tSlopetemperature slope (dimensionless), (default: 0).
BackgroundFileName(optional) path and name of file describing the background as function of frequency (in MHz), (default:None).
nH_flag: is deprecatedN_H(optional) Hydrogen column density (in cm\(^{-2}\)), (default:None).
beta_dust(optional) spectral index for dust (dimensionless), (default:None).
kappa_1300(optional) kappa at 1.3 mm (cm\(^2 \,\)g\(^{-1}\)), (default:None).
DustFileName(optional) path and name of file describing dust contribution (default:"").
Te_ff: (optional) electronic temperature (in K) used for free-free continuum, (default:None).
EM_ff: (optional) emission measure (in pc cm\(^{-6}\)) used for free-free continuum, (default:None).
kappa_sync: (optional) kappa of energy spectrum of electrons for synchrotron continuum (electrons m\(^{-3}\) GeV\(^{-1}\)), (default:None).
B_sync: (optional) magnetic field (in Gauss) used for synchrotron continuum, (default:None).
p_sync: (optional) energy spectral index (dimensionless) used for synchrotron continuum, (default:None).
l_sync: (optional) thickness of slab (in AU) used for synchrotron continuum (default:None).
ContPhenFuncID: (optional) describes which phenomenological function is used to describe the continuum, (default:None).
ContPhenFuncParam1: (optional) first parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam2: (optional) second parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam3: (optional) third parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam4: (optional) fourth parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam5: (optional) fifth parameter for phenomenological function describing the continuum, (default:None).
iso_flag: use isotopologues (True/False). Ifiso_flagis set toTrueisotopologues defined in the iso ratio file are used, (default:False).
IsoTableFileName(has to be given only ifiso_flagis set toTrue): path and name of an ASCII file including the iso ratios between certain molecules. A detailed description of the iso ratio file is given in Sect. “The iso ratio file”. If no file name is given (default), the so-called iso-flag is set toFalse.NOTE, if the given path is a relative path, i.e. the path does not start with
/, this path has to be defined relative to the current working directory!
CollisionFileName: (has to be defined only for Non-LTE calculation): path and name of file describing collision partners, (default:"").
DBFileName: (optional) path and name of a database file, which is used instead of the default filecdms_sqlite.db, (default:""").
dbFilename: (optional) alias ofDBFileName, (default:""").NumModelPixelXX: (optional) used for sub-beam modeling, describes the number of pixels used in x-direction, (default: 100).
NumModelPixelYY: (optional) used for sub-beam modeling, describes the number of pixels used in y-direction, (default: 100).
LocalOverlapFlag(True/False): (optional) take local overlap into account (default:False), see Sect. “Local-overlap of neighboring lines”.
NoSubBeamFlag(True/False): (optional) do not use sub-beam description (default:False), see Sect. “Sub-beam description”.
BibTeXFlag(True/False): flag indicating that a BibTeX file calledreferences_db.bibis generated within the job directory describing all papers used to generate the molecular parameters of transitions within the selected frequency ranges.
Output parameters:
None
IterativeFitting
# import myXCLASSFit function module
from xclass.task_myXCLASSFit import IterativeFitting
iter_class = IterativeFitting(MolfitFile_Class = None,
MolfitFileName = "",
ObsXMLFile_Class = None,
ObsXMLFileName = "",
AlgXMLFile_Class = None,
AlgXMLFileName = "",
JobDir = "",
ParameterMapDir = "",
ClusterFileName = "",
NumberProcessors = 1,
NumberWorkers = 4,
ServerList = [],
ClusterFlag = False,
CheckMoleculeInDB = True,
AdditionMolfitFileContentFileName = "",
RangevelLowLimit = -50.0,
RangevelUpLimit = 50.0,
IncludeAllMolFlag = True,
FirstMoleculeOnlyFlag = False,
LimitRangesFlag = False,
dryrunFlag = False,
RedRangeFlag = False,
SetNewRangeFlag = False,
MinColumnDensityAbs = 1.e12,
MinColumnDensityEmis = 1.e12,
SourceLow = None, SourceUp = None,
TLow = None, TUp = None,
NLow = None, NUp = None,
VwidthLow = None, VwidthUp = None,
vOffLow = None, vOffUp = None,
DoFinalPlotFlag = True,
DoPlotPerIterationFlag = False,
UpdateFlag = False,
NumRep = 1,
PlotNum = 32,
backgroundFlag = False,
BackgroundDir = "",
EmsAbsFlag = True,
FastCompBackFlag = True,
FastFlag = False,
AddCompFileName = "",
PrintFlag = True,
IndentString = "\t")
The IterativeFitting function offers an iterative fitting procedure for single spectra and entire data cubes.
Input parameters:
MolfitFile_Class: (optional) molfit file class object (default:None).
MolfitFileName: (optional) path and name of molfit file (default:"").
ObsXMLFile_Class: (optional) obs. xml file class object (default:None).
ObsXMLFileName: (optional) path and name of obs. xml file (default:"").
AlgXMLFile_Class: (optional) alg. xml file class object (default:None).
AlgXMLFileName: (optional) path and name of algorithm xml file (default:"").
JobDir: (optional) current working directory (default:"").
ParameterMapDir: (optional) path of subdirectory containing FITS images for different molfit parameters, (default:None).
ClusterFileName: (optional) path and name of file describing cluster information (default:"").
NumberProcessors: (optional) total number of available cores (default: 1).
NumberWorkers: (optional) number of workers (default: 4).
ServerList: (optional) list of all available nodes (default: []).
ClusterFlag(True/False): (optional) flag indicating usage of cluster (default:False).
CheckMoleculeInDB(True/False): (optional) check, if molecules included in molfit file have transitions within given frequency range (default:True).
AdditionMolfitFileContentFileName: (optional) path and name of file describing additional information for molfit file (default:"").
RangevelLowLimit: (optional) lower velocity limit (default: -50.0).
RangevelUpLimit: (optional) upper velocity limit (default: 50.0).
IncludeAllMolFlag(True/False): (optional) flag to include all molecules in each single molecule fit (default:True).
FirstMoleculeOnlyFlag(True/False): (optional) flag for taking only first mol. in given molfit file (default:False).
LimitRangesFlag(True/False): (optional) flag for limit parameter ranges in molfit file (default:False).
dryrunFlag(True/False): (optional) flag for executing script without fitting (default:False).
RedRangeFlag(True/False): (optional) flag for reducing frequency ranges in obs. xml file for each molecule (default:False).
SetNewRangeFlag(True/False): (optional) flag for allowing modifications on parameter ranges in molfit file (default:False).
MinColumnDensityAbs: (optional) minimal column density for absorption components (default: 1.e12).
MinColumnDensityEmis: (optional) minimal column density for emission components (default: 1.e12).
SourceLow: (optional) relative lower limit for source size (default: 100.0).
SourceUp: (optional) relative upper limit for source size (default: 100.0).
TLow: (optional) relative lower limit for temperature (default: 50.0).
TUp: (optional) relative upper limit for temperature (default: 50.0).
NLow: (optional) relative lower limit for column density (default: 3.0).
NUp: (optional) relative upper limit for column density (default: 3.0).
VwidthLow: (optional) relative lower limit for velocity width (default: 5.0).
VwidthUp: (optional) relative upper limit for velocity width (default: 5.0).
vOffLow: (optional) relative lower limit for velocity offset (default: 5.0).
vOffUp: (optional) relative upper limit for velocity offset (default: 5.0).
DoFinalPlotFlag(True/False): (optional) flag for final plot with all optimized molecules (default:True).
DoPlotPerIterationFlag(True/False): (optional) flag for creating transition plots for each iteration (default:False).
UpdateFlag(True/False): (optional) flag for using new optimized molfit parameter for further molecules (only in serial mode, i.e. without clustering) (default:False).
NumRep: (optional) number of repetitions / iterations (default: 1).
PlotNum: (optional) number of subplots (default: 32).
backgroundFlag(True/False): (optional) background flag (default:True).
BackgroundDir: (optional) directory for background (default:"").
EmsAbsFlag(True/False): (optional) flag for computing emission / absorption function (default:True).
FastCompBackFlag(True/False): (optional) flag for fast computing of emission / absorption function at the beginning of each iteration for all molecules (default:True).
FastFlag(True/False): (optional) fast flag for usage of lite XCLASS routines (default:False).
AddCompFileName: (optional) path and name of file describing additional components for each molecule molfit file (default:"").
PrintFlag(True/False): (optional) flag for screen output (default:True).
IndentString: (optional) string describing indent for screen output (default: “t”).
Output parameters:
None
DatabaseQuery
# import DatabaseQuery function
from xclass.task_DatabaseQuery import DatabaseQuery
DatabaseQuery(QueryString,
DBFileName = "",
dbFileName = "")
This function sends a given SQL query string to the XCLASS database.
Input parameters:
QueryString: the query string (default “ “).
DBFileName: (optional) path and name of a database file, which is used instead of the default filecdms_sqlite.db, (default:""").
dbFilename: (optional) alias ofDBFileName, (default:""").
Output parameters:
Contents: contains the screen output of the query command.Please note, the user is free to define a different name for the output parameter.
UpdateDatabase
# import UpdateDatabase function
from xclass.task_UpdateDatabase import UpdateDatabase
UpdateDatabase(DBUpdateNew,
DBFileName = None,
ForceRRLAdd = False,
RRLMaxN = 900,
RRLMaxDeltaN = 6)
This function updates the XCLASS database file located in the
directory ~/.xclass/db/ from spectroscopic databases connected to the
Virtual Atomic and Molecular Datacentre, VAMDC (www.vamdc.eu).
Currently, data access to the CDMS database and the JPL catalog, is
supported.
Input parameters:
DBUpdateNew: indicates, if, a complete SQLite3 database file is downloaded from the CDMS web server (“new”) or if the existing database file is updated (“update”) via the VAMDC infrastructure, i.e. new entries are added to the tables and existing entries are overwritten with new data. The complete SQLite3 database file is regularly updated via the VAMDC infrastructure and contains data from CDMS and JPL (defaultnew).
DBFileName(optional): path and name of database file (default:None).
ForceRRLAdd(optional,True/False): flag for forcing the addition and calculation of RRLs (default:False)
RRLMaxN(optional): max. main quantum number n for which RRLs are computed, (default:900)
RRLMaxDeltaN(optional): max. Delta n for which RRLs are computed, (default:6)
Output parameters:
None
ListDatabase
# import ListDatabase function
from xclass.task_ListDatabase import ListDatabaseCore
ListDatabaseCore(FreqMin = 1.0,
FreqMax = 1.e8,
ElowMin = 0.0,
ElowMax = 1.e6,
SelectMolecule = [],
DBFileName = "",
dbFileName = "",
OutputDevice = "")
This function reads in entries from the table transitions located
in the SQLite3 database file and prints out the contents to the screen
or file (defined by the input parameter OutputDevice). The user
can limit the output by defining a minimum and maximum for the
frequency (or for the lower energy) for the transitions.
Input parameters:
FreqMin: minimum frequency (in MHz) for transitions in thetransitionstable (default: 0)
FreqMax: maximum frequency (in MHz) for transitions in thetransitionstable (default: \(10^8\))
ElowMin: minimum for lower energy (in K) (default: 0)
ElowMax: maximum for lower energy (in K) (default: \(10^6\))
SelectMolecule: a (python) list containing the names of all molecules which should be considered or a string defining the path and the name of an ASCII file including the molecules which should be selected.Note, if the parameter defines a relative path, this path has to be defined relative to the current working directory!
DBFileName: (optional) path and name of a database file, which is used instead of the default filecdms_sqlite.db, (default:""").
dbFilename: (optional) alias ofDBFileName, (default:""").
OutputDevice: path and name of the file where the output is written to. If no file name is defined, the content of the database is written to the screen. If this parameter is set toquietno information is printed to screen.Note, if the parameter defines a relative path, this path has to be defined relative to the current working directory!
Output parameters:
Contents: content of the database (as a python array, e.g.Contents[0]contains the entries for the first molecule within the frequency/energy range).Note, the user is free to define a different name for the output parameter.
GetTransitions
# import GetTransitions function
from xclass.task_GetTransitions import GetTransitionsCore
GetTransitionsCore(obsdata = None,
FreqMin = None,
FreqMax = None,
SelectMolecule = [],
FrequencyWidth = 5,
ElowMin = 0,
ElowMax = 1.e9,
modeldata = None,
DBFileName = "",
dbFileName = "")
Similar to the ListDatabasefunction, the GetTransitions
function displays information from the embedded database about
transitions within a user-defined frequency range. Here, the frequency
range selection is done in an interactive way. Therefore, the
GetTransitions function starts a graphical user interface (GUI)
which describes observational data (defined by parameter obsdata)
within a frequency range defined by the parameters FreqMin and
FreqMax. Additionally, the function offers the possibility to plot a
synthetic spectrum (defined by parameter modeldata) of a previous
myXCLASS or myXCLASSFit function call as well.
By using the mouse in combination with the right mouse button, the user
defines the central frequency (blue vertical line) of the frequency
range (grey vertical bar) which is used for the database query. The
initial half width of this range is defined by parameter
FrequencyWidth. Note, the width of the range can be enlarged
(shrinked) by pressing the “+” (“-”) key. In the following,
XCLASS prints out information of all transitions within the selected
frequency range to the screen, starting with the transition with the
lowest transition frequency. In order to sort the transitions by lower
energies, the user has to press the “1” key before selecting the
central frequency. Using the “2” (“3”) key, sorts the
transitions by the Einstein A coefficients (by a products of upper state
degeneracies and Einstein A coefficients (gA). Pressing the “0” key
restores the initial order, i.e. by transition frequencies.
Additionally, the user can reduce the number of transitions by
defining a lower and a upper limit for the lower energies using
parameters ElowMin and ElowMax, respectively. Furthermore, the
screen output can be limited to a few molecules by using parameter
SelectMolecule.
Note, the GUI has to be closed to stop the function!
Input parameters:
obsdata: 2D numpy array containing the observational data.Note, the data has to be given as a function of frequency (in MHz).
modeldata: 2D numpy array containing the synthetic spectrum.Note, the data has to be given as a function of frequency (in MHz).
FreqMin: minimum frequency (in MHz) for the observational data (default:None).
FreqMax: maximum frequency (in MHz) for the observational data (default:None).
SelectMolecule: a (python) list containing the names of all molecules which should be considered or a file name including the molecules which should be selected.Note, if the parameter defines a relative path, this path has to be defined relative to the current working directory!
FrequencyWidth: defines the width of the frequency interval in MHz (selected central frequency \(\pm\)FrequencyWidthwhere the line information are printed out to the screen. (default: 5)
ElowMin: minimum for lower energy (in K) (default: 0) in tabletransitions.
ElowMax: maximum for lower energy (in K) (default: \(10^6\)) in tabletransitions.
modeldata2D numpy array containing the modeled data, (default:None)Note, the data has to be given as a function of frequency (in MHz).
DBFileName: (optional) path and name of a database file, which is used instead of the default filecdms_sqlite.db, (default:""").
dbFilename: (optional) alias ofDBFileName, (default:""").
Output parameters:
None
AddPrivateEntries
# import AddPrivateEntries function
from xclass.task_AddPrivateEntries import AddPrivateEntries
AddPrivateEntries(DBFileName = "",
dbFileName = "",
downloadFlag = False,
FinalPath = "",
MoleculeName = "",
TransFileName = "",
PFFileName = "",
EGYFileName = "",
OldFormatEGYFileFlag = False,
EinsteinAFlag = True,
logpfFlag = True,
plotpfFlag = False,
myxclassFlag = False)
The AddPrivateEntries function offers the possibility to add private entries for specific species to a given XCLASS database file.
Input parameters:
DBFileName: (optional) path and name of used database file (default:"").
dbFileName: (optional) alias of DBFileName, (default:"").
downloadFlag(True/False): (optional) flag for downloading a new database file from the CDMS server (default:False)
FinalPath: (optional) path and name of final file or path (default:"").
Molecule: name of molecule which is added (default:"").
TransFileName: path and name of ASCII file containing parameters for each transition (default:""). The format of this so called cat file is described in https://cdms.astro.uni-koeln.de/classic/general/#format_of_the_cat_file. XCLASS uses the format “(A13, 2A8, A2, A10, A3, A7, A4, 12A2)”, with parameters: ‘TransFreq’ (transition frequency in MHz), ‘UncerTransFreq’ (uncertainty of transition frequency), ‘EinsteinA’ (Einstein A coefficient on log10 scale), ‘Int’ (intensity in units of nm\(^2\) MHz on log10 scale), ‘Elow’ (lower energy in cm\(^{-1}\)), ‘gup’ (upper state degeneracy), ‘QNFMT’ (identifies the format of the quantum numbers, described in https://cdms.astro.uni-koeln.de/classic/general/#format_of_quantum_numbers), ‘QN’ (Quantum numbers for the state).
PFFileName: path and name of ASCII file containing parameters for partition function (default:""). XCLASS expects two columns, where the first column describes the temperature and the second column the corresponding partition function. Comments are indicated by the%sign.
EGYFileName: path and name of ASCII file describing energies (default:""). The format of the egy file is described in https://cdms.astro.uni-koeln.de/classic/general/#format_of_the_egy_file. XCLASS uses the following formats (indicated by input parameterOldFormatEGYFileFlag) for the egy file:old format: “(I5, I5, 3F18.6, 6A3)” (Fortran format), with parameters: ‘IBLK’ (Internal Hamiltonian block number), ‘INDX’ (Internal index Hamiltonian block), ‘EGY’ (Energy in cm\(^{-1}\)), ‘ERR’ (Expected error of the energy in cm\(^{-1}\)), ‘PMIX’ (mixing coefficient), ‘QN’ (Quantum numbers for the state).
new format: “(F30.10, I33, 1X, 2I3, 3A3)” (Fortran format), with parameters: ‘EGY’ (Energy in cm\(^{-1}\)), ‘gup’ (upper state degeneracy), ‘dum1’ (dummy parameter), ‘dum2’ (dummy parameter), ‘QN’ (Quantum numbers for the state).
OldFormatEGYFileFlag(True/False): (optional) flag for using old format for egy file (default:False).
EinsteinAFlag(True/False): (optional) flag for reading Einstein A coefficient (default:True). If this flag is set toFalse, XCLASS tries to compute the Einstein A coefficient from the given intensity using the Eq. (5) described in http://www.astro.uni-koeln.de/cdms/catalog#equations.
logpfFlag: (optional, only used if ASCII file describing the partition function is specified) flag indicating that the partition function is described on log-scale (default:True)
plotpfFlag(True/False): (optional) flag for plotting partition function (default:True)
myxclassFlag(True/False): (optional) flag for using myXCLASS function to check new entries (default:True), i.e. the myXCLASS function is used in combination with a default molfit file to check the new entries.
Output parameters:
None
LineID
# import LineID function
from xclass.task_LineIdentification import LineIdentificationCore
LineIdentificationCore(JobDir = None,
DefaultMolfitFileName = "",
DefaultMolfitFile = "",
SourceName = "",
MaxOverestimationHeight = 10.0,
Tolerance = 10.0,
SelectedMolecules = [],
StrongMoleculeList = [],
MinColumnDensityEmis = 1.0,
MinColumnDensityAbs = 1.0,
NumberIteration = 50,
FastFitFlag = False,
NumberProcessors = 4,
ChannelIntegrationFlag = True,
ClusterFileName = "",
clusterdef = "",
ParallelizationMethod = "process",
NumberWorkers = None,
AlgXMLFile_Class = None,
AlgorithmXMLFileName = "",
AlgorithmXMLFileISO = "",
AlgorithmXMLFileSMF = "",
AlgorithmXMLFileOverAll = "",
regionFileName = "",
FITSImageMaskFileName = "",
PlotSMFFlag = False,
EstimateParamFlag = False,
InternalParameterList = None,
PrintFlag = True,
ObsXMLFile_Class = None,
ObsXMLFileName = "",
ObsDataFileName = "",
experimentalData = "",
Noise = 2.0,
TelescopeSize = None,
BMIN = None,
BMAJ = None,
BPA = None,
Inter_Flag = False,
Redshift = 0.0,
vLSR = 0.0,
RestFreq = 0.0,
FreqMin = None,
FreqMax = None,
t_back_flag = True,
tBack = 0.0,
tSlope = 0.0,
BackgroundFileName = "",
nH_flag = None,
N_H = 0.0,
beta_dust = 0.0,
kappa_1300 = 0.0,
DustFileName = "",
Te_ff = None,
EM_ff = None,
kappa_sync = None,
B_sync = None,
p_sync = None,
l_sync = None,
ContPhenFuncID = None,
ContPhenFuncParam1 = None,
ContPhenFuncParam2 = None,
ContPhenFuncParam3 = None,
ContPhenFuncParam4 = None,
ContPhenFuncParam5 = None,
NoSubBeamFlag = None,
LocalOverlapFlag = None,
NumModelPixelXX = 100,
NumModelPixelYY = 100,
iso_flag = False,
IsoTableFileName = "",
CollisionFileName = None,
EmAbsPATH = None,
BibTeXFlag = False,
DBFileName = "",
dbFileName = "")
This function identifies molecules in a given spectrum and returns a quantitative description of the data.
In general, the LineIdentification function takes all molecules into account which have at least one transition within the user-defined frequency range(s) [1]. These molecules are stored into a so-called molecule list. In order to determine the contribution of each molecule, the LineIdentification function performs so-called single molecule fits for each molecule. If a molecule covers a defined fraction of the spectrum the molecule is “for now” identified and the optimized molfit file is append to a so-called overall molfit file which describes the contribution of all identified molecules. After all single molecule fits are done, the LineIdentification function performs a final fit, using the overall molfit file created before.
Input parameters:
JobDir: describes path of the current job directory, i.e.OutDict["JobDir"]is the path and name of the current job directory.
DefaultMolfitFileName: This parameter defines the path and name of a so-called default molfit file which is used to fit the contribution of each molecule, (default:"").The default molfit file defines for one molecule the number of components, the ranges and initial values for each parameter.
During the line identification process the name of the (first) molecule defined in the default molfit file is replaced by the name of the current molecule.
NOTE, if the parameter
DefaultMolfitFileNamedefines a relative path, the path has to be defined relative to the current working directory!
DefaultMolfitFile: (optional) alias ofDefaultMolfitFileName,(default:
"").
SourceName: This parameter defines the path and the name of a source (template) molfit file which is used for the single molecule fits. Whenever the LineIdentification function fits the contribution of one of the molecules included in the source molfit file the definitions herein are used instead of the definitions in the default molfit file, see description below. This offers the possibility to fit the contribution(s) of one or more molecules with a different number of components, ranges, initial values etc. as defined in the default molfit file. In order to perform the single molecule fits for molecules which are not described by the source molfit file, the LineIdentification function uses the definitions in the default molfit file, see below.
Noise: noise level (in K), (default: 2.0 K). All parts of the spectrum with intensities lower than the noise level are ignored.The usage of an observational xml file offers the possibility to specify a noise level for each frequency range. For that purpose, a tag called
NoiseLevel, see below, has to be defined for each frequency range for all spectra defined in the xml file. If the definition is not given for all frequency ranges, the value defined by this parameterNoiseis used instead for all frequency ranges of all spectra.
MaxOverestimationHeight: defines the overestimation factor of the modeled single molecule spectrum (in %), (default: 10 %).Please note, the input parameter
MaxOverestimationHeighthas to be given in % and does not define the final overestimation limit directly which is given byMaxOverestimationHeight+ 100 %.If the modeled spectrum does not overestimate the observed spectrum, the corresponding molecule is “for now” identified and the fitted molfit file is added to the overall molfit file which is used at the end of the line identification process to determine the final contribution of each molecule.
Tolerance: defines the max. fraction (in %) of overestimated lines. If the fraction of overestimated lines in the result of a single molecule fit is lower than the given threshold, the corresponding molecule is “for now” identified and the optimized molfit file (and the corresponding iso ratios) are considered in the final overall fit (default: 10.0).
SelectedMolecules: A (python) list defining molecules which should be excluded from or considered only by the line identification process.In general, the LineIdentification function considers all molecules which are included in the database within the defined frequency range(s).
If the list contains the command words “
known molecules”, the LineIdentification function adds all known molecules [1] to the list as well.In order to exclude a molecule from the line identification process, the name of the molecule has to start with “
--”, e.g. to exclude the moleculeCH3SH;v=0;the user has to defineSelectedMolecules = ["--CH3SH;v=0;"].If the molecule names do not start with “
--”, the LineIdentification function consider only these molecules.
StrongMoleculeList: A (python) list including the strong (highly abundant) molecules which should be fitted at the beginning (default: “[]”).
MinColumnDensityEmis: min. column density (for core lines) of an optimized component of a single molecule fit to be included in the overall molfit file (default: 0).
MinColumnDensityAbs: min. column density (for foreground lines) of an optimized component of a single molecule fit to be included in the overall molfit file (default: 0).
FastFitFlag(True/False): (optional) flag indicating usage of fast-fitting method, (default:True).
ChannelIntegrationFlag(True/False): (optional, forfull-fitting method only) flag indicating channel integration, (default:False).
NumberIteration: max. number of iterations (default: 50). This parameter is used only if no MAGIX xml files are defined by the parametersAlgorithmXMLFileNameandAlgorithmXMLFileOverAll.
NumberProcessors: (optional, only used if no alg. xml file is defined) number of cores, (default: 2).
NumberWorkers: (optional, only used if no cluster file is defined) number of workers (default: 4)
AlgorithmXMLFileName: path and name of an algorithm xml-file, (default:"").. (A relative path has to be defined relative to the current working directory!)
AlgorithmXMLFileSMF(optional) for each single molecule fit): path and name of a MAGIX xml-file defining settings for an algorithm or algorithm tree has to be given. (A relative path has to be defined relative to the current working directory!)NOTE, if the user specify a xml file, the number of iterations given by the parameter
NumberIterationis ignored. The number of iteration is then given by the xml file itself. In order to use the implemented fit algorithm (Levenberg-Marquardt) clear theAlgorithmXMLFileSMFparameter, i.e.AlgorithmXMLFileSMF= “”, and define the max. number of iterations by using parameterNumberIteration.
AlgorithmXMLFileOverAll: only necessary, if the user wants to use another fit algorithm (than Levenberg-Marquardt) for the overall fit. Therefore, the path and name of a MAGIX xml-file defining settings for an algorithm or algorithm tree has to be given. (A relative path has to be defined relative to the current working directory!)For more information see description of parameter
AlgorithmXMLFileSMF.
AlgXMLFile_Class: (optional) alg. xml file class object, (default:None).
ClusterFileName: (optional, only used for parallelization methods other thanprocess) path and name of a file describing a cluster, (default:"").
clusterdef: (optional) alias ofClusterFileName, (default:"")
ParallelizationMethod: (optional) method used for parallelization, (default:process). Methods other thanprocesswill be available in one of the next versions.
regionFileName: (optional) path and name of a so-called ds9 [2] region file, (default:"").
FITSImageMaskFileName: (optional) path and name of a FITS image used for masking. Here, each pixel with a real numerical value (except 0), i.e. which is not NaN or 0, is fitted, (default:""). The value of the selecting pixel is used to weight the contribution of the corresponding pixel to the total \(\chi^2\) value.
PlotSMFFlag(True/False): (optional, only used if single spectra are analyzed) flag for creating twf plots for each smf (default:False)
EstimateParamFlag(True/False): (optional) will be available in one of the next releases defines, if the calculation parameters for a single molecule fit are estimated directly from the given observational data (default:False). Parameter estimation is used for all single molecule fits, which use the default molfit file.
InternalParameterList: (optional, only used ifEstimateParamFlagis set toTrue) dictionary describing internal parameters (default:None)
printFlag(True/False): (optional) flag for printing messages to the screen (default:True).
ObsXMLFile_Class: (optional) obs. xml file class object, (default:None).
ObsXMLFileName: (optional) path and name of an obs. xml file (recommended procedure), (default:""), see Sect. “Observational xml file” for a detailed description.Note, the name of the xml file has to end with
".xml".
ObsDataFileName: (optional) describes the observational data, if no obs. xml file is given, (default:None). There are two options for transferring the data to XCLASS:the parameter describes path and name of an ASCII file called observational data file, where the first column describe the frequency (in MHz) and the second column the brightness temperature (in Kelvin)
Please note, without using an obs. xml file, it is not possible to fit more than one data cube, simultaneously.
the parameter describes a 2d numpy array, where the first column describes the frequency (in MHz) and the second column the brightness temperature (in Kelvin).
experimentalData: (optional) alias ofObsDataFileName, (default:None)
The following parameters are needed, if parameter ObsDataFileName is
used instead of parameter ObsXMLFileName:
TelescopeSize: for single dish observations (Inter_Flag = False):TelescopeSizedescribes the size of telescope (in m), (default: 1); for interferometric observations (Inter_Flag = True):TelescopeSizedescribes the interferometric beam FWHM size (in arcsec), (default: 1).
BMIN: (optional) beam minor axis length (in arcsec), (default:None).Note, in the header of FITS files, the beam minor and major axis lengths are often given in degrees and not in arcsec, i.e. the values have to be divided by 3600.
BMAJ: (optional) beam major axis length (in arcsec), (default:None).
BPA: (optional) beam position angle (in degree), (default:None).
Redshift: (optional) red shift (dimensionless), (default:None).
Inter_Flag(True/False): defines, if single dish (False) or interferometric observations (True) are described, (default:False).
t_back_flag(True/False): defines, if the user defined background temperature \(T_{\rm bg}\) and temperature slope \(T_{\rm slope}\) given by the input parameterstBackandtSlopedescribe the continuum contribution completely (t_back_flag = True) or not (t_back_flag = False) (default:True).
tBackbackground temperature (in K), (default: 0.0).
tSlopetemperature slope (dimensionless), (default: 0.0).
BackgroundFileName(optional) path and name of file describing the background as function of frequency (in MHz), (default:None).
nH_flag: is deprecatedN_H(optional) Hydrogen column density (in cm\(^{-2}\)), (default: 0.0).
beta_dust(optional) spectral index for dust (dimensionless), (default: 0.0).
kappa_1300(optional) kappa at 1.3 mm (cm\(^2 \,\)g\(^{-1}\)), (default: 0.0).
DustFileName(optional) path and name of file describing dust contribution, (default:None).
Te_ff: (optional) electronic temperature (in K) used for free-free continuum, (default:None).
EM_ff: (optional) emission measure (in pc cm\(^{-6}\)) used for free-free continuum, (default:None).
kappa_sync: (optional) kappa of energy spectrum of electrons for synchrotron continuum (electrons m\(^{-3}\) GeV\(^{-1}\)), (default:None).
B_sync: (optional) magnetic field (in Gauss) used for synchrotron continuum, (default:None).
p_sync: (optional) energy spectral index (dimensionless) used for synchrotron continuum, (default:None).
l_sync: (optional) thickness of slab (in AU) used for synchrotron continuum, (default:None).
ContPhenFuncID: (optional) describes which phenomenological function is used to describe the continuum, (default:None).
ContPhenFuncParam1: (optional) first parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam2: (optional) second parameter for phenomenological function describing the continuum,(default:None).
ContPhenFuncParam3: (optional) third parameter for phenomenological function describing the continuum, (default:None).
ContPhenFuncParam4: (optional) fourth parameter for phenomenological function describing the continuum,(default:None).
ContPhenFuncParam5: (optional) fifth parameter for phenomenological function describing the continuum, (default:None).
iso_flag: use isotopologues (True/False). Ifiso_flagis set toTrueisotopologues defined in the iso ratio file are used, (default:False).
IsoTableFileName(has to be given only ifiso_flagis set toTrue): path and name of an ASCII file including the iso ratios between certain molecules. A detailed description of the iso ratio file is given in Sect. “The iso ratio file”. If no file name is given (default), the so-called iso-flag is set toFalse.NOTE, if the given path is a relative path, i.e. the path does not start with
/, this path has to be defined relative to the current working directory!
CollisionFileName: (has to be defined only for Non-LTE calculation): path and name of file describing collision partners, (default:None).
NumModelPixelXX: (optional) used for sub-beam modeling, Sect. “Sub-beam description”, describes the number of pixels used in x-direction, (default: 100).
NumModelPixelYY: (optional) used for sub-beam modeling, Sect. “Sub-beam description”, describes the number of pixels used in y-direction, (default: 100).
LocalOverlapFlag(True/False): (optional) flag indicates if local-overlap description is used or not, see Sect. “Local-overlap of neighboring lines”, (default:False).
NoSubBeamFlag(True/False): (optional) do not use sub-beam description (default:False), see Sect. “Sub-beam description”, i.e. XCLASS uses the old source description, see Sect. “Simple 1-d description”.
DBFileName: (optional) path and name of a database file, which is used instead of the default filecdms_sqlite.db, (default:""").
dbFilename: (optional) alias ofDBFileName, (default:""").
Output parameters:
IdentifiedLines: contains the optimized molfit file including the fitted parameters of all identified molecules.
JobDir: absolute path of the job directory created for the current run.
CubeFit
# import CubeFit function
from xclass.task_CubeFit import CubeFitCore
CubeFitCore(MolfitFileName = None,
MolfitsFileName = None,
ObsXMLFile_Class = None,
ObsXMLFileName = "",
ObsDataFileName = "",
experimentalData = "",
AlgXMLFile_Class = None,
AlgorithmXMLFileName = "",
AlgorithmXMLFile = "",
NumberIteration = 10,
NumberProcessors = 4,
FastFitFlag = False,
ChannelIntegrationFlag = False,
OpacityFlag = False,
Parallelization_Class = None,
daskFlag = False,
ClusterFileName = "",
clusterdef = "",
ParallelizationMethod = "process",
regionFileName = "",
FITSImageMaskFileName = "",
dbgFlag = False,
CleanupFlag = False,
JobDir = None,
PrintFlag = True,
CubeFitObj = None,
CubeFitParamDict = {},
CubeFitParamDic = None)
This function offers the possibility to directly fit existing data cubes to physical source models, such as rotating disks, or infalling cores etc. It allows a fast exploration of the parameter space, and optimization of the model parameters.
Input parameters:
JobDir: describes path of the current job directory, i.e.OutDict["JobDir"]is the path and name of the current job directory.
MolfitFileName: is deprecated
MolfitsFileName: is deprecated
ObsXMLFile_Class: (optional) obs. xml file class object, (default:None).
ObsXMLFileName: the parameter defines the path and name of an observational xml-file.Please note, if more than one data cube is specified in the xml file, the data cubes must describe the same map, i.e. the right ascension and the declination has to be identical! The different data cubes must be differ only in the frequency/velocity axis. Note, it is possible to specify different units for the frequency axis for different data cubes.
ObsDataFileName: is deprecated.
experimentalData: is deprecated.
AlgXMLFile_Class: (optional) alg. xml file class object, (default:None).
AlgorithmXMLFileName: path and name of an algorithm xml-file describing settings for the applied optimization algorithm(s).
AlgorithmXMLFile: (optional) alias ofAlgorithmXMLFileName.
NumberIteration: max. number of iterations (default: 50). This parameter is used only if no MAGIX xml files are defined by the parametersAlgorithmXMLFileNameandAlgorithmXMLFileOverAll.
NumberProcessors: (optional, only used if no alg. xml file is defined) number of cores, (default: 2).
FastFitFlag(True/False): (optional) flag indicating usage of fast-fitting method, (default:True).
ChannelIntegrationFlag(True/False): (optional, forfull-fitting method only) flag indicating channel integration, (default:False).
OpacityFlag(True/False): (optional, only forFastFitFlagis set toFalse) flag indicating storage of opacities and intensities for each distance (ifLocalOverlapFlagis set toTrue) or for each component (ifLocalOverlapFlagis set toFalse) (default:False). If the OpacityFlag is set toTrue, the job directory contains an additional subdirectory calledintensities_and_opacities, which includes files describing the intensities and opacities per component.
ParallelizationMethod: (optional) method used for parallelization, (default:process). Methods other thanprocesswill be available in one of the next versions.
ClusterFileName: (optional, only used for parallelization methods other thanprocess) path and name of a file describing a cluster, (default:"").
clusterdef: (optional) alias ofClusterFileName, (default:"")
daskFlag: will be used in later releases.
Parallelization_Class: (optional) parallel class object, (default:None).
regionFileName: path and name of a ds9 [11] region file, see Sect. “myXCLASSMapFit”.
FITSImageMaskFileName: path and name of a FITS image, see Sect. “myXCLASSMapFit”, whose pixels that are not NaN or 0 indicate which pixels are considered by the CubeFit function. The value of the selecting pixel is used to weight the contribution of the corresponding pixel to the total \(\chi^2\) value.Note that the value \(\sigma_i^\mathsf{mask}\) of each pixel in the FITS image is used to weight the corresponding pixel in the calculation of the total \(\chi^2\) value according to the following formula
\[\chi^2 = \sum_{i=1}^N \sum_{j=1}^M \left[ \left(y_{i,j}^\mathsf{obs} - y_{i,j}^\mathsf{fit} \right)^2 \cdot \frac{1}{\left(\sigma_i^\mathsf{mask}\right)^2} \right],\]where the first sum runs over the all \(N\) selected pixels, and the second sum over each spectral channel, respectively. Here, \(y_{i,j}^\mathsf{obs}\) represents the value of the observational data at pixel \(i\) and channel \(j\), and \(y_{i,j}^\mathsf{fit}\) the corresponding value of the fit function. Additionally, \(\sigma_i^{\rm mask}\) represents the value of the FITS image at pixel \(i\).
FITS image and ds9 regions can be used together.
clusterdef: is deprecated.
FastFitFlag: is deprecated.
dbgFlag: is deprecated.
CleanupFlag: is deprecated.
CubeFitObj: objective function to be called as generator for input file(s).
CubeFitParamDict: python dictionary describing all model parameters used by the objective function, see description above.
Output parameters:
OutDict: python dictionary, e.g.OutDict, containing the following items"JobDir": describes path of the current job directory, i.e.OutDict["JobDir"]is the path and name of the current job directory.
"CubeFitParamDict": complete python dictionary describing all model parameters
"CubeFitFitParam": list of fitted model parameters.
Footnotes
References
Virtanen, Gommers, Oliphant, Haberland, Reddy, Cournapeau, Burovski, Peterson, Weckesser, Bright, van der Walt, Brett, Wilson, Millman, Mayorov, Nelson, Jones, Kern, Larson, Carey, Polat, Feng, Moore, VanderPlas, Laxalde, Perktold, Cimrman, Henriksen, Quintero, Harris, Archibald, Ribeiro, Pedregosa, van Mulbregt, & SciPy 1.0 Contributors
Sanchez-Monge, A., Schilke, P., Ginsburg, A., Cesaroni, R., & Schmiedeke, A. 2018, Astronomy & Astrophysics, 609, A101
Rohlfs, K. & Wilson, T.~L. 1996, Tools of Radio Astronomy
Greisen, E.~W. & Calabretta, M.~R. 2002, Astronomy & Astrophysics, 395, 1061
Calabretta, M.~R. & Greisen, E.~W. 2002, Astronomy & Astrophysics, 395, 1077
Greisen, E., Calabretta, M., Valdes, F., & Allen, S. 2006, Astronomy & Astrophysics, 446, 747