8 DESCRIPTION OF DONLP2's COMMON BLOCKS

The common blocks used by DONLP2 are listed here: (may be copied from O8COMM.INC and its further included *.INC files).

       IMPLICIT NONE
       INCLUDE 'O8PARA.INC'
c accinf contains information about intermediate results useful for a
c postmortem analysis .
C****
C    ACCINF A C C U M U L A T E D I N F O R M A T I O N
C    ON ITERATION SEQUENCE
C    1: STEP-NR
C    2: F(X-K)    CURRENT VALUE OF OBJECTIVE (ZERO IN FEASIBILITY IMPROVEMENT
C                        PHASE (-1) )
C    3: SCF         INTERNAL SCALING OF OBJECTIVE (ZERO IN PHASE -1)
C    4: PSI          THE WEIGHTED PENALTY-TERM
C    5: UPSI       THE UNWEIGHTED PENALTY-TERM (L1-NORM OF CONSTRAINT VECTOR)
C    6: DEL.K.1  BOUND FOR CURRENTLY ACTIVE CONSTRAINTS
C    7: B2N0       L2-NORM OF PROJECTED GRADIENT, BASED ON CONSTRAINTS IN LEVEL DELMIN
C                        AND BELOW
C    8: B2N         L2-NORM OF PROJECTED GRADIENT BASED ON DEL.K.1
C    9: NR           NUMBER OF BINDING CONSTRAINTS
C  10: SING       IF 1, THE BINDING CONSTRAINTS DON'T SATISFY THE REGULARITY CONDITION
C  11: UMIN      INFINITY NORM OF NEGATIVE PART OF MULTIPLIER
C  12: -------------
C  13: COND.R  CONDITION NUMBER OF DIAGONAL PART OF QR-DECOMP. OF NORMALIZED
C                         GRADIENTS OF BINDING CONSTRAINTS
C  14: COND.H  CONDITION NUMBER OF DIAGONAL OF CHOLESKY-FACTOR OF UPDATED FULL HESSIAN
C  15: SCF0        THE RELATIVE DAMPING OF TANGENTIAL COMPONENT IF UPSI?TAU0/2
C  16: XNORM   L2-NORM OF X
C  17: DNORM   L2-NORM OF D (CORRECTION FROM QP -SUBPROBLEM, UNSCALED)
C  18: PHASE      -1 : INFEASIBILITY IMPROVEMENT PHASE, 0: INITIAL OPTIMIZATION
C                         1: BINDING CONSTRAINTS UNCHANGED , 2: D SMALL, MARATOS CORRECTION
C                          IN USE
C  19: C.K           NUMBER OF DECREASES OF PENALTY WEIGHTS
C  20: WMAX     INFINITY NORM OF WEIGHTS
C  21: SIG.K        STEPSIZE FROM UNIDIMENSIONAL MINIMIZATION (BACKTRACKING)
C  22: CFINCR    NUMBER OF OBJECTIVE EVALUATIONS FOR STEPSIZE-ALGORITHM
C  23: DIRDER    DIRECTIONAL DERIVATIVE OF PENALTY-FUNCTION ALONG D (SCALED)
C  24: DSCAL     SCALING FACTOR FOR D
C  25: COSPHI   COS OF ARC BETWEEN D AND D.PREVIOUS. IF LARGER THETA , SIG LARGER THAN
C                         ONE (UP TO SIGLA) IS TRIED
C  26: VIOLIS(0) NUMBER OF CONSTRAINTS NOT BINDING AT X BUT HIT DURING LINE SEARCH
C  27:                   TYPE OF UPDATE FOR H: 1 NORMAL P&M-BFGS-UPDATE, 0 UPDATE SUPPRESSED,
C                        -1 RESTART WITH SCALED UNIT MATRIX , 2 STANDARD BFGS, 3 BFGS MODIFIED
C                         BY POWELLS DEVICE
C  28: NY.K/TK MODIFICATION FACTOR FOR DAMPING THE PROJECTOR IN BFGS/PANTOJA
C                        MAYNE-TERM RESPECTIVELY
C  29: 1-MY.K/XSIK MODIFICATION FACTOR FOR DAMPING THE QUASI-NEWTON-RELATION IN BFGS 
C                        FOR UNMODIFIED BFGS NY.K SHOULD BE LARGER THAN UPDMY0 (NEAR ONE) 
C                        AND 1-MY.K EQUAL ONE./PANTOJA-MAYNE TERM RESPECTIVELYC**** 
C  30: QPTERM 0, IF SING=-1, TERMINATION INDICATOR OF QP-SOLVER OTHERWISE 
C                        1: SUCCESSFUL, -1: TAU BECOMES LARGER THAN TAUQP WITHOUT SLACK
C                        VARIABLES BECOMING SUFFICIENTLY SMALL . 
C                        -2: INFEASIBLE QP-PROBLEM (THEORETICALLY IMPOSSIBLE) 
C  31: TAUQP: WEIGHT OF SLACK-VARIABLES IN QP-SOLVER 
C  32: INFEAS  L1-NORM OF SLACK-VARIABLES IN QP-SOLVER 

