source: JEM-EUSO/esaf_cc_at_lal/packages/common/atmosphere/src/LowtranAtmosphere.cc @ 114

Last change on this file since 114 was 114, checked in by moretto, 11 years ago

actual version of ESAF at CCin2p3

File size: 8.6 KB
Line 
1// ESAF : Euso Simulation and Analysis Framework
2// $Id: LowtranAtmosphere.cc 2749 2006-07-04 13:16:23Z moreggia $
3// Corinne Berat created Apr,  2 2004
4
5/*****************************************************************************
6 * ESAF: Euso Simulation and Analysis Framework                              *
7 *                                                                           *
8 *  Id: LowtranAtmosphere                                                    *
9 *  Package: atmosphere                                                      *
10 *  Coordinator: S. Moreggia, C. Berat                                       *
11 *                                                                           *
12 *****************************************************************************/
13
14//_____________________________________________________________________________
15//
16// LowtranAtmosphere
17//
18// <extensive class description>
19//
20//   Config file parameters
21//   ======================
22//
23//   <parameter name>: <parameter description>
24//   -Valid options: <available options>
25//
26
27#include "LowtranAtmosphere.hh"
28#include "Config.hh"
29#include "LowtranAtmosphereData.hh"
30#include <TH1F.h>
31#include <TROOT.h>
32#include <math.h>
33#include <TSpline.h>
34
35using namespace sou;
36using namespace TMath;
37
38ClassImp(LowtranAtmosphere)
39
40//_____________________________________________________________________________
41LowtranAtmosphere::LowtranAtmosphere(string type) : Atmosphere(),  fData(0) {
42    //
43    // ctor
44    //
45    Msg(EsafMsg::Info) << "Building with the configuration:"<< MsgDispatch;
46    fType = type;
47    Build();
48}
49
50//_____________________________________________________________________________
51LowtranAtmosphere::~LowtranAtmosphere() {
52    //
53    // ctor
54    //
55    SafeDelete(fData);
56}
57
58//___________________________________________________________________________________
59void LowtranAtmosphere::CreateInstance() {
60    //
61    // Create the single atmosphere instance as being a LowtranAtmosphere object
62    //
63    if(!Atmosphere::fChild) Atmosphere::fChild = new LowtranAtmosphere("lowtran");
64    else Atmosphere::fChild->Msg(EsafMsg::Warning) << "Trying to create two atmospheres !!"<< MsgDispatch;
65}
66
67//_____________________________________________________________________________
68void LowtranAtmosphere::Build() {
69    //
70    // Build atmosphere according to config parameters
71    //
72    fData =  new LowtranAtmosphereData();
73    if(!fData)  Msg(EsafMsg::Panic) << "No new LowtranData, memory pb" << MsgDispatch ;
74   
75    const Double_t* data = 0;
76    // interpolators initialization
77    data =  fData->GetTemperatureTable();
78    TArrayD t_arr(fData->NumberOfElements(),data);
79    fTempInterpol = new TSpline3("fTempInterpol",0,fData->NumberOfElements()*km,t_arr.GetArray(),fData->NumberOfElements());   
80   
81    data =  fData->GetAir_DensityTable();
82    TArrayD rho_arr(fData->NumberOfElements(),data);
83    fDensInterpol = new TSpline3("fDensInterpol",0,fData->NumberOfElements()*km,rho_arr.GetArray(),fData->NumberOfElements());   
84   
85    data =  fData->GetPressureTable();
86    TArrayD press_arr(fData->NumberOfElements(),data);
87    fPressInterpol = new TSpline3("fPressInterpol",0,fData->NumberOfElements()*km,press_arr.GetArray(),fData->NumberOfElements());   
88}
89
90//_____________________________________________________________________________
91Double_t LowtranAtmosphere::Interpolate(string& val,string& mod,Double_t h) const {
92    //
93    // Calculate simple interpolation for secondary data
94    //
95    const Double_t* data = 0;
96
97    if(val == "Odens")            data = fData->GetO_DensityTable();
98    else if(val == "O2dens")      data = fData->GetO2_DensityTable();
99    else if(val == "H2Odens")     data = fData->GetH2O_DensityTable();
100    else if(val == "N2dens")      data = fData->GetN2_DensityTable();
101    else if(val == "O3dens")      data = fData->GetO3_DensityTable();
102    else if(val == "O3densPPMV")  data = fData->GetO3_DensityPPMVTable();
103    else   Msg(EsafMsg::Panic)<< "Wrong val argument in LOWTRAN_00Atmosphere::Interpolate" << MsgDispatch ;
104   
105   
106    if(h <0) return data[0];
107    Int_t nb = fData->NumberOfElements();
108    Double_t step = 1*km;
109    if(h > step*(nb-1)) return data[nb-1];
110    Double_t val1, val2;
111    Double_t h1;
112    Int_t i = 0;
113    i = Int_t(floor(h/step));
114    h1 = i * step;
115    val1 = data[i];
116    val2 = data[i+1];
117   
118    if(mod == "linear" || val1 == 0 || val2 == 0)  return (val2 - val1)*(h - h1)/step + val1;
119    else if(mod == "expon")                        return val1 * exp(-(h - h1) * log(val1/val2)/step);
120    else Msg(EsafMsg::Panic)<<"Wrong mod argument in LowtranAtmosphere::Interpolate" << MsgDispatch ;
121   
122    return 0;
123}
124
125//_____________________________________________________________________________
126Double_t LowtranAtmosphere::Pressure(Double_t h) const {
127    //
128    // Pressure at a given altitude
129    //
130    if(h > fTOA)  return 0.;
131    // may occur with interfaces with external codes
132    if(IsNaN(h)) Msg(EsafMsg::Panic)<<"<Pressure> Given altitude is NaN" << MsgDispatch ;
133    return fPressInterpol->Eval(h);
134}
135
136//_____________________________________________________________________________
137Double_t LowtranAtmosphere::Temperature(Double_t h) const {
138    //
139    // Temperature at a given altitude
140    //
141    if(h > fTOA)  return 0.;
142    // may occur with interfaces with external codes
143    if(IsNaN(h)) Msg(EsafMsg::Panic)<<"<Pressure> Given altitude is NaN" << MsgDispatch ;
144    return fTempInterpol->Eval(h);
145}
146
147//_____________________________________________________________________________
148Double_t LowtranAtmosphere::AbsoluteHumidity(Double_t h) const {
149    //
150    // Absolute humidity at a given altitude (water vapor density)
151    //
152   
153    if(h > fTOA)  return 0.;
154    string val = "H2Odens";
155    string mod = "linear";
156   
157    return Interpolate(val,mod,h);
158}
159
160//_____________________________________________________________________________
161Double_t LowtranAtmosphere::Air_Density(Double_t h) const {
162    //
163    // Atmospheric air density at a given altitude
164    //
165    if(h > fTOA)  return 0.;
166    // may occur with interfaces with external codes
167    if(IsNaN(h)) Msg(EsafMsg::Panic)<<"<Pressure> Given altitude is NaN" << MsgDispatch ;
168    return fDensInterpol->Eval(h);
169}
170
171//_____________________________________________________________________________
172Double_t LowtranAtmosphere::O_Density(Double_t h) const {
173    //
174    // Oxygen density at a given altitude
175    //
176    if(h > fTOA) return 0.;
177
178    // exponential interpolation between table values
179    string val = "Odens";
180    string mod = "expon";
181   
182    return Interpolate(val,mod,h);
183}
184
185//_____________________________________________________________________________
186Double_t LowtranAtmosphere::O2_Density(Double_t h) const {
187    //
188    // Di-oxygen density at a given altitude
189    //
190    if(h > fTOA) return 0.;
191
192        // exponential interpolation between table values
193        string val = "O2dens";
194        string mod = "expon";
195        return Interpolate(val,mod,h);
196}
197
198//_____________________________________________________________________________
199Double_t LowtranAtmosphere::O3_Density(Double_t h) const {
200    //
201    // Ozone density at a given altitude
202    //
203    if(h > fTOA) return 0.;
204
205    // linear interpolation between table values
206    string val = "O3dens";
207    string mod = "linear";
208   
209    return Interpolate(val,mod,h);
210}
211
212//_____________________________________________________________________________
213Double_t LowtranAtmosphere::O3_DensityPPMV(Double_t h) const {
214    //
215    // Ozone PPMV density at a given altitude
216    //
217    if(h > fTOA) return 0.;
218
219    // linear interpolation between table values
220    string val = "O3densPPMV";
221    string mod = "linear";
222   
223    return Interpolate(val,mod,h);
224}
225
226//_____________________________________________________________________________
227Double_t LowtranAtmosphere::N2_Density(Double_t h) const {
228    //
229    // Di-nitrogen density at a given altitude
230    //
231    if(h > fTOA) return 0.;
232
233    // exponential interpolation between table values
234    string val = "N2dens";
235    string mod = "expon";
236   
237    return Interpolate(val,mod,h);
238}
239
240//_____________________________________________________________________________
241Double_t LowtranAtmosphere::GetLatitude() const {
242    //
243    // used to identify the model
244    //
245   
246    return fData->GetLatitude();
247}
248
249//_____________________________________________________________________________
250Double_t LowtranAtmosphere::GetLongitude() const {
251    //
252    // dummy
253    //
254   
255    return fData->GetLongitude();
256}
257
258//_____________________________________________________________________________
259Double_t LowtranAtmosphere::GetDate() const {
260    //
261    // dummy
262    //
263   
264    return fData->GetDate();
265}
266
267
268
269
Note: See TracBrowser for help on using the repository browser.