mpqldec.f90 33.7 KB
Newer Older
1
2
3
4
5
6
!> \file
!! QL decompostion.
!!
!! \author Claus Kleinwort, DESY, 2015 (Claus.Kleinwort@desy.de)
!!
!! \copyright
7
!! Copyright (c) 2015-2021 Deutsches Elektronen-Synchroton,
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
!! Member of the Helmholtz Association, (DESY), HAMBURG, GERMANY \n\n
!! This library is free software; you can redistribute it and/or modify
!! it under the terms of the GNU Library General Public License as
!! published by the Free Software Foundation; either version 2 of the
!! License, or (at your option) any later version. \n\n
!! This library is distributed in the hope that it will be useful,
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!! GNU Library General Public License for more details. \n\n
!! You should have received a copy of the GNU Library General Public
!! License along with this program (see the file COPYING.LIB for more
!! details); if not, write to the Free Software Foundation, Inc.,
!! 675 Mass Ave, Cambridge, MA 02139, USA.
!!
!! QL decomposition of constraints matrix by Householder transformations
23
!! for solution by elimination. Optionally split into disjoint blocks.
24
25
26
27
28
29
30
31
32
!!

!> QL data.
MODULE mpqldec
    USE mpdef
    IMPLICIT NONE

    INTEGER(mpi) :: npar   !< number of parameters
    INTEGER(mpi) :: ncon   !< number of constraints
33
    INTEGER(mpi) :: nblock !< number of blocks
34
    INTEGER(mpl) :: matsize !< size of contraints matrix 
35
    INTEGER(mpi) :: iblock !< active block
36
    INTEGER(mpi) :: monpg  !< flag for progress monitoring
37
38
39
40
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: matV  !< unit normals (v_i) of Householder reflectors
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: vecVk !< secondary diagonal of matV (including last element)
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: matL  !< lower diagonal matrix L
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: vecN  !< normal vector
41
42
    INTEGER(mpi), DIMENSION(:), ALLOCATABLE :: nparBlock !< number of parameters in block
    INTEGER(mpi), DIMENSION(:), ALLOCATABLE :: ioffBlock !< block offset (1. constraint -1)
