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
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
|
"""Masked arrays add-ons.
A collection of utilities for maskedarray
:author: Pierre Gerard-Marchant
:contact: pierregm_at_uga_dot_edu
:version: $Id: extras.py 3473 2007-10-29 15:18:13Z jarrod.millman $
"""
__author__ = "Pierre GF Gerard-Marchant ($Author: jarrod.millman $)"
__version__ = '1.0'
__revision__ = "$Revision: 3473 $"
__date__ = '$Date: 2007-10-29 17:18:13 +0200 (Mon, 29 Oct 2007) $'
__all__ = ['apply_along_axis', 'atleast_1d', 'atleast_2d', 'atleast_3d',
'average',
'column_stack','compress_cols','compress_rowcols', 'compress_rows',
'count_masked', 'corrcoef', 'cov',
'dot','dstack',
'ediff1d','expand_dims',
'flatnotmasked_contiguous','flatnotmasked_edges',
'hsplit','hstack',
'mask_cols','mask_rowcols','mask_rows','masked_all','masked_all_like',
'median','mr_',
'notmasked_contiguous','notmasked_edges',
'polyfit',
'row_stack',
'vander','vstack',
]
from itertools import groupby
import warnings
import core as ma
from core import MaskedArray, MAError, add, array, asarray, concatenate, count,\
filled, getmask, getmaskarray, masked, masked_array, mask_or, nomask, ones,\
sort, zeros
#from core import *
import numpy as np
from numpy import ndarray, array as nxarray
import numpy.core.umath as umath
from numpy.lib.index_tricks import AxisConcatenator
from numpy.lib.polynomial import _lstsq
#...............................................................................
def issequence(seq):
"""Is seq a sequence (ndarray, list or tuple)?"""
if isinstance(seq, ndarray):
return True
elif isinstance(seq, tuple):
return True
elif isinstance(seq, list):
return True
return False
def count_masked(arr, axis=None):
"""Count the number of masked elements along the given axis.
Parameters
----------
axis : int, optional
Axis along which to count.
If None (default), a flattened version of the array is used.
"""
m = getmaskarray(arr)
return m.sum(axis)
def masked_all(shape, dtype=float):
"""Return an empty masked array of the given shape and dtype,
where all the data are masked.
Parameters
----------
dtype : dtype, optional
Data type of the output.
"""
a = masked_array(np.empty(shape, dtype),
mask=np.ones(shape, bool))
return a
def masked_all_like(arr):
"""Return an empty masked array of the same shape and dtype as
the array `a`, where all the data are masked.
"""
a = masked_array(np.empty_like(arr),
mask=np.ones(arr.shape, bool))
return a
#####--------------------------------------------------------------------------
#---- --- Standard functions ---
#####--------------------------------------------------------------------------
class _fromnxfunction:
"""Defines a wrapper to adapt numpy functions to masked arrays."""
def __init__(self, funcname):
self._function = funcname
self.__doc__ = self.getdoc()
def getdoc(self):
"Retrieves the __doc__ string from the function."
return getattr(np, self._function).__doc__ +\
"*Notes*:\n (The function is applied to both the _data and the _mask, if any.)"
def __call__(self, *args, **params):
func = getattr(np, self._function)
if len(args)==1:
x = args[0]
if isinstance(x, ndarray):
_d = func(np.asarray(x), **params)
_m = func(getmaskarray(x), **params)
return masked_array(_d, mask=_m)
elif isinstance(x, tuple) or isinstance(x, list):
_d = func(tuple([np.asarray(a) for a in x]), **params)
_m = func(tuple([getmaskarray(a) for a in x]), **params)
return masked_array(_d, mask=_m)
else:
arrays = []
args = list(args)
while len(args)>0 and issequence(args[0]):
arrays.append(args.pop(0))
res = []
for x in arrays:
_d = func(np.asarray(x), *args, **params)
_m = func(getmaskarray(x), *args, **params)
res.append(masked_array(_d, mask=_m))
return res
atleast_1d = _fromnxfunction('atleast_1d')
atleast_2d = _fromnxfunction('atleast_2d')
atleast_3d = _fromnxfunction('atleast_3d')
vstack = row_stack = _fromnxfunction('vstack')
hstack = _fromnxfunction('hstack')
column_stack = _fromnxfunction('column_stack')
dstack = _fromnxfunction('dstack')
hsplit = _fromnxfunction('hsplit')
def expand_dims(a, axis):
"""Expands the shape of a by including newaxis before axis.
"""
if not isinstance(a, MaskedArray):
return np.expand_dims(a, axis)
elif getmask(a) is nomask:
return np.expand_dims(a, axis).view(MaskedArray)
m = getmaskarray(a)
return masked_array(np.expand_dims(a, axis),
mask=np.expand_dims(m, axis))
#####--------------------------------------------------------------------------
#----
#####--------------------------------------------------------------------------
def flatten_inplace(seq):
"""Flatten a sequence in place."""
k = 0
while (k != len(seq)):
while hasattr(seq[k],'__iter__'):
seq[k:(k+1)] = seq[k]
k += 1
return seq
def apply_along_axis(func1d, axis, arr, *args, **kwargs):
"""Execute func1d(arr[i],*args) where func1d takes 1-D arrays and
arr is an N-d array. i varies so as to apply the function along
the given axis for each 1-d subarray in arr.
"""
arr = array(arr, copy=False, subok=True)
nd = arr.ndim
if axis < 0:
axis += nd
if (axis >= nd):
raise ValueError("axis must be less than arr.ndim; axis=%d, rank=%d."
% (axis,nd))
ind = [0]*(nd-1)
i = np.zeros(nd,'O')
indlist = range(nd)
indlist.remove(axis)
i[axis] = slice(None,None)
outshape = np.asarray(arr.shape).take(indlist)
i.put(indlist, ind)
j = i.copy()
res = func1d(arr[tuple(i.tolist())], *args, **kwargs)
# if res is a number, then we have a smaller output array
asscalar = np.isscalar(res)
if not asscalar:
try:
len(res)
except TypeError:
asscalar = True
# Note: we shouldn't set the dtype of the output from the first result...
#...so we force the type to object, and build a list of dtypes
#...we'll just take the largest, to avoid some downcasting
dtypes = []
if asscalar:
dtypes.append(np.asarray(res).dtype)
outarr = zeros(outshape, object)
outarr[tuple(ind)] = res
Ntot = np.product(outshape)
k = 1
while k < Ntot:
# increment the index
ind[-1] += 1
n = -1
while (ind[n] >= outshape[n]) and (n > (1-nd)):
ind[n-1] += 1
ind[n] = 0
n -= 1
i.put(indlist, ind)
res = func1d(arr[tuple(i.tolist())], *args, **kwargs)
outarr[tuple(ind)] = res
dtypes.append(asarray(res).dtype)
k += 1
else:
res = array(res, copy=False, subok=True)
j = i.copy()
j[axis] = ([slice(None, None)] * res.ndim)
j.put(indlist, ind)
Ntot = np.product(outshape)
holdshape = outshape
outshape = list(arr.shape)
outshape[axis] = res.shape
dtypes.append(asarray(res).dtype)
outshape = flatten_inplace(outshape)
outarr = zeros(outshape, object)
outarr[tuple(flatten_inplace(j.tolist()))] = res
k = 1
while k < Ntot:
# increment the index
ind[-1] += 1
n = -1
while (ind[n] >= holdshape[n]) and (n > (1-nd)):
ind[n-1] += 1
ind[n] = 0
n -= 1
i.put(indlist, ind)
j.put(indlist, ind)
res = func1d(arr[tuple(i.tolist())], *args, **kwargs)
outarr[tuple(flatten_inplace(j.tolist()))] = res
dtypes.append(asarray(res).dtype)
k += 1
max_dtypes = np.dtype(np.asarray(dtypes).max())
if not hasattr(arr, '_mask'):
result = np.asarray(outarr, dtype=max_dtypes)
else:
result = asarray(outarr, dtype=max_dtypes)
result.fill_value = ma.default_fill_value(result)
return result
def average(a, axis=None, weights=None, returned=False):
"""Average the array over the given axis.
Parameters
----------
axis : {None,int}, optional
Axis along which to perform the operation.
If None, applies to a flattened version of the array.
weights : {None, sequence}, optional
Sequence of weights.
The weights must have the shape of a, or be 1D with length
the size of a along the given axis.
If no weights are given, weights are assumed to be 1.
returned : {False, True}, optional
Flag indicating whether a tuple (result, sum of weights/counts)
should be returned as output (True), or just the result (False).
"""
a = asarray(a)
mask = a.mask
ash = a.shape
if ash == ():
ash = (1,)
if axis is None:
if mask is nomask:
if weights is None:
n = a.sum(axis=None)
d = float(a.size)
else:
w = filled(weights, 0.0).ravel()
n = umath.add.reduce(a._data.ravel() * w)
d = umath.add.reduce(w)
del w
else:
if weights is None:
n = a.filled(0).sum(axis=None)
d = umath.add.reduce((-mask).ravel().astype(int))
else:
w = array(filled(weights, 0.0), float, mask=mask).ravel()
n = add.reduce(a.ravel() * w)
d = add.reduce(w)
del w
else:
if mask is nomask:
if weights is None:
d = ash[axis] * 1.0
n = add.reduce(a._data, axis, dtype=float)
else:
w = filled(weights, 0.0)
wsh = w.shape
if wsh == ():
wsh = (1,)
if wsh == ash:
w = np.array(w, float, copy=0)
n = add.reduce(a*w, axis)
d = add.reduce(w, axis)
del w
elif wsh == (ash[axis],):
ni = ash[axis]
r = [None]*len(ash)
r[axis] = slice(None, None, 1)
w = eval ("w["+ repr(tuple(r)) + "] * ones(ash, float)")
n = add.reduce(a*w, axis, dtype=float)
d = add.reduce(w, axis, dtype=float)
del w, r
else:
raise ValueError, 'average: weights wrong shape.'
else:
if weights is None:
n = add.reduce(a, axis, dtype=float)
d = umath.add.reduce((-mask), axis=axis, dtype=float)
else:
w = filled(weights, 0.0)
wsh = w.shape
if wsh == ():
wsh = (1,)
if wsh == ash:
w = array(w, dtype=float, mask=mask, copy=0)
n = add.reduce(a*w, axis, dtype=float)
d = add.reduce(w, axis, dtype=float)
elif wsh == (ash[axis],):
ni = ash[axis]
r = [None]*len(ash)
r[axis] = slice(None, None, 1)
w = eval ("w["+ repr(tuple(r)) + \
"] * masked_array(ones(ash, float), mask)")
n = add.reduce(a*w, axis, dtype=float)
d = add.reduce(w, axis, dtype=float)
else:
raise ValueError, 'average: weights wrong shape.'
del w
if n is masked or d is masked:
return masked
result = n/d
del n
if isinstance(result, MaskedArray):
if ((axis is None) or (axis==0 and a.ndim == 1)) and \
(result.mask is nomask):
result = result._data
if returned:
if not isinstance(d, MaskedArray):
d = masked_array(d)
if isinstance(d, ndarray) and (not d.shape == result.shape):
d = ones(result.shape, dtype=float) * d
if returned:
return result, d
else:
return result
def median(a, axis=None, out=None, overwrite_input=False):
"""Compute the median along the specified axis.
Returns the median of the array elements.
Parameters
----------
a : array-like
Input array or object that can be converted to an array
axis : {None, int}, optional
Axis along which the medians are computed. The default (axis=None) is to
compute the median along a flattened version of the array.
out : {None, ndarray}, optional
Alternative output array in which to place the result. It must
have the same shape and buffer length as the expected output
but the type will be cast if necessary.
overwrite_input : {False, True}, optional
If True, then allow use of memory of input array (a) for
calculations. The input array will be modified by the call to
median. This will save memory when you do not need to preserve
the contents of the input array. Treat the input as undefined,
but it will probably be fully or partially sorted. Default is
False. Note that, if overwrite_input is true, and the input
is not already an ndarray, an error will be raised.
Returns
-------
median : ndarray.
A new array holding the result is returned unless out is
specified, in which case a reference to out is returned.
Return datatype is float64 for ints and floats smaller than
float64, or the input datatype otherwise.
See Also
-------
mean
Notes
-----
Given a vector V with N non masked values, the median of V is the middle
value of a sorted copy of V (Vs) - i.e. Vs[(N-1)/2], when N is odd, or
{Vs[N/2 - 1] + Vs[N/2]}/2. when N is even.
"""
def _median1D(data):
counts = filled(count(data, axis),0)
(idx, rmd) = divmod(counts, 2)
if rmd:
choice = slice(idx, idx+1)
else:
choice = slice(idx-1, idx+1)
return data[choice].mean(0)
#
if overwrite_input:
if axis is None:
asorted = a.ravel()
asorted.sort()
else:
a.sort(axis=axis)
asorted = a
else:
asorted = sort(a, axis=axis)
if axis is None:
result = _median1D(asorted)
else:
result = apply_along_axis(_median1D, axis, asorted)
if out is not None:
out = result
return result
#..............................................................................
def compress_rowcols(x, axis=None):
"""Suppress the rows and/or columns of a 2D array that contains
masked values.
The suppression behavior is selected with the `axis`parameter.
- If axis is None, rows and columns are suppressed.
- If axis is 0, only rows are suppressed.
- If axis is 1 or -1, only columns are suppressed.
Parameters
----------
axis : int, optional
Axis along which to perform the operation.
If None, applies to a flattened version of the array.
Returns
-------
compressed_array : an ndarray.
"""
x = asarray(x)
if x.ndim != 2:
raise NotImplementedError, "compress2d works for 2D arrays only."
m = getmask(x)
# Nothing is masked: return x
if m is nomask or not m.any():
return x._data
# All is masked: return empty
if m.all():
return nxarray([])
# Builds a list of rows/columns indices
(idxr, idxc) = (range(len(x)), range(x.shape[1]))
masked = m.nonzero()
if not axis:
for i in np.unique(masked[0]):
idxr.remove(i)
if axis in [None, 1, -1]:
for j in np.unique(masked[1]):
idxc.remove(j)
return x._data[idxr][:,idxc]
def compress_rows(a):
"""Suppress whole rows of a 2D array that contain masked values.
"""
return compress_rowcols(a, 0)
def compress_cols(a):
"""Suppress whole columnss of a 2D array that contain masked values.
"""
return compress_rowcols(a, 1)
def mask_rowcols(a, axis=None):
"""Mask whole rows and/or columns of a 2D array that contain
masked values. The masking behavior is selected with the
`axis`parameter.
- If axis is None, rows and columns are masked.
- If axis is 0, only rows are masked.
- If axis is 1 or -1, only columns are masked.
Parameters
----------
axis : int, optional
Axis along which to perform the operation.
If None, applies to a flattened version of the array.
Returns
-------
a *pure* ndarray.
"""
a = asarray(a)
if a.ndim != 2:
raise NotImplementedError, "compress2d works for 2D arrays only."
m = getmask(a)
# Nothing is masked: return a
if m is nomask or not m.any():
return a
maskedval = m.nonzero()
a._mask = a._mask.copy()
if not axis:
a[np.unique(maskedval[0])] = masked
if axis in [None, 1, -1]:
a[:,np.unique(maskedval[1])] = masked
return a
def mask_rows(a, axis=None):
"""Mask whole rows of a 2D array that contain masked values.
Parameters
----------
axis : int, optional
Axis along which to perform the operation.
If None, applies to a flattened version of the array.
"""
return mask_rowcols(a, 0)
def mask_cols(a, axis=None):
"""Mask whole columns of a 2D array that contain masked values.
Parameters
----------
axis : int, optional
Axis along which to perform the operation.
If None, applies to a flattened version of the array.
"""
return mask_rowcols(a, 1)
def dot(a,b, strict=False):
"""Return the dot product of two 2D masked arrays a and b.
Like the generic numpy equivalent, the product sum is over the last
dimension of a and the second-to-last dimension of b. If strict is True,
masked values are propagated: if a masked value appears in a row or column,
the whole row or column is considered masked.
Parameters
----------
strict : {boolean}
Whether masked data are propagated (True) or set to 0 for the computation.
Notes
-----
The first argument is not conjugated.
"""
#!!!: Works only with 2D arrays. There should be a way to get it to run with higher dimension
if strict and (a.ndim == 2) and (b.ndim == 2):
a = mask_rows(a)
b = mask_cols(b)
#
d = np.dot(filled(a, 0), filled(b, 0))
#
am = (~getmaskarray(a))
bm = (~getmaskarray(b))
m = ~np.dot(am, bm)
return masked_array(d, mask=m)
#...............................................................................
def ediff1d(array, to_end=None, to_begin=None):
"""Return the differences between consecutive elements of an
array, possibly with prefixed and/or appended values.
Parameters
----------
array : {array}
Input array, will be flattened before the difference is taken.
to_end : {number}, optional
If provided, this number will be tacked onto the end of the returned
differences.
to_begin : {number}, optional
If provided, this number will be taked onto the beginning of the
returned differences.
Returns
-------
ed : {array}
The differences. Loosely, this will be (ary[1:] - ary[:-1]).
"""
a = masked_array(array, copy=True)
if a.ndim > 1:
a.reshape((a.size,))
(d, m, n) = (a._data, a._mask, a.size-1)
dd = d[1:]-d[:-1]
if m is nomask:
dm = nomask
else:
dm = m[1:]-m[:-1]
#
if to_end is not None:
to_end = asarray(to_end)
nend = to_end.size
if to_begin is not None:
to_begin = asarray(to_begin)
nbegin = to_begin.size
r_data = np.empty((n+nend+nbegin,), dtype=a.dtype)
r_mask = np.zeros((n+nend+nbegin,), dtype=bool)
r_data[:nbegin] = to_begin._data
r_mask[:nbegin] = to_begin._mask
r_data[nbegin:-nend] = dd
r_mask[nbegin:-nend] = dm
else:
r_data = np.empty((n+nend,), dtype=a.dtype)
r_mask = np.zeros((n+nend,), dtype=bool)
r_data[:-nend] = dd
r_mask[:-nend] = dm
r_data[-nend:] = to_end._data
r_mask[-nend:] = to_end._mask
#
elif to_begin is not None:
to_begin = asarray(to_begin)
nbegin = to_begin.size
r_data = np.empty((n+nbegin,), dtype=a.dtype)
r_mask = np.zeros((n+nbegin,), dtype=bool)
r_data[:nbegin] = to_begin._data
r_mask[:nbegin] = to_begin._mask
r_data[nbegin:] = dd
r_mask[nbegin:] = dm
#
else:
r_data = dd
r_mask = dm
return masked_array(r_data, mask=r_mask)
def cov(x, y=None, rowvar=True, bias=False, allow_masked=True):
"""Estimates the covariance matrix.
Normalization is by (N-1) where N is the number of observations (unbiased
estimate). If bias is True then normalization is by N.
By default, masked values are recognized as such. If x and y have the same
shape, a common mask is allocated: if x[i,j] is masked, then y[i,j] will also
be masked.
Setting `allow_masked` to False will raise an exception if values are missing
in either of the input arrays.
Parameters
----------
x : array-like
Input data.
If x is a 1D array, returns the variance.
If x is a 2D array, returns the covariance matrix.
y : {None, array-like}, optional
Optional set of variables.
rowvar : {False, True} optional
If rowvar is true, then each row is a variable with observations in columns.
If rowvar is False, each column is a variable and the observations are in
the rows.
bias : {False, True} optional
Whether to use a biased (True) or unbiased (False) estimate of the covariance.
If bias is True, then the normalization is by N, the number of observations.
Otherwise, the normalization is by (N-1).
allow_masked : {True, False} optional
If True, masked values are propagated pair-wise: if a value is masked in x,
the corresponding value is masked in y.
If False, raises a ValueError exception when some values are missing.
Raises
------
ValueError:
Raised if some values are missing and allow_masked is False.
"""
x = array(x, ndmin=2, copy=True, dtype=float)
xmask = getmaskarray(x)
# Quick exit if we can't process masked data
if not allow_masked and xmask.any():
raise ValueError("Cannot process masked data...")
#
if x.shape[0] == 1:
rowvar = True
# Make sure that rowvar is either 0 or 1
rowvar = int(bool(rowvar))
axis = 1-rowvar
if rowvar:
tup = (slice(None), None)
else:
tup = (None, slice(None))
#
if y is None:
xnotmask = np.logical_not(xmask).astype(int)
else:
y = array(y, copy=False, ndmin=2, dtype=float)
ymask = getmaskarray(y)
if not allow_masked and ymask.any():
raise ValueError("Cannot process masked data...")
if xmask.any() or ymask.any():
if y.shape == x.shape:
# Define some common mask
common_mask = np.logical_or(xmask, ymask)
if common_mask is not nomask:
x.unshare_mask()
y.unshare_mask()
xmask = x._mask = y._mask = ymask = common_mask
x = concatenate((x,y),axis)
xnotmask = np.logical_not(np.concatenate((xmask, ymask), axis)).astype(int)
x -= x.mean(axis=rowvar)[tup]
if not rowvar:
fact = dot(xnotmask.T, xnotmask)*1. - (1 - bool(bias))
result = (dot(x.T, x.conj(), strict=False) / fact).squeeze()
else:
fact = dot(xnotmask, xnotmask.T)*1. - (1 - bool(bias))
result = (dot(x, x.T.conj(), strict=False) / fact).squeeze()
return result
def corrcoef(x, y=None, rowvar=True, bias=False, allow_masked=True):
"""The correlation coefficients formed from the array x, where the
rows are the observations, and the columns are variables.
corrcoef(x,y) where x and y are 1d arrays is the same as
corrcoef(transpose([x,y]))
Parameters
----------
x : ndarray
Input data. If x is a 1D array, returns the variance.
If x is a 2D array, returns the covariance matrix.
y : {None, ndarray} optional
Optional set of variables.
rowvar : {False, True} optional
If True, then each row is a variable with observations in columns.
If False, each column is a variable and the observations are in the rows.
bias : {False, True} optional
Whether to use a biased (True) or unbiased (False) estimate of the
covariance.
If True, then the normalization is by N, the number of non-missing
observations.
Otherwise, the normalization is by (N-1).
allow_masked : {True, False} optional
If True, masked values are propagated pair-wise: if a value is masked
in x, the corresponding value is masked in y.
If False, raises an exception.
See Also
--------
cov
"""
c = cov(x, y, rowvar, bias, allow_masked)
try:
d = ma.diag(c)
except ValueError: # scalar covariance
return 1
return c/ma.sqrt(ma.multiply.outer(d,d))
#####--------------------------------------------------------------------------
#---- --- Concatenation helpers ---
#####--------------------------------------------------------------------------
class MAxisConcatenator(AxisConcatenator):
"""Translate slice objects to concatenation along an axis.
"""
def __init__(self, axis=0):
AxisConcatenator.__init__(self, axis, matrix=False)
def __getitem__(self,key):
if isinstance(key, str):
raise MAError, "Unavailable for masked array."
if type(key) is not tuple:
key = (key,)
objs = []
scalars = []
final_dtypedescr = None
for k in range(len(key)):
scalar = False
if type(key[k]) is slice:
step = key[k].step
start = key[k].start
stop = key[k].stop
if start is None:
start = 0
if step is None:
step = 1
if type(step) is type(1j):
size = int(abs(step))
newobj = np.linspace(start, stop, num=size)
else:
newobj = np.arange(start, stop, step)
elif type(key[k]) is str:
if (key[k] in 'rc'):
self.matrix = True
self.col = (key[k] == 'c')
continue
try:
self.axis = int(key[k])
continue
except (ValueError, TypeError):
raise ValueError, "Unknown special directive"
elif type(key[k]) in np.ScalarType:
newobj = asarray([key[k]])
scalars.append(k)
scalar = True
else:
newobj = key[k]
objs.append(newobj)
if isinstance(newobj, ndarray) and not scalar:
if final_dtypedescr is None:
final_dtypedescr = newobj.dtype
elif newobj.dtype > final_dtypedescr:
final_dtypedescr = newobj.dtype
if final_dtypedescr is not None:
for k in scalars:
objs[k] = objs[k].astype(final_dtypedescr)
res = concatenate(tuple(objs),axis=self.axis)
return self._retval(res)
class mr_class(MAxisConcatenator):
"""Translate slice objects to concatenation along the first axis.
For example:
>>> np.ma.mr_[np.ma.array([1,2,3]), 0, 0, np.ma.array([4,5,6])]
array([1, 2, 3, 0, 0, 4, 5, 6])
"""
def __init__(self):
MAxisConcatenator.__init__(self, 0)
mr_ = mr_class()
#####--------------------------------------------------------------------------
#---- Find unmasked data ---
#####--------------------------------------------------------------------------
def flatnotmasked_edges(a):
"""Find the indices of the first and last not masked values in a
1D masked array. If all values are masked, returns None.
"""
m = getmask(a)
if m is nomask or not np.any(m):
return [0,-1]
unmasked = np.flatnonzero(~m)
if len(unmasked) > 0:
return unmasked[[0,-1]]
else:
return None
def notmasked_edges(a, axis=None):
"""Find the indices of the first and last not masked values along
the given axis in a masked array.
If all values are masked, return None. Otherwise, return a list
of 2 tuples, corresponding to the indices of the first and last
unmasked values respectively.
Parameters
----------
axis : int, optional
Axis along which to perform the operation.
If None, applies to a flattened version of the array.
"""
a = asarray(a)
if axis is None or a.ndim == 1:
return flatnotmasked_edges(a)
m = getmask(a)
idx = array(np.indices(a.shape), mask=np.asarray([m]*a.ndim))
return [tuple([idx[i].min(axis).compressed() for i in range(a.ndim)]),
tuple([idx[i].max(axis).compressed() for i in range(a.ndim)]),]
def flatnotmasked_contiguous(a):
"""Find contiguous unmasked data in a flattened masked array.
Return a sorted sequence of slices (start index, end index).
"""
m = getmask(a)
if m is nomask:
return (a.size, [0,-1])
unmasked = np.flatnonzero(~m)
if len(unmasked) == 0:
return None
result = []
for k, group in groupby(enumerate(unmasked), lambda (i,x):i-x):
tmp = np.array([g[1] for g in group], int)
# result.append((tmp.size, tuple(tmp[[0,-1]])))
result.append( slice(tmp[0], tmp[-1]) )
result.sort()
return result
def notmasked_contiguous(a, axis=None):
"""Find contiguous unmasked data in a masked array along the given
axis.
Parameters
----------
axis : int, optional
Axis along which to perform the operation.
If None, applies to a flattened version of the array.
Returns
-------
A sorted sequence of slices (start index, end index).
Notes
-----
Only accepts 2D arrays at most.
"""
a = asarray(a)
nd = a.ndim
if nd > 2:
raise NotImplementedError,"Currently limited to atmost 2D array."
if axis is None or nd == 1:
return flatnotmasked_contiguous(a)
#
result = []
#
other = (axis+1)%2
idx = [0,0]
idx[axis] = slice(None, None)
#
for i in range(a.shape[other]):
idx[other] = i
result.append( flatnotmasked_contiguous(a[idx]) )
return result
#####--------------------------------------------------------------------------
#---- Polynomial fit ---
#####--------------------------------------------------------------------------
def vander(x, n=None):
"""%s
Notes
-----
Masked values in x will result in rows of zeros.
"""
_vander = np.vander(x, n)
m = getmask(x)
if m is not nomask:
_vander[m] = 0
return _vander
def polyfit(x, y, deg, rcond=None, full=False):
"""
Least squares polynomial fit.
Do a best fit polynomial of degree 'deg' of 'x' to 'y'. Return value is a
vector of polynomial coefficients [pk ... p1 p0]. Eg, for ``deg = 2``::
p2*x0^2 + p1*x0 + p0 = y1
p2*x1^2 + p1*x1 + p0 = y1
p2*x2^2 + p1*x2 + p0 = y2
.....
p2*xk^2 + p1*xk + p0 = yk
Parameters
----------
x : array_like
1D vector of sample points.
y : array_like
1D vector or 2D array of values to fit. The values should run down the
columns in the 2D case.
deg : integer
Degree of the fitting polynomial
rcond: {None, float}, optional
Relative condition number of the fit. Singular values smaller than this
relative to the largest singular value will be ignored. The defaul value
is len(x)*eps, where eps is the relative precision of the float type,
about 2e-16 in most cases.
full : {False, boolean}, optional
Switch determining nature of return value. When it is False just the
coefficients are returned, when True diagnostic information from the
singular value decomposition is also returned.
Returns
-------
coefficients, [residuals, rank, singular_values, rcond] : variable
When full=False, only the coefficients are returned, running down the
appropriate colume when y is a 2D array. When full=True, the rank of the
scaled Vandermonde matrix, its effective rank in light of the rcond
value, its singular values, and the specified value of rcond are also
returned.
Warns
-----
RankWarning : if rank is reduced and not full output
The warnings can be turned off by:
>>> import warnings
>>> warnings.simplefilter('ignore',np.RankWarning)
See Also
--------
polyval : computes polynomial values.
Notes
-----
If X is a the Vandermonde Matrix computed from x (see
http://mathworld.wolfram.com/VandermondeMatrix.html), then the
polynomial least squares solution is given by the 'p' in
X*p = y
where X.shape is a matrix of dimensions (len(x), deg + 1), p is a vector of
dimensions (deg + 1, 1), and y is a vector of dimensions (len(x), 1).
This equation can be solved as
p = (XT*X)^-1 * XT * y
where XT is the transpose of X and -1 denotes the inverse. However, this
method is susceptible to rounding errors and generally the singular value
decomposition of the matrix X is preferred and that is what is done here.
The singular value method takes a paramenter, 'rcond', which sets a limit on
the relative size of the smallest singular value to be used in solving the
equation. This may result in lowering the rank of the Vandermonde matrix, in
which case a RankWarning is issued. If polyfit issues a RankWarning, try a
fit of lower degree or replace x by x - x.mean(), both of which will
generally improve the condition number. The routine already normalizes the
vector x by its maximum absolute value to help in this regard. The rcond
parameter can be set to a value smaller than its default, but the resulting
fit may be spurious. The current default value of rcond is len(x)*eps, where
eps is the relative precision of the floating type being used, generally
around 1e-7 and 2e-16 for IEEE single and double precision respectively.
This value of rcond is fairly conservative but works pretty well when x -
x.mean() is used in place of x.
DISCLAIMER: Power series fits are full of pitfalls for the unwary once the
degree of the fit becomes large or the interval of sample points is badly
centered. The problem is that the powers x**n are generally a poor basis for
the polynomial functions on the sample interval, resulting in a Vandermonde
matrix is ill conditioned and coefficients sensitive to rounding erros. The
computation of the polynomial values will also sensitive to rounding errors.
Consequently, the quality of the polynomial fit should be checked against
the data whenever the condition number is large. The quality of polynomial
fits *can not* be taken for granted. If all you want to do is draw a smooth
curve through the y values and polyfit is not doing the job, try centering
the sample range or look into scipy.interpolate, which includes some nice
spline fitting functions that may be of use.
For more info, see
http://mathworld.wolfram.com/LeastSquaresFittingPolynomial.html,
but note that the k's and n's in the superscripts and subscripts
on that page. The linear algebra is correct, however.
Notes
-----
Any masked values in x is propagated in y, and vice-versa.
"""
order = int(deg) + 1
x = asarray(x)
mx = getmask(x)
y = asarray(y)
if y.ndim == 1:
m = mask_or(mx, getmask(y))
elif y.ndim == 2:
y = mask_rows(y)
my = getmask(y)
if my is not nomask:
m = mask_or(mx, my[:,0])
else:
m = mx
else:
raise TypeError,"Expected a 1D or 2D array for y!"
if m is not nomask:
x[m] = y[m] = masked
# Set rcond
if rcond is None :
rcond = len(x)*np.finfo(x.dtype).eps
# Scale x to improve condition number
scale = abs(x).max()
if scale != 0 :
x = x / scale
# solve least squares equation for powers of x
v = vander(x, order)
c, resids, rank, s = _lstsq(v, y.filled(0), rcond)
# warn on rank reduction, which indicates an ill conditioned matrix
if rank != order and not full:
warnings.warn("Polyfit may be poorly conditioned", np.RankWarning)
# scale returned coefficients
if scale != 0 :
if c.ndim == 1 :
c /= np.vander([scale], order)[0]
else :
c /= np.vander([scale], order).T
if full :
return c, resids, rank, s, rcond
else :
return c
################################################################################
|