Salome HOME
c04fd1d4af14548eacd8d9cca931138f019c402d
[modules/adao.git] / src / daComposant / daAlgorithms / UnscentedKalmanFilter.py
1 # -*- coding: utf-8 -*-
2 #
3 # Copyright (C) 2008-2017 EDF R&D
4 #
5 # This library is free software; you can redistribute it and/or
6 # modify it under the terms of the GNU Lesser General Public
7 # License as published by the Free Software Foundation; either
8 # version 2.1 of the License.
9 #
10 # This library is distributed in the hope that it will be useful,
11 # but WITHOUT ANY WARRANTY; without even the implied warranty of
12 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
13 # Lesser General Public License for more details.
14 #
15 # You should have received a copy of the GNU Lesser General Public
16 # License along with this library; if not, write to the Free Software
17 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
18 #
19 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 #
21 # Author: Jean-Philippe Argaud, jean-philippe.argaud@edf.fr, EDF R&D
22
23 import logging
24 from daCore import BasicObjects
25 import numpy, math
26
27 # ==============================================================================
28 class ElementaryAlgorithm(BasicObjects.Algorithm):
29     def __init__(self):
30         BasicObjects.Algorithm.__init__(self, "UNSCENTEDKALMANFILTER")
31         self.defineRequiredParameter(
32             name     = "ConstrainedBy",
33             default  = "EstimateProjection",
34             typecast = str,
35             message  = "Prise en compte des contraintes",
36             listval  = ["EstimateProjection"],
37             )
38         self.defineRequiredParameter(
39             name     = "EstimationOf",
40             default  = "State",
41             typecast = str,
42             message  = "Estimation d'etat ou de parametres",
43             listval  = ["State", "Parameters"],
44             )
45         self.defineRequiredParameter(
46             name     = "Alpha",
47             default  = 1.,
48             typecast = float,
49             message  = "",
50             minval   = 1.e-4,
51             maxval   = 1.,
52             )
53         self.defineRequiredParameter(
54             name     = "Beta",
55             default  = 2,
56             typecast = float,
57             message  = "",
58             )
59         self.defineRequiredParameter(
60             name     = "Kappa",
61             default  = 0,
62             typecast = int,
63             message  = "",
64             maxval   = 2,
65             )
66         self.defineRequiredParameter(
67             name     = "Reconditioner",
68             default  = 1.,
69             typecast = float,
70             message  = "",
71             minval   = 1.e-3,
72             maxval   = 1.e+1,
73             )
74         self.defineRequiredParameter(
75             name     = "StoreInternalVariables",
76             default  = False,
77             typecast = bool,
78             message  = "Stockage des variables internes ou intermédiaires du calcul",
79             )
80         self.defineRequiredParameter(
81             name     = "StoreSupplementaryCalculations",
82             default  = [],
83             typecast = tuple,
84             message  = "Liste de calculs supplémentaires à stocker et/ou effectuer",
85             listval  = ["APosterioriCorrelations", "APosterioriCovariance", "APosterioriStandardDeviations", "APosterioriVariances", "BMA", "CurrentState", "CostFunctionJ", "CostFunctionJb", "CostFunctionJo", "Innovation"]
86             )
87         self.defineRequiredParameter( # Pas de type
88             name     = "Bounds",
89             message  = "Liste des valeurs de bornes",
90             )
91
92     def run(self, Xb=None, Y=None, U=None, HO=None, EM=None, CM=None, R=None, B=None, Q=None, Parameters=None):
93         self._pre_run(Parameters)
94         #
95         if self._parameters["EstimationOf"] == "Parameters":
96             self._parameters["StoreInternalVariables"] = True
97         #
98         L     = Xb.size
99         Alpha = self._parameters["Alpha"]
100         Beta  = self._parameters["Beta"]
101         if self._parameters["Kappa"] == 0:
102             if self._parameters["EstimationOf"] == "State":
103                 Kappa = 0
104             elif self._parameters["EstimationOf"] == "Parameters":
105                 Kappa = 3 - L
106         else:
107             Kappa = self._parameters["Kappa"]
108         Lambda = float( Alpha**2 ) * ( L + Kappa ) - L
109         Gamma  = math.sqrt( L + Lambda )
110         #
111         Ww = []
112         Ww.append( 0. )
113         for i in range(2*L):
114             Ww.append( 1. / (2.*(L + Lambda)) )
115         #
116         Wm = numpy.array( Ww )
117         Wm[0] = Lambda / (L + Lambda)
118         Wc = numpy.array( Ww )
119         Wc[0] = Lambda / (L + Lambda) + (1. - Alpha**2 + Beta)
120         #
121         # Opérateurs
122         # ----------
123         if B is None:
124             raise ValueError("Background error covariance matrix has to be properly defined!")
125         if R is None:
126             raise ValueError("Observation error covariance matrix has to be properly defined!")
127         #
128         H = HO["Direct"].appliedControledFormTo
129         #
130         if self._parameters["EstimationOf"] == "State":
131             M = EM["Direct"].appliedControledFormTo
132         #
133         if CM is not None and "Tangent" in CM and U is not None:
134             Cm = CM["Tangent"].asMatrix(Xb)
135         else:
136             Cm = None
137         #
138         # Nombre de pas identique au nombre de pas d'observations
139         # -------------------------------------------------------
140         if hasattr(Y,"stepnumber"):
141             duration = Y.stepnumber()
142         else:
143             duration = 2
144         #
145         # Précalcul des inversions de B et R
146         # ----------------------------------
147         if self._parameters["StoreInternalVariables"]:
148             BI = B.getI()
149             RI = R.getI()
150         #
151         # Initialisation
152         # --------------
153         Xn = Xb
154         if hasattr(B,"asfullmatrix"):
155             Pn = B.asfullmatrix(Xn.size)
156         else:
157             Pn = B
158         #
159         self.StoredVariables["Analysis"].store( Xn.A1 )
160         if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
161             self.StoredVariables["APosterioriCovariance"].store( Pn )
162             covarianceXa = Pn
163         Xa               = Xn
164         previousJMinimum = numpy.finfo(float).max
165         #
166         for step in range(duration-1):
167             if hasattr(Y,"store"):
168                 Ynpu = numpy.asmatrix(numpy.ravel( Y[step+1] )).T
169             else:
170                 Ynpu = numpy.asmatrix(numpy.ravel( Y )).T
171             #
172             if U is not None:
173                 if hasattr(U,"store") and len(U)>1:
174                     Un = numpy.asmatrix(numpy.ravel( U[step] )).T
175                 elif hasattr(U,"store") and len(U)==1:
176                     Un = numpy.asmatrix(numpy.ravel( U[0] )).T
177                 else:
178                     Un = numpy.asmatrix(numpy.ravel( U )).T
179             else:
180                 Un = None
181             #
182             Pndemi = numpy.linalg.cholesky(Pn)
183             Xnp = numpy.hstack([Xn, Xn+Gamma*Pndemi, Xn-Gamma*Pndemi])
184             nbSpts = 2*Xn.size+1
185             #
186             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
187                 for point in range(nbSpts):
188                     Xnp[:,point] = numpy.max(numpy.hstack((Xnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
189                     Xnp[:,point] = numpy.min(numpy.hstack((Xnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
190             #
191             XEtnnp = []
192             for point in range(nbSpts):
193                 if self._parameters["EstimationOf"] == "State":
194                     XEtnnpi = numpy.asmatrix(numpy.ravel( M( (Xnp[:,point], Un) ) )).T
195                     if Cm is not None and Un is not None: # Attention : si Cm est aussi dans M, doublon !
196                         Cm = Cm.reshape(Xn.size,Un.size) # ADAO & check shape
197                         XEtnnpi = XEtnnpi + Cm * Un
198                     if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
199                         XEtnnpi = numpy.max(numpy.hstack((XEtnnpi,numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
200                         XEtnnpi = numpy.min(numpy.hstack((XEtnnpi,numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
201                 elif self._parameters["EstimationOf"] == "Parameters":
202                     # --- > Par principe, M = Id, Q = 0
203                     XEtnnpi = Xnp[:,point]
204                 XEtnnp.append( XEtnnpi )
205             XEtnnp = numpy.hstack( XEtnnp )
206             #
207             Xncm = numpy.matrix( XEtnnp.getA()*numpy.array(Wm) ).sum(axis=1)
208             #
209             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
210                 Xncm = numpy.max(numpy.hstack((Xncm,numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
211                 Xncm = numpy.min(numpy.hstack((Xncm,numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
212             #
213             if self._parameters["EstimationOf"] == "State":        Pnm = Q
214             elif self._parameters["EstimationOf"] == "Parameters": Pnm = 0.
215             for point in range(nbSpts):
216                 Pnm += Wc[i] * (XEtnnp[:,point]-Xncm) * (XEtnnp[:,point]-Xncm).T
217             #
218             if self._parameters["EstimationOf"] == "Parameters" and self._parameters["Bounds"] is not None:
219                 Pnmdemi = self._parameters["Reconditioner"] * numpy.linalg.cholesky(Pnm)
220             else:
221                 Pnmdemi = numpy.linalg.cholesky(Pnm)
222             #
223             Xnnp = numpy.hstack([Xncm, Xncm+Gamma*Pnmdemi, Xncm-Gamma*Pnmdemi])
224             #
225             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
226                 for point in range(nbSpts):
227                     Xnnp[:,point] = numpy.max(numpy.hstack((Xnnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
228                     Xnnp[:,point] = numpy.min(numpy.hstack((Xnnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
229             #
230             Ynnp = []
231             for point in range(nbSpts):
232                 if self._parameters["EstimationOf"] == "State":
233                     Ynnpi = numpy.asmatrix(numpy.ravel( H( (Xnnp[:,point], None) ) )).T
234                 elif self._parameters["EstimationOf"] == "Parameters":
235                     Ynnpi = numpy.asmatrix(numpy.ravel( H( (Xnnp[:,point], Un) ) )).T
236                 Ynnp.append( Ynnpi )
237             Ynnp = numpy.hstack( Ynnp )
238             #
239             Yncm = numpy.matrix( Ynnp.getA()*numpy.array(Wm) ).sum(axis=1)
240             #
241             Pyyn = R
242             Pxyn = 0.
243             for point in range(nbSpts):
244                 Pyyn += Wc[i] * (Ynnp[:,point]-Yncm) * (Ynnp[:,point]-Yncm).T
245                 Pxyn += Wc[i] * (Xnnp[:,point]-Xncm) * (Ynnp[:,point]-Yncm).T
246             #
247             d  = Ynpu - Yncm
248             if self._parameters["EstimationOf"] == "Parameters":
249                 if Cm is not None and Un is not None: # Attention : si Cm est aussi dans H, doublon !
250                     d = d - Cm * Un
251             #
252             Kn = Pxyn * Pyyn.I
253             Xn = Xncm + Kn * d
254             Pn = Pnm - Kn * Pyyn * Kn.T
255             #
256             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
257                 Xn = numpy.max(numpy.hstack((Xn,numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
258                 Xn = numpy.min(numpy.hstack((Xn,numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
259             #
260             self.StoredVariables["Analysis"].store( Xn.A1 )
261             if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
262                 self.StoredVariables["APosterioriCovariance"].store( Pn )
263             if "Innovation" in self._parameters["StoreSupplementaryCalculations"]:
264                 self.StoredVariables["Innovation"].store( numpy.ravel( d.A1 ) )
265             if self._parameters["StoreInternalVariables"]:
266                 Jb  = 0.5 * (Xn - Xb).T * BI * (Xn - Xb)
267                 Jo  = 0.5 * d.T * RI * d
268                 J   = float( Jb ) + float( Jo )
269                 if self._parameters["StoreInternalVariables"] or "CurrentState" in self._parameters["StoreSupplementaryCalculations"]:
270                     self.StoredVariables["CurrentState"].store( Xn )
271                 self.StoredVariables["CostFunctionJb"].store( Jb )
272                 self.StoredVariables["CostFunctionJo"].store( Jo )
273                 self.StoredVariables["CostFunctionJ" ].store( J )
274                 if J < previousJMinimum:
275                     previousJMinimum  = J
276                     Xa                = Xn
277                     if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
278                         covarianceXa  = Pn
279             else:
280                 Xa = Xn
281             #
282         #
283         # Stockage supplementaire de l'optimum en estimation de parametres
284         # ----------------------------------------------------------------
285         if self._parameters["EstimationOf"] == "Parameters":
286             self.StoredVariables["Analysis"].store( Xa.A1 )
287             if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
288                 self.StoredVariables["APosterioriCovariance"].store( covarianceXa )
289         #
290         if "BMA" in self._parameters["StoreSupplementaryCalculations"]:
291             self.StoredVariables["BMA"].store( numpy.ravel(Xb) - numpy.ravel(Xa) )
292         #
293         self._post_run(HO)
294         return 0
295
296 # ==============================================================================
297 if __name__ == "__main__":
298     print('\n AUTODIAGNOSTIC \n')