c optite is the termination parameter. 
c optite ranges from -10 to 7, a negative value indicating some irregular 
c event. optite+11 corresponds one the following messages: 


c         'CONSTRAINT EVALUATION RETURNS ERROR WITH CURRENT POINT', 
c         'OBJECTIVE EVALUATION RETURNS ERROR WITH CURRENT POINT', 
c         'QPSOLVER: WORKING SET SINGULAR IN DUAL EXTENDED QP ', 
c         'QPSOLVER: EXTENDED QP-PROBLEM SEEMINGLY INFEASIBLE ', 
c         'QPSOLVER: NO DESCENT FOR INFEAS FROM QP FOR TAU=TAU.MAX', 
c         'QPSOLVER: ON EXIT CORRECTION SMALL, INFEASIBLE POINT', 
c         'STEPSIZESELECTION: COMPUTED D FROM QP NOT A DIR. OF DESCENT', 
c         'MORE THAN MAXIT ITERATION STEPS', 
c         'STEPSIZESELECTION: NO ACCEPTABLE STEPSIZE IN [SIGSM,SIGLA]', 
c         'SMALL CORRECTION FROM QP, INFEASIBLE POINT', 
c         'KT-CONDITIONS SATISFIED, NO FURTHER CORRECTION COMPUTED', 
c         'COMPUTED CORRECTION SMALL, REGULAR CASE ', 
c         'STEPSIZESELECTION: X ALMOST FEASIBLE, DIR. DERIV. VERY SMALL', 
c         'KT-CONDITIONS (RELAXED) SATISFIED, SINGULAR POINT', 
c         'VERY SLOW PRIMAL PROGRESS, SINGULAR OR ILLCONDITONED PROBLEM', 
c         'MORE THAN NRESET SMALL CORRECTIONS IN X ', 
c         'CORRECTION FROM QP VERY SMALL, ALMOST FEASIBLE, SINGULAR ', 
c          'NUMSM SMALL DIFFERENCES IN PENALTY FUNCTION,TERMINATE' 

c    a "good" termination is case 8,9,10.(value 0,1,2) 
c    In the cases 12-15 the final 
c    precision might be quite poor, since in the singular case DONLP2 
c    reduces its requirements automatically by multiplying epsx by 100. 

c runtim :  the processes cpu-time 
c phase :      PHASE -1 : INFEASIBILITY IMPROVEMENT PHASE, 0: INITIAL OPTIMIZATION 
c                 1: BINDING CONSTRAINTS UNCHANGED , 2: D SMALL, MARATOS CORRECTION 
c                 IN USE 
c
         REAL RUNTIM 
         DOUBLE PRECISION ACCINF 
         REAL OPTITE 
         INTEGER ITSTEP,PHASE 
         COMMON/O8ITIN/OPTITE,ITSTEP,PHASE,RUNTIM, 
       F                  ACCINF(0:MAXIT,32) 
c  component "0" introduced to overcome errors of some optimizing compilers 
c  o8xdat contains information pertaining the value of current x's 
c  x = current point, x0= previous point, x1 = tentative new point 
c  xst = starting point
c  same convention for names upsi , psi , fx , sig, -norm 