43
44
45
    INTEGER(mpl), DIMENSION(:), ALLOCATABLE :: ioffRow !< row offsets in matV (for constrint block)
    INTEGER(mpi), DIMENSION(:), ALLOCATABLE :: ioffPar !< parameter number offsets for matV ( " )
    INTEGER(mpi), DIMENSION(:,:), ALLOCATABLE :: irangeParNZ !< range for non zero part (except vecVk)
46

47
48
49
50
51
52
END MODULE mpqldec

!> Initialize QL decomposition.
!!
!! \param [in]     n  number of rows (parameters)
!! \param [in]     m  number of columns (constraints)
53
!! \param [in]     l  number of disjoint blocks
54
!! \param [in]     s  size of constrints matrix (compressed or n*m)
55
!! \param [in]     k  flag for progress monitoring
56
!!
57
SUBROUTINE qlini(n,m,l,s,k)
58
59
60
61
62
63
64
65
    USE mpqldec
    USE mpdalc

    IMPLICIT NONE
    INTEGER(mpl) :: length

    INTEGER(mpi), INTENT(IN)          :: n
    INTEGER(mpi), INTENT(IN)          :: m
66
    INTEGER(mpi), INTENT(IN)          :: l
67
    INTEGER(mpl), INTENT(IN)          :: s
68
    INTEGER(mpi), INTENT(IN)          :: k
69
70
71

    npar=n
    ncon=m
72
    nblock=l
73
    matsize=s
74
    iblock=1
75
    monpg=k
76
    ! allocate 
77
78
    length=matsize
    !print *, ' full length (V)', length
79
    CALL mpalloc(matV,length,'QLDEC: V')
80
    matV=0.
81
    length=INT(ncon,mpl)*INT(ncon,mpl)
82
    CALL mpalloc(matL,length,'QLDEC: L')
83
    matL=0.
84
    length=npar
85
    CALL mpalloc(vecN,length,'QLDEC: v') 
86
87
88
89
90
91
92
93
94
    length=ncon
    CALL mpalloc(vecVk,length,'QLDEC: sec. diag(V)')
    vecVk=0. 
    CALL mpalloc(ioffPar,length,'QLDEC: parameter offsets (V)')
    ioffPar=0
    CALL mpalloc(irangeParNZ,2_mpl,length,'QLDEC: parameter non zero range (V)')
    length=ncon+1
    CALL mpalloc(ioffRow,length,'QLDEC: row offsets (V)') 
    ioffRow=0   
95
96
97
98
99
100
    length=nblock   
    CALL mpalloc(nparBlock,length,'QLDEC: npar in block')
    nparBlock=0
    length=nblock+1   
    CALL mpalloc(ioffBlock,length,'QLDEC: ioff for block')
    ioffBlock=0   
101
102
103
END SUBROUTINE qlini

!                                                 141217 C. Kleinwort, DESY-FH1
104
!> QL decomposition (as single block).
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
!!
!! QL decomposition with Householder transformations.
!! Decompose N-By-M matrix A into orthogonal N-by-N matrix Q and a
!! N-by-M matrix containing zeros except for a lower triangular
!! M-by-M matrix L (at the bottom):
!!
!!              | 0 |
!!      A = Q * |   |
!!              | L |
!!
!! The decomposition is stored in a N-by-M matrix matV containing the unit
!! normal vectors v_i of the hyperplanes (Householder reflectors) defining Q.
!! The lower triangular matrix L is stored in the M-by-M matrix matL.
!!
!! \param [in]     a  Npar-by-Ncon matrix
!!
SUBROUTINE qldec(a)
    USE mpqldec
    USE mpdalc

    ! cost[dot ops] ~= Npar*Ncon*Ncon

    IMPLICIT NONE
    INTEGER(mpi) :: i
    INTEGER(mpl) :: ioff1
    INTEGER(mpl) :: ioff2
    INTEGER(mpl) :: ioff3
    INTEGER(mpi) :: k
    INTEGER(mpi) :: kn
    INTEGER(mpl) :: length
    REAL(mpd) :: nrm
    REAL(mpd) :: sp

    REAL(mpd), INTENT(IN)             :: a(*)

140
141
    ! prepare
    vecVk=0._mpd
142
    length=INT(npar,mpl)*INT(ncon,mpl)
143
144
    matV=a(1:length)
    matL=0.0_mpd
145
146
147
    ! implemented as single block
    nblock=1
    nparBlock(1)=npar
148
149
150
151
    ioffBlock(2)=ncon
    DO k=1,ncon
        ioffRow(k+1)=ioffRow(k)+npar
    END DO   
152
153
154

    ! Householder procedure
    DO k=ncon,1,-1
155
156
        ! monitoring
        IF(monpg>0) CALL monpgs(ncon+1-k)
157
158
        kn=npar+k-ncon
        ! column offset
159
        ioff1=INT(k-1,mpl)*INT(npar,mpl)
160
        ! get column
161
162
        vecN(1:kn)=matV(ioff1+1:ioff1+kn)
        nrm = SQRT(dot_product(vecN(1:kn),vecN(1:kn)))
163
164
        IF (nrm == 0.0_mpd) CYCLE
        !
165
166
        IF (vecN(kn) >= 0.0_mpd) THEN
            vecN(kn)=vecN(kn)+nrm
167
        ELSE
168
            vecN(kn)=vecN(kn)-nrm
169
170
        END IF
        ! create normal vector
171
172
        nrm = SQRT(dot_product(vecN(1:kn),vecN(1:kn)))
        vecN(1:kn)=vecN(1:kn)/nrm
173
174
175
        ! transformation
        ioff2=0
        DO i=1,k
176
177
            sp=dot_product(vecN(1:kn),matV(ioff2+1:ioff2+kn))
            matV(ioff2+1:ioff2+kn)=matV(ioff2+1:ioff2+kn)-2.0_mpd*vecN(1:kn)*sp
178
179
180
            ioff2=ioff2+npar
        END DO
        ! store column of L
181
        ioff3=INT(k-1,mpl)*INT(ncon,mpl)
182
183
        matL(ioff3+k:ioff3+ncon)=matV(ioff1+kn:ioff1+npar)
        ! store normal vector
184
185
186
187
188
189
        matV(ioff1+1:ioff1+kn-1)=vecN(1:kn-1)
        matV(ioff1+kn:ioff1+npar)=0.0_mpd
        irangeParNZ(1,k)=1
        irangeParNZ(2,k)=kn-1
        ! store secondary diagonal 
        vecVk(k)=vecN(kn)
190
191
192
193
    END DO

END SUBROUTINE qldec

194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
!                                                 190312 C. Kleinwort, DESY-BELLE
!> QL decomposition (for disjoint block matrix).
!!
!! QL decomposition with Householder transformations.
!! Decompose N-By-M matrix A into orthogonal N-by-N matrix Q and a
!! N-by-M matrix containing zeros except for a lower triangular
!! M-by-M matrix L (at the bottom):
!!
!!              | 0 |
!!      A = Q * |   |
!!              | L |
!!
!! The decomposition is stored in a N-by-M matrix matV containing the unit
!! normal vectors v_i of the hyperplanes (Householder reflectors) defining Q.
!! The lower triangular matrix L is stored in the M-by-M matrix matL.
!!
210
!! \param [in]     a     block compressed Npar-by-Ncon matrix
211
!! \param [in]     bpar  2-by-NparBlock+1 matrix (with parameter block definition)
212
!! \param [in]     bcon  3-by-NconBlock+1 matrix (with constraint block definition)
213
!! \param [in]     cons  3-by-Ncons+1     matrix (with constraint definition)
214
!!
215
SUBROUTINE qldecb(a,bpar,bcon,cons)
216
217
218
219
220
221
222
    USE mpqldec
    USE mpdalc

    ! cost[dot ops] ~= Npar*Ncon*Ncon

    IMPLICIT NONE
    INTEGER(mpi) :: i
223
    INTEGER(mpi) :: ibcon
224
225
    INTEGER(mpi) :: iblast
    INTEGER(mpi) :: iboff
226
    INTEGER(mpi) :: ibpar
227
228
    INTEGER(mpi) :: ifirst
    INTEGER(mpi) :: ilast
229
    INTEGER(mpi) :: in
230
231
232
    INTEGER(mpl) :: ioff1
    INTEGER(mpl) :: ioff2
    INTEGER(mpl) :: ioff3
233
234
235
236
    INTEGER(mpi) :: iclast
    INTEGER(mpi) :: icoff
    INTEGER(mpi) :: iplast
    INTEGER(mpi) :: ipoff
237
238
239
240
    INTEGER(mpi) :: k
    INTEGER(mpi) :: k1
    INTEGER(mpi) :: kn
    INTEGER(mpi) :: ncb
241
    INTEGER(mpi) :: ncol
242
243
244
245
246
    INTEGER(mpi) :: npb
    REAL(mpd) :: nrm
    REAL(mpd) :: sp

    REAL(mpd), INTENT(IN)             :: a(*)
247
248
    INTEGER(mpi), INTENT(IN)          :: bpar(2,*)
    INTEGER(mpi), INTENT(IN)          :: bcon(3,*)
249
    INTEGER(mpi), INTENT(IN)          :: cons(3,*)
250
251

    !$POMP INST BEGIN(qldecb)
252
    ! prepare 
253
254
    vecVk=0.0_mpd
    matV=a(1:matsize)
255
    matL=0.0_mpd
256

257
    ! prepare offsets
258
259
    icoff=0
    DO ibpar=1,nblock ! parameter block
260
        iclast=icoff
261
262
263
264
265
        DO ibcon=bpar(2,ibpar)+1, bpar(2,ibpar+1)! constraint block
            ncb=bcon(1,ibcon+1)-bcon(1,ibcon) ! number of constraints in constraint block
            npb=bcon(3,ibcon)+1-bcon(2,ibcon) ! number of parameters in constraint block
            ifirst=bcon(2,ibcon)
            ilast=bcon(3,ibcon)
266
267
268
269
270
            DO i=bcon(1,ibcon),bcon(1,ibcon+1)-1
                ioffRow(i+1)=ioffRow(i)+npb  ! row offset
                ioffPar(i)=bcon(2,ibcon)-1   ! parameter offset
                irangeParNZ(1,i)=cons(1,i)
                irangeParNZ(2,i)=cons(2,i)
271
            END DO
272
            iclast=iclast+ncb
273
        END DO
274
275
276
277
278
279
280
281
282
283
284
285
286
287
        ! set up matL
        iblast=bpar(1,ibpar+1) ! last parameter in parameter block
        DO k=icoff+1,iclast
            kn=iblast+k-iclast
            ioff1=ioffRow(k)
            npb=INT(ioffRow(k+1)-ioff1,mpi)
            ! size of overlap
            ncol=ioffPar(k)+npb-kn
            IF (ncol >= 0) THEN
                ioff3=INT(k-1,mpl)*INT(ncon,mpl)
                matL(ioff3+iclast-ncol-icoff:ioff3+iclast-icoff)=matV(ioff1+npb-ncol:ioff1+npb)
            END IF
        END DO    
        icoff=iclast
288
289
        nparBlock(ibpar)=bpar(1,ibpar+1)-bpar(1,ibpar)
        ioffBlock(ibpar+1)=icoff
290
291
    END DO
 
292
    DO ibpar=1,nblock ! parameter block
293
294
        iboff=bpar(1,ibpar)    ! last offset in parameter block
        iblast=bpar(1,ibpar+1) ! last parameter in parameter block
295
296
297
298
299
300
        icoff=ioffBlock(ibpar) ! constraint offset in parameter block
        iclast=ioffBlock(ibpar+1) ! last constraint in parameter block
        ibcon=bpar(2,ibpar+1) ! start with last constraint block
        k1=bcon(1,ibcon) ! first constraint in block
        ! Householder procedure
        DO k=iclast,icoff+1,-1
301
302
            ! monitoring
            IF(monpg>0) CALL monpgs(ncon+1-k)
303
            kn=iblast+k-iclast
304
305
306
307
308
            ! different constraint block?
            IF (k < k1) THEN
                ibcon=ibcon-1
                k1=bcon(1,ibcon)
            END IF
309
310
311
312
            ! parameter offset
            ipoff=ioffPar(k)
            ! index if first non-zero parameter
            ifirst=ipoff+1
313
            IF (ifirst > kn) CYCLE
314
315
316
317
318
319
320
321
322
323
            ! column offset
            ioff1=ioffRow(k)
            ! number of parameter
            npb=INT(ioffRow(k+1)-ioff1,mpi)
            ! index of last parameter
            iplast=ioffPar(k)+npb
            ! index of last used parameter
            ilast=min(iplast,kn)
            ! number of used columns
            ncol=ilast-ipoff
324
            ! get column
325
326
            in=iblast+k1-iclast
            vecN(in:kn)=0.0_mpd
327
            vecN(ifirst:ilast)=matV(ioff1+1:ioff1+ncol)
328
            nrm = SQRT(dot_product(vecN(ifirst:ilast),vecN(ifirst:ilast)))
329
330
331
332
333
334
335
            IF (nrm == 0.0_mpd) CYCLE
            !
            IF (vecN(kn) >= 0.0_mpd) THEN
                vecN(kn)=vecN(kn)+nrm
            ELSE
                vecN(kn)=vecN(kn)-nrm
            END IF
336
            ! create normal vector        
337
338
339
340
341
342
343
            IF (ilast < kn) THEN
                nrm = SQRT(dot_product(vecN(ifirst:ilast),vecN(ifirst:ilast))+vecN(kn)*vecN(kn))
                vecN(ifirst:ilast)=vecN(ifirst:ilast)/nrm
                vecN(kn)=vecN(kn)/nrm
            ELSE
                nrm = SQRT(dot_product(vecN(ifirst:ilast),vecN(ifirst:ilast)))
                vecN(ifirst:ilast)=vecN(ifirst:ilast)/nrm
344
            END IF
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
            ! update L too
            ioff3=INT(k1-1,mpl)*INT(ncon,mpl)
            ! transformation
            DO i=k1,k
                ioff2=ioffRow(i)
                sp=dot_product(vecN(ifirst:ilast),matV(ioff2+1:ioff2+ncol))
                IF (sp /= 0.0_mpd) THEN
                    ! update matV
                    matV(ioff2+1:ioff2+ncol)=matV(ioff2+1:ioff2+ncol)-2.0_mpd*vecN(ifirst:ilast)*sp
                    ! update matL
                    in=iblast+i-iclast
                    matL(ioff3+i-icoff:ioff3+k-icoff)=matL(ioff3+i-icoff:ioff3+k-icoff)-2.0_mpd*vecN(in:kn)*sp
                    ! update non zero range
                    irangeParNZ(1,i)=min(irangeParNZ(1,i),irangeParNZ(1,k))
                    irangeParNZ(2,i)=max(irangeParNZ(2,i),irangeParNZ(2,k))
                END IF
                ioff3=ioff3+ncon
            END DO  
            ! store secondary diagonal
            vecVk(icoff+k)=vecN(kn)         
            ! store normal vector (non zero part)
            ifirst=irangeParNZ(1,k)
            ilast=min(irangeParNZ(2,k),kn-1)
            ncol=ilast-ifirst+1
            matV(ioff1+1:ioff1+ncol)=vecN(ifirst:ilast)
            matV(ioff1+ncol+1:ioff1+npb)=0.0_mpd
            ! local to parameter block
            irangeParNZ(1,k)=ifirst-iboff
            irangeParNZ(2,k)=ilast-iboff
374
        END DO
375
    END DO
376
377
    !$POMP INST END(qldecb)
    
378
379
380
    
END SUBROUTINE qldecb

381

382
!> Multiply left by Q(t) (per block).
383
384
385
!!
!! Multiply left by Q(t) from QL decomposition.
!!
386
!! \param [in,out] x    NparBlock-by-M matrix, overwritten with Q*X (t=false) or Q^t*X (t=true)
387
388
389
390
391
392
393
394
395
396
!! \param [in]     m    number of columns
!! \param [in]     t    use transposed of Q
!!
SUBROUTINE qlmlq(x,m,t)
    USE mpqldec

    ! cost[dot ops] ~= N*M*Nhr

    IMPLICIT NONE
    INTEGER(mpi) :: i
397
398
    INTEGER(mpi) :: icoff
    INTEGER(mpi) :: iclast
399
400
    INTEGER(mpi) :: ifirst
    INTEGER(mpi) :: ilast
401
402
403
    INTEGER(mpl) :: ioff2
    INTEGER(mpi) :: j
    INTEGER(mpi) :: k
404
    INTEGER(mpi) :: l
405
    INTEGER(mpi) :: kn
406
407
    INTEGER(mpi) :: nconb
    INTEGER(mpi) :: nparb
408
409
410
411
412
413
    REAL(mpd) :: sp

    REAL(mpd), INTENT(IN OUT)         :: x(*)
    INTEGER(mpi), INTENT(IN)          :: m
    LOGICAL, INTENT(IN)               :: t

414
415
416
417
418
    icoff=ioffBlock(iblock) ! constraint offset in parameter block
    iclast=ioffBlock(iblock+1) ! last constraint in parameter block
    nconb=iclast-icoff ! number of constraints in block
    nparb=nparBlock(iblock) ! number of parameters in block
    DO j=1,nconb
419
        k=j
420
421
        IF (t) k=nconb+1-j
        kn=nparb+k-nconb
422
423
424
425
426
427
428
429
        ! expand row 'l' of matV into vecN
        l=k+icoff
        ! non-zero range (excluding 'kn')
        ifirst=irangeParNZ(1,l) 
        ilast=irangeParNZ(2,l)
        vecN(1:nparb)=0._mpd
        vecN(ifirst:ilast)=matV(ioffRow(l)+1:ioffRow(l)+1+ilast-ifirst) 
        vecN(kn)=vecVk(l)
430
431
432
        ! transformation
        ioff2=0
        DO i=1,m
433
434
435
            sp=dot_product(vecN(ifirst:ilast),x(ioff2+ifirst:ioff2+ilast))+vecN(kn)*x(ioff2+kn)
            x(ioff2+ifirst:ioff2+ilast)=x(ioff2+ifirst:ioff2+ilast)-2.0_mpd*vecN(ifirst:ilast)*sp
            x(ioff2+kn)=x(ioff2+kn)-2.0_mpd*vecN(kn)*sp
436
            ioff2=ioff2+nparb
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
        END DO
    END DO

END SUBROUTINE qlmlq


!> Multiply right by Q(t).
!!
!! Multiply right by Q(t) from QL decomposition.
!!
!! \param [in,out] x    M-by-Npar matrix, overwritten with X*Q (t=false) or X*Q^t (t=true)
!! \param [in]     m    number of rows
!! \param [in]     t    use transposed of Q
!!
SUBROUTINE qlmrq(x,m,t)
    USE mpqldec

    ! cost[dot ops] ~= N*M*Nhr

    IMPLICIT NONE
    INTEGER(mpi) :: i
    INTEGER(mpl) :: iend
459
460
    INTEGER(mpi) :: ifirst
    INTEGER(mpi) :: ilast
461
462
463
464
465
466
467
468
469
470
471
472
473
    INTEGER(mpi) :: j
    INTEGER(mpi) :: k
    INTEGER(mpi) :: kn
    REAL(mpd) :: sp

    REAL(mpd), INTENT(IN OUT)         :: x(*)
    INTEGER(mpi), INTENT(IN)          :: m
    LOGICAL, INTENT(IN)               :: t

    DO j=1,ncon
        k=j
        IF (.not.t) k=ncon+1-j
        kn=npar+k-ncon
474
475
476
477
478
479
480
        ! expand row 'k' of matV into vecN
        ! non-zero range (excluding 'kn')
        ifirst=irangeParNZ(1,k) 
        ilast=irangeParNZ(2,k)
        vecN=0._mpd
        vecN(ifirst:ilast)=matV(ioffRow(k)+1:ioffRow(k)+1+ilast-ifirst) 
        vecN(kn)=vecVk(k)
481
482
483
        ! transformation
        iend=m*kn
        DO i=1,npar
484
485
            sp=dot_product(vecN(1:kn),x(i:iend:m))
            x(i:iend:m)=x(i:iend:m)-2.0_mpd*vecN(1:kn)*sp
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
        END DO
    END DO

END SUBROUTINE qlmrq


!> Similarity transformation by Q(t).
!!
!! Similarity transformation by Q from QL decomposition.
!!
!! \param [in,out] x    Npar-by-Npar matrix, overwritten with Q*X*Q^t (t=false) or Q^t*X*Q (t=true)
!! \param [in]     t    use transposed of Q
!!
SUBROUTINE qlsmq(x,t)
    USE mpqldec

    ! cost[dot ops] ~= N*N*Nhr

    IMPLICIT NONE
    INTEGER(mpi) :: i
    INTEGER(mpl) :: ioff2
    INTEGER(mpl) :: iend
508
509
    INTEGER(mpi) :: ifirst
    INTEGER(mpi) :: ilast
510
511
512
513
514
515
516
517
518
    INTEGER(mpi) :: j
    INTEGER(mpi) :: k
    INTEGER(mpi) :: kn
    REAL(mpd) :: sp

    REAL(mpd), INTENT(IN OUT)         :: x(*)
    LOGICAL, INTENT(IN)               :: t

    DO j=1,ncon
519
520
        ! monitoring
        IF(monpg>0) CALL monpgs(j)
521
522
523
        k=j
        IF (t) k=ncon+1-j
        kn=npar+k-ncon
524
525
526
527
528
529
530
        ! expand row 'k' of matV into vecN
        ! non-zero range (excluding 'kn')
        ifirst=irangeParNZ(1,k) 
        ilast=irangeParNZ(2,k)
        vecN=0._mpd
        vecN(ifirst:ilast)=matV(ioffRow(k)+1:ioffRow(k)+1+ilast-ifirst) 
        vecN(kn)=vecVk(k)
531
        ! transformation
532
        iend=INT(npar,mpl)*INT(kn,mpl)
533
        DO i=1,npar
534
535
            sp=dot_product(vecN(1:kn),x(i:iend:npar))
            x(i:iend:npar)=x(i:iend:npar)-2.0_mpd*vecN(1:kn)*sp
536
537
538
        END DO
        ioff2=0
        DO i=1,npar
539
540
            sp=dot_product(vecN(1:kn),x(ioff2+1:ioff2+kn))
            x(ioff2+1:ioff2+kn)=x(ioff2+1:ioff2+kn)-2.0_mpd*vecN(1:kn)*sp
541
542
543
544
545
546
547
548
549
550
551
552
            ioff2=ioff2+npar
        END DO
    END DO

END SUBROUTINE qlsmq


!> Similarity transformation by Q(t).
!!
!! Similarity transformation for symmetric matrix by Q from QL decomposition.
!!
!! \param [in]     aprod    external procedure to calculate A*v
553
!! \param [in,out] A        symmetric Npar-by-Npar matrix A in symmetric or unpacked storage mode
554
!!                          overwritten with Q*A*Q^t (t=false) or Q^t*A*Q (t=true)
555
!! \param [in]     roff     row offsets for A
556
557
!! \param [in]     t        use transposed of Q
!!
558
SUBROUTINE qlssq(aprod,A,roff,t)
559
560
561
562
563
564
565
    USE mpqldec
    USE mpdalc

    ! cost[dot ops] ~= N*N*Nhr

    IMPLICIT NONE
    INTEGER(mpi) :: i
566
567
568
    INTEGER(mpi) :: ibpar
    INTEGER(mpi) :: icoff
    INTEGER(mpi) :: iclast
569
570
571
    INTEGER(mpi) :: ifirst
    INTEGER(mpi) :: ilast
    INTEGER(mpi) :: ilasti
572
    INTEGER(mpl) :: ioff2
573
    INTEGER(mpi) :: ioffp
574
575
    INTEGER(mpi) :: j
    INTEGER(mpi) :: k
576
    INTEGER(mpi) :: l
577
578
    INTEGER(mpi) :: kn
    INTEGER(mpl) :: length
579
580
    INTEGER(mpi) :: nconb
    INTEGER(mpi) :: nparb
581
582
583
584
    REAL(mpd) :: vtAv
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: Av

    REAL(mpd), INTENT(IN OUT)         :: A(*)
585
    INTEGER(mpl), INTENT(IN)          :: roff(*)
586
587
588
    LOGICAL, INTENT(IN)               :: t

    INTERFACE
589
        SUBROUTINE aprod(n,l,x,is,ie,y) ! y=A*x
590
591
            USE mpdef
            INTEGER(mpi), INTENT(in) :: n
592
            INTEGER(mpl), INTENT(in) :: l
593
            REAL(mpd), INTENT(IN)    :: x(n)
594
595
            INTEGER(mpi), INTENT(in) :: is
            INTEGER(mpi), INTENT(in) :: ie
596
597
598
            REAL(mpd), INTENT(OUT)   :: y(n)
        END SUBROUTINE aprod
    END INTERFACE
599
600
    !$POMP INST BEGIN(qlssq)
    
601
602
603
    length=npar
    CALL mpalloc(Av,length,'qlssq: A*v')

604
    ioffp=0 ! parameter offset for block
605
606
607
608
609
610
611
    DO ibpar=1,nblock ! parameter block
        icoff=ioffBlock(ibpar) ! constraint offset in parameter block
        iclast=ioffBlock(ibpar+1) ! last constraint in parameter block
        nconb=iclast-icoff ! number of constraints in block
        nparb=nparBlock(ibpar) ! number of parameters in block
        DO j=1,nconb
            k=j
612
613
            ! monitoring
            IF(monpg>0) CALL monpgs(icoff+k)
614
615
            IF (t) k=nconb+1-j
            kn=nparb+k-nconb
616
617
618
619
620
621
622
623
            ! expand row 'l' of matV into vecN
            l=k+icoff
            ! non-zero range (excluding 'kn')
            ifirst=irangeParNZ(1,l) 
            ilast=irangeParNZ(2,l)
            vecN(1:nparb)=0._mpd
            vecN(ifirst:ilast)=matV(ioffRow(l)+1:ioffRow(l)+1+ilast-ifirst) 
            vecN(kn)=vecVk(k)
624
            ! A*v
625
626
627
            Av(1:nparb)=0._mpd
            CALL aprod(nparb,INT(ioffp,mpl),vecN(1:nparb),ifirst,ilast,Av(1:nparb))
            CALL aprod(nparb,INT(ioffp,mpl),vecN(1:nparb),kn,kn,Av(1:nparb))
628
629
630
            ! transformation
            ! diagonal block
            ! v^t*A*v
631
            vtAv=dot_product(vecN(ifirst:ilast),Av(ifirst:ilast))+vecN(kn)*Av(kn)
632
            ! update
633
634
635
            ! parallelize row loop
            ! slot of 8 'I' for next idle thread 
            !$OMP PARALLEL DO &
636
            !$OMP PRIVATE(IOFF2,ILASTI) &
637
            !$OMP SCHEDULE(DYNAMIC,8)
638
            DO i=1,kn
639
                ioff2=roff(i+ioffp)+ioffp
640
641
642
643
                ilasti=min(ilast,i)
                ! correct with  2*(2v*vtAv*v^t - Av*v^t)
                A(ioff2+ifirst:ioff2+ilasti)=A(ioff2+ifirst:ioff2+ilasti)+2.0_mpd* &
                    ((2.0_mpd*vecN(i)*vtAv-Av(i))*vecN(ifirst:ilasti))
644
            END DO
645
            !$OMP END PARALLEL DO
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
            
            ! parallelize row loop
            ! slot of 8 'I' for next idle thread 
            !$OMP PARALLEL DO &
            !$OMP PRIVATE(IOFF2) &
            !$OMP SCHEDULE(DYNAMIC,8)
            DO i=ifirst,ilast
                ioff2=roff(i+ioffp)+ioffp
                ! correct with  -2(Av*v^t)^t)
                A(ioff2+1:ioff2+i)=A(ioff2+1:ioff2+i)-2.0_mpd*Av(1:i)*vecN(i)
            END DO            
            !$OMP END PARALLEL DO
            ! i=kn, add secondary diagonal element
            ioff2=roff(kn+ioffp)+ioffp
            A(ioff2+kn)=A(ioff2+kn)+2.0_mpd*((2.0_mpd*vecVk(l)*vtAv-Av(kn))*vecVk(l)-Av(kn)*vecVk(l))
