domino.plsr
1import numpy as np 2import xarray as xr 3import sklearn.cross_decomposition as skcc 4 5class PLSR_Reduction(object): 6 """ 7 Wraps around the scikit-learn partial-least-squares-regression algorithm, supporting prediction-focused dimensionality reduction. 8 9 **Arguments** 10 11 *mode_num* 12 13 An integer number of PLSR modes to retain. Defaults to the second dimension of the predictor variables *X* passed to *PLSR_Reduction.fit* 14 15 """ 16 def __init__(self,mode_num): 17 """Initialise a PLSR_Reduction object""" 18 self.mode_num=mode_num 19 """@private""" 20 return 21 22 def __repr__(self): 23 return 'A PLSR_Reduction object' 24 25 def __str__(self): 26 return self.__repr__() 27 28 def _validate_shape(self,da): 29 30 #Expand if necessary 31 if da.ndim==1: 32 da=da.expand_dims('feature') 33 34 assert da.ndim==2 35 assert self.sample_dim in da.dims 36 37 d0,d1=da.dims 38 #Transpose if necessary 39 if d0!=self.sample_dim: 40 da=da.T 41 return da 42 43 44 def fit(self,X,y,sample_dim='time'): 45 """ 46 Performs the partial-least-squares regression of *X* against *y*, returning the transformed PLSR modes. 47 48 **Arguments** 49 50 *X* 51 52 The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by *sample_dim*. 53 54 *y* 55 56 The target variables to be predicted. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by *sample_dim*. 57 58 **Optional arguments** 59 60 *sample_dim* 61 62 A string specifying the sample dimension which must be shared between *X* and *y*. Defaults to *time*. 63 64 **Returns** 65 66 *plsr_X* 67 68 A DataArray containing PLSR modes, with a sample coordinate given by the intersection of the sample coordinates of *X* and *y* 69 70 """ 71 X,y=xr.align(X.dropna(sample_dim),y.dropna(sample_dim)) 72 assert len(y>1) 73 self.sample_dim=sample_dim 74 self.X=self._validate_shape(X) 75 self.y=self._validate_shape(y) 76 77 d0,d1=self.X.dims 78 assert d0==self.sample_dim 79 self.feature_dim=d1 80 81 model=skcc.PLSRegression(n_components=self.mode_num) 82 model.fit(self.X.values,self.y.values) 83 rotations=model.x_rotations_ 84 coords={d1:self.X[d1].values,'PLSR_mode':range(1,self.mode_num+1)} 85 rotation_da=xr.DataArray(data=rotations,coords=coords,dims=[d1,'PLSR_mode']) 86 87 plsr_X=self.X.values@rotations 88 plsr_X=xr.DataArray(data=plsr_X,coords={sample_dim:self.X[sample_dim],'PLSR_mode':range(1,self.mode_num+1)}) 89 90 self.rotations=rotations 91 self.rotation_da=rotation_da 92 self.PLSR_X=plsr_X 93 94 return plsr_X 95 96 def project_quasimodes(self,X): 97 98 """ 99 Use precomputed rotation matrix from calling *PLSR_Reduction.fit* to project *X* into the space of PLSR modes. 100 This can be used to extend the PLSR modes to times when the target events are not defined. 101 102 **Arguments** 103 104 105 *X* 106 107 The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by the *sample_dim* used to fit the PLSR regression. 108 109 **Returns** 110 111 *proj_plsr_X* 112 113 A DataArray containing projected PLSR modes, with the same sample coordinate as *X*. 114 115 """ 116 117 #PLSR_projection 118 119 X=self._validate_shape(X) 120 121 122 d0,d1=X.dims 123 proj_plsr_X=X.values@self.rotations 124 proj_plsr_X=xr.DataArray(data=proj_plsr_X,coords={self.sample_dim:X[self.sample_dim],'PLSR_mode':range(1,self.mode_num+1)}) 125 return proj_plsr_X 126 127 def project_pattern(self,X): 128 """ 129 Compute the patterns corresponding to each PLSR mode, by computing a weighted sum of the spatial patterns corresponding to the predictor indices used to fit the partial least squares regression. 130 131 **Arguments** 132 133 134 *X* 135 136 The pattern to transform. An xarray.Datarray of arbitrary shape, and with an identical feature coordinate to the predictor variables passed into *PLSR_Reduction.fit*. 137 138 **Returns** 139 140 *pattern* 141 142 A DataArray containing the resulting weighted sums, with the same coordinates as *X*, except for the feature coordinate is replaced with the PLSR_mode coordinate. 143 """ 144 assert self.feature_dim in X.dims 145 pattern=(self.rotation_da*X).sum(self.feature_dim) 146 return pattern 147 148
6class PLSR_Reduction(object): 7 """ 8 Wraps around the scikit-learn partial-least-squares-regression algorithm, supporting prediction-focused dimensionality reduction. 9 10 **Arguments** 11 12 *mode_num* 13 14 An integer number of PLSR modes to retain. Defaults to the second dimension of the predictor variables *X* passed to *PLSR_Reduction.fit* 15 16 """ 17 def __init__(self,mode_num): 18 """Initialise a PLSR_Reduction object""" 19 self.mode_num=mode_num 20 """@private""" 21 return 22 23 def __repr__(self): 24 return 'A PLSR_Reduction object' 25 26 def __str__(self): 27 return self.__repr__() 28 29 def _validate_shape(self,da): 30 31 #Expand if necessary 32 if da.ndim==1: 33 da=da.expand_dims('feature') 34 35 assert da.ndim==2 36 assert self.sample_dim in da.dims 37 38 d0,d1=da.dims 39 #Transpose if necessary 40 if d0!=self.sample_dim: 41 da=da.T 42 return da 43 44 45 def fit(self,X,y,sample_dim='time'): 46 """ 47 Performs the partial-least-squares regression of *X* against *y*, returning the transformed PLSR modes. 48 49 **Arguments** 50 51 *X* 52 53 The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by *sample_dim*. 54 55 *y* 56 57 The target variables to be predicted. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by *sample_dim*. 58 59 **Optional arguments** 60 61 *sample_dim* 62 63 A string specifying the sample dimension which must be shared between *X* and *y*. Defaults to *time*. 64 65 **Returns** 66 67 *plsr_X* 68 69 A DataArray containing PLSR modes, with a sample coordinate given by the intersection of the sample coordinates of *X* and *y* 70 71 """ 72 X,y=xr.align(X.dropna(sample_dim),y.dropna(sample_dim)) 73 assert len(y>1) 74 self.sample_dim=sample_dim 75 self.X=self._validate_shape(X) 76 self.y=self._validate_shape(y) 77 78 d0,d1=self.X.dims 79 assert d0==self.sample_dim 80 self.feature_dim=d1 81 82 model=skcc.PLSRegression(n_components=self.mode_num) 83 model.fit(self.X.values,self.y.values) 84 rotations=model.x_rotations_ 85 coords={d1:self.X[d1].values,'PLSR_mode':range(1,self.mode_num+1)} 86 rotation_da=xr.DataArray(data=rotations,coords=coords,dims=[d1,'PLSR_mode']) 87 88 plsr_X=self.X.values@rotations 89 plsr_X=xr.DataArray(data=plsr_X,coords={sample_dim:self.X[sample_dim],'PLSR_mode':range(1,self.mode_num+1)}) 90 91 self.rotations=rotations 92 self.rotation_da=rotation_da 93 self.PLSR_X=plsr_X 94 95 return plsr_X 96 97 def project_quasimodes(self,X): 98 99 """ 100 Use precomputed rotation matrix from calling *PLSR_Reduction.fit* to project *X* into the space of PLSR modes. 101 This can be used to extend the PLSR modes to times when the target events are not defined. 102 103 **Arguments** 104 105 106 *X* 107 108 The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by the *sample_dim* used to fit the PLSR regression. 109 110 **Returns** 111 112 *proj_plsr_X* 113 114 A DataArray containing projected PLSR modes, with the same sample coordinate as *X*. 115 116 """ 117 118 #PLSR_projection 119 120 X=self._validate_shape(X) 121 122 123 d0,d1=X.dims 124 proj_plsr_X=X.values@self.rotations 125 proj_plsr_X=xr.DataArray(data=proj_plsr_X,coords={self.sample_dim:X[self.sample_dim],'PLSR_mode':range(1,self.mode_num+1)}) 126 return proj_plsr_X 127 128 def project_pattern(self,X): 129 """ 130 Compute the patterns corresponding to each PLSR mode, by computing a weighted sum of the spatial patterns corresponding to the predictor indices used to fit the partial least squares regression. 131 132 **Arguments** 133 134 135 *X* 136 137 The pattern to transform. An xarray.Datarray of arbitrary shape, and with an identical feature coordinate to the predictor variables passed into *PLSR_Reduction.fit*. 138 139 **Returns** 140 141 *pattern* 142 143 A DataArray containing the resulting weighted sums, with the same coordinates as *X*, except for the feature coordinate is replaced with the PLSR_mode coordinate. 144 """ 145 assert self.feature_dim in X.dims 146 pattern=(self.rotation_da*X).sum(self.feature_dim) 147 return pattern
Wraps around the scikit-learn partial-least-squares-regression algorithm, supporting prediction-focused dimensionality reduction.
Arguments
mode_num
An integer number of PLSR modes to retain. Defaults to the second dimension of the predictor variables X passed to PLSR_Reduction.fit
17 def __init__(self,mode_num): 18 """Initialise a PLSR_Reduction object""" 19 self.mode_num=mode_num 20 """@private""" 21 return
Initialise a PLSR_Reduction object
45 def fit(self,X,y,sample_dim='time'): 46 """ 47 Performs the partial-least-squares regression of *X* against *y*, returning the transformed PLSR modes. 48 49 **Arguments** 50 51 *X* 52 53 The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by *sample_dim*. 54 55 *y* 56 57 The target variables to be predicted. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by *sample_dim*. 58 59 **Optional arguments** 60 61 *sample_dim* 62 63 A string specifying the sample dimension which must be shared between *X* and *y*. Defaults to *time*. 64 65 **Returns** 66 67 *plsr_X* 68 69 A DataArray containing PLSR modes, with a sample coordinate given by the intersection of the sample coordinates of *X* and *y* 70 71 """ 72 X,y=xr.align(X.dropna(sample_dim),y.dropna(sample_dim)) 73 assert len(y>1) 74 self.sample_dim=sample_dim 75 self.X=self._validate_shape(X) 76 self.y=self._validate_shape(y) 77 78 d0,d1=self.X.dims 79 assert d0==self.sample_dim 80 self.feature_dim=d1 81 82 model=skcc.PLSRegression(n_components=self.mode_num) 83 model.fit(self.X.values,self.y.values) 84 rotations=model.x_rotations_ 85 coords={d1:self.X[d1].values,'PLSR_mode':range(1,self.mode_num+1)} 86 rotation_da=xr.DataArray(data=rotations,coords=coords,dims=[d1,'PLSR_mode']) 87 88 plsr_X=self.X.values@rotations 89 plsr_X=xr.DataArray(data=plsr_X,coords={sample_dim:self.X[sample_dim],'PLSR_mode':range(1,self.mode_num+1)}) 90 91 self.rotations=rotations 92 self.rotation_da=rotation_da 93 self.PLSR_X=plsr_X 94 95 return plsr_X
Performs the partial-least-squares regression of X against y, returning the transformed PLSR modes.
Arguments
X
The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by sample_dim.
y
The target variables to be predicted. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by sample_dim.
Optional arguments
sample_dim
A string specifying the sample dimension which must be shared between X and y. Defaults to time.
Returns
plsr_X
A DataArray containing PLSR modes, with a sample coordinate given by the intersection of the sample coordinates of X and y
97 def project_quasimodes(self,X): 98 99 """ 100 Use precomputed rotation matrix from calling *PLSR_Reduction.fit* to project *X* into the space of PLSR modes. 101 This can be used to extend the PLSR modes to times when the target events are not defined. 102 103 **Arguments** 104 105 106 *X* 107 108 The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by the *sample_dim* used to fit the PLSR regression. 109 110 **Returns** 111 112 *proj_plsr_X* 113 114 A DataArray containing projected PLSR modes, with the same sample coordinate as *X*. 115 116 """ 117 118 #PLSR_projection 119 120 X=self._validate_shape(X) 121 122 123 d0,d1=X.dims 124 proj_plsr_X=X.values@self.rotations 125 proj_plsr_X=xr.DataArray(data=proj_plsr_X,coords={self.sample_dim:X[self.sample_dim],'PLSR_mode':range(1,self.mode_num+1)}) 126 return proj_plsr_X
Use precomputed rotation matrix from calling PLSR_Reduction.fit to project X into the space of PLSR modes. This can be used to extend the PLSR modes to times when the target events are not defined.
Arguments
X
The predictor variables to transform. An xarray.Datarray with two coordinates representing the sample and feature dimensions as specified by the sample_dim used to fit the PLSR regression.
Returns
proj_plsr_X
A DataArray containing projected PLSR modes, with the same sample coordinate as X.
128 def project_pattern(self,X): 129 """ 130 Compute the patterns corresponding to each PLSR mode, by computing a weighted sum of the spatial patterns corresponding to the predictor indices used to fit the partial least squares regression. 131 132 **Arguments** 133 134 135 *X* 136 137 The pattern to transform. An xarray.Datarray of arbitrary shape, and with an identical feature coordinate to the predictor variables passed into *PLSR_Reduction.fit*. 138 139 **Returns** 140 141 *pattern* 142 143 A DataArray containing the resulting weighted sums, with the same coordinates as *X*, except for the feature coordinate is replaced with the PLSR_mode coordinate. 144 """ 145 assert self.feature_dim in X.dims 146 pattern=(self.rotation_da*X).sum(self.feature_dim) 147 return pattern
Compute the patterns corresponding to each PLSR mode, by computing a weighted sum of the spatial patterns corresponding to the predictor indices used to fit the partial least squares regression.
Arguments
X
The pattern to transform. An xarray.Datarray of arbitrary shape, and with an identical feature coordinate to the predictor variables passed into PLSR_Reduction.fit.
Returns
pattern
A DataArray containing the resulting weighted sums, with the same coordinates as X, except for the feature coordinate is replaced with the PLSR_mode coordinate.