Salome HOME
Code improvements for warning on iteration control
[modules/adao.git] / src / daComposant / daCore / NumericObjects.py
1 # -*- coding: utf-8 -*-
2 #
3 # Copyright (C) 2008-2022 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 __doc__ = """
24     Définit les objets numériques génériques.
25 """
26 __author__ = "Jean-Philippe ARGAUD"
27
28 import os, copy, types, sys, logging, numpy
29 from daCore.BasicObjects import Operator, Covariance, PartialAlgorithm
30 from daCore.PlatformInfo import PlatformInfo
31 mpr = PlatformInfo().MachinePrecision()
32 mfp = PlatformInfo().MaximumPrecision()
33 # logging.getLogger().setLevel(logging.DEBUG)
34
35 # ==============================================================================
36 def ExecuteFunction( triplet ):
37     assert len(triplet) == 3, "Incorrect number of arguments"
38     X, xArgs, funcrepr = triplet
39     __X = numpy.ravel( X ).reshape((-1,1))
40     __sys_path_tmp = sys.path ; sys.path.insert(0,funcrepr["__userFunction__path"])
41     __module = __import__(funcrepr["__userFunction__modl"], globals(), locals(), [])
42     __fonction = getattr(__module,funcrepr["__userFunction__name"])
43     sys.path = __sys_path_tmp ; del __sys_path_tmp
44     if isinstance(xArgs, dict):
45         __HX  = __fonction( __X, **xArgs )
46     else:
47         __HX  = __fonction( __X )
48     return numpy.ravel( __HX )
49
50 # ==============================================================================
51 class FDApproximation(object):
52     """
53     Cette classe sert d'interface pour définir les opérateurs approximés. A la
54     création d'un objet, en fournissant une fonction "Function", on obtient un
55     objet qui dispose de 3 méthodes "DirectOperator", "TangentOperator" et
56     "AdjointOperator". On contrôle l'approximation DF avec l'incrément
57     multiplicatif "increment" valant par défaut 1%, ou avec l'incrément fixe
58     "dX" qui sera multiplié par "increment" (donc en %), et on effectue de DF
59     centrées si le booléen "centeredDF" est vrai.
60     """
61     def __init__(self,
62             name                  = "FDApproximation",
63             Function              = None,
64             centeredDF            = False,
65             increment             = 0.01,
66             dX                    = None,
67             extraArguments        = None,
68             reducingMemoryUse     = False,
69             avoidingRedundancy    = True,
70             toleranceInRedundancy = 1.e-18,
71             lenghtOfRedundancy    = -1,
72             mpEnabled             = False,
73             mpWorkers             = None,
74             mfEnabled             = False,
75             ):
76         self.__name = str(name)
77         self.__extraArgs = extraArguments
78         #
79         if mpEnabled:
80             try:
81                 import multiprocessing
82                 self.__mpEnabled = True
83             except ImportError:
84                 self.__mpEnabled = False
85         else:
86             self.__mpEnabled = False
87         self.__mpWorkers = mpWorkers
88         if self.__mpWorkers is not None and self.__mpWorkers < 1:
89             self.__mpWorkers = None
90         logging.debug("FDA Calculs en multiprocessing : %s (nombre de processus : %s)"%(self.__mpEnabled,self.__mpWorkers))
91         #
92         self.__mfEnabled = bool(mfEnabled)
93         logging.debug("FDA Calculs en multifonctions : %s"%(self.__mfEnabled,))
94         #
95         self.__rmEnabled = bool(reducingMemoryUse)
96         logging.debug("FDA Calculs avec réduction mémoire : %s"%(self.__rmEnabled,))
97         #
98         if avoidingRedundancy:
99             self.__avoidRC = True
100             self.__tolerBP = float(toleranceInRedundancy)
101             self.__lenghtRJ = int(lenghtOfRedundancy)
102             self.__listJPCP = [] # Jacobian Previous Calculated Points
103             self.__listJPCI = [] # Jacobian Previous Calculated Increment
104             self.__listJPCR = [] # Jacobian Previous Calculated Results
105             self.__listJPPN = [] # Jacobian Previous Calculated Point Norms
106             self.__listJPIN = [] # Jacobian Previous Calculated Increment Norms
107         else:
108             self.__avoidRC = False
109         logging.debug("FDA Calculs avec réduction des doublons : %s"%self.__avoidRC)
110         if self.__avoidRC:
111             logging.debug("FDA Tolérance de détermination des doublons : %.2e"%self.__tolerBP)
112         #
113         if self.__mpEnabled:
114             if isinstance(Function,types.FunctionType):
115                 logging.debug("FDA Calculs en multiprocessing : FunctionType")
116                 self.__userFunction__name = Function.__name__
117                 try:
118                     mod = os.path.join(Function.__globals__['filepath'],Function.__globals__['filename'])
119                 except:
120                     mod = os.path.abspath(Function.__globals__['__file__'])
121                 if not os.path.isfile(mod):
122                     raise ImportError("No user defined function or method found with the name %s"%(mod,))
123                 self.__userFunction__modl = os.path.basename(mod).replace('.pyc','').replace('.pyo','').replace('.py','')
124                 self.__userFunction__path = os.path.dirname(mod)
125                 del mod
126                 self.__userOperator = Operator( name = self.__name, fromMethod = Function, avoidingRedundancy = self.__avoidRC, inputAsMultiFunction = self.__mfEnabled, extraArguments = self.__extraArgs )
127                 self.__userFunction = self.__userOperator.appliedTo # Pour le calcul Direct
128             elif isinstance(Function,types.MethodType):
129                 logging.debug("FDA Calculs en multiprocessing : MethodType")
130                 self.__userFunction__name = Function.__name__
131                 try:
132                     mod = os.path.join(Function.__globals__['filepath'],Function.__globals__['filename'])
133                 except:
134                     mod = os.path.abspath(Function.__func__.__globals__['__file__'])
135                 if not os.path.isfile(mod):
136                     raise ImportError("No user defined function or method found with the name %s"%(mod,))
137                 self.__userFunction__modl = os.path.basename(mod).replace('.pyc','').replace('.pyo','').replace('.py','')
138                 self.__userFunction__path = os.path.dirname(mod)
139                 del mod
140                 self.__userOperator = Operator( name = self.__name, fromMethod = Function, avoidingRedundancy = self.__avoidRC, inputAsMultiFunction = self.__mfEnabled, extraArguments = self.__extraArgs )
141                 self.__userFunction = self.__userOperator.appliedTo # Pour le calcul Direct
142             else:
143                 raise TypeError("User defined function or method has to be provided for finite differences approximation.")
144         else:
145             self.__userOperator = Operator( name = self.__name, fromMethod = Function, avoidingRedundancy = self.__avoidRC, inputAsMultiFunction = self.__mfEnabled, extraArguments = self.__extraArgs )
146             self.__userFunction = self.__userOperator.appliedTo
147         #
148         self.__centeredDF = bool(centeredDF)
149         if abs(float(increment)) > 1.e-15:
150             self.__increment  = float(increment)
151         else:
152             self.__increment  = 0.01
153         if dX is None:
154             self.__dX     = None
155         else:
156             self.__dX     = numpy.ravel( dX )
157
158     # ---------------------------------------------------------
159     def __doublon__(self, e, l, n, v=None):
160         __ac, __iac = False, -1
161         for i in range(len(l)-1,-1,-1):
162             if numpy.linalg.norm(e - l[i]) < self.__tolerBP * n[i]:
163                 __ac, __iac = True, i
164                 if v is not None: logging.debug("FDA Cas%s déja calculé, récupération du doublon %i"%(v,__iac))
165                 break
166         return __ac, __iac
167
168     # ---------------------------------------------------------
169     def __listdotwith__(self, __LMatrix, __dotWith = None, __dotTWith = None):
170         "Produit incrémental d'une matrice liste de colonnes avec un vecteur"
171         if not isinstance(__LMatrix, (list,tuple)):
172             raise TypeError("Columnwise list matrix has not the proper type: %s"%type(__LMatrix))
173         if __dotWith is not None:
174             __Idwx = numpy.ravel( __dotWith )
175             assert len(__LMatrix) == __Idwx.size, "Incorrect size of elements"
176             __Produit = numpy.zeros(__LMatrix[0].size)
177             for i, col in enumerate(__LMatrix):
178                 __Produit += float(__Idwx[i]) * col
179             return __Produit
180         elif __dotTWith is not None:
181             _Idwy = numpy.ravel( __dotTWith ).T
182             assert __LMatrix[0].size == _Idwy.size, "Incorrect size of elements"
183             __Produit = numpy.zeros(len(__LMatrix))
184             for i, col in enumerate(__LMatrix):
185                 __Produit[i] = float( _Idwy @ col)
186             return __Produit
187         else:
188             __Produit = None
189         return __Produit
190
191     # ---------------------------------------------------------
192     def DirectOperator(self, X, **extraArgs ):
193         """
194         Calcul du direct à l'aide de la fonction fournie.
195
196         NB : les extraArgs sont là pour assurer la compatibilité d'appel, mais
197         ne doivent pas être données ici à la fonction utilisateur.
198         """
199         logging.debug("FDA Calcul DirectOperator (explicite)")
200         if self.__mfEnabled:
201             _HX = self.__userFunction( X, argsAsSerie = True )
202         else:
203             _HX = numpy.ravel(self.__userFunction( numpy.ravel(X) ))
204         #
205         return _HX
206
207     # ---------------------------------------------------------
208     def TangentMatrix(self, X, dotWith = None, dotTWith = None ):
209         """
210         Calcul de l'opérateur tangent comme la Jacobienne par différences finies,
211         c'est-à-dire le gradient de H en X. On utilise des différences finies
212         directionnelles autour du point X. X est un numpy.ndarray.
213
214         Différences finies centrées (approximation d'ordre 2):
215         1/ Pour chaque composante i de X, on ajoute et on enlève la perturbation
216            dX[i] à la  composante X[i], pour composer X_plus_dXi et X_moins_dXi, et
217            on calcule les réponses HX_plus_dXi = H( X_plus_dXi ) et HX_moins_dXi =
218            H( X_moins_dXi )
219         2/ On effectue les différences (HX_plus_dXi-HX_moins_dXi) et on divise par
220            le pas 2*dXi
221         3/ Chaque résultat, par composante, devient une colonne de la Jacobienne
222
223         Différences finies non centrées (approximation d'ordre 1):
224         1/ Pour chaque composante i de X, on ajoute la perturbation dX[i] à la
225            composante X[i] pour composer X_plus_dXi, et on calcule la réponse
226            HX_plus_dXi = H( X_plus_dXi )
227         2/ On calcule la valeur centrale HX = H(X)
228         3/ On effectue les différences (HX_plus_dXi-HX) et on divise par
229            le pas dXi
230         4/ Chaque résultat, par composante, devient une colonne de la Jacobienne
231
232         """
233         logging.debug("FDA Début du calcul de la Jacobienne")
234         logging.debug("FDA   Incrément de............: %s*X"%float(self.__increment))
235         logging.debug("FDA   Approximation centrée...: %s"%(self.__centeredDF))
236         #
237         if X is None or len(X)==0:
238             raise ValueError("Nominal point X for approximate derivatives can not be None or void (given X: %s)."%(str(X),))
239         #
240         _X = numpy.ravel( X )
241         #
242         if self.__dX is None:
243             _dX  = self.__increment * _X
244         else:
245             _dX = numpy.ravel( self.__dX )
246         assert len(_X) == len(_dX), "Inconsistent dX increment length with respect to the X one"
247         assert _X.size == _dX.size, "Inconsistent dX increment size with respect to the X one"
248         #
249         if (_dX == 0.).any():
250             moyenne = _dX.mean()
251             if moyenne == 0.:
252                 _dX = numpy.where( _dX == 0., float(self.__increment), _dX )
253             else:
254                 _dX = numpy.where( _dX == 0., moyenne, _dX )
255         #
256         __alreadyCalculated  = False
257         if self.__avoidRC:
258             __bidon, __alreadyCalculatedP = self.__doublon__(_X,  self.__listJPCP, self.__listJPPN, None)
259             __bidon, __alreadyCalculatedI = self.__doublon__(_dX, self.__listJPCI, self.__listJPIN, None)
260             if __alreadyCalculatedP == __alreadyCalculatedI > -1:
261                 __alreadyCalculated, __i = True, __alreadyCalculatedP
262                 logging.debug("FDA Cas J déjà calculé, récupération du doublon %i"%__i)
263         #
264         if __alreadyCalculated:
265             logging.debug("FDA   Calcul Jacobienne (par récupération du doublon %i)"%__i)
266             _Jacobienne = self.__listJPCR[__i]
267             logging.debug("FDA Fin du calcul de la Jacobienne")
268             if dotWith is not None:
269                 return numpy.dot(_Jacobienne,   numpy.ravel( dotWith ))
270             elif dotTWith is not None:
271                 return numpy.dot(_Jacobienne.T, numpy.ravel( dotTWith ))
272         else:
273             logging.debug("FDA   Calcul Jacobienne (explicite)")
274             if self.__centeredDF:
275                 #
276                 if self.__mpEnabled and not self.__mfEnabled:
277                     funcrepr = {
278                         "__userFunction__path" : self.__userFunction__path,
279                         "__userFunction__modl" : self.__userFunction__modl,
280                         "__userFunction__name" : self.__userFunction__name,
281                     }
282                     _jobs = []
283                     for i in range( len(_dX) ):
284                         _dXi            = _dX[i]
285                         _X_plus_dXi     = numpy.array( _X, dtype=float )
286                         _X_plus_dXi[i]  = _X[i] + _dXi
287                         _X_moins_dXi    = numpy.array( _X, dtype=float )
288                         _X_moins_dXi[i] = _X[i] - _dXi
289                         #
290                         _jobs.append( (_X_plus_dXi,  self.__extraArgs, funcrepr) )
291                         _jobs.append( (_X_moins_dXi, self.__extraArgs, funcrepr) )
292                     #
293                     import multiprocessing
294                     self.__pool = multiprocessing.Pool(self.__mpWorkers)
295                     _HX_plusmoins_dX = self.__pool.map( ExecuteFunction, _jobs )
296                     self.__pool.close()
297                     self.__pool.join()
298                     #
299                     _Jacobienne  = []
300                     for i in range( len(_dX) ):
301                         _Jacobienne.append( numpy.ravel( _HX_plusmoins_dX[2*i] - _HX_plusmoins_dX[2*i+1] ) / (2.*_dX[i]) )
302                     #
303                 elif self.__mfEnabled:
304                     _xserie = []
305                     for i in range( len(_dX) ):
306                         _dXi            = _dX[i]
307                         _X_plus_dXi     = numpy.array( _X, dtype=float )
308                         _X_plus_dXi[i]  = _X[i] + _dXi
309                         _X_moins_dXi    = numpy.array( _X, dtype=float )
310                         _X_moins_dXi[i] = _X[i] - _dXi
311                         #
312                         _xserie.append( _X_plus_dXi )
313                         _xserie.append( _X_moins_dXi )
314                     #
315                     _HX_plusmoins_dX = self.DirectOperator( _xserie )
316                      #
317                     _Jacobienne  = []
318                     for i in range( len(_dX) ):
319                         _Jacobienne.append( numpy.ravel( _HX_plusmoins_dX[2*i] - _HX_plusmoins_dX[2*i+1] ) / (2.*_dX[i]) )
320                     #
321                 else:
322                     _Jacobienne  = []
323                     for i in range( _dX.size ):
324                         _dXi            = _dX[i]
325                         _X_plus_dXi     = numpy.array( _X, dtype=float )
326                         _X_plus_dXi[i]  = _X[i] + _dXi
327                         _X_moins_dXi    = numpy.array( _X, dtype=float )
328                         _X_moins_dXi[i] = _X[i] - _dXi
329                         #
330                         _HX_plus_dXi    = self.DirectOperator( _X_plus_dXi )
331                         _HX_moins_dXi   = self.DirectOperator( _X_moins_dXi )
332                         #
333                         _Jacobienne.append( numpy.ravel( _HX_plus_dXi - _HX_moins_dXi ) / (2.*_dXi) )
334                 #
335             else:
336                 #
337                 if self.__mpEnabled and not self.__mfEnabled:
338                     funcrepr = {
339                         "__userFunction__path" : self.__userFunction__path,
340                         "__userFunction__modl" : self.__userFunction__modl,
341                         "__userFunction__name" : self.__userFunction__name,
342                     }
343                     _jobs = []
344                     _jobs.append( (_X, self.__extraArgs, funcrepr) )
345                     for i in range( len(_dX) ):
346                         _X_plus_dXi    = numpy.array( _X, dtype=float )
347                         _X_plus_dXi[i] = _X[i] + _dX[i]
348                         #
349                         _jobs.append( (_X_plus_dXi, self.__extraArgs, funcrepr) )
350                     #
351                     import multiprocessing
352                     self.__pool = multiprocessing.Pool(self.__mpWorkers)
353                     _HX_plus_dX = self.__pool.map( ExecuteFunction, _jobs )
354                     self.__pool.close()
355                     self.__pool.join()
356                     #
357                     _HX = _HX_plus_dX.pop(0)
358                     #
359                     _Jacobienne = []
360                     for i in range( len(_dX) ):
361                         _Jacobienne.append( numpy.ravel(( _HX_plus_dX[i] - _HX ) / _dX[i]) )
362                     #
363                 elif self.__mfEnabled:
364                     _xserie = []
365                     _xserie.append( _X )
366                     for i in range( len(_dX) ):
367                         _X_plus_dXi    = numpy.array( _X, dtype=float )
368                         _X_plus_dXi[i] = _X[i] + _dX[i]
369                         #
370                         _xserie.append( _X_plus_dXi )
371                     #
372                     _HX_plus_dX = self.DirectOperator( _xserie )
373                     #
374                     _HX = _HX_plus_dX.pop(0)
375                     #
376                     _Jacobienne = []
377                     for i in range( len(_dX) ):
378                         _Jacobienne.append( numpy.ravel(( _HX_plus_dX[i] - _HX ) / _dX[i]) )
379                    #
380                 else:
381                     _Jacobienne  = []
382                     _HX = self.DirectOperator( _X )
383                     for i in range( _dX.size ):
384                         _dXi            = _dX[i]
385                         _X_plus_dXi     = numpy.array( _X, dtype=float )
386                         _X_plus_dXi[i]  = _X[i] + _dXi
387                         #
388                         _HX_plus_dXi = self.DirectOperator( _X_plus_dXi )
389                         #
390                         _Jacobienne.append( numpy.ravel(( _HX_plus_dXi - _HX ) / _dXi) )
391             #
392             if (dotWith is not None) or (dotTWith is not None):
393                 __Produit = self.__listdotwith__(_Jacobienne, dotWith, dotTWith)
394             else:
395                 __Produit = None
396             if __Produit is None or self.__avoidRC:
397                 _Jacobienne = numpy.transpose( numpy.vstack( _Jacobienne ) )
398                 if self.__avoidRC:
399                     if self.__lenghtRJ < 0: self.__lenghtRJ = 2 * _X.size
400                     while len(self.__listJPCP) > self.__lenghtRJ:
401                         self.__listJPCP.pop(0)
402                         self.__listJPCI.pop(0)
403                         self.__listJPCR.pop(0)
404                         self.__listJPPN.pop(0)
405                         self.__listJPIN.pop(0)
406                     self.__listJPCP.append( copy.copy(_X) )
407                     self.__listJPCI.append( copy.copy(_dX) )
408                     self.__listJPCR.append( copy.copy(_Jacobienne) )
409                     self.__listJPPN.append( numpy.linalg.norm(_X) )
410                     self.__listJPIN.append( numpy.linalg.norm(_Jacobienne) )
411             logging.debug("FDA Fin du calcul de la Jacobienne")
412             if __Produit is not None:
413                 return __Produit
414         #
415         return _Jacobienne
416
417     # ---------------------------------------------------------
418     def TangentOperator(self, paire, **extraArgs ):
419         """
420         Calcul du tangent à l'aide de la Jacobienne.
421
422         NB : les extraArgs sont là pour assurer la compatibilité d'appel, mais
423         ne doivent pas être données ici à la fonction utilisateur.
424         """
425         if self.__mfEnabled:
426             assert len(paire) == 1, "Incorrect length of arguments"
427             _paire = paire[0]
428             assert len(_paire) == 2, "Incorrect number of arguments"
429         else:
430             assert len(paire) == 2, "Incorrect number of arguments"
431             _paire = paire
432         X, dX = _paire
433         if dX is None or len(dX) == 0:
434             #
435             # Calcul de la forme matricielle si le second argument est None
436             # -------------------------------------------------------------
437             _Jacobienne = self.TangentMatrix( X )
438             if self.__mfEnabled: return [_Jacobienne,]
439             else:                return _Jacobienne
440         else:
441             #
442             # Calcul de la valeur linéarisée de H en X appliqué à dX
443             # ------------------------------------------------------
444             _HtX = self.TangentMatrix( X, dotWith = dX )
445             if self.__mfEnabled: return [_HtX,]
446             else:                return _HtX
447
448     # ---------------------------------------------------------
449     def AdjointOperator(self, paire, **extraArgs ):
450         """
451         Calcul de l'adjoint à l'aide de la Jacobienne.
452
453         NB : les extraArgs sont là pour assurer la compatibilité d'appel, mais
454         ne doivent pas être données ici à la fonction utilisateur.
455         """
456         if self.__mfEnabled:
457             assert len(paire) == 1, "Incorrect length of arguments"
458             _paire = paire[0]
459             assert len(_paire) == 2, "Incorrect number of arguments"
460         else:
461             assert len(paire) == 2, "Incorrect number of arguments"
462             _paire = paire
463         X, Y = _paire
464         if Y is None or len(Y) == 0:
465             #
466             # Calcul de la forme matricielle si le second argument est None
467             # -------------------------------------------------------------
468             _JacobienneT = self.TangentMatrix( X ).T
469             if self.__mfEnabled: return [_JacobienneT,]
470             else:                return _JacobienneT
471         else:
472             #
473             # Calcul de la valeur de l'adjoint en X appliqué à Y
474             # --------------------------------------------------
475             _HaY = self.TangentMatrix( X, dotTWith = Y )
476             if self.__mfEnabled: return [_HaY,]
477             else:                return _HaY
478
479 # ==============================================================================
480 def SetInitialDirection( __Direction = [], __Amplitude = 1., __Position = None ):
481     "Établit ou élabore une direction avec une amplitude"
482     #
483     if len(__Direction) == 0 and __Position is None:
484         raise ValueError("If initial direction is void, current position has to be given")
485     if abs(float(__Amplitude)) < mpr:
486         raise ValueError("Amplitude of perturbation can not be zero")
487     #
488     if len(__Direction) > 0:
489         __dX0 = numpy.asarray(__Direction)
490     else:
491         __dX0 = []
492         __X0 = numpy.ravel(numpy.asarray(__Position))
493         __mX0 = numpy.mean( __X0, dtype=mfp )
494         if abs(__mX0) < 2*mpr: __mX0 = 1. # Évite le problème de position nulle
495         for v in __X0:
496             if abs(v) > 1.e-8:
497                 __dX0.append( numpy.random.normal(0.,abs(v)) )
498             else:
499                 __dX0.append( numpy.random.normal(0.,__mX0) )
500     #
501     __dX0 = numpy.asarray(__dX0,float) # Évite le problème d'array de taille 1
502     __dX0 = numpy.ravel( __dX0 )       # Redresse les vecteurs
503     __dX0 = float(__Amplitude) * __dX0
504     #
505     return __dX0
506
507 # ==============================================================================
508 def EnsembleOfCenteredPerturbations( __bgCenter, __bgCovariance, __nbMembers ):
509     "Génération d'un ensemble de taille __nbMembers-1 d'états aléatoires centrés"
510     #
511     __bgCenter = numpy.ravel(__bgCenter)[:,None]
512     if __nbMembers < 1:
513         raise ValueError("Number of members has to be strictly more than 1 (given number: %s)."%(str(__nbMembers),))
514     #
515     if __bgCovariance is None:
516         _Perturbations = numpy.tile( __bgCenter, __nbMembers)
517     else:
518         _Z = numpy.random.multivariate_normal(numpy.zeros(__bgCenter.size), __bgCovariance, size=__nbMembers).T
519         _Perturbations = numpy.tile( __bgCenter, __nbMembers) + _Z
520     #
521     return _Perturbations
522
523 # ==============================================================================
524 def EnsembleOfBackgroundPerturbations( __bgCenter, __bgCovariance, __nbMembers, __withSVD = True):
525     "Génération d'un ensemble de taille __nbMembers-1 d'états aléatoires centrés"
526     def __CenteredRandomAnomalies(Zr, N):
527         """
528         Génère une matrice de N anomalies aléatoires centrées sur Zr selon les
529         notes manuscrites de MB et conforme au code de PS avec eps = -1
530         """
531         eps = -1
532         Q = numpy.identity(N-1)-numpy.ones((N-1,N-1))/numpy.sqrt(N)/(numpy.sqrt(N)-eps)
533         Q = numpy.concatenate((Q, [eps*numpy.ones(N-1)/numpy.sqrt(N)]), axis=0)
534         R, _ = numpy.linalg.qr(numpy.random.normal(size = (N-1,N-1)))
535         Q = numpy.dot(Q,R)
536         Zr = numpy.dot(Q,Zr)
537         return Zr.T
538     #
539     __bgCenter = numpy.ravel(__bgCenter).reshape((-1,1))
540     if __nbMembers < 1:
541         raise ValueError("Number of members has to be strictly more than 1 (given number: %s)."%(str(__nbMembers),))
542     if __bgCovariance is None:
543         _Perturbations = numpy.tile( __bgCenter, __nbMembers)
544     else:
545         if __withSVD:
546             _U, _s, _V = numpy.linalg.svd(__bgCovariance, full_matrices=False)
547             _nbctl = __bgCenter.size
548             if __nbMembers > _nbctl:
549                 _Z = numpy.concatenate((numpy.dot(
550                     numpy.diag(numpy.sqrt(_s[:_nbctl])), _V[:_nbctl]),
551                     numpy.random.multivariate_normal(numpy.zeros(_nbctl),__bgCovariance,__nbMembers-1-_nbctl)), axis = 0)
552             else:
553                 _Z = numpy.dot(numpy.diag(numpy.sqrt(_s[:__nbMembers-1])), _V[:__nbMembers-1])
554             _Zca = __CenteredRandomAnomalies(_Z, __nbMembers)
555             _Perturbations = __bgCenter + _Zca
556         else:
557             if max(abs(__bgCovariance.flatten())) > 0:
558                 _nbctl = __bgCenter.size
559                 _Z = numpy.random.multivariate_normal(numpy.zeros(_nbctl),__bgCovariance,__nbMembers-1)
560                 _Zca = __CenteredRandomAnomalies(_Z, __nbMembers)
561                 _Perturbations = __bgCenter + _Zca
562             else:
563                 _Perturbations = numpy.tile( __bgCenter, __nbMembers)
564     #
565     return _Perturbations
566
567 # ==============================================================================
568 def EnsembleMean( __Ensemble ):
569     "Renvoie la moyenne empirique d'un ensemble"
570     return numpy.asarray(__Ensemble).mean(axis=1, dtype=mfp).astype('float').reshape((-1,1))
571
572 # ==============================================================================
573 def EnsembleOfAnomalies( __Ensemble, __OptMean = None, __Normalisation = 1.):
574     "Renvoie les anomalies centrées à partir d'un ensemble"
575     if __OptMean is None:
576         __Em = EnsembleMean( __Ensemble )
577     else:
578         __Em = numpy.ravel( __OptMean ).reshape((-1,1))
579     #
580     return __Normalisation * (numpy.asarray( __Ensemble ) - __Em)
581
582 # ==============================================================================
583 def EnsembleErrorCovariance( __Ensemble, __Quick = False ):
584     "Renvoie l'estimation empirique de la covariance d'ensemble"
585     if __Quick:
586         # Covariance rapide mais rarement définie positive
587         __Covariance = numpy.cov( __Ensemble )
588     else:
589         # Résultat souvent identique à numpy.cov, mais plus robuste
590         __n, __m = numpy.asarray( __Ensemble ).shape
591         __Anomalies = EnsembleOfAnomalies( __Ensemble )
592         # Estimation empirique
593         __Covariance = ( __Anomalies @ __Anomalies.T ) / (__m-1)
594         # Assure la symétrie
595         __Covariance = ( __Covariance + __Covariance.T ) * 0.5
596         # Assure la positivité
597         __epsilon    = mpr*numpy.trace( __Covariance )
598         __Covariance = __Covariance + __epsilon * numpy.identity(__n)
599     #
600     return __Covariance
601
602 # ==============================================================================
603 def EnsemblePerturbationWithGivenCovariance(
604         __Ensemble,
605         __Covariance,
606         __Seed = None,
607         ):
608     "Ajout d'une perturbation à chaque membre d'un ensemble selon une covariance prescrite"
609     if hasattr(__Covariance,"assparsematrix"):
610         if (abs(__Ensemble).mean() > mpr) and (abs(__Covariance.assparsematrix())/abs(__Ensemble).mean() < mpr).all():
611             # Traitement d'une covariance nulle ou presque
612             return __Ensemble
613         if (abs(__Ensemble).mean() <= mpr) and (abs(__Covariance.assparsematrix()) < mpr).all():
614             # Traitement d'une covariance nulle ou presque
615             return __Ensemble
616     else:
617         if (abs(__Ensemble).mean() > mpr) and (abs(__Covariance)/abs(__Ensemble).mean() < mpr).all():
618             # Traitement d'une covariance nulle ou presque
619             return __Ensemble
620         if (abs(__Ensemble).mean() <= mpr) and (abs(__Covariance) < mpr).all():
621             # Traitement d'une covariance nulle ou presque
622             return __Ensemble
623     #
624     __n, __m = __Ensemble.shape
625     if __Seed is not None: numpy.random.seed(__Seed)
626     #
627     if hasattr(__Covariance,"isscalar") and __Covariance.isscalar():
628         # Traitement d'une covariance multiple de l'identité
629         __zero = 0.
630         __std  = numpy.sqrt(__Covariance.assparsematrix())
631         __Ensemble += numpy.random.normal(__zero, __std, size=(__m,__n)).T
632     #
633     elif hasattr(__Covariance,"isvector") and __Covariance.isvector():
634         # Traitement d'une covariance diagonale avec variances non identiques
635         __zero = numpy.zeros(__n)
636         __std  = numpy.sqrt(__Covariance.assparsematrix())
637         __Ensemble += numpy.asarray([numpy.random.normal(__zero, __std) for i in range(__m)]).T
638     #
639     elif hasattr(__Covariance,"ismatrix") and __Covariance.ismatrix():
640         # Traitement d'une covariance pleine
641         __Ensemble += numpy.random.multivariate_normal(numpy.zeros(__n), __Covariance.asfullmatrix(__n), size=__m).T
642     #
643     elif isinstance(__Covariance, numpy.ndarray):
644         # Traitement d'une covariance numpy pleine, sachant qu'on arrive ici en dernier
645         __Ensemble += numpy.random.multivariate_normal(numpy.zeros(__n), __Covariance, size=__m).T
646     #
647     else:
648         raise ValueError("Error in ensemble perturbation with inadequate covariance specification")
649     #
650     return __Ensemble
651
652 # ==============================================================================
653 def CovarianceInflation(
654         __InputCovOrEns,
655         __InflationType   = None,
656         __InflationFactor = None,
657         __BackgroundCov   = None,
658         ):
659     """
660     Inflation applicable soit sur Pb ou Pa, soit sur les ensembles EXb ou EXa
661
662     Synthèse : Hunt 2007, section 2.3.5
663     """
664     if __InflationFactor is None:
665         return __InputCovOrEns
666     else:
667         __InflationFactor = float(__InflationFactor)
668     #
669     __InputCovOrEns = numpy.asarray(__InputCovOrEns)
670     if __InputCovOrEns.size == 0: return __InputCovOrEns
671     #
672     if __InflationType in ["MultiplicativeOnAnalysisCovariance", "MultiplicativeOnBackgroundCovariance"]:
673         if __InflationFactor < 1.:
674             raise ValueError("Inflation factor for multiplicative inflation has to be greater or equal than 1.")
675         if __InflationFactor < 1.+mpr: # No inflation = 1
676             return __InputCovOrEns
677         __OutputCovOrEns = __InflationFactor**2 * __InputCovOrEns
678     #
679     elif __InflationType in ["MultiplicativeOnAnalysisAnomalies", "MultiplicativeOnBackgroundAnomalies"]:
680         if __InflationFactor < 1.:
681             raise ValueError("Inflation factor for multiplicative inflation has to be greater or equal than 1.")
682         if __InflationFactor < 1.+mpr: # No inflation = 1
683             return __InputCovOrEns
684         __InputCovOrEnsMean = __InputCovOrEns.mean(axis=1, dtype=mfp).astype('float')
685         __OutputCovOrEns = __InputCovOrEnsMean[:,numpy.newaxis] \
686             + __InflationFactor * (__InputCovOrEns - __InputCovOrEnsMean[:,numpy.newaxis])
687     #
688     elif __InflationType in ["AdditiveOnAnalysisCovariance", "AdditiveOnBackgroundCovariance"]:
689         if __InflationFactor < 0.:
690             raise ValueError("Inflation factor for additive inflation has to be greater or equal than 0.")
691         if __InflationFactor < mpr: # No inflation = 0
692             return __InputCovOrEns
693         __n, __m = __InputCovOrEns.shape
694         if __n != __m:
695             raise ValueError("Additive inflation can only be applied to squared (covariance) matrix.")
696         __tr = __InputCovOrEns.trace()/__n
697         if __InflationFactor > __tr:
698             raise ValueError("Inflation factor for additive inflation has to be small over %.0e."%__tr)
699         __OutputCovOrEns = (1. - __InflationFactor)*__InputCovOrEns + __InflationFactor * numpy.identity(__n)
700     #
701     elif __InflationType == "HybridOnBackgroundCovariance":
702         if __InflationFactor < 0.:
703             raise ValueError("Inflation factor for hybrid inflation has to be greater or equal than 0.")
704         if __InflationFactor < mpr: # No inflation = 0
705             return __InputCovOrEns
706         __n, __m = __InputCovOrEns.shape
707         if __n != __m:
708             raise ValueError("Additive inflation can only be applied to squared (covariance) matrix.")
709         if __BackgroundCov is None:
710             raise ValueError("Background covariance matrix B has to be given for hybrid inflation.")
711         if __InputCovOrEns.shape != __BackgroundCov.shape:
712             raise ValueError("Ensemble covariance matrix has to be of same size than background covariance matrix B.")
713         __OutputCovOrEns = (1. - __InflationFactor) * __InputCovOrEns + __InflationFactor * __BackgroundCov
714     #
715     elif __InflationType == "Relaxation":
716         raise NotImplementedError("Relaxation inflation type not implemented")
717     #
718     else:
719         raise ValueError("Error in inflation type, '%s' is not a valid keyword."%__InflationType)
720     #
721     return __OutputCovOrEns
722
723 # ==============================================================================
724 def HessienneEstimation(__selfA, __nb, __HaM, __HtM, __BI, __RI):
725     "Estimation de la Hessienne"
726     #
727     __HessienneI = []
728     for i in range(int(__nb)):
729         __ee    = numpy.zeros((__nb,1))
730         __ee[i] = 1.
731         __HtEE  = numpy.dot(__HtM,__ee).reshape((-1,1))
732         __HessienneI.append( numpy.ravel( __BI * __ee + __HaM * (__RI * __HtEE) ) )
733     #
734     __A = numpy.linalg.inv(numpy.array( __HessienneI ))
735     __A = (__A + __A.T) * 0.5 # Symétrie
736     __A = __A + mpr*numpy.trace( __A ) * numpy.identity(__nb) # Positivité
737     #
738     if min(__A.shape) != max(__A.shape):
739         raise ValueError("The %s a posteriori covariance matrix A is of shape %s, despites it has to be a squared matrix. There is an error in the observation operator, please check it."%(__selfA._name,str(__A.shape)))
740     if (numpy.diag(__A) < 0).any():
741         raise ValueError("The %s a posteriori covariance matrix A has at least one negative value on its diagonal. There is an error in the observation operator, please check it."%(__selfA._name,))
742     if logging.getLogger().level < logging.WARNING: # La vérification n'a lieu qu'en debug
743         try:
744             numpy.linalg.cholesky( __A )
745         except:
746             raise ValueError("The %s a posteriori covariance matrix A is not symmetric positive-definite. Please check your a priori covariances and your observation operator."%(__selfA._name,))
747     #
748     return __A
749
750 # ==============================================================================
751 def QuantilesEstimations(selfA, A, Xa, HXa = None, Hm = None, HtM = None):
752     "Estimation des quantiles a posteriori à partir de A>0 (selfA est modifié)"
753     nbsamples = selfA._parameters["NumberOfSamplesForQuantiles"]
754     #
755     # Traitement des bornes
756     if "StateBoundsForQuantiles" in selfA._parameters:
757         LBounds = selfA._parameters["StateBoundsForQuantiles"] # Prioritaire
758     elif "Bounds" in selfA._parameters:
759         LBounds = selfA._parameters["Bounds"]  # Défaut raisonnable
760     else:
761         LBounds = None
762     if LBounds is not None:
763         LBounds = ForceNumericBounds( LBounds )
764     __Xa = numpy.ravel(Xa)
765     #
766     # Échantillonnage des états
767     YfQ  = None
768     EXr  = None
769     for i in range(nbsamples):
770         if selfA._parameters["SimulationForQuantiles"] == "Linear" and HtM is not None and HXa is not None:
771             dXr = (numpy.random.multivariate_normal(__Xa,A) - __Xa).reshape((-1,1))
772             if LBounds is not None: # "EstimateProjection" par défaut
773                 dXr = numpy.max(numpy.hstack((dXr,LBounds[:,0].reshape((-1,1))) - __Xa.reshape((-1,1))),axis=1)
774                 dXr = numpy.min(numpy.hstack((dXr,LBounds[:,1].reshape((-1,1))) - __Xa.reshape((-1,1))),axis=1)
775             dYr = HtM @ dXr
776             Yr = HXa.reshape((-1,1)) + dYr
777             if selfA._toStore("SampledStateForQuantiles"): Xr = __Xa + numpy.ravel(dXr)
778         elif selfA._parameters["SimulationForQuantiles"] == "NonLinear" and Hm is not None:
779             Xr = numpy.random.multivariate_normal(__Xa,A)
780             if LBounds is not None: # "EstimateProjection" par défaut
781                 Xr = numpy.max(numpy.hstack((Xr.reshape((-1,1)),LBounds[:,0].reshape((-1,1)))),axis=1)
782                 Xr = numpy.min(numpy.hstack((Xr.reshape((-1,1)),LBounds[:,1].reshape((-1,1)))),axis=1)
783             Yr = numpy.asarray(Hm( Xr ))
784         else:
785             raise ValueError("Quantile simulations has only to be Linear or NonLinear.")
786         #
787         if YfQ is None:
788             YfQ = Yr.reshape((-1,1))
789             if selfA._toStore("SampledStateForQuantiles"): EXr = Xr.reshape((-1,1))
790         else:
791             YfQ = numpy.hstack((YfQ,Yr.reshape((-1,1))))
792             if selfA._toStore("SampledStateForQuantiles"): EXr = numpy.hstack((EXr,Xr.reshape((-1,1))))
793     #
794     # Extraction des quantiles
795     YfQ.sort(axis=-1)
796     YQ = None
797     for quantile in selfA._parameters["Quantiles"]:
798         if not (0. <= float(quantile) <= 1.): continue
799         indice = int(nbsamples * float(quantile) - 1./nbsamples)
800         if YQ is None: YQ = YfQ[:,indice].reshape((-1,1))
801         else:          YQ = numpy.hstack((YQ,YfQ[:,indice].reshape((-1,1))))
802     if YQ is not None: # Liste non vide de quantiles
803         selfA.StoredVariables["SimulationQuantiles"].store( YQ )
804     if selfA._toStore("SampledStateForQuantiles"):
805         selfA.StoredVariables["SampledStateForQuantiles"].store( EXr )
806     #
807     return 0
808
809 # ==============================================================================
810 def ForceNumericBounds( __Bounds ):
811     "Force les bornes à être des valeurs numériques, sauf si globalement None"
812     # Conserve une valeur par défaut à None s'il n'y a pas de bornes
813     if __Bounds is None: return None
814     # Converti toutes les bornes individuelles None à +/- l'infini
815     __Bounds = numpy.asarray( __Bounds, dtype=float )
816     if len(__Bounds.shape) != 2 or min(__Bounds.shape) <= 0 or __Bounds.shape[1] != 2:
817         raise ValueError("Incorrectly shaped bounds data")
818     __Bounds[numpy.isnan(__Bounds[:,0]),0] = -sys.float_info.max
819     __Bounds[numpy.isnan(__Bounds[:,1]),1] =  sys.float_info.max
820     return __Bounds
821
822 # ==============================================================================
823 def RecentredBounds( __Bounds, __Center, __Scale = None):
824     "Recentre les bornes autour de 0, sauf si globalement None"
825     # Conserve une valeur par défaut à None s'il n'y a pas de bornes
826     if __Bounds is None: return None
827     if __Scale is None:
828         # Recentre les valeurs numériques de bornes
829         return ForceNumericBounds( __Bounds ) - numpy.ravel( __Center ).reshape((-1,1))
830     else:
831         # Recentre les valeurs numériques de bornes et change l'échelle par une matrice
832         return __Scale @ (ForceNumericBounds( __Bounds ) - numpy.ravel( __Center ).reshape((-1,1)))
833
834 # ==============================================================================
835 def ApplyBounds( __Vector, __Bounds, __newClip = True):
836     "Applique des bornes numériques à un point"
837     # Conserve une valeur par défaut s'il n'y a pas de bornes
838     if __Bounds is None: return __Vector
839     #
840     if not isinstance(__Vector, numpy.ndarray): # Is an array
841         raise ValueError("Incorrect array definition of vector data")
842     if not isinstance(__Bounds, numpy.ndarray): # Is an array
843         raise ValueError("Incorrect array definition of bounds data")
844     if 2*__Vector.size != __Bounds.size: # Is a 2 column array of vector lenght
845         raise ValueError("Incorrect bounds number (%i) to be applied for this vector (of size %i)"%(__Bounds.size,__Vector.size))
846     if len(__Bounds.shape) != 2 or min(__Bounds.shape) <= 0 or __Bounds.shape[1] != 2:
847         raise ValueError("Incorrectly shaped bounds data")
848     #
849     if __newClip:
850         __Vector = __Vector.clip(
851             __Bounds[:,0].reshape(__Vector.shape),
852             __Bounds[:,1].reshape(__Vector.shape),
853             )
854     else:
855         __Vector = numpy.max(numpy.hstack((__Vector.reshape((-1,1)),numpy.asmatrix(__Bounds)[:,0])),axis=1)
856         __Vector = numpy.min(numpy.hstack((__Vector.reshape((-1,1)),numpy.asmatrix(__Bounds)[:,1])),axis=1)
857         __Vector = numpy.asarray(__Vector)
858     #
859     return __Vector
860
861 # ==============================================================================
862 def Apply3DVarRecentringOnEnsemble(__EnXn, __EnXf, __Ynpu, __HO, __R, __B, __Betaf):
863     "Recentre l'ensemble Xn autour de l'analyse 3DVAR"
864     #
865     Xf = EnsembleMean( __EnXf )
866     Pf = Covariance( asCovariance=EnsembleErrorCovariance(__EnXf) )
867     Pf = (1 - __Betaf) * __B.asfullmatrix(Xf.size) + __Betaf * Pf
868     #
869     selfB = PartialAlgorithm("3DVAR")
870     selfB._parameters["Minimizer"] = "LBFGSB"
871     selfB._parameters["MaximumNumberOfIterations"] = 15000
872     selfB._parameters["CostDecrementTolerance"] = 1.e-7
873     selfB._parameters["ProjectedGradientTolerance"] = -1
874     selfB._parameters["GradientNormTolerance"] = 1.e-05
875     selfB._parameters["StoreInternalVariables"] = False
876     selfB._parameters["optiprint"] = -1
877     selfB._parameters["optdisp"] = 0
878     selfB._parameters["Bounds"] = None
879     selfB._parameters["InitializationPoint"] = Xf
880     from daAlgorithms.Atoms import std3dvar
881     std3dvar.std3dvar(selfB, Xf, __Ynpu, None, __HO, None, __R, Pf)
882     Xa = selfB.get("Analysis")[-1].reshape((-1,1))
883     del selfB
884     #
885     return Xa + EnsembleOfAnomalies( __EnXn )
886
887 # ==============================================================================
888 def multiXOsteps(selfA, Xb, Y, U, HO, EM, CM, R, B, Q, oneCycle,
889         __CovForecast = False, __LinEvolution = False,
890         ):
891     """
892     Prévision multi-pas avec une correction par pas (multi-méthodes)
893     """
894     #
895     # Initialisation
896     # --------------
897     if selfA._parameters["EstimationOf"] == "State":
898         if len(selfA.StoredVariables["Analysis"])==0 or not selfA._parameters["nextStep"]:
899             Xn = numpy.asarray(Xb)
900             if __CovForecast: Pn = B
901             selfA.StoredVariables["Analysis"].store( Xn )
902             if selfA._toStore("APosterioriCovariance"):
903                 if hasattr(B,"asfullmatrix"):
904                     selfA.StoredVariables["APosterioriCovariance"].store( B.asfullmatrix(Xn.size) )
905                 else:
906                     selfA.StoredVariables["APosterioriCovariance"].store( B )
907             selfA._setInternalState("seed", numpy.random.get_state())
908         elif selfA._parameters["nextStep"]:
909             Xn = selfA._getInternalState("Xn")
910             if __CovForecast: Pn = selfA._getInternalState("Pn")
911     else:
912         Xn = numpy.asarray(Xb)
913         if __CovForecast: Pn = B
914     #
915     if hasattr(Y,"stepnumber"):
916         duration = Y.stepnumber()
917     else:
918         duration = 2
919     #
920     # Multi-steps
921     # -----------
922     for step in range(duration-1):
923         selfA.StoredVariables["CurrentStepNumber"].store( len(selfA.StoredVariables["Analysis"]) )
924         #
925         if hasattr(Y,"store"):
926             Ynpu = numpy.asarray( Y[step+1] ).reshape((-1,1))
927         else:
928             Ynpu = numpy.asarray( Y ).reshape((-1,1))
929         #
930         if U is not None:
931             if hasattr(U,"store") and len(U)>1:
932                 Un = numpy.asarray( U[step] ).reshape((-1,1))
933             elif hasattr(U,"store") and len(U)==1:
934                 Un = numpy.asarray( U[0] ).reshape((-1,1))
935             else:
936                 Un = numpy.asarray( U ).reshape((-1,1))
937         else:
938             Un = None
939         #
940         # Predict (Time Update)
941         # ---------------------
942         if selfA._parameters["EstimationOf"] == "State":
943             if __CovForecast or __LinEvolution:
944                 Mt = EM["Tangent"].asMatrix(Xn)
945                 Mt = Mt.reshape(Xn.size,Xn.size) # ADAO & check shape
946             if __CovForecast:
947                 Ma = EM["Adjoint"].asMatrix(Xn)
948                 Ma = Ma.reshape(Xn.size,Xn.size) # ADAO & check shape
949                 Pn_predicted = Q + Mt @ (Pn @ Ma)
950             if __LinEvolution:
951                 Xn_predicted = Mt @ Xn
952             else:
953                 M  = EM["Direct"].appliedControledFormTo
954                 Xn_predicted = M( (Xn, Un) )
955             if CM is not None and "Tangent" in CM and Un is not None: # Attention : si Cm est aussi dans M, doublon !
956                 Cm = CM["Tangent"].asMatrix(Xn_predicted)
957                 Cm = Cm.reshape(Xn.size,Un.size) # ADAO & check shape
958                 Xn_predicted = Xn_predicted + (Cm @ Un).reshape((-1,1))
959         elif selfA._parameters["EstimationOf"] == "Parameters": # No forecast
960             # --- > Par principe, M = Id, Q = 0
961             Xn_predicted = Xn
962             if __CovForecast: Pn_predicted = Pn
963         Xn_predicted = numpy.asarray(Xn_predicted).reshape((-1,1))
964         if selfA._toStore("ForecastState"):
965             selfA.StoredVariables["ForecastState"].store( Xn_predicted )
966         if __CovForecast:
967             if hasattr(Pn_predicted,"asfullmatrix"):
968                 Pn_predicted = Pn_predicted.asfullmatrix(Xn.size)
969             else:
970                 Pn_predicted = numpy.asarray(Pn_predicted).reshape((Xn.size,Xn.size))
971             if selfA._toStore("ForecastCovariance"):
972                 selfA.StoredVariables["ForecastCovariance"].store( Pn_predicted )
973         #
974         # Correct (Measurement Update)
975         # ----------------------------
976         if __CovForecast:
977             oneCycle(selfA, Xn_predicted, Ynpu, Un, HO, CM, R, Pn_predicted, True)
978         else:
979             oneCycle(selfA, Xn_predicted, Ynpu, Un, HO, CM, R, B, True)
980         #
981         #--------------------------
982         Xn = selfA._getInternalState("Xn")
983         if __CovForecast: Pn = selfA._getInternalState("Pn")
984     #
985     return 0
986
987 # ==============================================================================
988 if __name__ == "__main__":
989     print('\n AUTODIAGNOSTIC\n')