661
662
            ! off diagonal block
            DO i=kn+1,nparb
663
                ioff2=roff(i+ioffp)+ioffp
664
                ! correct with -2Av*v^t
665
666
                A(ioff2+ifirst:ioff2+ilast)=A(ioff2+ifirst:ioff2+ilast)-2.0_mpd*vecN(ifirst:ilast)*Av(i)
                A(ioff2+kn)=A(ioff2+kn)-2.0_mpd*vecVk(l)*Av(i)
667
668
            END DO
        END DO
669
670
        ! update parameter offset
        ioffp=ioffp+nparb
671
672
673
    END DO

    CALL mpdealloc(Av)
674
    !$POMP INST END(qlssq)
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696

END SUBROUTINE qlssq

!> Partial similarity transformation by Q(t).
!!
!! Partial similarity transformation for symmetric matrix by Q from QL decomposition.
!! Calculate corrections to band part of matrix.
!!
!! \param [in]     aprod    external procedure to calculate A*v
!! \param [in,out] B        band part of symmetric Npar-by-Npar matrix A in symmetric storage mode,
!!                          overwritten with band part of Q^t*A*Q (t=false) or Q^t*A*Q (t=true)
!! \param [in]     m        band width (including diagonal)
!! \param [in]     t        use transposed of Q
!!
SUBROUTINE qlpssq(aprod,B,m,t)
    USE mpqldec
    USE mpdalc

    ! cost[dot ops] ~= N*N*Nhr

    IMPLICIT NONE
    INTEGER(mpi) :: i
