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