icemc
icemodel.hh
Go to the documentation of this file.
1 #ifndef ICEMODEL_HH_
2 #define ICEMODEL_HH_
3 
4 #include <cmath>
5 #include <iostream>
6 #include <fstream>
7 #include <vector>
8 #include "TVector3.h"
9 
10 #include "earthmodel.hh"
11 #include "source.hh"
12 #include "TH2.h"
13 
14 class Balloon;
15 class Position;
16 class Vector;
17 class EarthModel;
18 class Settings;
19 class Interaction;
20 class Ray;
21 class TFile;
22 
23 //Constants relating to all ice models
24 const double FIRNDEPTH=-150.; // depth of the firn, in meters: currently a constant over all ice
25 using std::ofstream;
26 using std::vector;
27 using std::string;
28 
29 
30 
31 //Variables for conversion between BEDMAP polar stereographic coordinates and lat/lon. Conversion equations from ftp://164.214.2.65/pub/gig/tm8358.2/TM8358_2.pdf
32 const double scale_factor=0.97276901289; //scale factor at pole corresponding to 71 deg S latitude of true scale (used in BEDMAP)
33 const double ellipsoid_inv_f = 298.257223563; //of Earth
35 const double eccentricity = sqrt((1/ellipsoid_inv_f)*(2-(1/ellipsoid_inv_f)));
36 const double bedmap_a_bar = pow(eccentricity,2)/2 + 5*pow(eccentricity,4)/24 + pow(eccentricity,6)/12 + 13*pow(eccentricity,8)/360;
37 const double bedmap_b_bar = 7*pow(eccentricity,4)/48 + 29*pow(eccentricity,6)/240 + 811*pow(eccentricity,8)/11520;
38 const double bedmap_c_bar = 7*pow(eccentricity,6)/120 + 81*pow(eccentricity,8)/1120;
39 const double bedmap_d_bar = 4279*pow(eccentricity,8)/161280;
40 const double bedmap_c_0 = (2*EarthModel::R_EARTH / sqrt(1-pow(eccentricity,2))) * pow(( (1-eccentricity) / (1+eccentricity) ),eccentricity/2);
41 
42 #ifdef ICEMODEL_DEBUG_TREE
43 struct icemodel_debug
44 {
45  TVector3 balloon;
46  TVector3 ref;
47  int nintersect;
48  TVector3 int1;
49  TVector3 int2;
50  TVector3 pint;
51  TVector3 nudir;
52  bool rightdir1;
53  bool rightdir2;
54  bool noway;
55  bool leave_err;
56  TVector3 nupos;
57  TVector3 enterice;
58  TVector3 exitice;
59  TVector3 exitearth;
60  double x,y;
61 
62  void reset()
63  {
64  balloon.SetXYZ(0,0,0);
65  ref.SetXYZ(0,0,0);
66  int1.SetXYZ(0,0,0);
67  int2.SetXYZ(0,0,0);
68  pint.SetXYZ(0,0,0);
69  nudir.SetXYZ(0,0,0);
70  nupos.SetXYZ(0,0,0);
71  enterice.SetXYZ(0,0,0);
72  exitice.SetXYZ(0,0,0);
73  exitearth.SetXYZ(0,0,0);
74  nintersect = 0;
75  rightdir1 = false;
76  rightdir2 = false;
77  noway = false;
78  x = 0;
79  y = 0;
80  leave_err = false;
81  }
82 };
83 #endif
84 
85 
86 
88 class IceModel : public EarthModel {
89 
90 
91 public:
92 
93  //BEDMAP data
97 
98 
99 double bedmap_R; //varies with latitude, defined here for 71 deg S latitude
100 double bedmap_nu;
101 
102 
103 //Parameters of the BEDMAP ice model. (See http://www.antarctica.ac.uk/aedc/bedmap/download/)
104 int nCols_ice; //number of columns in data, set by header file (should be 1200)
105 int nRows_ice; //number of rows in data, set by header file (should be 1000)
106 int cellSize; //in meters, set by header file (should be 5000) - same for both files
117 int NODATA;
118 
119 
120  void IceENtoLonLat(int e,
121  int n,
122 
123  double& lon,
124  double& lat);
125 
126  void GroundENtoLonLat(int e,
127  int n,
128  double& lon,
129  double& lat);
130  void WaterENtoLonLat(int e,
131  int n,
132  double& lon,
133  double& lat);
134 
135 
136 
137  vector<double> volume_inhorizon; // volume of ice within horizon for each balloon phi position
138  IceModel(int model=0,int earth_mode=0,int WEIGHTABSORPTION_SETTING=1);
139  virtual ~IceModel();
140  double IceThickness(double lon,double lat);
141  double IceThickness(const Position& pos) ;
142  double Surface(double lon,double lat) ;
143  double Surface(const Position& pos) ;
144  double SurfaceAboveGeoid(double lon,double lat);
145  double SurfaceAboveGeoid(const Position& pos) ;
146  double WaterDepth(double lon,double lat);
147  double WaterDepth(const Position& pos);
148  Position PickInteractionLocation(int ibnposition, Settings *settings1, const Position &rbn, Interaction *interaction1);
150  void GetMAXHORIZON(Balloon *bn1); // get upper limit on the horizon wrt the balloon.
151  int RossIceShelf(const Position &position);
152  int IceOnWater(const Position &postition);
153  int RossExcept(const Position &position);
154  int RonneIceShelf(const Position &position);
155  int WestLand(const Position &pos);
156  int AcceptableRfexit(const Vector &nsurf_rfexit,const Position &rfexit,const Vector &n_exit2rx);
157  double GetBalloonPositionWeight(int ibnpos);
158  int OutsideAntarctica(const Position &pos);
159  int OutsideAntarctica(double lat);
160  int WhereDoesItEnterIce(const Position &posnu,
161  const Vector &nnu,
162  double stepsize,
163  Position &r_enterice);
164 
165  int WhereDoesItExitIce(const Position &posnu,
166  const Vector &nnu,
167  double stepsize,
168  Position &r_enterice);
169  int WhereDoesItExitIceForward(const Position &posnu,
170  const Vector &nnu,
171  double stepsize,
172  Position &r_enterice);
173 
174 
175  //resolution is in cartesian meters
176  void CreateCartesianTopAndBottom(int resolution, bool force_new =false);
177  const TH2 * GetCartesianTop() const { return cart_ice_top; }
178  const TH2 * GetCartesianBottom() const { return cart_ice_bot; }
179  bool CartesianIsInIce(double x, double y, double z);
180 
181  int GetIceIntersectionsCartesian(const Position &posnu, const Vector &nnu,
182  std::vector<std::pair<double,double> > & intersections, double initial_step_size = 50, int map_resolution=1000);
183 
184 
185  void CreateHorizons(Settings *settings1,Balloon *bn1,double theta_bn,double phi_bn,double altitude_bn,ofstream &foutput);
186  Vector GetSurfaceNormal(const Position &r_out); //overloaded from EarthModel to include procedures for new ice models.
187  double GetN(double depth);
188  double GetN(const Position &pos);
189  double EffectiveAttenuationLength(Settings *settings1,const Position &pos, const int &whichray);
190 
191  void FillArraysforTree(double lon_ground[1068][869],double lat_ground[1068][869],double lon_ice[1200][1000],double lat_ice[1200][1000],double lon_water[1200][1000],double lat_water[1200][1000]);
192  int PickUnbiased(Interaction *interaction1, double len_int_kgm2, double & position_weight, double chord_step, Vector * force_dir = 0);
193 
195  const Position * balloon_position, double max_ps_distance, double chord_step,
196  double len_int_kgm2, const Vector * force_dir = 0);
197 
198  double getSampleX() const { return sample_x; }
199  double getSampleY() const { return sample_y; }
200 
201  void LonLattoEN(double lon,
202  double lat,
203  double& E,
204  double& N);
205 
206 
207 protected:
210 
211 
212  //Information on horizons - what ice the balloon can see at each position along its path.
213 
214 
215  double volume_inhorizon_average; // average volume of ice seen by balloon
216  vector< vector<int> > ilon_inhorizon; // indices in lon and lat for bins in horizon for NPHI balloon positions along 80 deg latitude line.
217  vector< vector<int> >ilat_inhorizon;
218  vector< vector<int> >easting_inhorizon; //indicies in easting and northing for bins in horizon for NPHI balloon positions along 80 deg latitude line.
219  vector< vector<int> >northing_inhorizon;
220  vector<double> maxvol_inhorizon; // maximum volume of ice for a bin
221 
222  double cart_max_z;
223  double cart_min_z;
224  double sample_x, sample_y;
225 
228  TFile *cart_ice_file;
230 
231  //BEDMAP utility methods
232  double Area(double latitude);
233 
234  void ENtoLonLat(int e_coord,
235  int n_coord,
236  double xLowerLeft,
237  double yLowerLeft,
238 
239  double& lon,
240  double& lat) ;
241 
242 
243 
244 
245  //BEDMAP data input methods
246  void ReadIceThickness();
247  void ReadGroundBed();
248  void ReadWaterDepth();
249 
250 private:
251 
252  const static int N_sheetup=2810;
254  const static int N_shelfup=420;
256  const static int N_westlandup=420;
258 
259  const static int N_sheetdown=2810;
261  const static int N_shelfdown=420;
263  const static int N_westlanddown=420;
265 
266 
267 }; //class IceModel
268 
269 #endif
void LonLattoEN(double lon, double lat, double &E, double &N)
Definition: icemodel.cc:1262
static constexpr double R_EARTH
Definition: earthmodel.hh:46
const double bedmap_c_0
Definition: icemodel.hh:40
void ReadIceThickness()
Definition: icemodel.cc:1646
int GetIceIntersectionsCartesian(const Position &posnu, const Vector &nnu, std::vector< std::pair< double, double > > &intersections, double initial_step_size=50, int map_resolution=1000)
Definition: icemodel.cc:2013
int WhereDoesItExitIceForward(const Position &posnu, const Vector &nnu, double stepsize, Position &r_enterice)
Definition: icemodel.cc:902
double EffectiveAttenuationLength(Settings *settings1, const Position &pos, const int &whichray)
Definition: icemodel.cc:1199
int yLowerLeft_ice
Definition: icemodel.hh:108
TH2 * cart_ice_top
Definition: icemodel.hh:226
Position PickInteractionLocation(int ibnposition, Settings *settings1, const Position &rbn, Interaction *interaction1)
Definition: icemodel.cc:260
double l_westlandup[N_westlandup]
Definition: icemodel.hh:257
double d_westlanddown[N_westlanddown]
Definition: icemodel.hh:264
vector< double > volume_inhorizon
Definition: icemodel.hh:137
virtual ~IceModel()
Definition: icemodel.cc:245
static const int N_shelfup
Definition: icemodel.hh:254
double sample_x
Definition: icemodel.hh:224
const TH2 * GetCartesianBottom() const
Definition: icemodel.hh:178
const double ellipsoid_b
Definition: icemodel.hh:34
const double eccentricity
Definition: icemodel.hh:35
static const int N_sheetup
Definition: icemodel.hh:252
double l_sheetdown[N_sheetdown]
Definition: icemodel.hh:260
static const int N_sheetdown
Definition: icemodel.hh:259
int nCols_ground
Definition: icemodel.hh:109
This class represents a three-vector. Operators are overloaded to provide for the familiar operations...
Definition: vector.hh:27
int nRows_water
Definition: icemodel.hh:114
int yLowerLeft_ground
Definition: icemodel.hh:112
Vector GetSurfaceNormal(const Position &r_out)
Definition: icemodel.cc:620
vector< double > maxvol_inhorizon
Definition: icemodel.hh:220
int RonneIceShelf(const Position &position)
Definition: icemodel.cc:1117
int whichray
Definition: icemc.cc:162
TFile * cart_ice_file
Definition: icemodel.hh:228
int AcceptableRfexit(const Vector &nsurf_rfexit, const Position &rfexit, const Vector &n_exit2rx)
Definition: icemodel.cc:1157
int cellSize
Definition: icemodel.hh:106
int WhereDoesItEnterIce(const Position &posnu, const Vector &nnu, double stepsize, Position &r_enterice)
Definition: icemodel.cc:722
static const int N_westlanddown
Definition: icemodel.hh:263
void GroundENtoLonLat(int e, int n, double &lon, double &lat)
Definition: icemodel.cc:1331
bool CartesianIsInIce(double x, double y, double z)
Definition: icemodel.cc:1997
static const int N_westlandup
Definition: icemodel.hh:256
int IceOnWater(const Position &postition)
Definition: icemodel.cc:1087
void IceENtoLonLat(int e, int n, double &lon, double &lat)
Definition: icemodel.cc:1326
vector< vector< int > > ilon_inhorizon
Definition: icemodel.hh:216
double d_shelfup[N_shelfup]
Definition: icemodel.hh:255
void ReadWaterDepth()
Definition: icemodel.cc:1789
const double bedmap_c_bar
Definition: icemodel.hh:38
TH2 * cart_ice_bot
Definition: icemodel.hh:227
Ice thicknesses and water depth.
Definition: icemodel.hh:88
double volume_inhorizon_average
Definition: icemodel.hh:215
double Surface(double lon, double lat)
Definition: icemodel.cc:1053
const double bedmap_a_bar
Definition: icemodel.hh:36
int OutsideAntarctica(const Position &pos)
Definition: icemodel.cc:1149
void FillArraysforTree(double lon_ground[1068][869], double lat_ground[1068][869], double lon_ice[1200][1000], double lat_ice[1200][1000], double lon_water[1200][1000], double lat_water[1200][1000])
double d_westlandup[N_westlandup]
Definition: icemodel.hh:257
This class is a 3-vector that represents a position on the Earth&#39;s surface.
Definition: position.hh:26
double d_sheetup[N_sheetup]
Definition: icemodel.hh:253
double d_sheetdown[N_sheetup]
Definition: icemodel.hh:260
double l_shelfdown[N_shelfdown]
Definition: icemodel.hh:262
const double FIRNDEPTH
Definition: icemodel.hh:24
Position PickBalloonPosition()
Definition: icemodel.cc:254
double SurfaceAboveGeoid(double lon, double lat)
Definition: icemodel.cc:1022
int nCols_ice
Definition: icemodel.hh:104
int ice_model
Definition: icemodel.hh:208
int WestLand(const Position &pos)
Definition: icemodel.cc:1131
double getSampleX() const
Definition: icemodel.hh:198
int nRows_ground
Definition: icemodel.hh:110
double l_shelfup[N_shelfup]
Definition: icemodel.hh:255
TH2D h_water_depth
Definition: icemodel.hh:95
void ReadGroundBed()
Definition: icemodel.cc:1724
int xLowerLeft_ground
Definition: icemodel.hh:111
int xLowerLeft_ice
Definition: icemodel.hh:107
vector< vector< int > > northing_inhorizon
Definition: icemodel.hh:219
Ray tracing.
Definition: ray.hh:20
Handles everything related to balloon positions, payload orientation over the course of a flight...
Definition: balloon.hh:30
Reads in and stores input settings for the run.
Definition: Settings.h:35
const TH2 * GetCartesianTop() const
Definition: icemodel.hh:177
const double bedmap_d_bar
Definition: icemodel.hh:39
double Area(double latitude)
Definition: icemodel.cc:1255
TH2D h_ice_thickness
Definition: icemodel.hh:94
double sample_y
Definition: icemodel.hh:224
void WaterENtoLonLat(int e, int n, double &lon, double &lat)
Definition: icemodel.cc:1335
void GetMAXHORIZON(Balloon *bn1)
Definition: icemodel.cc:1342
int xLowerLeft_water
Definition: icemodel.hh:115
void ENtoLonLat(int e_coord, int n_coord, double xLowerLeft, double yLowerLeft, double &lon, double &lat)
Definition: icemodel.cc:1277
const double scale_factor
Definition: icemodel.hh:32
IceModel(int model=0, int earth_mode=0, int WEIGHTABSORPTION_SETTING=1)
Definition: icemodel.cc:104
int PickUnbiased(Interaction *interaction1, double len_int_kgm2, double &position_weight, double chord_step, Vector *force_dir=0)
Definition: icemodel.cc:542
double cart_min_z
Definition: icemodel.hh:223
vector< vector< int > > ilat_inhorizon
Definition: icemodel.hh:217
TH2D h_ground_elevation
Definition: icemodel.hh:95
int RossIceShelf(const Position &position)
Definition: icemodel.cc:1093
int nCols_water
Definition: icemodel.hh:113
double l_westlanddown[N_westlanddown]
Definition: icemodel.hh:264
double getSampleY() const
Definition: icemodel.hh:199
double GetN(double depth)
Definition: icemodel.cc:1175
int nRows_ice
Definition: icemodel.hh:105
int RossExcept(const Position &position)
Definition: icemodel.cc:1107
double cart_max_z
Definition: icemodel.hh:222
double WaterDepth(double lon, double lat)
Definition: icemodel.cc:1061
int PickUnbiasedPointSourceNearBalloon(Interaction *interaction1, const Position *balloon_position, double max_ps_distance, double chord_step, double len_int_kgm2, const Vector *force_dir=0)
Definition: icemodel.cc:421
int NODATA
Definition: icemodel.hh:117
const double bedmap_b_bar
Definition: icemodel.hh:37
int WhereDoesItExitIce(const Position &posnu, const Vector &nnu, double stepsize, Position &r_enterice)
Definition: icemodel.cc:811
double bedmap_R
Definition: icemodel.hh:96
int yLowerLeft_water
Definition: icemodel.hh:116
double GetBalloonPositionWeight(int ibnpos)
Definition: icemodel.cc:1140
double bedmap_nu
Definition: icemodel.hh:100
double l_sheetup[N_sheetup]
Definition: icemodel.hh:253
double IceThickness(double lon, double lat)
Definition: icemodel.cc:991
void CreateCartesianTopAndBottom(int resolution, bool force_new=false)
Definition: icemodel.cc:1859
vector< vector< int > > easting_inhorizon
Definition: icemodel.hh:218
void CreateHorizons(Settings *settings1, Balloon *bn1, double theta_bn, double phi_bn, double altitude_bn, ofstream &foutput)
Definition: icemodel.cc:1353
int cart_resolution
Definition: icemodel.hh:229
Shape of the earth, ice thicknesses, profiles of earth layers, densities, neutrino absorption...
Definition: earthmodel.hh:40
static const int N_shelfdown
Definition: icemodel.hh:261
Stores everything about a particular neutrino interaction. Interaction.
Definition: Primaries.h:136
const double ellipsoid_inv_f
Definition: icemodel.hh:33
int DEPTH_DEPENDENT_N
Definition: icemodel.hh:209
double d_shelfdown[N_shelfdown]
Definition: icemodel.hh:262