adagsl_335d13f0/examples/statistics/src/statistics.adb

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
with Ada.Text_Io; use Ada.Text_IO ;
with Ada.Long_Float_Text_IO ; use Ada.Long_Float_Text_IO ;
with Ada.Integer_Text_Io; use Ada.Integer_Text_Io;
with Interfaces.C; use Interfaces.C ;
with Interfaces.C.Strings;
with gsl ;
with gsl.rng ;
with gsl.randist ;
with gsl.statistics ;
with gsl.rstat ;
with gsl.sort_double ;

procedure statistics is

-- rand(Float64,64)

   random_values : constant gsl.Long_Float_array :=
   (0.6839381853605606, 0.9969779599687887, 0.6494799416396935, 0.6738563793810768, 
    0.4466766434069719, 0.5468592763191233, 0.5005775319913994, 0.24104970185527097, 
    0.721870773818154, 0.24420878862651374, 0.3210445173795178, 0.8172867136455346, 
    0.39704768557558456, 0.10609804690837332, 0.3945131035759507, 0.5839216839010354, 
    0.9104335733318845, 0.31634769072179325, 0.5895070365220522, 0.27132384222531336, 
    0.6988626948612783, 0.18059903196885874, 0.08336299160749083, 0.32020647353421194, 
    0.40297612787898773, 0.7235828674700184, 0.8347223450852815, 0.29785908898899704, 
    0.1181012910541206, 0.5461417721567989, 0.634215657292513, 0.2878742359492805, 
    0.9977499192830578, 0.809157181797072, 0.593448576909427, 0.04982609704658092, 
    0.7342102740376639, 0.8609372491744139, 0.7498670011359, 0.2184094951012312, 
    0.061617757003315066, 0.03667014848542027, 0.04614381563591252, 0.7195117811325236, 
    0.18691528314341987, 0.482672855180976, 0.6921117159031446, 0.8998497561326896, 
    0.11872098633853478, 0.5785075636309084, 0.5321582005125715, 0.6658189583772859, 
    0.9527903919192431, 0.5568841478976171, 0.6443321406578448, 0.9598618957298377, 
    0.5886627917010728, 0.8577310882535446, 0.4337618388364225, 0.043704531558826365, 
    0.04944233237964346, 0.21555069639611824, 0.14851620562404466, 0.4600302846819925);

-- julia> maximum(r)
--0.9977499192830578

--julia> minimum(r)
--0.03667014848542027

--julia> mean(r)
--0.4919858846187607

--julia> std(r)
--0.2840111098108367

--julia> var(r)
--0.08066231049598314

--julia> skewness(r)
-- -0.034201534655284256