c  fx = function value , psi =scaled penalty-term 
c  upsi = unscaled penalty term , phi=scf*f+psi 
c  xmin = during stepsize increase current best point 
c  (with values sigmin, upsim, psimin, fmin, phimin, resmin ) 
c  sig : current step-size 
c  dirder = directional derivative of phi along d 
c  cosphi=cos(arc(d,d0)) 
c  dscal: scaling of d (if too long, resulting in a too large change of x) 
c  d = direction , dd = second order correction 
c  difx = x.new-x.old 
c  b2n = the projected gradient (in a system transformed by the 
c  cholesky-factor of H) , b2n0 = same based on a reduced minimal 
c  working set ,
         DOUBLE PRECISION UPSI,UPSI0,UPSI1,UPSIST,PSI,PSI0, 
       F         PSI1,PSIST,PSIMIN, F PHI,PHI0,PHI1,PHIMIN,FX,FX0,FX1, 
       F         FXST,FMIN,B2N,B2N0,XNORM,X0NORM,SIG0,DSCAL,DNORM,D0NORM
         DOUBLE PRECISION SIG,SIGMIN,DIRDER,COSPHI,UPSIM 
         DOUBLE PRECISION X,X0,X1,XMIN,D,D0,DD,DIFX,RESMIN
         COMMON/O8XDAT/X(NX),X0(NX),X1(NX),XMIN(NX),RESMIN(NRESM), 
       F                                  D(NX),D0(NX),DD(NX),DIFX(NX), 
       F                                  XNORM,X0NORM,DNORM,D0NORM,SIG,SIG0,SIGMIN,DSCAL, 
       F                                  UPSI,UPSI0,UPSI1,UPSIST,UPSIM, 
       F                                  PSI,PSI0,PSI1,PSIST,PSIMIN, 
       F                                  PHI,PHI0,PHI1,PHIMIN, 
       F                                  FX,FX0,FX1,FXST,FMIN, 
       F                                  B2N,B2N0,DIRDER,COSPHI 
c  gradf= gradient(f), gfn= norm of gradient of f , 
c  qgf = internal transform of gradf 
c  gphi1 and gphi0 : gradient of lagrangian at new and old point 
c  gres : gradients of constraints, 
c  gresn(i)=max(1,norm of gradient(constraint.i)))
         DOUBLE PRECISION GRADF,GFN,QGF,GPHI0,GPHI1,GRES,GRESN
         COMMON/O8GRD/GRADF(NX),GFN,QGF(NX),GRES(NX,NRESM), 
      F                       GRESN(NRESM),GPHI0(NX),GPHI1(NX) 
C** 
c  information pertaining the qr-decomposition of the working set's 
c  gradients 
c  perm = current valid row permutation making it continuous 
c  perm1 = new permutation 
c  colno = column permutation using column pivoting 
c  rank = rank of working set 
c  qr : matrix containing information of the householder-decomposition 
c  diag : diagonal of the r-part 
c  betaq : factor for reflection normals 
c  cscal: the column scaling 
c  colle: column lenght for the qp.problem (extended columns)
         INTEGER PERM,PERM1,COLNO,RANK 
         DOUBLE PRECISION QR,BETAQ,DIAG,CSCAL,COLLE
         COMMON/O8RDAT/QR(NX,NRESM),BETAQ(NRESM),DIAG(NRESM),CSCAL(NRESM), 
      F                                 COLLE(NRESM), 
      F                                 COLNO(2*NRESM),PERM(NX),
      F                                 PERM1(NX),RANK 
C  COLNO ALSO USED IN O8QPSO WITH DOUBLE LENGTH ! 
c   information pertaining the constraints 
c   val(i)=true, if gradient of constraint i has been evaluated at 
c   the new point 
c   i=0 corresponds to f 
c   gconst(i)=.true. if constraint i is linear. 
c   such a constraint gradient is evaluated once only 
c   gunit : see descrition above 
c   llow and lup are computed by -"tt DONLP2"" using gunit 
c   llow(i)=.true. if there is a finite lower bound for x(i) 
c   lup(i)=.true. if there is a finite upper bound for x(i) 
c   sort: internal order of evaluation of constraints
        LOGICAL VAL,GCONST,LLOW,LUP 
        INTEGER GUNIT 
        COMMON/O8GRI/VAL(0:NRESM),GCONST(0:NRESM),GUNIT(3,0:NRESM), 
      F                                LLOW(NX),LUP(NX) 