697
698
699
700
    INTEGER(mpi) :: ifirst
    INTEGER(mpi) :: ilast
    INTEGER(mpi) :: ifirst2
    INTEGER(mpi) :: ilast2
701
702
    INTEGER(mpl) :: ioff1
    INTEGER(mpl) :: ioff2
703
    INTEGER(mpi) :: istat(3)
704
705
706
707
708
    INTEGER(mpi) :: j
    INTEGER(mpi) :: j2
    INTEGER(mpi) :: k
    INTEGER(mpi) :: k2
    INTEGER(mpi) :: kn
709
    INTEGER(mpi) :: kn2
710
    INTEGER(mpi) :: l
711
712
    INTEGER(mpi) :: l1
    INTEGER(mpi) :: l2
713
714
    INTEGER(mpl) :: length
    INTEGER(mpi) :: mbnd
715
    REAL(mpd) :: v2kn
716
717
718
    REAL(mpd) :: vtAv
    REAL(mpd) :: vtAvp
    REAL(mpd) :: vtvp
719
720
721
722
723
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: vecAv    ! A*v 
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: matvtvp  ! v^t*v' 
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: matvtAvp ! v^t*(A*v')
    REAL(mpd), DIMENSION(:), ALLOCATABLE :: matCoeff ! coefficients (d(A*v)=sum(c_i*v_i)) 
    INTEGER(mpi), DIMENSION(:,:), ALLOCATABLE :: irangeCoeff !< range for non zero part
724
725
726
727
728
729

    REAL(mpd), INTENT(IN OUT)         :: B(*)
    INTEGER(mpi), INTENT(IN)          :: m
    LOGICAL, INTENT(IN)               :: t

    INTERFACE
730
        SUBROUTINE aprod(n,l,x,is,ie,y) ! y=A*x
731
732
            USE mpdef
            INTEGER(mpi), INTENT(in) :: n
733
            INTEGER(mpl), INTENT(in) :: l
734
            REAL(mpd), INTENT(IN)    :: x(n)
735
736
            INTEGER(mpi), INTENT(in) :: is
            INTEGER(mpi), INTENT(in) :: ie
737
738
739
            REAL(mpd), INTENT(OUT)   :: y(n)
        END SUBROUTINE aprod
    END INTERFACE
740
    !$POMP INST BEGIN(qlpssq)
741

742
743
744
745
746
747
748
749
750
751
752
    length=npar
    CALL mpalloc(vecAv,length,'qlpssq: A*v')    
    length=INT(ncon,mpl)*INT(ncon,mpl)
    CALL mpalloc(matvtvp,length,"qlpssq: v^t*v'")
    matvtvp=0._mpd
    CALL mpalloc(matvtAvp,length,"qlpssq: v^t*(A*v')")
    matvtAvp=0._mpd
    CALL mpalloc(matCoeff,length,'qlpssq: coefficients')
    matCoeff=0._mpd
    length=ncon
    CALL mpalloc(irangeCoeff,2_mpl,length,'qlpssq: non zero coefficient range')
753
754

    mbnd=max(0,m-1) ! band width without diagonal
755
 
756
757
    DO j=1,ncon
        k=j
758
759
        ! monitoring
        IF(monpg>0) CALL monpgs(k)
760
761
        IF (t) k=ncon+1-j
        kn=npar+k-ncon
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
        ioff1=INT(k-1,mpl)*INT(ncon,mpl)
        irangeCoeff(1,k)=ncon
        irangeCoeff(2,k)=1
        ! expand row 'k' of matV into vecN
        ! non-zero range (excluding 'kn')
        ifirst=irangeParNZ(1,k) 
        ilast=irangeParNZ(2,k)
        vecN=0._mpd
        vecN(ifirst:ilast)=matV(ioffRow(k)+1:ioffRow(k)+1+ilast-ifirst) 
        vecN(kn)=vecVk(k)  
        ! transformation A*v
        vecAv(1:npar)=0._mpd
        CALL aprod(npar,0_mpl,vecN(1:npar),ifirst,ilast,vecAv(1:npar))
        CALL aprod(npar,0_mpl,vecN(1:npar),kn,kn,vecAv(1:npar))
        ! products v^t*v'               
        DO j2=j+1,ncon
            k2=j2
            IF (t) k2=ncon+1-j2
            kn2=npar+k2-ncon
            ioff2=INT(k2-1,mpl)*INT(ncon,mpl)
            ! non-zero range (excluding 'kn')
            ifirst2=irangeParNZ(1,k2) 
            ilast2=irangeParNZ(2,k2)
            v2kn=0._mpd
            IF (kn >= ifirst2.AND.kn <= ilast2) v2kn=matV(ioffRow(k2)+1+kn-ifirst2)
            ! overlap regions
            l1=max(ifirst,ifirst2)
            l2=min(ilast,ilast2)
            vtvp=vecN(kn2)*vecVk(k2)+vecN(kn)*v2kn ! v^t*v'
            IF (l1 <= l2) vtvp=vtvp+dot_product(vecN(l1:l2), &
                matV(ioffRow(k2)+1+l1-ifirst2:ioffRow(k2)+1+l2-ifirst2))
            ! significant term?
            IF (abs(vtvp) > 16.0_mpd*epsilon(vtvp)) THEN
                matvtvp(ioff1+k2)=vtvp 
                matvtvp(ioff2+k)=vtvp 
            END IF    
        END DO
        matvtvp(ioff1+k)=1.0_mpd
        ! products v^t*(A*v')                 
        DO j2=1,j
            k2=j2
            IF (t) k2=ncon+1-j2
            kn2=npar+k2-ncon
            ! non-zero range (excluding 'kn')
            ifirst2=irangeParNZ(1,k2) 
            ilast2=irangeParNZ(2,k2)
            ! non-zero regions
            matvtAvp(ioff1+k2)=vecVk(k2)*vecAv(kn2)+dot_product(vecAv(ifirst2:ilast2), &
                matV(ioffRow(k2)+1:ioffRow(k2)+1+ilast2-ifirst2)) ! v'^t*(A*v)
        END DO
        ! update with (initial) A*v
        ioff2=0
        vtAv=matvtAvp(ioff1+k)
        DO i=1,kn
            ! correct with  2*(2v*vtAv*v^t - Av*v^t - (Av*v^t)^t)
            DO l=max(1,i-mbnd),i
                ioff2=ioff2+1
                B(ioff2)=B(ioff2)+2.0_mpd*((2.0_mpd*vecN(i)*vtAv-vecAv(i))*vecN(l)-vecAv(l)*vecN(i))
            END DO
        END DO
        ! off diagonal block
        DO i=kn+1,npar
            ! correct with -2Av*v^t
            DO l=max(1,i-mbnd),i
                ioff2=ioff2+1
                B(ioff2)=B(ioff2)-2.0_mpd*vecAv(i)*vecN(l)
            END DO
        END DO
    END DO

    ! corrections for A*v (as linear combination of v's)
    DO j=1,ncon
        k=j
        IF (t) k=ncon+1-j
        kn=npar+k-ncon
        ioff1=INT(k-1,mpl)*INT(ncon,mpl)
        ! expand row 'k' of matV into vecN
        ! non-zero range (excluding 'kn')
        ifirst=irangeParNZ(1,k) 
        ilast=irangeParNZ(2,k)
        vecN=0._mpd
        vecN(ifirst:ilast)=matV(ioffRow(k)+1:ioffRow(k)+1+ilast-ifirst) 
        vecN(kn)=vecVk(k)
845
        ! transformation (diagonal block)
846
847
        l1=irangeCoeff(1,k)
        l2=irangeCoeff(2,k)
848
849
        ! diagonal block
        ! v^t*A*v
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
        vtAv=matvtAvp(ioff1+k)+dot_product(matCoeff(ioff1+l1:ioff1+l2),matvtvp(ioff1+l1:ioff1+l2))
        ! expand correction to initial A*v
        vecAv(1:npar)=0._mpd
        istat=0
        DO k2=l1,l2
            IF (matCoeff(ioff1+k2) == 0._mpd) CYCLE
            if (istat(1)==0) istat(1)=k2
            istat(2)=k2
            istat(3)=istat(3)+1
            kn2=npar+k2-ncon
            ! expand row 'k2' of matV directly into vecAv
            ! non-zero range (excluding 'kn')
            ifirst2=irangeParNZ(1,k2) 
            ilast2=irangeParNZ(2,k2)          
            vecAv(ifirst2:ilast2)=vecAv(ifirst2:ilast2)+matCoeff(ioff1+k2)* &
                matV(ioffRow(k2)+1:ioffRow(k2)+1+ilast2-ifirst2)            
            vecAv(kn2)=vecAv(kn2)+matCoeff(ioff1+k2)*vecVk(k2)
        END DO
868
869
870
871
872
873
        ! update
        ioff2=0
        DO i=1,kn
            ! correct with  2*(2v*vtAv*v^t - Av*v^t - (Av*v^t)^t)
            DO l=max(1,i-mbnd),i
                ioff2=ioff2+1
874
                B(ioff2)=B(ioff2)+2.0_mpd*((2.0_mpd*vecN(i)*vtAv-vecAv(i))*vecN(l)-vecAv(l)*vecN(i))
875
876
877
878
879
880
881
            END DO
        END DO
        ! off diagonal block
        DO i=kn+1,npar
            ! correct with -2Av*v^t
            DO l=max(1,i-mbnd),i
                ioff2=ioff2+1
882
                B(ioff2)=B(ioff2)-2.0_mpd*vecAv(i)*vecN(l)
883
884
885
886
887
888
            END DO
        END DO
        ! correct A*v for the remainung v
        DO j2=j+1,ncon
            k2=j2
            IF (t) k2=ncon+1-j2
889
890
891
892
893
894
895
896
897
898
899
            kn2=npar+k2-ncon
            ioff2=INT(k2-1,mpl)*INT(ncon,mpl)
            vtvp=matvtvp(ioff1+k2) ! v^t*v'
            ! non-zero regions
            l1=irangeCoeff(1,k2)
            l2=irangeCoeff(2,k2)
            vtAvp=matvtAvp(ioff2+k)
            IF (l1 <= l2) vtAvp=vtAvp+dot_product(matCoeff(ioff2+l1:ioff2+l2),matvtvp(ioff1+l1:ioff1+l2)) ! v^t*(A*v')
            l1=min(l1,k)
            l2=max(l2,k)
            matCoeff(ioff2+k)=matCoeff(ioff2+k)+2.0_mpd*(2.0_mpd*vtAv*vtvp-vtAvp)
900
            IF (vtvp /= 0._mpd) THEN
901
902
903
                l1=min(l1,irangeCoeff(1,k))
                l2=max(l2,irangeCoeff(2,k))
                matCoeff(ioff2+l1:ioff2+l2)=matCoeff(ioff2+l1:ioff2+l2)-2.0_mpd*matCoeff(ioff1+l1:ioff1+l2)*vtvp
904
            END IF
905
906
            irangeCoeff(1,k2)=l1
            irangeCoeff(2,k2)=l2
907
908
909
        END DO
    END DO

910
911
912
913
914
    CALL mpdealloc(irangeCoeff)
    CALL mpdealloc(matCoeff)
    CALL mpdealloc(matvtAvp)
    CALL mpdealloc(matvtvp)
    CALL mpdealloc(vecAv)
915
    !$POMP INST END(qlpssq)
916
917
918
919
920
921

END SUBROUTINE qlpssq


!> Get eigenvalues.
!!
922
!! Get smallest and largest |eigenvalue| of L.
923
!!
924
925
!! \param [out]    emin  eigenvalue with smallest absolute value
!! \param [out]    emax  eigenvalue with largest absolute value
926
927
928
929
930
931
!!
SUBROUTINE qlgete(emin,emax)
    USE mpqldec

    IMPLICIT NONE
    INTEGER(mpi) :: i
932
933
934
    INTEGER(mpi) :: ibpar
    INTEGER(mpi) :: icoff
    INTEGER(mpi) :: iclast
935
936
937
938
939
940
941
    INTEGER(mpl) :: idiag

    REAL(mpd), INTENT(OUT)         :: emin
    REAL(mpd), INTENT(OUT)         :: emax

    emax=matL(1)
    emin=emax
942
943
944
    DO ibpar=1,nblock ! parameter block
        icoff=ioffBlock(ibpar) ! constraint offset in parameter block
        iclast=ioffBlock(ibpar+1) ! last constraint in parameter block
945
        idiag=INT(ncon,mpl)*INT(icoff,mpl)+1
946
947
948
949
950
        DO i=icoff+1,iclast
            IF (ABS(emax) < ABS(matL(idiag))) emax=matL(idiag)
            IF (ABS(emin) > ABS(matL(idiag))) emin=matL(idiag)
            idiag=idiag+ncon+1
        END DO
951
952
953
954
955
    END DO

END SUBROUTINE qlgete


956
!> Backward substitution (per block).
957
958
959
960
961
962
963
964
965
966
!!
!! Get y from L^t*y=d.
!!
!! \param [in]     d  Ncon vector, resdiduals
!! \param [out]    y  Ncon vector, solution
!!
SUBROUTINE qlbsub(d,y)
    USE mpqldec

    IMPLICIT NONE
967
968
    INTEGER(mpi) :: icoff
    INTEGER(mpi) :: iclast
969
970
    INTEGER(mpl) :: idiag
    INTEGER(mpi) :: k
971
    INTEGER(mpi) :: nconb
972

973
974
    REAL(mpd), INTENT(IN)         :: d(*)
    REAL(mpd), INTENT(OUT)        :: y(*)
975
976

    ! solve L*y=d by forward substitution
977
978
979
    icoff=ioffBlock(iblock) ! constraint offset in parameter block
    iclast=ioffBlock(iblock+1) ! last constraint in parameter block
    nconb=iclast-icoff ! number of constraints in block
980
    idiag=INT(ncon,mpl)*INT(iclast-1,mpl)+nconb
981
982
    DO k=nconb,1,-1
        y(k)=(d(k)-dot_product(matL(idiag+1:idiag+nconb-k),y(k+1:nconb)))/matL(idiag)
983
984
985
986
        idiag=idiag-ncon-1
    END DO

END SUBROUTINE qlbsub
987
988
989
990
991
992
993
994
995
996
997
998
999
1000

!> Set block
!!
SUBROUTINE qlsetb(ib)
    USE mpqldec

    IMPLICIT NONE
    INTEGER(mpi), INTENT(IN)      :: ib

    iblock=ib

END SUBROUTINE qlsetb

!> Print statistics