--julia> kurtosis(r)
-- -1.1146180363311742

   procedure Test is
      vals : gsl.double_array(1..5);
      vals2 : gsl.double_array := (17.2, 18.1, 16.5, 18.3, 12.6) ;
   begin 

      for i in vals'range
      loop
         vals(i) := double(i) ;
      end loop ;
      Put( "Mean "); Put( Long_Float(gsl.statistics.mean(vals))  ); New_Line;
      Put( "Variance "); Put( Long_Float(gsl.statistics.variance(vals))  ); New_Line;
      Put( "Sd "); Put( Long_Float(gsl.statistics.sd(vals))  ); New_Line;
      Put("Predefined values");New_Line;
      Put( "Mean "); Put( Long_Float(gsl.statistics.mean(vals2))  , aft => 2 , exp => 0 ); New_Line;
      Put( "Variance "); Put( Long_Float(gsl.statistics.variance(vals2)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Sd "); Put( Long_Float(gsl.statistics.sd(vals2)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Skew "); Put( Long_Float(gsl.statistics.skew(vals2)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Kurtosis "); Put( Long_Float(gsl.statistics.kurtosis(vals2)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Tss "); Put( Long_Float(gsl.statistics.tss(vals2)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Max "); Put( Long_Float(gsl.statistics.max(vals2)  ) , aft => 2 , exp => 0 ); New_Line;
      Put( "Min "); Put( Long_Float(gsl.statistics.min(vals2)) , aft => 2 , exp => 0 ); New_Line;

      Put("Random Numbers from Julia");New_Line;
      Put( "Max "); Put( Long_Float(gsl.statistics.max(random_values)  ) , aft => 2 , exp => 0 ); New_Line;
      Put( "Min "); Put( Long_Float(gsl.statistics.min(random_values)) , aft => 2 , exp => 0 ); New_Line;

      Put( "Mean "); Put( Long_Float(gsl.statistics.mean(random_values))  , aft => 2 , exp => 0 ); New_Line;
      Put( "Variance "); Put( Long_Float(gsl.statistics.variance(random_values)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Sd "); Put( Long_Float(gsl.statistics.sd(random_values)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Skew "); Put( Long_Float(gsl.statistics.skew(random_values)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Kurtosis "); Put( Long_Float(gsl.statistics.kurtosis(random_values)) , aft => 2 , exp => 0  ); New_Line;
      Put( "Tss "); Put( Long_Float(gsl.statistics.tss(random_values)) , aft => 2 , exp => 0  ); New_Line;
   
   end Test ;
   procedure TestRunning is
      ws : access gsl.rstat.gsl_rstat_workspace := gsl.rstat.alloc ;
      status : Int ;
   begin
      for i in random_values'Range
      loop
         Status := gsl.rstat.add(random_values(i),ws);
      end loop;

      Put( "Count "); Put(Integer(gsl.rstat.n(ws))); New_Line;
      Put( "Max "); Put( Long_Float(gsl.rstat.max(ws))  ); New_Line;
      Put( "Min "); Put( Long_Float(gsl.rstat.min(ws))  ); New_Line;

      Put( "Mean "); Put( Long_Float(gsl.rstat.mean(ws))  ); New_Line;
      Put( "Variance "); Put( Long_Float(gsl.rstat.variance(ws))  ); New_Line;
      Put( "Sd "); Put( Long_Float(gsl.rstat.sd(ws))   ); New_Line;
      Put( "Skew "); Put( Long_Float(gsl.rstat.skew(ws))  ); New_Line;
      Put( "Kurtosis "); Put( Long_Float(gsl.rstat.kurtosis(ws))  ); New_Line;
      Put( "RMS "); Put(Long_Float(gsl.rstat.rms(ws))); New_Line;
      gsl.rstat.free(ws);
   end TestRunning;

   procedure TestQuantile is
      def : access constant gsl.rng.gsl_rng_type := gsl.rng.default ;
      rng : access gsl.rng.gsl_rng := gsl.rng.alloc(def) ;

      ws_25 : access gsl.rstat.gsl_rstat_quantile_workspace;
      ws_5 : access gsl.rstat.gsl_rstat_quantile_workspace;
      ws_75 : access gsl.rstat.gsl_rstat_quantile_workspace;
      ran : double ;
      N : constant := 10_000;
      ranvals : gsl.double_array(1..N);
      Status : Int ;
      qval : double ;
   begin
      Ada.Long_Float_Text_IO.Default_aft := 4;
      Put_Line("Quantile estimate test");
      ws_25 := gsl.rstat.quantile_alloc(0.25);
      ws_5 := gsl.rstat.quantile_alloc(0.5);
      ws_75 := gsl.rstat.quantile_alloc(0.75);
      for i in 1..N
      loop
         ran := gsl.randist.rayleigh(rng,1.0);
         Status := gsl.rstat.quantile_add(ran,ws_25);
         Status := gsl.rstat.quantile_add(ran,ws_5);
         Status := gsl.rstat.quantile_add(ran,ws_75);
         ranvals(i) := ran;
      end loop;
      qval := gsl.rstat.quantile_get(ws_25);
      Put("Quantile 0.25 estimate "); Put(Long_Float(qval)); New_Line;

      qval := gsl.rstat.quantile_get(ws_5);
      Put("Quantile 0.5 estimate "); Put(Long_Float(qval)); New_Line;

      qval := gsl.rstat.quantile_get(ws_75);
      Put("Quantile 0.75 estimate "); Put(Long_Float(qval)); New_Line;

      gsl.rstat.quantile_free(ws_25);
      gsl.rstat.quantile_free(ws_5);
      gsl.rstat.quantile_free(ws_75);

      gsl.sort_double.sort(ranvals) ;
      Put_Line("From sorted data");
      Put("Median "); Put(Long_Float(gsl.statistics.median_from_sorted_data(ranvals))); New_Line;
      Put("Quantile 0.25 "); Put(Long_Float(gsl.statistics.quantile_from_sorted_data(ranvals,0.25))); New_Line;
      Put("Quantile 0.5 "); Put(Long_Float(gsl.statistics.quantile_from_sorted_data(ranvals,0.5))); New_Line;
      Put("Quantile 0.75 "); Put(Long_Float(gsl.statistics.quantile_from_sorted_data(ranvals,0.75))); New_Line;

      Put_Line("Estimators");
      declare
         ws_13 : double := gsl.statistics.quantile_from_sorted_data(ranvals,0.3333) ;
         ws_5 : double := gsl.statistics.quantile_from_sorted_data(ranvals,0.5) ;
         ws_23 : double := gsl.statistics.quantile_from_sorted_data(ranvals,0.6666) ;
      begin
         Put("Computed Gastwirth estimate ");
         Put(Long_Float(0.3 * ws_13 + 0.4 * ws_5 + 0.3 * ws_23));
         New_Line ;
         Put("Library result ");
         Put(Long_Float(gsl.statistics.gastwirth_from_sorted_data(ranvals)));
         New_Line;
      end ;

   end TestQuantile;

begin
   Ada.Long_Float_Text_IO.Default_exp := 0;
   Ada.Long_Float_Text_IO.Default_aft := 2;
   Test;
   TestRunning;
   TestQuantile;
end statistics ;