-
Notifications
You must be signed in to change notification settings - Fork 15
/
lmdif.f90
658 lines (597 loc) · 27.8 KB
/
lmdif.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
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
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
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
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
#include "alias.inc"
subroutine leasqr_lm (get_eig, NN_TABLE, EDFT, PWGHT, PINPT, PPRAM, PKPTS, PGEOM, fnorm)
use parameters
use mpi_setup
use print_io
implicit none
type (incar) :: PINPT
type (params) :: PPRAM
type (hopping), dimension(PINPT%nsystem) :: NN_TABLE
type (weight) , dimension(PINPT%nsystem) :: PWGHT
type (kpoints), dimension(PINPT%nsystem) :: PKPTS
type (poscar) , dimension(PINPT%nsystem) :: PGEOM
type (energy) , dimension(PINPT%nsystem) :: EDFT
integer*4 info, maxfev
integer*4 nparam_free, mpierr
integer*4 ldjac, imode ! ldjac > nparam_free
integer*4 i
real*8 epsfcn,factor, tol, xtol, ftol, gtol
external get_eig
logical flag_write_info, flag_fit_orbital
real*8 fnorm, fnorm_plain, fnorm_orb ! fnorm of last step
if(PPRAM%slater_koster_type .gt. 10) then
nparam_free = PPRAM%nparam_nrl_free ! total number of free parameters
! PPRAM%nparam total number of sk parameter set
! PPRAM%nparam_free total number of sk free parameter set
else
nparam_free = PPRAM%nparam_free ! total number of free parameters
endif
imode = 13
ldjac = 0
do i =1, PINPT%nsystem
ldjac = ldjac + PKPTS(i)%nkpoint * PGEOM(i)%nband * PINPT%nspin
enddo
flag_fit_orbital = PINPT%flag_fit_orbital
if( PINPT%ls_type == 'LMDIF' ) then
write(message,*)' Start: fitting procedures with ',PINPT%ls_type,' method.' ; write_msg
if ( nparam_free <= 0 ) then
write(message,*)' !!!! WARN !!!! nparam_free (total number of free parameters) <= 0' ; write_msg
write(message,*)' Check parameter settings. Exit program.' ; write_msg
kill_job
elseif( ldjac < nparam_free ) then
write(message,*)' !!!! WARN !!!! number of kpoints * number of eigenvalues <= number of free parameters' ; write_msg
write(message,*)' => Increase kpoints, and try again. Exit program.' ; write_msg
kill_job
endif
factor = 100.0D+00
maxfev = PINPT%miter * ( nparam_free + 1 )
!maxfev = 2000 * ( nparam_free + 1 )
ftol = PINPT%ftol ;xtol = PINPT%ptol ; gtol = 0.0D+00;epsfcn = 0.000D+00
if(PINPT%flag_python_module) then
if(iverbose .eq. 1) then
flag_write_info = .true.
else
flag_write_info = .false.
endif
else
flag_write_info = .true.
endif
call lmdif(get_eig, NN_TABLE, ldjac, imode, PINPT, PPRAM, PKPTS, PGEOM, EDFT, nparam_free, PWGHT, &
ftol, xtol, gtol, fnorm, fnorm_plain, fnorm_orb, maxfev, epsfcn, factor, info, flag_write_info, flag_fit_orbital)
if_main call infostamp(info,PINPT%ls_type)
!write(message,*)" End: fitting procedures" ; write_msg
elseif( PINPT%ls_type == 'GA' ) then
if ( nparam_free <= 0 ) then
write(message,*)' !!!! WARN !!!! nparam_free (number of free parameters) <= 0' ; write_msg
write(message,*)' Check parameter settings. Exit program.' ; write_msg
kill_job
elseif( ldjac < nparam_free ) then
write(message,*)' !!!! WARN !!!! number of kpoints * number of eigenvalues <= number of free parameters' ; write_msg
write(message,*)' => Increase kpoints, and try again. Exit program.' ; write_msg
kill_job
endif
factor = 100.0D+00
maxfev = PINPT%miter * ( nparam_free + 1 )
ftol = PINPT%ftol ;xtol = PINPT%ptol ; gtol = 0.0D+00;epsfcn = 0.0D+00
flag_write_info = .false.
call lmdif(get_eig, NN_TABLE, ldjac, imode, PINPT, PPRAM, PKPTS, PGEOM, EDFT, nparam_free, PWGHT, &
ftol, xtol, gtol, fnorm, fnorm_plain, maxfev, epsfcn, factor, info, flag_write_info, flag_fit_orbital)
endif
return
endsubroutine
subroutine lmdif(get_eig, NN_TABLE, ldjac, imode, PINPT, PPRAM, PKPTS, PGEOM, EDFT, nparam_free, PWGHT, &
ftol, xtol, gtol, fnorm, fnorm_plain, fnorm_orb, maxfev, epsfcn, factor, info, flag_write_info_, flag_fit_orbital)
!*****************************************************************************80
!
!! LMDIF minimizes M functions in N variables by the Levenberg-Marquardt method.
!
! Discussion:
!
! LMDIF minimizes the sum of the squares of M nonlinear functions in
! N variables by a modification of the Levenberg-Marquardt algorithm.
! The user must provide a subroutine which calculates the functions.
! The jacobian is then calculated by a forward-difference approximation.
!
! Licensing:
!
! This code is distributed under the GNU LGPL license.
!
! Modified:
!
! 06 April 2010
!
! Author:
!
! Original FORTRAN77 version by Jorge More, Burton Garbow, Kenneth Hillstrom.
! FORTRAN90 version by John Burkardt.
!
! Jorge More, Burton Garbow, Kenneth Hillstrom,
! User Guide for MINPACK-1,
! Technical Report ANL-80-74,
! Argonne National Laboratory, 1980.
!
! Modified: by Hyun-Jung Kim (KIAS, Infant@kias.re.kr)
! 14 December 2017
! Modified for the TBFIT purpose.
! The original source code can be found in: https://people.sc.fsu.edu/~jburkardt/f_src/minpack/minpack.html
! Reference:
!
use parameters
use cost_function
use mpi_setup
use reorder_band
use print_io
use projected_band
implicit none
type(incar) :: PINPT
type(params) :: PPRAM
type(weight) , dimension(PINPT%nsystem) :: PWGHT
type(poscar) , dimension(PINPT%nsystem) :: PGEOM
type(hopping) , dimension(PINPT%nsystem) :: NN_TABLE
type(energy) , dimension(PINPT%nsystem) :: ETBA_FIT, EDFT
type(kpoints) , dimension(PINPT%nsystem) :: PKPTS
integer*4 nparam_free ! number of free parameters (PPRAM%nparam_nrl_free if slater_koster_type >= 10, PPRAM%nparam_free if < 10)
integer*4 k, nsub
integer*4 i,j,l,iter,info,ipvt(nparam_free),maxfev,nfev, info_temp
integer*4 irange(PPRAM%param_nsub_max), irange_(PPRAM%param_nsub_max)
integer*4 ldjac ! leading demension of jacobean fjac
integer*4 imode
real*8 actred,delta,diag(nparam_free),dirder,enorm,epsfcn,epsmch,factor
real*8 fjac(ldjac, nparam_free)
real*8 fvec(ldjac)
real*8 fvec_plain(ldjac)
real*8 fvec_orb(ldjac)
real*8 wa4(ldjac)
real*8 fnorm,fnorm1,fnorm_,fnorm_plain,fnorm_orb, ftol,gnorm,gtol,par
real*8 pnorm,prered,qtf(nparam_free),ratio,sum2,temp,temp1,temp2,xnorm,xtol
real*8 wa1(nparam_free),wa2(nparam_free),wa3(nparam_free)
real*8 wa2_temp(nparam_free)
character*132 pfileoutnm_temp
logical flag_write_info_, flag_write_info
external get_eig
integer*4 i_dummy
character*132 gnu_command
logical flag_wait_plot
logical flag_order, flag_order_weight, flag_cost_history, flag_fit_orbital
integer*4 mpierr
! for the debug purpose.
! info_temp = info
flag_write_info = flag_write_info_
! if(myid .eq. 1) then
! flag_write_info = .TRUE.
! endif
!!!!!!!!!!!
if(allocated(PPRAM%cost_history)) then
flag_cost_history = .TRUE.
else
flag_cost_history = .FALSE.
endif
flag_order = PINPT%flag_get_band_order .and. (.not. PINPT%flag_get_band_order_print_only)
flag_order_weight = .false. ! experimental feature
! manually allocate instead of using allocate_ETBA routine.
do i = 1, PINPT%nsystem
allocate(ETBA_FIT(i)%E(PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
allocate(ETBA_FIT(i)%V(PGEOM(i)%neig*PINPT%ispin,PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
allocate(ETBA_FIT(i)%ORB(PINPT%lmmax,PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
allocate(ETBA_FIT(i)%SV(PGEOM(i)%neig*PINPT%ispin,PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
if(flag_order) then
allocate(ETBA_FIT(i)%E_ORD(PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
allocate(ETBA_FIT(i)%V_ORD(PGEOM(i)%neig*PINPT%ispin,PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
allocate(ETBA_FIT(i)%SV_ORD(PGEOM(i)%neig*PINPT%ispin,PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
allocate(ETBA_FIT(i)%IDX(PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint))
if(PINPT%flag_python_module) then
write(6,'(A)')' ! WARN ! LORDER tag is not supported with python module. Exit program...'
kill_job
endif
endif
if(PINPT%flag_fit_degeneracy) then
allocate(ETBA_FIT(i)%D(3, PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint) )
if(PINPT%flag_python_module) then
write(6,'(A)')' ! WARN ! DEGENW tag in "SET WEIGHT" is not supported with python module. Exit program...'
kill_job
endif
endif
enddo
fnorm_ = 0d0
i_dummy = 0
epsmch = epsilon ( epsmch )
info = 0 ; nfev = 0
fvec = 0d0; fvec_plain = 0d0 ; fvec_orb = 0d0
if (ftol < 0.0D+00 .or. xtol < 0.0D+00 .or. gtol < 0.0D+00 .or. maxfev <= 0) go to 300
! Evaluate degeneracy information for DFT band : in the beginning
do i = 1, PINPT%nsystem
call get_degeneracy(EDFT(i), PGEOM(i)%neig*PINPT%ispin, PKPTS(i)%nkpoint, PINPT)
enddo
call get_dE(fvec, fvec_plain, fvec_orb, ldjac, imode, PINPT, PPRAM, NN_TABLE, EDFT, ETBA_FIT, PWGHT, PGEOM, PKPTS, flag_fit_orbital)
nfev = 1
fnorm = enorm ( ldjac , fvec )
fnorm_plain = enorm(ldjac, fvec_plain)
fnorm_orb = enorm(ldjac, fvec_orb)
! Initialize Levenberg-Marquardt parameter and iteration counter.
iter = 1 ; par = 0.0D+00
if(flag_cost_history) PPRAM%cost_history(iter) = fnorm
if(flag_write_info) then
write(message,'(A)')' ' ; write_msg
write(message,'(A)')' # dE: (EDFT - ETBA), WT= weight ' ; write_msg
if(.not.flag_fit_orbital) then
write(message,'(A,I8, 2(A,F16.6))')' ITER=',iter,' ,rt(sum((dE*WT)^2)) = ',fnorm, &
' ,rt(sum(dE^2)) = ', fnorm_plain ; write_msg
elseif(flag_fit_orbital) then
write(message,'(A,I8, 3(A,F16.6))')' ITER=',iter,' ,rt(sum((dE*WT)^2)) = ',fnorm, &
' ,rt(sum(dE^2)) = ', fnorm_plain, &
' ,rt(sum(dORB^2)) =',fnorm_orb ; write_msg
endif
!', rt(sum(dE^2)) = ', enorm ( ldjac , fvec_plain ) ; write_msg
endif
30 continue ! Beginning of the outer loop.
! Calculate the jacobian matrix.
if(flag_order) then
do i = 1, PINPT%nsystem
call get_ordered_band(ETBA_FIT(i), PKPTS(i), PGEOM(i), PWGHT(i), PINPT, flag_order_weight, PPRAM%flag_use_overlap)
enddo
endif
call fdjac2 (get_eig,NN_TABLE,ldjac,imode,PINPT,PPRAM,PGEOM,fvec,ETBA_FIT,EDFT,nparam_free,PWGHT,fjac,epsfcn,PKPTS)
nfev = nfev + nparam_free
! Compute the QR factorization of the jacobian.
call qrfac ( ldjac, nparam_free, fjac, ipvt, wa1, wa2 )
! On the first iteration, scale according to the norms of the columns of the initial jacobian.
if ( iter == 1 ) then
do j = 1, PPRAM%nparam_free
if(PPRAM%slater_koster_type .gt. 10) then
nsub = PPRAM%param_nsub(PPRAM%iparam_free(j))
i = PPRAM%iparam_free_nrl(j)
diag(i:i+nsub-1) = wa2(i:i+nsub-1)
else
diag(j) = wa2(j)
endif
enddo
do j = 1, nparam_free
if ( wa2(j) == 0.0D+00 ) diag(j) = 1.0D+00
enddo
! On the first iteration, calculate the norm of the scaled X
! and initialize the step bound DELTA.
do j = 1, PPRAM%nparam_free
if(PPRAM%slater_koster_type .gt. 10) then
nsub = PPRAM%param_nsub(PPRAM%iparam_free(j))
i = PPRAM%iparam_free_nrl(j)
wa3(i:i+nsub-1) = diag(i:i+nsub-1)*PPRAM%param_nrl(1:nsub,PPRAM%iparam_free(j))
else
wa3(j) = diag(j)*PPRAM%param(PPRAM%iparam_free(j))
endif
enddo
xnorm = enorm ( nparam_free, wa3 )
delta = factor * xnorm
if ( delta == 0.0D+00 ) delta = factor
endif ! if iter = 1
! Form Q' * FVEC and store the first N components in QTF.
wa4(:) = fvec(:)
do j = 1, nparam_free
if ( fjac(j,j) /= 0.0D+00 ) then
sum2 = dot_product ( wa4(j:ldjac), fjac(j:ldjac,j) )
temp = -sum2 / fjac(j,j)
wa4(j:ldjac) = wa4(j:ldjac) + fjac(j:ldjac,j) * temp
endif
fjac(j,j) = wa1(j)
qtf(j) = wa4(j)
enddo
! Compute the norm of the scaled gradient.
gnorm = 0.0D+00
if ( fnorm /= 0.0D+00 ) then
do j = 1, nparam_free
l = ipvt(j)
if (wa2(l) /= 0.0D+00) then
sum2 = sum( fjac(:,j)*(qtf(:)/fnorm) )
gnorm = max( gnorm, abs(sum2/wa2(l)) )
endif
enddo
endif
! Test for convergence of the gradient norm.
if (gnorm <= gtol) then
info = 4
go to 300
endif
! Rescale if necessary.
diag(1:nparam_free) = max ( diag(1:nparam_free), wa2(1:nparam_free) )
200 continue ! Beginning of the inner loop.
! Determine the Levenberg-Marquardt parameter.
call lmpar ( nparam_free, fjac, ldjac, ipvt, diag, qtf, delta, par, wa1, wa2 )
! Store the direction P and X + P.
! Calculate the norm of P.
wa1(1:nparam_free) = -wa1(1:nparam_free)
! NOTE: paramter constraint for upper/lower bound is applied only for (PPRAM%slater_koster_type .le. 10) cases.
! KHJ, 19.Aug.2022
do j = 1, PPRAM%nparam_free
if(PPRAM%slater_koster_type .gt. 10) then
nsub = PPRAM%param_nsub(PPRAM%iparam_free(j))
i = PPRAM%iparam_free_nrl(j)
wa2(i:i+nsub-1) = PPRAM%param_nrl(1:nsub,PPRAM%iparam_free(j)) + wa1(i:i+nsub-1)
else
!wa2(j) = PPRAM%param(PPRAM%iparam_free(j)) + wa1(j)
if(PPRAM%param(PPRAM%iparam_free(j)) + wa1(j) .lt. PPRAM%param_const(3,PPRAM%iparam_free(j))) then ! lower bound
!wa2(j) = PPRAM%param(PPRAM%iparam_free(j))
!wa1(j) = 0.d0
wa2(j) = PPRAM%param_const(3,PPRAM%iparam_free(j))
wa1(j) = PPRAM%param_const(3,PPRAM%iparam_free(j)) - PPRAM%param(PPRAM%iparam_free(j))
elseif( PPRAM%param(PPRAM%iparam_free(j)) + wa1(j) .gt. PPRAM%param_const(2,PPRAM%iparam_free(j))) then ! upper bound
!wa2(j) = PPRAM%param(PPRAM%iparam_free(j))
!wa1(j) = 0.d0
wa2(j) = PPRAM%param_const(2,PPRAM%iparam_free(j))
wa1(j) = PPRAM%param_const(2,PPRAM%iparam_free(j)) - PPRAM%param(PPRAM%iparam_free(j))
else
wa2(j) = PPRAM%param(PPRAM%iparam_free(j)) + wa1(j)
endif
endif
enddo
wa3(1:nparam_free) = diag(1:nparam_free) * wa1(1:nparam_free)
pnorm = enorm (nparam_free, wa3)
! On the first iteration, adjust the initial step bound.
if ( iter == 1 ) delta = min (delta, pnorm)
! Evaluate the function at X + P and calculate its norm.
if(PPRAM%slater_koster_type .gt. 10) then
do j = 1, PPRAM%nparam_free
nsub = PPRAM%param_nsub(PPRAM%iparam_free(j)) ; i = PPRAM%iparam_free_nrl(j)
wa2_temp(i:i+nsub-1) = PPRAM%param_nrl(1:nsub,PPRAM%iparam_free(j)) ! store temp
PPRAM%param_nrl(1:nsub,PPRAM%iparam_free(j)) = wa2(i:i+nsub-1) ! update param
enddo
else
wa2_temp(1:nparam_free) = PPRAM%param(PPRAM%iparam_free) ! store temp
PPRAM%param(PPRAM%iparam_free) = wa2(1:nparam_free) ! update param
endif
call get_dE(wa4, fvec_plain, fvec_orb, ldjac, imode, PINPT, PPRAM, NN_TABLE, EDFT, ETBA_FIT, PWGHT, PGEOM, PKPTS, flag_fit_orbital)
nfev = nfev + 1
fnorm1 = enorm ( ldjac, wa4 )
if(PPRAM%slater_koster_type .gt. 10) then
do j = 1, PPRAM%nparam_free
nsub = PPRAM%param_nsub(PPRAM%iparam_free(j)) ; i = PPRAM%iparam_free_nrl(j)
PPRAM%param_nrl(1:nsub,PPRAM%iparam_free(j)) = wa2_temp(i:i+nsub-1) ! restore param
enddo
else
PPRAM%param(PPRAM%iparam_free) = wa2_temp(1:nparam_free) ! restore param
endif
! Compute the scaled actual reduction.
if ( 0.1D+00 * fnorm1 < fnorm ) then
actred = 1.0D+00 - (fnorm1/fnorm)**2
else
actred = -1.0D+00
endif
! Compute the scaled predicted reduction and the scaled directional derivative.
do j = 1, nparam_free
wa3(j) = 0.0D+00
l = ipvt(j)
temp = wa1(l)
wa3(1:j) = wa3(1:j) + fjac(1:j,j) * temp
enddo
temp1 = enorm ( nparam_free, wa3 ) / fnorm
temp2 = ( sqrt ( par ) * pnorm ) / fnorm
prered = temp1 ** 2 + temp2 ** 2 / 0.5D+00
dirder = - ( temp1 ** 2 + temp2 ** 2 )
! Compute the ratio of the actual to the predicted reduction.
ratio = 0.0D+00
if ( prered /= 0.0D+00 ) ratio = actred / prered
! Update the step bound.
if ( ratio <= 0.25D+00 ) then
if (actred >= 0.0D+00) temp = 0.5D+00
if (actred < 0.0D+00) temp = 0.5D+00*dirder/(dirder+0.5D+00*actred)
if (0.1D+00*fnorm1 >= fnorm .or. temp < 0.1D+00) temp = 0.1D+00
delta = temp * min ( delta, pnorm / 0.1D+00 )
par = par / temp
else
if ( par == 0.0D+00 .or. ratio >= 0.75D+00 ) then
delta = 2.0D+00 * pnorm
par = 0.5D+00 * par
endif
endif
! Test for successful iteration.
! Successful iteration. update X, FVEC, and their norms.
if ( 0.0001D+00 <= ratio ) then
do j = 1, PPRAM%nparam_free
if(PPRAM%slater_koster_type .gt. 10) then
nsub = PPRAM%param_nsub(PPRAM%iparam_free(j)) ; i = PPRAM%iparam_free_nrl(j)
PPRAM%param_nrl(1:nsub,PPRAM%iparam_free(j)) = wa2(i:i+nsub-1) ! save
wa2(i:i+nsub-1) = diag(i:i+nsub-1)*PPRAM%param_nrl(1:nsub,PPRAM%iparam_free(j)) ! update
else
PPRAM%param(PPRAM%iparam_free(j)) = wa2(j)
wa2(j) = diag(j) * PPRAM%param(PPRAM%iparam_free(j))
endif
enddo
fvec(1:ldjac) = wa4(1:ldjac)
xnorm = enorm ( nparam_free, wa2 )
fnorm = fnorm1
fnorm_plain = enorm ( ldjac , fvec_plain )
if(flag_fit_orbital) then
fnorm_orb = enorm ( ldjac , fvec_orb )
endif
iter = iter + 1
if(iter .le. PINPT%miter .and. flag_cost_history) then
PPRAM%cost_history(iter) = fnorm
PPRAM%niter = iter
endif
if(flag_write_info) then
write(message,'(A)')' ' ; write_msg
if(.not. flag_fit_orbital) then
write(message,'(A,I8, 2(A,F16.6))')' ITER=',iter,' ,rt(sum((dE*WT)^2)) = ',fnorm, &
' ,rt(sum(dE^2)) = ', fnorm_plain ; write_msg
elseif(flag_fit_orbital) then
write(message,'(A,I8, 3(A,F16.6))')' ITER=',iter,' ,rt(sum((dE*WT)^2)) = ',fnorm, &
' ,rt(sum(dE^2)) = ', fnorm_plain, &
' ,rt(sum(dORB^2)) =',fnorm_orb ; write_msg
endif
if_main write(pfileoutnm_temp,'(A,A)')trim(PPRAM%pfileoutnm),'_temp'
if_main call print_param(PINPT,PPRAM,PWGHT(1),pfileoutnm_temp,.TRUE.) ! only main system will be printed..
fnorm_ = fnorm ! fnorm of previous step
endif
endif
!if_main write(6,*)"ZZZZ "
!kill_job
! Tests for convergence.
if (abs(actred) <= ftol .and. prered <= ftol .and. 0.5D+00*ratio <= 1.0D+00) info=1
if (delta <= xtol*xnorm ) info=2
if (abs(actred) <= ftol .and. prered <= ftol .and. 0.5D+00*ratio <= 1.0D+00 .and. info == 2) info=3
if ( info /= 0 ) then
go to 300
endif
! Tests for termination and stringent tolerances.
if ( maxfev <= nfev ) info=5
if (abs(actred) <= epsmch .and. prered <= epsmch .and. 0.5D+00*ratio <= 1.0D+00) info=6
if (delta <= epsmch*xnorm) info=7
if (gnorm <= epsmch) info=4
if (info /= 0) then
go to 300
endif
! End of the inner loop. Repeat if iteration unsuccessful.
if ( ratio < 0.0001D+00 ) go to 200
! check kill
call check_kill_tbfit(PINPT,PPRAM, PWGHT)
! End of the outer loop.
go to 30
300 continue
! Termination, either normal or user imposed.
if(flag_write_info) then
if_main_then
if(info .eq. 1) then
write(message,*)' ' ; write_msg
write(message,101)" Termination INFO=",info,' , condition: |actred|,prered <= ftol, ratio <= 2 :',' actred=',actred, &
' prered=',prered, &
' ratio=',ratio ; write_msg
elseif(info .eq. 2) then
write(message,*)' ' ; write_msg
write(message,102)" Termination INFO=",info,' , condition: delta <= xtol*xnorm :',' delta=',delta, &
' xtol*xnorm=',xtol*xnorm ; write_msg
elseif(info .eq. 3) then
write(message,*)' ' ; write_msg
write(message,103)" Termination INFO=",info,' , condition: |actred|,prered <= ftol, delta<=xtol*xnorm :',' actred=',actred, &
' prered=',prered, &
' delta=',delta, &
' xtol*xnorm=',xtol*xnorm ; write_msg
elseif(info .eq. 5) then
write(message,*)' ' ; write_msg
write(message,104)" Termination INFO=",info,' , condition: maxfev(=miter*(nparam_free+1)) <= nfev :',' miter=',PINPT%miter, &
' nparam_free=',nparam_free, &
' nfev=',nfev ; write_msg
elseif(info .eq. 6) then
write(message,*)' ' ; write_msg
write(message,105)" Termination INFO=",info,' , condition: |actred|,prered <= epsmch, ratio <= 2 :',' actred=',actred, &
' prered=',prered, &
' epsmch=',epsmch, &
' ratio=',ratio ; write_msg
elseif(info .eq. 7) then
write(message,*)' ' ; write_msg
write(message,106)" Termination INFO=",info,' , condition: delta <= epsmch*xnorm :',' delta=',delta, &
' epsmch*xnorm=',epsmch*xnorm; write_msg
elseif(info .eq. 4) then
write(message,*)' ' ; write_msg
write(message,107)" Termination INFO=",info,' , condition: gnorm <= epsmch :',' gnorm=',gnorm, &
' epsmch=',epsmch; write_msg
else
write(message,*)' ' ; write_msg
write(message,108)" Termination INFO=",info ; write_msg
endif
if_main_end
endif
101 format(A,I2,A,3(A,E14.7))
102 format(A,I2,A,2(A,E14.7))
103 format(A,I2,A,4(A,E14.7))
104 format(A,I2,A,3(A,I6 ))
105 format(A,I2,A,4(A,E14.7))
106 format(A,I2,A,2(A,E14.7))
107 format(A,I2,A,2(A,E14.7))
108 format(A,I2)
if( PPRAM%slater_koster_type .gt. 10) then
call update_param_nrl( PPRAM )
else
call update_param( PPRAM )
endif
return
endsubroutine
subroutine fdjac2 (get_eig,NN_TABLE,ldjac,imode,PINPT,PPRAM,PGEOM,fvec,ETBA_FIT,EDFT,nparam_free,PWGHT,fjac,epsfcn,PKPTS)
!*****************************************************************************80
!
!! FDJAC2 estimates an M by N jacobian matrix using forward differences.
!
! Discussion:
!
! This subroutine computes a forward-difference approximation
! to the M by N jacobian matrix associated with a specified
! problem of M functions in N variables.
!
! Licensing:
!
! This code is distributed under the GNU LGPL license.
!
! Modified:
!
! 06 April 2010
! Author:
!
! Original FORTRAN77 version by Jorge More, Burton Garbow, Kenneth Hillstrom.
! FORTRAN90 version by John Burkardt.
!
! Modified: by Hyun-Jung Kim (KIAS, Infant@kias.re.kr)
! 14 December 2017
! Modified for the TBFIT purpose.
! The original source code can be found in: https://people.sc.fsu.edu/~jburkardt/f_src/minpack/minpack.html
! Reference:
!
! Jorge More, Burton Garbow, Kenneth Hillstrom,
! User Guide for MINPACK-1,
! Technical Report ANL-80-74,
! Argonne National Laboratory, 1980.
!
use parameters
use cost_function
use mpi_setup
use reorder_band
use projected_band
implicit none
type (incar ) :: PINPT
type (params ) :: PPRAM
type (hopping), dimension(PINPT%nsystem) :: NN_TABLE
type (weight) , dimension(PINPT%nsystem) :: PWGHT
type (energy) , dimension(PINPT%nsystem) :: ETBA_FIT, EDFT
type (kpoints), dimension(PINPT%nsystem) :: PKPTS
type (poscar) , dimension(PINPT%nsystem) :: PGEOM
integer*4 i,j,ii, nparam_free
integer*4 ldjac, imode
real*8 eps,epsfcn,epsmch,h,temp,fjac(ldjac,nparam_free)
real*8 wa(ldjac), fvec(ldjac), fvec_plain(ldjac), fvec_orb(ldjac)
logical flag_order_weight, flag_order, flag_fit_orbital
external get_eig
character*20, external :: int2str
flag_order = PINPT%flag_get_band_order .and. (.not. PINPT%flag_get_band_order_print_only)
flag_order_weight = .false.
flag_fit_orbital = .false. ! NOTE: even if the PINPT%flag_fit_orbital is .TRUE. in obtainding Jacobian, we will not use orbital norm
!flag_fit_orbital = .true. ! NOTE: even if the PINPT%flag_fit_orbital is .TRUE. in obtainding Jacobian, we will not use orbital norm
! since it is not well defined (not smooth) in derivation.
! Normally, in PSO routine, selecting Global or Personal best particle, the orbital norm can be used.
! See pso.f90 routine for the detail.. H.-J. Kim (11. Mar. 2021)
if(PPRAM%slater_koster_type .gt. 10) ii = 0
epsmch = epsilon(epsmch)
eps = sqrt( max(epsfcn,epsmch) )
if(PPRAM%slater_koster_type .gt. 10) then
do j = 1, PPRAM%nparam_free
do i = 1, PPRAM%param_nsub(PPRAM%iparam_free(j))
! update param with h added, original param will be saved to temp, and calculate wa
temp = PPRAM%param_nrl(i,PPRAM%iparam_free(j))
h = eps*abs(temp)
if (h == 0.0D+00 ) h=eps
PPRAM%param_nrl(i,PPRAM%iparam_free(j)) = temp+h
call get_dE(wa, fvec_plain, fvec_orb, ldjac, imode, PINPT, PPRAM, NN_TABLE, EDFT, ETBA_FIT, PWGHT, PGEOM, PKPTS, flag_fit_orbital)
! restore param from temp, and calculate derivation fjac from wa and fvec
PPRAM%param_nrl(i,PPRAM%iparam_free(j)) = temp
ii = ii + 1
fjac(:,ii) = ( wa(:) - fvec(:) ) / h
enddo
enddo
else
do j = 1, PPRAM%nparam_free
temp = PPRAM%param(PPRAM%iparam_free(j))
h = eps*abs(temp)
if (h == 0.0D+00 ) h=eps
PPRAM%param(PPRAM%iparam_free(j)) = temp+h
call get_dE(wa, fvec_plain, fvec_orb, ldjac, imode, PINPT, PPRAM, NN_TABLE, EDFT, ETBA_FIT, PWGHT, PGEOM, PKPTS, flag_fit_orbital)
PPRAM%param(PPRAM%iparam_free(j)) = temp
fjac(:,j) = ( wa(:) - fvec(:) ) / h
enddo
endif
return
endsubroutine