mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
phoenixSpectrum.hpp
1/** \file phoenixSpectra.hpp
2 * \author Jared R. Males
3 * \brief Utilities for working with spectra from phoenix code.
4 * \ingroup astrophot
5 *
6 */
7
8#ifndef phoenixSpectrum_hpp
9#define phoenixSpectrum_hpp
10
11#include "../math/vectorUtils.hpp"
12#include "../ioutils/fileUtils.hpp"
13
14namespace mx
15{
16namespace astro
17{
18
19/// A spectrum from the Phoenix model, for use with ioutils::astro::astroSpectrum.
20/** \ingroup astrophot_spectra
21 */
22template <typename _units>
24{
25 typedef _units units;
26 typedef typename units::realT realT;
27
28 static const bool freq = false;
29
30 /// Convert from A to SI m
31 static constexpr realT wavelengthUnits = static_cast<realT>( 1e10 );
32
33 /// Convert from erg s-1 cm-2 A-1 to SI W m-3
34 static constexpr realT fluxUnits =
35 static_cast<realT>( 1e7 ) / ( static_cast<realT>( 1e4 ) * static_cast<realT>( 1e10 ) );
36
37 static constexpr const char *dataDirEnvVar = "PHOENIX_DATADIR";
38
39 typedef std::string paramsT; ///< The parameter is a string name.
40
41 static std::string fileName( const std::string &name )
42 {
43 return name;
44 }
45
46 static int
47 readSpectrum( std::vector<realT> &rawLambda,
48 std::vector<realT> &rawSpectrum,
49 const std::string &path,
50 const paramsT &params ///< [in] the parameters are passed in case needed to construct the spectrum
51 )
52 {
53 if( ioutils::readColumns( path, rawSpectrum ) < 0 )
54 return -1;
55
56 if( ioutils::readColumns( ioutils::parentPath( path ) + "/wavelength.dat", rawLambda ) < 0 )
57 return -1;
58
59 return 0;
60 }
61
62 static void scaleSpectrum( std::vector<realT> &spectrum, realT radius, realT distance )
63 {
64 for( int i = 0; i < spectrum.size(); ++i )
65 {
66 spectrum[i] *=
67 pow( radius / distance * ( constants::radJupiter<units>() / constants::parsec<units>() ), 2 );
68 }
69 }
70};
71
72/// Read in, crop, scale, and re-write a Phoenix spectrum data file.
73/** For working with spectra obtained from: http://perso.ens-lyon.fr/france.allard/
74 *
75 * The Phoenix spectra contain many columns of line data which are not often used, are
76 * not necessarily sorted by wavelength, and are sometimes formated so that points with leading
77 * minus signs are not in a separate column. We also have to change the fortan D to e.
78 *
79 * We also want to apply the dilution factor and take the power of 10.
80 *
81 * This function deals with these issues, producing a two-column space-delimited
82 * file in a specified wavelength range.
83 *
84 * References:
85 * - https://phoenix.ens-lyon.fr/Grids/BT-Settl/README
86 * - https://phoenix.ens-lyon.fr/Grids/FORMAT
87 *
88 * NOTE: This overwrites the input file!
89 *
90 *
91 * \tparam floatT the floating point type in which to work
92 *
93 * \ingroup astrophot_spectra
94 */
95template <typename floatT>
97 const std::string &filename, ///< [in.out] complete name of the file to rewrite
98 floatT lmin, ///< [in] minimum wavelength to rewrite [microns]
99 floatT lmax, ///< [in] maximum wavelemngth to rewrite [microns]
100 int sepWavelength =
101 0, ///< [in] [optional] controls how wavelength is handled. 0=> wavelength is included in file as first column.
102 ///< 1=> wavelength is written to a separate 'wavelength.dat' file. -1=> means don't write wavelength.
103 floatT DF = -8.0 ///< [in] [optional] the dilution factor. See references.
104)
105{
106 float lambda, flambda;
107
108 std::ifstream fin;
109
110 fin.open( filename );
111
112 if( !fin.good() )
113 {
114 std::cerr << "Error opening file: " << filename << "\n";
115 return;
116 }
117
118 std::vector<floatT> lambdas;
119 std::vector<floatT> flambdas;
120
121 int lineSize = 1024;
122 char line[lineSize];
123 std::string sline;
124
125 fin.getline( line, lineSize );
126 sline = line;
127
128 if( sline.length() < 25 )
129 {
130 std::cerr << "Error reading file: " << filename << "\n";
131 }
132
133 lambda = ioutils::convertFromString<floatT>( sline.substr( 1, 12 ) );
134
135 // First have to diagnose if this has the column problem
136 int nst;
137 if( sline[13] == '-' )
138 nst = 13;
139 else
140 nst = 14;
141
142 // convert the D to e
143 size_t dpos = sline.find( 'D', nst );
144 sline[dpos] = 'e';
145
146 flambda = ioutils::convertFromString<floatT>( sline.substr( nst, 12 ) );
147
148 if( lambda >= lmin * 1e4 && lambda <= lmax * 1e4 )
149 {
150 lambdas.push_back( lambda );
151 flambdas.push_back( flambda );
152 }
153
154 fin.getline( line, lineSize );
155 sline = line;
156
157 while( fin.good() )
158 {
159 if( sline.length() < 25 )
160 continue;
161
162 dpos = sline.find( 'D', nst );
163 sline[dpos] = 'e';
164
165 lambda = ioutils::convertFromString<floatT>( sline.substr( 1, 12 ) );
166 flambda = ioutils::convertFromString<floatT>( sline.substr( nst, 12 ) );
167
168 if( lambda >= lmin * 1e4 && lambda <= lmax * 1e4 )
169 {
170 lambdas.push_back( lambda );
171 flambdas.push_back( flambda );
172 }
173
174 fin.getline( line, lineSize );
175 sline = line;
176 }
177
178 fin.close();
179
180 std::vector<size_t> idx = math::vectorSortOrder( lambdas );
181
182 std::ofstream fout;
183
184 std::string fname = filename; // + ".t";
185 fout.open( filename );
186 fout.precision( 10 );
187
188 for( int i = 0; i < lambdas.size(); ++i )
189 {
190 if( sepWavelength == 0 )
191 {
192 fout << lambdas[idx[i]] << " ";
193 }
194 fout << pow( 10, flambdas[idx[i]] + DF ) << "\n";
195 }
196 fout.close();
197
198 if( sepWavelength == 1 )
199 {
200 fname = ioutils::parentPath( filename ) + "/wavelength.dat";
201
202 fout.open( fname );
203 fout.precision( 10 );
204 for( int i = 0; i < lambdas.size(); ++i )
205 {
206 fout << lambdas[idx[i]] << "\n";
207 }
208 fout.close();
209 }
210
211} // rewritePhoenixSpectrum
212
213/// Call rewritePhoenixSpectrum for all files in a directory.
214/** \ingroup astropho_spectra
215 */
216template <typename floatT>
217void rewritePhoenixSpectrumBatch( const std::string &dir, floatT lmin, floatT lmax, floatT DF = -8.0 )
218{
219 std::vector<std::string> flist = ioutils::getFileNames( dir, "lte", "", ".7" );
220
221 int sepWavelength = 1;
222 for( int i = 0; i < flist.size(); ++i )
223 {
224 rewritePhoenixSpectrum( flist[i], lmin, lmax, sepWavelength, DF );
225 if( i == 0 )
226 sepWavelength = -1;
227 }
228}
229
230} // namespace astro
231} // namespace mx
232
233#endif // phoenixSpectrum_hpp
int readColumns(const std::string &fname, arrTs &...arrays)
Read in columns from a text file.
void rewritePhoenixSpectrum(const std::string &filename, floatT lmin, floatT lmax, int sepWavelength=0, floatT DF=-8.0)
Read in, crop, scale, and re-write a Phoenix spectrum data file.
std::vector< std::string > getFileNames(const std::string &directory, const std::string &prefix, const std::string &substr, const std::string &extension)
Get a list of file names from the specified directory, specifying a prefix, a substring to match,...
std::string parentPath(const std::string &fname)
Get the parent path from a filename.
Definition fileUtils.cpp:77
std::vector< size_t > vectorSortOrder(std::vector< memberT > const &values)
Return the indices of the vector in sorted order, without altering the vector itself.
The mxlib c++ namespace.
Definition mxError.hpp:106
A spectrum from the Phoenix model, for use with ioutils::astro::astroSpectrum.
static constexpr realT fluxUnits
Convert from erg s-1 cm-2 A-1 to SI W m-3.
static int readSpectrum(std::vector< realT > &rawLambda, std::vector< realT > &rawSpectrum, const std::string &path, const paramsT &params)
static constexpr realT wavelengthUnits
Convert from A to SI m.
std::string paramsT
The parameter is a string name.