Salome HOME
Improving warning message for available DFO minimizer choice
[modules/adao.git] / src / daComposant / daAlgorithms / UnscentedKalmanFilter.py
1 # -*- coding: utf-8 -*-
2 #
3 # Copyright (C) 2008-2018 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         self.requireInputArguments(
92             mandatory= ("Xb", "Y", "HO", "R", "B" ),
93             optional = ("U", "EM", "CM", "Q"),
94             )
95
96     def run(self, Xb=None, Y=None, U=None, HO=None, EM=None, CM=None, R=None, B=None, Q=None, Parameters=None):
97         self._pre_run(Parameters, Xb, Y, R, B, Q)
98         #
99         if self._parameters["EstimationOf"] == "Parameters":
100             self._parameters["StoreInternalVariables"] = True
101         #
102         L     = Xb.size
103         Alpha = self._parameters["Alpha"]
104         Beta  = self._parameters["Beta"]
105         if self._parameters["Kappa"] == 0:
106             if self._parameters["EstimationOf"] == "State":
107                 Kappa = 0
108             elif self._parameters["EstimationOf"] == "Parameters":
109                 Kappa = 3 - L
110         else:
111             Kappa = self._parameters["Kappa"]
112         Lambda = float( Alpha**2 ) * ( L + Kappa ) - L
113         Gamma  = math.sqrt( L + Lambda )
114         #
115         Ww = []
116         Ww.append( 0. )
117         for i in range(2*L):
118             Ww.append( 1. / (2.*(L + Lambda)) )
119         #
120         Wm = numpy.array( Ww )
121         Wm[0] = Lambda / (L + Lambda)
122         Wc = numpy.array( Ww )
123         Wc[0] = Lambda / (L + Lambda) + (1. - Alpha**2 + Beta)
124         #
125         # Opérateurs
126         # ----------
127         H = HO["Direct"].appliedControledFormTo
128         #
129         if self._parameters["EstimationOf"] == "State":
130             M = EM["Direct"].appliedControledFormTo
131         #
132         if CM is not None and "Tangent" in CM and U is not None:
133             Cm = CM["Tangent"].asMatrix(Xb)
134         else:
135             Cm = None
136         #
137         # Nombre de pas identique au nombre de pas d'observations
138         # -------------------------------------------------------
139         if hasattr(Y,"stepnumber"):
140             duration = Y.stepnumber()
141         else:
142             duration = 2
143         #
144         # Précalcul des inversions de B et R
145         # ----------------------------------
146         if self._parameters["StoreInternalVariables"]:
147             BI = B.getI()
148             RI = R.getI()
149         #
150         # Initialisation
151         # --------------
152         Xn = Xb
153         if hasattr(B,"asfullmatrix"):
154             Pn = B.asfullmatrix(Xn.size)
155         else:
156             Pn = B
157         #
158         self.StoredVariables["Analysis"].store( Xn.A1 )
159         if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
160             self.StoredVariables["APosterioriCovariance"].store( Pn )
161             covarianceXa = Pn
162         Xa               = Xn
163         previousJMinimum = numpy.finfo(float).max
164         #
165         for step in range(duration-1):
166             if hasattr(Y,"store"):
167                 Ynpu = numpy.asmatrix(numpy.ravel( Y[step+1] )).T
168             else:
169                 Ynpu = numpy.asmatrix(numpy.ravel( Y )).T
170             #
171             if U is not None:
172                 if hasattr(U,"store") and len(U)>1:
173                     Un = numpy.asmatrix(numpy.ravel( U[step] )).T
174                 elif hasattr(U,"store") and len(U)==1:
175                     Un = numpy.asmatrix(numpy.ravel( U[0] )).T
176                 else:
177                     Un = numpy.asmatrix(numpy.ravel( U )).T
178             else:
179                 Un = None
180             #
181             Pndemi = numpy.linalg.cholesky(Pn)
182             Xnp = numpy.hstack([Xn, Xn+Gamma*Pndemi, Xn-Gamma*Pndemi])
183             nbSpts = 2*Xn.size+1
184             #
185             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
186                 for point in range(nbSpts):
187                     Xnp[:,point] = numpy.max(numpy.hstack((Xnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
188                     Xnp[:,point] = numpy.min(numpy.hstack((Xnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
189             #
190             XEtnnp = []
191             for point in range(nbSpts):
192                 if self._parameters["EstimationOf"] == "State":
193                     XEtnnpi = numpy.asmatrix(numpy.ravel( M( (Xnp[:,point], Un) ) )).T
194                     if Cm is not None and Un is not None: # Attention : si Cm est aussi dans M, doublon !
195                         Cm = Cm.reshape(Xn.size,Un.size) # ADAO & check shape
196                         XEtnnpi = XEtnnpi + Cm * Un
197                     if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
198                         XEtnnpi = numpy.max(numpy.hstack((XEtnnpi,numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
199                         XEtnnpi = numpy.min(numpy.hstack((XEtnnpi,numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
200                 elif self._parameters["EstimationOf"] == "Parameters":
201                     # --- > Par principe, M = Id, Q = 0
202                     XEtnnpi = Xnp[:,point]
203                 XEtnnp.append( XEtnnpi )
204             XEtnnp = numpy.hstack( XEtnnp )
205             #
206             Xncm = numpy.matrix( XEtnnp.getA()*numpy.array(Wm) ).sum(axis=1)
207             #
208             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
209                 Xncm = numpy.max(numpy.hstack((Xncm,numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
210                 Xncm = numpy.min(numpy.hstack((Xncm,numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
211             #
212             if self._parameters["EstimationOf"] == "State":        Pnm = Q
213             elif self._parameters["EstimationOf"] == "Parameters": Pnm = 0.
214             for point in range(nbSpts):
215                 Pnm += Wc[i] * (XEtnnp[:,point]-Xncm) * (XEtnnp[:,point]-Xncm).T
216             #
217             if self._parameters["EstimationOf"] == "Parameters" and self._parameters["Bounds"] is not None:
218                 Pnmdemi = self._parameters["Reconditioner"] * numpy.linalg.cholesky(Pnm)
219             else:
220                 Pnmdemi = numpy.linalg.cholesky(Pnm)
221             #
222             Xnnp = numpy.hstack([Xncm, Xncm+Gamma*Pnmdemi, Xncm-Gamma*Pnmdemi])
223             #
224             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
225                 for point in range(nbSpts):
226                     Xnnp[:,point] = numpy.max(numpy.hstack((Xnnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
227                     Xnnp[:,point] = numpy.min(numpy.hstack((Xnnp[:,point],numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
228             #
229             Ynnp = []
230             for point in range(nbSpts):
231                 if self._parameters["EstimationOf"] == "State":
232                     Ynnpi = numpy.asmatrix(numpy.ravel( H( (Xnnp[:,point], None) ) )).T
233                 elif self._parameters["EstimationOf"] == "Parameters":
234                     Ynnpi = numpy.asmatrix(numpy.ravel( H( (Xnnp[:,point], Un) ) )).T
235                 Ynnp.append( Ynnpi )
236             Ynnp = numpy.hstack( Ynnp )
237             #
238             Yncm = numpy.matrix( Ynnp.getA()*numpy.array(Wm) ).sum(axis=1)
239             #
240             Pyyn = R
241             Pxyn = 0.
242             for point in range(nbSpts):
243                 Pyyn += Wc[i] * (Ynnp[:,point]-Yncm) * (Ynnp[:,point]-Yncm).T
244                 Pxyn += Wc[i] * (Xnnp[:,point]-Xncm) * (Ynnp[:,point]-Yncm).T
245             #
246             d  = Ynpu - Yncm
247             if self._parameters["EstimationOf"] == "Parameters":
248                 if Cm is not None and Un is not None: # Attention : si Cm est aussi dans H, doublon !
249                     d = d - Cm * Un
250             #
251             Kn = Pxyn * Pyyn.I
252             Xn = Xncm + Kn * d
253             Pn = Pnm - Kn * Pyyn * Kn.T
254             #
255             if self._parameters["Bounds"] is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
256                 Xn = numpy.max(numpy.hstack((Xn,numpy.asmatrix(self._parameters["Bounds"])[:,0])),axis=1)
257                 Xn = numpy.min(numpy.hstack((Xn,numpy.asmatrix(self._parameters["Bounds"])[:,1])),axis=1)
258             #
259             self.StoredVariables["Analysis"].store( Xn.A1 )
260             if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
261                 self.StoredVariables["APosterioriCovariance"].store( Pn )
262             if "Innovation" in self._parameters["StoreSupplementaryCalculations"]:
263                 self.StoredVariables["Innovation"].store( numpy.ravel( d.A1 ) )
264             if self._parameters["StoreInternalVariables"]:
265                 Jb  = 0.5 * (Xn - Xb).T * BI * (Xn - Xb)
266                 Jo  = 0.5 * d.T * RI * d
267                 J   = float( Jb ) + float( Jo )
268                 if self._parameters["StoreInternalVariables"] or "CurrentState" in self._parameters["StoreSupplementaryCalculations"]:
269                     self.StoredVariables["CurrentState"].store( Xn )
270                 self.StoredVariables["CostFunctionJb"].store( Jb )
271                 self.StoredVariables["CostFunctionJo"].store( Jo )
272                 self.StoredVariables["CostFunctionJ" ].store( J )
273                 if J < previousJMinimum:
274                     previousJMinimum  = J
275                     Xa                = Xn
276                     if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
277                         covarianceXa  = Pn
278             else:
279                 Xa = Xn
280             #
281         #
282         # Stockage supplementaire de l'optimum en estimation de parametres
283         # ----------------------------------------------------------------
284         if self._parameters["EstimationOf"] == "Parameters":
285             self.StoredVariables["Analysis"].store( Xa.A1 )
286             if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
287                 self.StoredVariables["APosterioriCovariance"].store( covarianceXa )
288         #
289         if "BMA" in self._parameters["StoreSupplementaryCalculations"]:
290             self.StoredVariables["BMA"].store( numpy.ravel(Xb) - numpy.ravel(Xa) )
291         #
292         self._post_run(HO)
293         return 0
294
295 # ==============================================================================
296 if __name__ == "__main__":
297     print('\n AUTODIAGNOSTIC \n')