1 #-*-coding:iso-8859-1-*-
3 # Copyright (C) 2008-2014 EDF R&D
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.
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.
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
19 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 # Author: Jean-Philippe Argaud, jean-philippe.argaud@edf.fr, EDF R&D
24 from daCore import BasicObjects, PlatformInfo
25 m = PlatformInfo.SystemUsage()
28 # ==============================================================================
29 class ElementaryAlgorithm(BasicObjects.Algorithm):
31 BasicObjects.Algorithm.__init__(self, "UNSCENTEDKALMANFILTER")
32 self.defineRequiredParameter(
33 name = "ConstrainedBy",
34 default = "EstimateProjection",
36 message = "Prise en compte des contraintes",
37 listval = ["EstimateProjection"],
39 self.defineRequiredParameter(
40 name = "EstimationOf",
43 message = "Estimation d'etat ou de parametres",
44 listval = ["State", "Parameters"],
46 self.defineRequiredParameter(
54 self.defineRequiredParameter(
60 self.defineRequiredParameter(
67 self.defineRequiredParameter(
68 name = "Reconditioner",
75 self.defineRequiredParameter(
76 name = "StoreInternalVariables",
79 message = "Stockage des variables internes ou intermédiaires du calcul",
81 self.defineRequiredParameter(
82 name = "StoreSupplementaryCalculations",
85 message = "Liste de calculs supplémentaires à stocker et/ou effectuer",
86 listval = ["APosterioriCovariance", "BMA", "Innovation"]
89 def run(self, Xb=None, Y=None, U=None, HO=None, EM=None, CM=None, R=None, B=None, Q=None, Parameters=None):
90 logging.debug("%s Lancement"%self._name)
91 logging.debug("%s Taille mémoire utilisée de %.1f Mo"%(self._name, m.getUsedMemory("M")))
93 # Paramètres de pilotage
94 # ----------------------
95 self.setParameters(Parameters)
97 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):
98 Bounds = self._parameters["Bounds"]
99 logging.debug("%s Prise en compte des bornes effectuee"%(self._name,))
102 if self._parameters["EstimationOf"] == "Parameters":
103 self._parameters["StoreInternalVariables"] = True
106 Alpha = self._parameters["Alpha"]
107 Beta = self._parameters["Beta"]
108 if self._parameters["Kappa"] == 0:
109 if self._parameters["EstimationOf"] == "State":
111 elif self._parameters["EstimationOf"] == "Parameters":
114 Kappa = self._parameters["Kappa"]
115 Lambda = float( Alpha**2 ) * ( L + Kappa ) - L
116 Gamma = math.sqrt( L + Lambda )
121 Ww.append( 1. / (2.*(L + Lambda)) )
123 Wm = numpy.array( Ww )
124 Wm[0] = Lambda / (L + Lambda)
125 Wc = numpy.array( Ww )
126 Wc[0] = Lambda / (L + Lambda) + (1. - Alpha**2 + Beta)
131 raise ValueError("Background error covariance matrix has to be properly defined!")
133 raise ValueError("Observation error covariance matrix has to be properly defined!")
135 H = HO["Direct"].appliedControledFormTo
137 if self._parameters["EstimationOf"] == "State":
138 M = EM["Direct"].appliedControledFormTo
140 if CM is not None and CM.has_key("Tangent") and U is not None:
141 Cm = CM["Tangent"].asMatrix(Xb)
145 # Nombre de pas du Kalman identique au nombre de pas d'observations
146 # -----------------------------------------------------------------
147 if hasattr(Y,"stepnumber"):
148 duration = Y.stepnumber()
152 # Précalcul des inversions de B et R
153 # ----------------------------------
154 if self._parameters["StoreInternalVariables"]:
161 if hasattr(B,"asfullmatrix"):
162 Pn = B.asfullmatrix(Xn.size)
166 self.StoredVariables["Analysis"].store( Xn.A1 )
167 if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
168 self.StoredVariables["APosterioriCovariance"].store( Pn )
171 previousJMinimum = numpy.finfo(float).max
173 for step in range(duration-1):
174 if hasattr(Y,"store"):
175 Ynpu = numpy.asmatrix(numpy.ravel( Y[step+1] )).T
177 Ynpu = numpy.asmatrix(numpy.ravel( Y )).T
180 if hasattr(U,"store") and len(U)>1:
181 Un = numpy.asmatrix(numpy.ravel( U[step] )).T
182 elif hasattr(U,"store") and len(U)==1:
183 Un = numpy.asmatrix(numpy.ravel( U[0] )).T
185 Un = numpy.asmatrix(numpy.ravel( U )).T
189 Pndemi = numpy.linalg.cholesky(Pn)
190 Xnp = numpy.hstack([Xn, Xn+Gamma*Pndemi, Xn-Gamma*Pndemi])
193 if Bounds is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
194 for point in range(nbSpts):
195 Xnp[:,point] = numpy.max(numpy.hstack((Xnp[:,point],numpy.asmatrix(Bounds)[:,0])),axis=1)
196 Xnp[:,point] = numpy.min(numpy.hstack((Xnp[:,point],numpy.asmatrix(Bounds)[:,1])),axis=1)
199 for point in range(nbSpts):
200 if self._parameters["EstimationOf"] == "State":
201 XEtnnpi = numpy.asmatrix(numpy.ravel( M( (Xnp[:,point], Un) ) )).T
202 if Cm is not None and Un is not None: # Attention : si Cm est aussi dans M, doublon !
203 Cm = Cm.reshape(Xn.size,Un.size) # ADAO & check shape
204 XEtnnpi = XEtnnpi + Cm * Un
205 if Bounds is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
206 XEtnnpi = numpy.max(numpy.hstack((XEtnnpi,numpy.asmatrix(Bounds)[:,0])),axis=1)
207 XEtnnpi = numpy.min(numpy.hstack((XEtnnpi,numpy.asmatrix(Bounds)[:,1])),axis=1)
208 elif self._parameters["EstimationOf"] == "Parameters":
209 # --- > Par principe, M = Id, Q = 0
210 XEtnnpi = Xnp[:,point]
211 XEtnnp.append( XEtnnpi )
212 XEtnnp = numpy.hstack( XEtnnp )
214 Xncm = numpy.matrix( XEtnnp.getA()*numpy.array(Wm) ).sum(axis=1)
216 if Bounds is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
217 Xncm = numpy.max(numpy.hstack((Xncm,numpy.asmatrix(Bounds)[:,0])),axis=1)
218 Xncm = numpy.min(numpy.hstack((Xncm,numpy.asmatrix(Bounds)[:,1])),axis=1)
220 if self._parameters["EstimationOf"] == "State": Pnm = Q
221 elif self._parameters["EstimationOf"] == "Parameters": Pnm = 0.
222 for point in range(nbSpts):
223 Pnm += Wc[i] * (XEtnnp[:,point]-Xncm) * (XEtnnp[:,point]-Xncm).T
225 if self._parameters["EstimationOf"] == "Parameters" and Bounds is not None:
226 Pnmdemi = self._parameters["Reconditioner"] * numpy.linalg.cholesky(Pnm)
228 Pnmdemi = numpy.linalg.cholesky(Pnm)
230 Xnnp = numpy.hstack([Xncm, Xncm+Gamma*Pnmdemi, Xncm-Gamma*Pnmdemi])
232 if Bounds is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
233 for point in range(nbSpts):
234 Xnnp[:,point] = numpy.max(numpy.hstack((Xnnp[:,point],numpy.asmatrix(Bounds)[:,0])),axis=1)
235 Xnnp[:,point] = numpy.min(numpy.hstack((Xnnp[:,point],numpy.asmatrix(Bounds)[:,1])),axis=1)
238 for point in range(nbSpts):
239 if self._parameters["EstimationOf"] == "State":
240 Ynnpi = numpy.asmatrix(numpy.ravel( H( (Xnnp[:,point], None) ) )).T
241 elif self._parameters["EstimationOf"] == "Parameters":
242 Ynnpi = numpy.asmatrix(numpy.ravel( H( (Xnnp[:,point], Un) ) )).T
244 Ynnp = numpy.hstack( Ynnp )
246 Yncm = numpy.matrix( Ynnp.getA()*numpy.array(Wm) ).sum(axis=1)
250 for point in range(nbSpts):
251 Pyyn += Wc[i] * (Ynnp[:,point]-Yncm) * (Ynnp[:,point]-Yncm).T
252 Pxyn += Wc[i] * (Xnnp[:,point]-Xncm) * (Ynnp[:,point]-Yncm).T
255 if self._parameters["EstimationOf"] == "Parameters":
256 if Cm is not None and Un is not None: # Attention : si Cm est aussi dans H, doublon !
261 Pn = Pnm - Kn * Pyyn * Kn.T
263 if Bounds is not None and self._parameters["ConstrainedBy"] == "EstimateProjection":
264 Xn = numpy.max(numpy.hstack((Xn,numpy.asmatrix(Bounds)[:,0])),axis=1)
265 Xn = numpy.min(numpy.hstack((Xn,numpy.asmatrix(Bounds)[:,1])),axis=1)
267 self.StoredVariables["Analysis"].store( Xn.A1 )
268 if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
269 self.StoredVariables["APosterioriCovariance"].store( Pn )
270 if "Innovation" in self._parameters["StoreSupplementaryCalculations"]:
271 self.StoredVariables["Innovation"].store( numpy.ravel( d.A1 ) )
272 if self._parameters["StoreInternalVariables"]:
273 Jb = 0.5 * (Xn - Xb).T * BI * (Xn - Xb)
274 Jo = 0.5 * d.T * RI * d
275 J = float( Jb ) + float( Jo )
276 self.StoredVariables["CurrentState"].store( Xn )
277 self.StoredVariables["CostFunctionJb"].store( Jb )
278 self.StoredVariables["CostFunctionJo"].store( Jo )
279 self.StoredVariables["CostFunctionJ" ].store( J )
280 if J < previousJMinimum:
283 if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
289 # Stockage supplementaire de l'optimum en estimation de parametres
290 # ----------------------------------------------------------------
291 if self._parameters["EstimationOf"] == "Parameters":
292 self.StoredVariables["Analysis"].store( Xa.A1 )
293 if "APosterioriCovariance" in self._parameters["StoreSupplementaryCalculations"]:
294 self.StoredVariables["APosterioriCovariance"].store( covarianceXa )
296 if "BMA" in self._parameters["StoreSupplementaryCalculations"]:
297 self.StoredVariables["BMA"].store( numpy.ravel(Xb) - numpy.ravel(Xa) )
299 logging.debug("%s Nombre d'évaluation(s) de l'opérateur d'observation direct/tangent/adjoint : %i/%i/%i"%(self._name, HO["Direct"].nbcalls(0),HO["Tangent"].nbcalls(0),HO["Adjoint"].nbcalls(0)))
300 logging.debug("%s Taille mémoire utilisée de %.1f Mo"%(self._name, m.getUsedMemory("M")))
301 logging.debug("%s Terminé"%self._name)
305 # ==============================================================================
306 if __name__ == "__main__":
307 print '\n AUTODIAGNOSTIC \n'