c   intakt= .true. : give output to prou-channel and to console 
c   inx= presently not used 
c   std=.true. use standard initialization 
c   this may be superseded by statements in setup 
c   te0 = .true. : give one-line-output for every step on console 
c   te1 = .true. : give post-mortem-dump of accumulated information in accinf 
c   te2 = .true. : give detailed protocol of iteration 
c   te3 = .true. : also print the gradients and the approximated hessian 
c   singul = .true. : the working set is singular 
c   ident = .true. : we try a modified step with the same point 
c   eqres = .true. : the working set stays fixed 
c   silent= .true. : give no output at all (user must evaluate information 
c   in the common-blocks himself) 
c   analyt = .true. : assume full precision in the gradients 
c   otherwise {\tt DONLP2} assumes epsdif relative precision in the gradients and 
c   automatically reduces its precision requirements 
c   cold = .true. : initialize weights, hessian and so on, otherwise use them 
c   unchanged in the first step.
            LOGICAL INTAKT,INX,STD,TE0,TE1,TE2,TE3,SINGUL 
            LOGICAL IDENT,EQRES,SILENT,ANALYT,COLD 
            COMMON/O8STPA/INTAKT,INX,STD,TE0,TE1,TE2,TE3,SINGUL, 
          F                      IDENT,EQRES,SILENT,ANALYT,COLD 
c   information pertaining the quasi newton update: 
c   A :the quasi-newton-update in cholesky-factors: upper triangle :the new one 
c   lower triangle + diag0 : the previous one 
c   matsc : scaling for identity to be used in case of restart 
c   (computed internally)
           DOUBLE PRECISION A,DIAG0,MATSC 
           COMMON/O8QN/A(NX,NX),DIAG0(NX), 
         F                             ,MATSC 
c   information pertaining the working set: 
c   violis: list of constraints hit during current linesearch 
c   alist: indices of the working set 
c   bind : constraints binding at x (violated or at zero) bind(i)=1 
c            (=0 otherwise) 
c   bind0 : the same for x0 
c   sort : internal evaluaution order for constraints
            INTEGER VIOLIS,ALIST,BIND,BIND0,SORT
            COMMON /O8RESI/ BIND(NRESM),BIND0(NRESM),VIOLIS(0:NSTEP*NRESM), 
          F                                   ALIST(0:NRESM),SORT(NRESM) 
c   u = mulipliers , u0 = old multipliers , w = weights , 
c   w1 = tentative new weights , res = constraint values , res0 = previous 
c   constraint values, res1 = same at tentative point , resst = 
c   constraint values at starting point , scf = current scaling of f 
c   scf0 = damping factor for tangential part of d if infeasibility ? tau0/2 
c   yu internal varable, stores multipliers for working set , 
c   slack = slack values in qp-problem , infeas = one-norm of slack , 
c   work = work array
            DOUBLE PRECISION U,U0,W,W1,RES,RES0,RES1,RESST,SCF,SCF0, 
          F                                    YU,SLACK,INFEAS,WORK
            COMMON/O8RESD/RES(NRESM),RES0(NRESM),RES1(NRESM),RESST(NRESM), 
          F                               U(NRESM),U0(NRESM),W(NRESM),W1(NRESM),WORK(0:NRESM), 
          F                               YU(NRESM),SLACK(NRESM),SCF,SCF0,INFEAS 
c   n=dim(x), nh=dim(h), ng=dim(g), nres=nh+ng, nr=dim(working set) 
c   (may be larger than n)
            INTEGER N,NR,NRES,NH,NG 
            COMMON/O8DIM/N,NH,NG,NR,NRES 
c   epsmac = machine precision, tolmac= smallest positive machine number 
c   (roughly approximated) both computed automatically by DONLP2 
c   deldif is a suitable difference stepsize for the automatic 
c   numerical differentiation routine supplied with this distribution
            DOUBLE PRECISION EPSMAC,TOLMAC,DELDIF 
            COMMON/O8MPAR/EPSMAC,TOLMAC,DELDIF 
