icemodel.hh
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
34 const double ellipsoid_b = EarthModel::R_EARTH*(1-(1/ellipsoid_inv_f));
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 
43 {
44  TVector3 balloon;
45  TVector3 ref;
46  int nintersect;
47  TVector3 int1;
48  TVector3 int2;
49  TVector3 pint;
50  TVector3 nudir;
51  bool rightdir1;
52  bool rightdir2;
53  bool noway;
54  bool leave_err;
55  TVector3 nupos;
56  TVector3 enterice;
57  TVector3 exitice;
58  TVector3 exitearth;
59  std::vector<double> abs_weights;
60  std::vector<double> weight_probs;
61  std::vector<double> lengths;
62  double sum_weight_probs;
63  double x,y;
64 
65  void reset()
66  {
67  balloon.SetXYZ(0,0,0);
68  ref.SetXYZ(0,0,0);
69  int1.SetXYZ(0,0,0);
70  int2.SetXYZ(0,0,0);
71  pint.SetXYZ(0,0,0);
72  nudir.SetXYZ(0,0,0);
73  nupos.SetXYZ(0,0,0);
74  enterice.SetXYZ(0,0,0);
75  exitice.SetXYZ(0,0,0);
76  exitearth.SetXYZ(0,0,0);
77  abs_weights.clear();
78  weight_probs.clear();
79  lengths.clear();
80  nintersect = 0;
81  rightdir1 = false;
82  rightdir2 = false;
83  noway = false;
84  x = 0;
85  y = 0;
86  leave_err = false;
87  sum_weight_probs = 0;
88  }
89 
91  {
92  reset();
93  }
94  virtual ~icemodel_debug()
95  {; }
96 
97  ClassDef(icemodel_debug,1);
98 };
99 
100 
101 
103 class IceModel : public EarthModel {
104 
105 
106 public:
107 
108  //BEDMAP data
109  TH2D h_ice_thickness;
110  TH2D h_ground_elevation;;
111  TH2D h_water_depth;;
112 
113 
114 double bedmap_R; //varies with latitude, defined here for 71 deg S latitude
115 double bedmap_nu;
116 
117 
118 //Parameters of the BEDMAP ice model. (See http://www.antarctica.ac.uk/aedc/bedmap/download/)
119 int nCols_ice; //number of columns in data, set by header file (should be 1200)
120 int nRows_ice; //number of rows in data, set by header file (should be 1000)
121 int cellSize; //in meters, set by header file (should be 5000) - same for both files
122 int xLowerLeft_ice;
123 int yLowerLeft_ice;
124 int nCols_ground;
125 int nRows_ground;
126 int xLowerLeft_ground;
127 int yLowerLeft_ground;
128 int nCols_water;
129 int nRows_water;
130 int xLowerLeft_water;
131 int yLowerLeft_water;
132 int NODATA;
133 
134 
135  void IceENtoLonLat(int e,
136  int n,
137 
138  double& lon,
139  double& lat);
140 
141  void GroundENtoLonLat(int e,
142  int n,
143  double& lon,
144  double& lat);
145  void WaterENtoLonLat(int e,
146  int n,
147  double& lon,
148  double& lat);
149 
150 
151 
152  vector<double> volume_inhorizon; // volume of ice within horizon for each balloon phi position
153  IceModel(int model=0,int earth_mode=0,int WEIGHTABSORPTION_SETTING=1);
154  virtual ~IceModel();
155  double IceThickness(double lon,double lat);
156  double IceThickness(const Position& pos) ;
157  double Surface(double lon,double lat) ;
158  double Surface(const Position& pos) ;
159  double SurfaceAboveGeoid(double lon,double lat);
160  double SurfaceAboveGeoid(const Position& pos) ;
161  double WaterDepth(double lon,double lat);
162  double WaterDepth(const Position& pos);
163  Position PickInteractionLocation(int ibnposition, Settings *settings1, const Position &rbn, Interaction *interaction1);
164  Position PickBalloonPosition();
165  void GetMAXHORIZON(Balloon *bn1); // get upper limit on the horizon wrt the balloon.
166  int RossIceShelf(const Position &position);
167  int IceOnWater(const Position &postition);
168  int RossExcept(const Position &position);
169  int RonneIceShelf(const Position &position);
170  int WestLand(const Position &pos);
171  int AcceptableRfexit(const Vector &nsurf_rfexit,const Position &rfexit,const Vector &n_exit2rx);
172  double GetBalloonPositionWeight(int ibnpos);
173  int OutsideAntarctica(const Position &pos);
174  int OutsideAntarctica(double lat);
175  int WhereDoesItEnterIce(const Position &posnu,
176  const Vector &nnu,
177  double stepsize,
178  Position &r_enterice);
179 
180  int WhereDoesItExitIce(const Position &posnu,
181  const Vector &nnu,
182  double stepsize,
183  Position &r_enterice);
184  int WhereDoesItExitIceForward(const Position &posnu,
185  const Vector &nnu,
186  double stepsize,
187  Position &r_enterice);
188 
189 
190  //resolution is in cartesian meters
191  void CreateCartesianTopAndBottom(int resolution, bool force_new =false);
192  const TH2 * GetCartesianTop() const { return cart_ice_top; }
193  const TH2 * GetCartesianBottom() const { return cart_ice_bot; }
194  bool CartesianIsInIce(double x, double y, double z);
195 
196  int GetIceIntersectionsCartesian(const Position &posnu, const Vector &nnu,
197  std::vector<std::pair<double,double> > & intersections, double initial_step_size = 50, int map_resolution=2000);
198 
199 
200  void CreateHorizons(Settings *settings1,Balloon *bn1,double theta_bn,double phi_bn,double altitude_bn,ofstream &foutput);
201  Vector GetSurfaceNormal(const Position &r_out); //overloaded from EarthModel to include procedures for new ice models.
202  double GetN(double depth);
203  double GetN(const Position &pos);
204  double EffectiveAttenuationLength(Settings *settings1,const Position &pos, const int &whichray);
205 
206  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]);
207  int PickUnbiased(Interaction *interaction1, double len_int_kgm2, double & position_weight, double chord_step, Vector * force_dir = 0);
208 
209  int PickUnbiasedPointSourceNearBalloon(Interaction *interaction1,
210  const Position * balloon_position, double max_ps_distance, double chord_step,
211  double len_int_kgm2, const Vector * force_dir = 0);
212 
213  double getSampleX() const { return sample_x; }
214  double getSampleY() const { return sample_y; }
215 
216  void LonLattoEN(double lon,
217  double lat,
218  double& E,
219  double& N);
220 
221 
222 protected:
223  int ice_model;
224  int DEPTH_DEPENDENT_N;
225 
226 
227  //Information on horizons - what ice the balloon can see at each position along its path.
228 
229 
230  double volume_inhorizon_average; // average volume of ice seen by balloon
231  vector< vector<int> > ilon_inhorizon; // indices in lon and lat for bins in horizon for NPHI balloon positions along 80 deg latitude line.
232  vector< vector<int> >ilat_inhorizon;
233  vector< vector<int> >easting_inhorizon; //indicies in easting and northing for bins in horizon for NPHI balloon positions along 80 deg latitude line.
234  vector< vector<int> >northing_inhorizon;
235  vector<double> maxvol_inhorizon; // maximum volume of ice for a bin
236 
237  double cart_max_z;
238  double cart_min_z;
239  double sample_x, sample_y;
240 
241  TH2 *cart_ice_top;
242  TH2 *cart_ice_bot;
243  TFile *cart_ice_file;
244  int cart_resolution;
245 
246  //BEDMAP utility methods
247  double Area(double latitude);
248 
249  void ENtoLonLat(int e_coord,
250  int n_coord,
251  double xLowerLeft,
252  double yLowerLeft,
253 
254  double& lon,
255  double& lat) ;
256 
257 
258 
259 
260  //BEDMAP data input methods
261  void ReadIceThickness();
262  void ReadGroundBed();
263  void ReadWaterDepth();
264 
265 private:
266 
267  const static int N_sheetup=2810;
268  double d_sheetup[N_sheetup], l_sheetup[N_sheetup];
269  const static int N_shelfup=420;
270  double d_shelfup[N_shelfup], l_shelfup[N_shelfup];
271  const static int N_westlandup=420;
272  double d_westlandup[N_westlandup],l_westlandup[N_westlandup];
273 
274  const static int N_sheetdown=2810;
275  double d_sheetdown[N_sheetup], l_sheetdown[N_sheetdown];
276  const static int N_shelfdown=420;
277  double d_shelfdown[N_shelfdown], l_shelfdown[N_shelfdown];
278  const static int N_westlanddown=420;
279  double d_westlanddown[N_westlanddown],l_westlanddown[N_westlanddown];
280 
281 
282 }; //class IceModel
283 
284 #endif
Reads in and stores input settings for the run.
Definition: Settings.h:37
This class is a 3-vector that represents a position on the Earth&#39;s surface.
Definition: position.hh:26
Stores everything about a particular neutrino interaction. Interaction.
Definition: Primaries.h:135
Shape of the earth, ice thicknesses, profiles of earth layers, densities, neutrino absorption...
Definition: earthmodel.hh:40
Handles everything related to balloon positions, payload orientation over the course of a flight...
Definition: balloon.hh:30
This class represents a three-vector. Operators are overloaded to provide for the familiar operations...
Definition: vector.hh:27
Ray tracing.
Definition: ray.hh:20
Ice thicknesses and water depth.
Definition: icemodel.hh:103