c   iterma= maximum number of steps allowed. (30*n a good value from 
c   experience) 
c   del = current delta, del0 = maximum of delta , del01=del0/10 , 
c   delmin: constraint are satisfied if |h_i(x)|  <= delmin , g_j(x) >= -delmin 
c   respectively 
c   tau0 : upper bound for the unscaled penalty-term, i.e. the one-norm of the 
c   constraint violation 
c   tau gives a weight between descent for f and infeasibility. 
c   the smaller tau, the more may f be scaled down. 
c   also used as additive increase for the penalty-weights 
c   ny increase factor for penalty-weights when computed from the multipliers 
c   smalld : use second order correction if ||d|| <= smalld and the working 
c   set stays fixed 
c   smallw: multipliers absolutely smaller than smallw are considered as zero 
c   rho : working set considered singular if condition number of r-part 
c   of its qr-decomposition (after scaling the gradients to one) is larger than 
c   1/rho. 
c   rho1: if condition number of quasi-newton-update of hessian becomes larger 
c   than 1/rho1, a restart with a scaled unit matrix is performed 
c   eta: level required for descent if weights are to be diminished, 
c   computed automatically 
c   epsx : termination parameter. successful termination is indicated if 
c   (roughly speaking) the kuhn-tucker-criteria are satisfied within an error of 
c   epsx 
c   epsphi: variable entering the termination criterion for small progress 
c   in the penalty function, see under termination criteria 
c   c1d: inactivation of an inequality contraint occurs if a multiplier*c1d 
c   is less than - a scaled projected gradient at the current point
c  scfmax: the internal automatic scaling of f is bounded by 1/scfmax from below 
c  and scfmax from above. this scaling aims at one multiplier at least being 
c  in the order of magnitude of one. 
c  tauqp: the penalty-parameter in the qp-problem 
c  taufac: increase factor for tauqp 
c  taumax: maximum allowable tauqp 
c  updmy0 : presently not used 
c
            INTEGER ITERMA,IFILL1 
            DOUBLE PRECISION DEL,DEL0,DEL01,DELMIN,TAU0,TAU,NY 
            DOUBLE PRECISION SMALLD,SMALLW,RHO,RHO1,ETA,EPSX,EPSPHI,C1D, 
          F                                      SCFMAX,UPDMY0,TAUQP,TAUFAC,TAUMAX
            COMMON/O8PAR/DEL0,DEL01,DEL,DELMIN,TAU0,TAU, 
          F                                 SMALLD,SMALLW,RHO,RHO1,ETA,NY,EPSX,EPSPHI, 
          F                                 C1D,SCFMAX,TAUQP,TAUFAC,TAUMAX,UPDMY0, 
          F                                 ITERMA,IFILL1 
c   the parameters of the stepsize-algorithm: 
c   alpha=smallest reduction factor, beta=feasible increase factor 
c   for change of x-norm adding the correction d. theta: if theta!=cos(arc(d,d0)) 
c   then stepsizes > 1 are tried for d, where d0=previous correction. 
c   sigsm= smallest stepsize allowed, sigla = largest stepsize allowed, 
c   delta = the armijo-parameter, stptrm = termination indicator 
c   for stepsize-algorithm, delta1= armijo-parameter for reduction 
c   of infeasibility, stmaxl= min(sigla, maximum stepsize for which the 
c   point on the ray x+sigma*d, projected to the bounds, changes)
           DOUBLE PRECISION ALPHA,BETA,THETA,SIGSM,SIGLA,DELTA,STPTRM 
           DOUBLE PRECISION DELTA1,STMAXL 
           COMMON/O8STEP/ALPHA,BETA,THETA,SIGSM,SIGLA,DELTA,STPTRM, 
         F                                  DELTA1,STMAXL 
c   evaluation counters cf= number of calls to f, cgf= number of calls to grad(f) 
c   cfincr= number of function values used in the current stepsize-selection 
c   cres(i)= number of calls to individual constraints (bounds not counted) 
c   cgres = same for the gradients. 
c   i=1,..,nh corresponds h and i=nh+1,...,nh+ng g.
           INTEGER CF,CGF,CFINCR,CRES,CGRES 
           COMMON/O8CNT/CF,CGF,CFINCR,CRES(NRESM),CGRES(NRESM) 
c   FFUERR error indicator for evaluating f 
c   CFUERR(i) error indicator for evaluating a constraining function 
c   i=1,..,nh corresponds to h and i=nh+1,...,nh+ng to g. 
c   these are all initialized by .FALSE.
           LOGICAL FFUERR,CFUERR(NRESM) 
           COMMON/O8ERR/FFUERR,CFUERR 
c   clow = number of times the weights have been dimished 
c   lastdw = last step number with weights diminished 
c   lastup = last step number with weights increased 
c   lastch = last step number with weights changed
           INTEGER CLOW,LASTDW,LASTUP,LASTCH 
           COMMON/O8WEI/CLOW,LASTDW,LASTUP,LASTCH 
c   problem identifier
           CHARACTER*40 NAME 
           COMMON/O8ID/NAME 
c   precision given by the u s e r, computing gradients numerically. 
c   not used if analyt=.true.
           DOUBLE PRECISION EPSDIF
           COMMON/O8DER/EPSDIF 
c   channel numbers for i/o of intermediate and final results
           INTEGER PROU,MEU 
           COMMON/O8IO/PROU,MEU 
c   ug, og and delfac (lower and upper bounds for x(i) and a special 
c   scaling weight for detection of binding bounds) are computed by 
c   DONLP2 itself, using GUNIT and evaluation of the bounds using x=0
           DOUBLE PRECISION UG(NX),OG(NX),DELFAC(NRESM) 
           COMMON/O8BD/UG,OG,DELFAC 
c   reset the quasi-newton-update after "nreset" small or unsuccessful steps 
c  numsm the number of allowed consecutive small changes in the penalty function
           INTEGER NRESET,NUMSM 
           COMMON/O8RST/NRESET,NUMSM 
c   here we store the starting point , for later printing only
           DOUBLE PRECISION XST(NX) 
           COMMON/O8STV/XST 
           INTRINSIC SQRT,EXP,LOG,MAX,MIN,ABS

The following common-area serves as an interface for transmitting and receiving information from a black-box evaluation routine which a user might wish to use:

C**** IF BLOC=.TRUE. THEN IT IS ASSUMED THAT FUNCTIONEVALUATION TAKES PLACE 
C**** AT ONCE IN AN EXTERNAL ROUTINE AND 
C**** THAT USER.EVAL HAS BEEN CALLED (WHICH IN TURN CALLS 
C**** EVAL.EXTERN(MODE) BEFORE CALLING FOR EVALUATION OF FUNCTIONS 
C**** THE LATTER THEN SIMPLY CONSISTS IN COPYING DATA FROM COMMON O8FUEXT 
C**** TO DONLP2'S OWN DATA-AREAS. 
C**** CORR IS SET TO TRUE BY DONLP2, IF THE INITIAL X DOES NOT SATISFY 
C**** THE BOUND CONSTRAINTS. X IS MODIFIED IN THIS CASE 
C**** DIFFTYPE=1,2,3 NUMERICAL DIFFERENTIATION BY THE ORDINARY FORWARD 
C**** DIFFERENCES, BY CENTRAL DIFFERENCES OR BY RICHARDSON-EXTRAPOLATION 
C**** OF ORDER 6, REQUIRING N, 2N , 6N ADDITIONAL FUNCTION EVALUATIONS 
C**** RESPECTIVELY 
C**** EPSFCN IS THE ASSUMED PRECISION OF THE FUNCTION EVALUATION, TO BE 
C**** SET BY THE USER 
C**** TAUBND: AMOUNT BY WHICH BOUND CONSTRAINTS MAY BE VIOLATED DURING 
C**** FINITE DIFFERENCING, SET BY THE USER
            LOGICAL SCALE,BLOC,VALID,CORR 
            INTEGER DIFFTYPE 
C**** XTR : THE CURRENT DESIGN X 
C**** FU : THE VECTOR OF FUNCTION VALUES IN THE ODER 0=OBJECTIVE FUNCTION, 
C**** 1,...,NH EQUALITY CONSTRAINTS, NH+1,NH+NG INEQUALITY CONSTRAINTS. 
C**** SIMILARLY FUGRAD(.,I) THE CORRESPONDING GRADIENTS 
C**** WITH BLOC=.TRUE. AND ANALYT=.FALSE. ONLY FU MUST BE COMPUTED, USING XTR 
C**** AS INPUT
            DOUBLE PRECISION XTR(NX),XSC(NX),FU(0:NRESM),FUGRAD(NX,0:NRESM), 
          F                                      FUD(0:NRESM,1:6),EPSFCN,TAUBND
            COMMON/O8FUEXT/XTR,XSC,FU,FUGRAD,FUD,EPSFCN,TAUBND
            COMMON/O8FUPAR/SCALE,BLOC,VALID,CORR,DIFFTYPE