-
Notifications
You must be signed in to change notification settings - Fork 15
/
cost_function.f90
580 lines (533 loc) · 25.9 KB
/
cost_function.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
#include "alias.inc"
module cost_function
use mpi_setup
contains
subroutine get_dE(fvec, fvec_plain, fvec_orb, ldjac, imode, PINPT, PPRAM, NN_TABLE, EDFT, ETBA_FIT, PWGHT, PGEOM, PKPTS, flag_fit_orbital)
use parameters, only: weight, incar, energy, params, poscar, hopping, kpoints
use projected_band
use reorder_band
implicit none
type(incar) :: PINPT
type(params) :: PPRAM
type(hopping),dimension(PINPT%nsystem) :: NN_TABLE
type(poscar) ,dimension(PINPT%nsystem) :: PGEOM
type(weight) ,dimension(PINPT%nsystem) :: PWGHT
type(energy) ,dimension(PINPT%nsystem) :: ETBA_FIT, EDFT
type(kpoints),dimension(PINPT%nsystem) :: PKPTS
integer*4 mpierr
integer*4 i
integer*4 ldjac, my_ldjac
integer*4 imode
logical flag_order
logical flag_get_orbital
logical flag_order_weight
logical flag_fit_orbital
real*8 fvec(ldjac)
real*8 fvec_plain(ldjac)
real*8 fvec_orb(ldjac)
integer*4 ildjac, fldjac
flag_order_weight = .false. ! experimental feature
flag_order = PINPT%flag_get_band_order .and. (.not. PINPT%flag_get_band_order_print_only)
ildjac = 0
fldjac = 0
do i = 1, PINPT%nsystem
flag_get_orbital = (PWGHT(i)%flag_weight_orb .or. flag_order .or. flag_fit_orbital)
if(imode .eq. 1 .or. imode .eq. 11 .or. imode .eq. 13) then
my_ldjac = PKPTS(i)%nkpoint * PGEOM(i)%nband * PINPT%nspin
elseif(imode .eq. 2 .or. imode .eq. 12) then
my_ldjac = PKPTS(i)%nkpoint
endif
ildjac = fldjac + 1
fldjac = ildjac + my_ldjac - 1
! Evaluate the function at the starting point and calculate its norm.
call get_eig(NN_TABLE(i), PKPTS(i)%kpoint, PKPTS(i)%nkpoint, PINPT, PPRAM, &
ETBA_FIT(i)%E, ETBA_FIT(i)%V, ETBA_FIT(i)%SV, &
PGEOM(i)%neig, PGEOM(i)%init_erange, PGEOM(i)%nband, &
flag_get_orbital, .false., .false., PINPT%flag_phase)
! Evaluate degeneracy information for TBA band : after calling get_eig
call get_degeneracy(ETBA_FIT(i), PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint, PINPT)
! Evaluate ordered band for TBA band
if(flag_order) then
call get_ordered_band(ETBA_FIT(i), PKPTS(i), PGEOM(i), PWGHT(i), PINPT, flag_order_weight, PPRAM%flag_use_overlap)
endif
! Evaluate orbital projection for TBA band
if(flag_fit_orbital) then
call get_orbital_projection(ETBA_FIT(i), PKPTS(i), PINPT, PGEOM(i), PPRAM%flag_use_overlap)
endif
call get_cost_function(fvec(ildjac:fldjac), fvec_plain(ildjac:fldjac), fvec_orb(ildjac:fldjac), &
ETBA_FIT(i), EDFT(i), PGEOM(i), PINPT, PKPTS(i), PWGHT(i), my_ldjac, imode, flag_order, flag_fit_orbital)
enddo
return
endsubroutine
subroutine get_dE12(PINPT, PPRAM, NN_TABLE, EDFT, ETBA_FIT, PWGHT, PGEOM, PKPTS, flag_fit_orbital)
use parameters, only: weight, incar, energy, params, poscar, hopping, kpoints
use projected_band
use reorder_band
implicit none
type(incar) :: PINPT
type(params) :: PPRAM
type(hopping),dimension(PINPT%nsystem) :: NN_TABLE
type(poscar) ,dimension(PINPT%nsystem) :: PGEOM
type(weight) ,dimension(PINPT%nsystem) :: PWGHT
type(energy) ,dimension(PINPT%nsystem) :: ETBA_FIT, EDFT
type(kpoints),dimension(PINPT%nsystem) :: PKPTS
integer*4 mpierr
integer*4 i
integer*4 ldjac, my_ldjac
integer*4 imode
logical flag_order
logical flag_get_orbital
logical flag_order_weight
logical flag_fit_orbital
real*8 fvec(1)
real*8 fvec_plain(1)
real*8 fvec_orb(1)
integer*4 ildjac, fldjac
flag_order_weight = .false. ! experimental feature
flag_order = PINPT%flag_get_band_order .and. (.not. PINPT%flag_get_band_order_print_only)
imode = 12
do i = 1, PINPT%nsystem
flag_get_orbital = (PWGHT(i)%flag_weight_orb .or. flag_order .or. flag_fit_orbital)
! Evaluate the function at the starting point and calculate its norm.
call get_eig(NN_TABLE(i), PKPTS(i)%kpoint, PKPTS(i)%nkpoint, PINPT, PPRAM, &
ETBA_FIT(i)%E, ETBA_FIT(i)%V, ETBA_FIT(i)%SV, &
PGEOM(i)%neig, PGEOM(i)%init_erange, PGEOM(i)%nband, &
flag_get_orbital, .false., .false., PINPT%flag_phase)
! Evaluate degeneracy information for TBA band : after calling get_eig
call get_degeneracy(ETBA_FIT(i), PGEOM(i)%nband*PINPT%nspin, PKPTS(i)%nkpoint, PINPT)
! Evaluate ordered band for TBA band
if(flag_order) then
call get_ordered_band(ETBA_FIT(i), PKPTS(i), PGEOM(i), PWGHT(i), PINPT, flag_order_weight, PPRAM%flag_use_overlap)
endif
! Evaluate orbital projection for TBA band
if(flag_fit_orbital) then
call get_orbital_projection(ETBA_FIT(i), PKPTS(i), PINPT, PGEOM(i), PPRAM%flag_use_overlap)
endif
call get_cost_function(fvec(1), fvec_plain(1), fvec_orb(1), ETBA_FIT(i), EDFT(i), PGEOM(i), PINPT, PKPTS(i), PWGHT(i), 1, imode, flag_order, flag_fit_orbital)
enddo
return
endsubroutine
subroutine get_cost_function(fvec, fvec_plain, fvec_orb, ETBA_FIT, EDFT, PGEOM, PINPT, PKPTS, PWGHT, ldjac, imode, flag_order, flag_fit_orbital)
use parameters, only: weight, incar, energy, pi2, poscar, kpoints
use do_math, only: fgauss
use mpi_setup
implicit none
type (weight) :: PWGHT
type (incar ) :: PINPT
type (poscar) :: PGEOM
type (kpoints) :: PKPTS
type (energy) :: ETBA_FIT, EDFT
integer*4 i, ik
integer*4 is, ie, ie_, iband
integer*4 je, je_
integer*4 my_ik, sizebuff
integer*4 mpierr
integer*4 ldjac, imode
real*8 fvec(ldjac)
real*8 fvec_(ldjac)
real*8 fvec_plain(ldjac)
real*8 fvec_plain_(ldjac)
real*8 fvec_orb(ldjac)
real*8 fvec_orb_(ldjac)
real*8 OW, dD(3)
real*8 dE, dE_plain, dE_smooth, dORB_smooth
real*8 cost_func(PGEOM%nband*PINPT%nspin)
real*8 cost_func_plain(PGEOM%nband*PINPT%nspin)
real*8 sigma, maxgauss, max_gauss
real*8 enorm
external enorm
logical flag_weight_orbital, flag_fit_degeneracy
logical flag_order
real*8 E_TBA(PGEOM%nband*PINPT%nspin,PKPTS%nkpoint)
real*8 dE_TBA(PGEOM%nband*PINPT%nspin,PKPTS%nkpoint)
real*8 dORB_TBA(PGEOM%nband*PINPT%nspin,PKPTS%nkpoint)
real*8 E_DFT(PGEOM%neig*PINPT%ispin,PKPTS%nkpoint)
complex*16, allocatable :: myV(:,:,:)
real*8, allocatable :: myORB_TBA(:,:,:), myORB_DFT(:,:,:)
real*8, allocatable :: orb_TBA(:), orb_DFT(:,:)
real*8, allocatable :: dE_gauss(:) , dE_ref(:), dORB(:), dORB_ref(:)
logical flag_fit_orbital
integer*4 ie_cutoff
!integer*4 ourjob(nprocs), ourjob_disp(0:nprocs-1), my_i
integer*4 ncpu, id
integer*4, allocatable :: ourjob(:), ourjob_disp(:), my_i
! imode : 1, if ldjac < nparam (total number of parameters) -> unusual cases. try to avoid.
! ldjac = nkpoint * nband * PINPT%nspin
! imode : 2, if ldjac > nparam -> usual cases
! ldjac = nkpoint
! Important NOTE:
! imode = 1 ldjac != 1 : return fvec, fvec_plain (ldjac = nkpoint * nband)
! imode = 2 ldjac != 1 : return fvec, fvec_plain (ldjac = nkpoint)
! imode = 11 with ldjac != 1 : return ETBA%dE, fvec, fvec_plain (ldjac = nkpoint * nband)
! imode = 12 with ldjac != 1 : return ETBA%dE, fvec, fvec_plain (ldjac = nkpoint)
! imode = 13 with ldjac = 1 : return ETBA%dE
! imode = 13 with ldjac != 1 : return fvec, fvec_plain (ldjac = nkpoint * nband)
flag_fit_degeneracy = PINPT%flag_fit_degeneracy
flag_weight_orbital = PWGHT%flag_weight_orb
dE = 0d0 ; sigma = PINPT%orbital_fit_smearing ! default
cost_func = 0d0 ; cost_func_plain = 0d0
if(ldjac .ne. 1) then
fvec = 0d0 ; fvec_plain = 0d0 ; fvec_orb = 0d0
fvec_ = 0d0 ; fvec_plain_ = 0d0 ; fvec_orb_ = 0d0
endif
maxgauss = 1d0/(sigma*sqrt(2d0*pi2)) *real(PGEOM%neig)
max_gauss = fgauss(sigma, 0d0)
iband = PGEOM%init_erange
ie_cutoff = PWGHT%ie_cutoff
if(flag_fit_degeneracy) dD= 0d0
if(flag_weight_orbital) OW= 0d0
OW = 0d0
if(imode .gt. 10) dE_TBA = 0d0
if(imode .gt. 10) dORB_TBA = 0d0
if(COMM_KOREA%flag_split) then
ncpu = COMM_KOREA%nprocs
id = COMM_KOREA%myid
else
ncpu = nprocs
id = myid
endif
allocate(ourjob(ncpu))
allocate( ourjob_disp(0:ncpu-1))
call mpi_job_distribution_chain(PKPTS%nkpoint, ncpu, ourjob, ourjob_disp)
if(flag_order) then
E_TBA = ETBA_FIT%E_ORD
E_DFT = EDFT%E_ORD
if(flag_weight_orbital) then
sizebuff = PGEOM%neig*PINPT%ispinor*PGEOM%nband*PINPT%nspin
allocate(myV(PGEOM%neig*PINPT%ispinor, PGEOM%nband*PINPT%nspin, ourjob(id+1)))
#ifdef MPI
! Distribute wavefunction data to neighboring nodes
if(COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%V_ORD, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_COMPLEX16, myV, ourjob(id+1)*sizebuff, &
MPI_COMPLEX16, 0, COMM_KOREA%mpi_comm, mpierr)
elseif(.not. COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%V_ORD, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_COMPLEX16, myV, ourjob(id+1)*sizebuff, &
MPI_COMPLEX16, 0, mpi_comm_earth, mpierr)
endif
#else
myV = ETBA_FIT%V_ORD
#endif
endif
if(flag_fit_orbital) then
sizebuff = PINPT%lmmax * PGEOM%nband*PINPT%nspin
allocate(dE_gauss(PGEOM%nband))
allocate(dE_ref( PGEOM%nband))
allocate(dORB_ref(PGEOM%nband))
allocate(myORB_TBA(PINPT%lmmax, PGEOM%nband*PINPT%nspin, ourjob(id+1)))
allocate(myORB_DFT(PINPT%lmmax, PGEOM%nband*PINPT%nspin, ourjob(id+1)))
allocate(orb_TBA(PINPT%lmmax))
allocate(orb_DFT(PINPT%lmmax, PGEOM%nband))
! note: the size of second column, nband*PINPT%nspin is differ from read_energy routine
! but should work anyway, since ERANGE or EWINDOW tag should not be applied if
! fitting is requested by TBFIT tag of INCAR
#ifdef MPI
! Distribute orbital info to neighboring nodes
if(COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%ORB, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_TBA, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, COMM_KOREA%mpi_comm, mpierr)
call MPI_SCATTERV(EDFT%ORB , ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_DFT, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, COMM_KOREA%mpi_comm, mpierr)
elseif(.not. COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%ORB, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_TBA, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, mpi_comm_earth, mpierr)
call MPI_SCATTERV(EDFT%ORB , ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_DFT, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, mpi_comm_earth, mpierr)
endif
#else
myORB_TBA = ETBA_FIT%ORB
myORB_DFT = EDFT%ORB
#endif
endif
elseif(.not. flag_order) then
E_TBA = ETBA_FIT%E
E_DFT = EDFT%E
if(flag_weight_orbital) then
sizebuff = PGEOM%neig*PINPT%ispinor*PGEOM%nband*PINPT%nspin
allocate(myV(PGEOM%neig*PINPT%ispinor, PGEOM%nband*PINPT%nspin, ourjob(id+1)))
#ifdef MPI
! Distribute wavefunction data to neighboring nodes
if(COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%V, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_COMPLEX16, myV, ourjob(id+1)*sizebuff, &
MPI_COMPLEX16, 0, COMM_KOREA%mpi_comm, mpierr)
elseif(.not. COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%V, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_COMPLEX16, myV, ourjob(id+1)*sizebuff, &
MPI_COMPLEX16, 0, mpi_comm_earth, mpierr)
endif
#else
myV = ETBA_FIT%V
#endif
endif
if(flag_fit_orbital) then
sizebuff = PINPT%lmmax * PGEOM%nband*PINPT%nspin
allocate(dE_gauss(PGEOM%nband))
allocate(dE_ref( PGEOM%nband))
allocate(dORB_ref(PGEOM%nband))
allocate(dORB(PGEOM%nband))
allocate(myORB_TBA(PINPT%lmmax, PGEOM%nband*PINPT%nspin, ourjob(id+1)))
allocate(myORB_DFT(PINPT%lmmax, PGEOM%nband*PINPT%nspin, ourjob(id+1)))
allocate(orb_TBA(PINPT%lmmax))
allocate(orb_DFT(PINPT%lmmax, PGEOM%nband))
! note: the size of second column, nband*PINPT%nspin is differ from read_energy routine
! but should work anyway, since ERANGE or EWINDOW tag should not be applied if
! fitting is requested by TBFIT tag of INCAR
#ifdef MPI
! Distribute orbital info to neighboring nodes
if(COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%ORB, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_TBA, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, COMM_KOREA%mpi_comm, mpierr)
call MPI_SCATTERV(EDFT%ORB , ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_DFT, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, COMM_KOREA%mpi_comm, mpierr)
elseif(.not. COMM_KOREA%flag_split) then
call MPI_SCATTERV(ETBA_FIT%ORB, ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_TBA, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, mpi_comm_earth, mpierr)
call MPI_SCATTERV(EDFT%ORB , ourjob*sizebuff, ourjob_disp*sizebuff, &
MPI_REAL8, myORB_DFT, ourjob(id+1)*sizebuff, &
MPI_REAL8, 0, mpi_comm_earth, mpierr)
endif
#else
myORB_TBA = ETBA_FIT%ORB
myORB_DFT = EDFT%ORB
#endif
endif
endif
!#### IMODE = 1 : if nkpoint is less than number of free parameters to be fitted, this routine works
if(imode .eq. 1 .or. imode .eq. 11) then
do ik = sum(ourjob(1:id))+1, sum(ourjob(1:id+1))
my_i = (ik - 1) * PGEOM%nband * PINPT%nspin
my_ik= ik - sum(ourjob(1:id))
do is = 1, PINPT%nspin
do ie = 1, PGEOM%nband
ie_ = ie + iband - 1 + (is-1)*PGEOM%neig
dE_plain = E_TBA(ie+(is-1)*PGEOM%nband,ik) - E_DFT(ie_,ik)
dE = dE_plain * PWGHT%WT(ie_,ik)
if(flag_fit_degeneracy) dD= (ETBA_FIT%D(:,ie+(is-1)*PGEOM%nband,ik) - EDFT%D(:,ie_,ik)) * PWGHT%DEGENERACY_WT(ie_,ik)
if(flag_weight_orbital) OW = sum( PWGHT%PENALTY_ORB(:,ie_,ik)*abs(myV(:,ie+(is-1)*PGEOM%nband,my_ik)) )
if(imode .eq. 11) then
dE_TBA(ie+(is-1)*PGEOM%nband,ik) = dE_plain
else
fvec(my_i + ie+(is-1)*PGEOM%nband) = abs(dE) + (sum(abs(dD(:)))) + OW
if(ie_cutoff .gt. 0) then
if(ie .le. ie_cutoff) then
fvec_plain(my_i + ie+(is-1)*PGEOM%nband) = abs(dE_plain)
endif
else
fvec_plain(my_i + ie+(is-1)*PGEOM%nband) = abs(dE_plain)
endif
endif
enddo ! ie
enddo ! is
enddo ! ik
!#### IMODE = 2 ; ldjac > nparam -> usual cases
elseif(imode .eq. 2 .or. imode .eq. 12) then
do ik = sum(ourjob(1:id))+1, sum(ourjob(1:id+1))
my_ik = ik - sum(ourjob(1:id))
do is = 1, PINPT%nspin
do ie = 1, PGEOM%nband
ie_ = ie + iband - 1 + (is-1) * PGEOM%neig
dE_plain = E_TBA(ie+(is-1)*PGEOM%nband,ik) - E_DFT(ie_,ik)
if(flag_fit_degeneracy) dD = (ETBA_FIT%D(:,ie+(is-1)*PGEOM%nband,ik) - EDFT%D(:,ie_,ik)) * PWGHT%DEGENERACY_WT(ie_,ik)
if(flag_weight_orbital) OW = sum( PWGHT%PENALTY_ORB(:,ie_,ik)*abs(myV(:,ie+(is-1)*PGEOM%nband,my_ik)) )
dE = abs(dE_plain * PWGHT%WT(ie_,ik))
if(flag_fit_orbital) then
orb_TBA = myORB_TBA(:,ie+(is-1)*PGEOM%nband,my_ik)
orb_DFT = myORB_DFT(:,1+(is-1)*PGEOM%nband:PGEOM%nband+(is-1)*PGEOM%nband,my_ik)
dORB = fdORB(orb_TBA, orb_DFT, PINPT%lmmax, PGEOM%nband)
dORB_ref = fdORB(orb_DFT(:,ie), orb_DFT, PINPT%lmmax, PGEOM%nband)
endif
if(imode .eq. 12) then
dE_TBA(ie+(is-1)*PGEOM%nband,ik) = dE_plain
if(flag_fit_orbital) dORB_TBA(ie+(is-1)*PGEOM%nband,ik)= sum(abs(dORB - dORB_ref))
else
cost_func(ie+(is-1)*PGEOM%nband) = dE + sum(abs(dD)) + OW
if(ie_cutoff .gt. 0) then
if(ie .le. ie_cutoff) then
cost_func_plain(ie+(is-1)*PGEOM%nband) = abs(dE_plain)
endif
else
cost_func_plain(ie+(is-1)*PGEOM%nband) = abs(dE_plain)
endif
endif
enddo
enddo ! is
if(imode .eq. 2) then
fvec(ik) = enorm(PGEOM%nband*PINPT%nspin, cost_func)
fvec_plain(ik) = enorm(PGEOM%nband*PINPT%nspin, cost_func_plain)
endif
enddo ! ik
elseif(imode .eq. 13) then
do ik = sum(ourjob(1:id))+1, sum(ourjob(1:id+1))
my_i = (ik - 1) * PGEOM%nband * PINPT%nspin
my_ik= ik - sum(ourjob(1:id))
do is = 1, PINPT%nspin
do ie = 1, PGEOM%nband
ie_ = ie + iband - 1 + (is-1) * PGEOM%neig
if(flag_weight_orbital) then
OW = sum( PWGHT%PENALTY_ORB(:,ie_,ik)*abs(myV(:,ie+(is-1)*PGEOM%nband,my_ik)) )
!if(ie_ .eq. 11 .and. ik .eq. 90) write(6,*)"ZZZ ", abs(myV(:,ie+(is-1)*PGEOM%nband,my_ik))
endif
dE_plain = E_TBA(ie+(is-1)*PGEOM%nband,ik) - E_DFT(ie_,ik)
dE = dE_plain * PWGHT%WT(ie_,ik) + OW
dE_TBA(ie+(is-1)*PGEOM%nband,ik) = dE_plain
if(flag_fit_orbital) then
! NOTE: what is the best way to make the orbital different to be a "smooth" function?
! orbital : s px py pz dz2 dx2 dxy dxz dyz"
orb_TBA = myORB_TBA(:,ie+(is-1)*PGEOM%nband,my_ik)
orb_DFT = myORB_DFT(:,1+(is-1)*PGEOM%nband:PGEOM%nband+(is-1)*PGEOM%nband,my_ik)
dORB = fdORB(orb_TBA, orb_DFT, PINPT%lmmax, PGEOM%nband)
dORB_ref = fdORB(orb_DFT(:,ie), orb_DFT, PINPT%lmmax, PGEOM%nband)
dORB_smooth= sum(abs(dORB - dORB_ref)) * PWGHT%WT(ie_,ik)
endif
if(ldjac .ne. 1) then
if(ie_cutoff .gt. 0) then
if(ie .le. ie_cutoff) then
fvec(my_i + ie+(is-1)*PGEOM%nband) = abs(dE)
fvec_plain(my_i + ie+(is-1)*PGEOM%nband) = abs(dE_plain)
if(flag_fit_orbital) then
fvec_orb(my_i + ie+(is-1)*PGEOM%nband) = dORB_smooth
endif
endif
else
fvec(my_i + ie+(is-1)*PGEOM%nband) = abs(dE)
fvec_plain(my_i + ie+(is-1)*PGEOM%nband) = abs(dE_plain)
if(flag_fit_orbital) then
fvec_orb(my_i + ie+(is-1)*PGEOM%nband) = dORB_smooth
endif
endif
endif
enddo ! ie
enddo ! is
enddo ! ik
endif ! IMODE ?
if(imode .lt. 10) then
#ifdef MPI
if(COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(fvec, fvec_, size(fvec), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
call MPI_ALLREDUCE(fvec_plain, fvec_plain_, size(fvec_plain), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
elseif(.not. COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(fvec, fvec_, size(fvec), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
call MPI_ALLREDUCE(fvec_plain, fvec_plain_, size(fvec_plain), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
endif
fvec = fvec_
fvec_plain = fvec_plain_
#endif
elseif(imode .gt. 10 .and. imode .le. 12) then
#ifdef MPI
if(COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(dE_TBA, ETBA_FIT%dE, size(dE_TBA), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
if(flag_fit_orbital) then
call MPI_ALLREDUCE(dORB_TBA, ETBA_FIT%dORB, size(dORB_TBA), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
endif
if(ldjac .ne. 1) then
call MPI_ALLREDUCE(fvec, fvec_, size(fvec), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
fvec = fvec_
call MPI_ALLREDUCE(fvec_plain, fvec_plain_, size(fvec_plain), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
fvec_plain = fvec_plain_
if(flag_fit_orbital) then
call MPI_ALLREDUCE(fvec_orb, fvec_orb_, size(fvec_orb), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
fvec_orb = fvec_orb_
endif
endif
elseif(.not. COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(dE_TBA, ETBA_FIT%dE, size(dE_TBA), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
if(flag_fit_orbital) then
call MPI_ALLREDUCE(dORB_TBA, ETBA_FIT%dORB, size(dORB_TBA), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
endif
if(ldjac .ne. 1) then
call MPI_ALLREDUCE(fvec, fvec_, size(fvec), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
fvec = fvec_
call MPI_ALLREDUCE(fvec_plain, fvec_plain_, size(fvec_plain), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
fvec_plain = fvec_plain_
if(flag_fit_orbital) then
call MPI_ALLREDUCE(fvec_orb, fvec_orb_, size(fvec_orb), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
fvec_orb = fvec_orb_
endif
endif
endif
#else
ETBA_FIT%dE = dE_TBA
if(flag_fit_orbital) then
ETBA_FIT%dORB = dORB_TBA
endif
#endif
elseif(imode .ge. 13) then
#ifdef MPI
if(COMM_KOREA%flag_split) then
if(ldjac .ne. 1) then
call MPI_ALLREDUCE(fvec, fvec_, size(fvec), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
fvec = fvec_
call MPI_ALLREDUCE(fvec_plain, fvec_plain_, size(fvec_plain), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
fvec_plain = fvec_plain_
if(flag_fit_orbital) then
call MPI_ALLREDUCE(fvec_orb, fvec_orb_, size(fvec_orb), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
fvec_orb = fvec_orb_
endif
elseif(ldjac .eq. 1) then
call MPI_ALLREDUCE(dE_TBA, ETBA_FIT%dE, size(dE_TBA), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
endif
elseif(.not. COMM_KOREA%flag_split) then
if(ldjac .ne. 1) then
call MPI_ALLREDUCE(fvec, fvec_, size(fvec), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
fvec = fvec_
call MPI_ALLREDUCE(fvec_plain, fvec_plain_, size(fvec_plain), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
fvec_plain = fvec_plain_
if(flag_fit_orbital) then
call MPI_ALLREDUCE(fvec_orb, fvec_orb_, size(fvec_orb), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
fvec_orb = fvec_orb_
endif
elseif(ldjac .eq. 1) then
call MPI_ALLREDUCE(dE_TBA, ETBA_FIT%dE, size(dE_TBA), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
endif
endif
#else
if(ldjac .eq. 1) then
ETBA_FIT%dE = dE_TBA
endif
#endif
endif
return
endsubroutine
function fdORB(myORB_TBA, myORB_DFT, lmmax, nband)
use do_math, only: fgauss
implicit none
integer*4 ie
integer*4 lmmax
integer*4 nband
real*8 myORB_TBA(lmmax)
real*8 myORB_DFT(lmmax, nband)
real*8 fdORB(nband)
external enorm
real*8 enorm
do ie = 1, nband
fdORB(ie) = dot_product(myORB_TBA, myORB_DFT(:,ie))/(sqrt(dot_product(myORB_TBA, myORB_TBA)) * sqrt(dot_product(myORB_DFT(:,ie),myORB_DFT(:,ie))))
enddo
return
endfunction
! function fdORB2(E_DFT, myORB_TBA, myORB_DFT, lmmax, nband)
! use do_math, only: fgauss
! implicit none
! integer*4 ie
! integer*4 lmmax
! integer*4 nband
! real*8 myORB_TBA(lmmax)
! real*8 myORB_DFT(lmmax, nband)
! real*8 fdORB(nband)
! real*8 E_DFT(nband)
! external enorm
! real*8 enorm
! do ie = 1, nband
! fdORB(ie) = dot_product(myORB_TBA, myORB_DFT(:,ie))/(sqrt(dot_product(myORB_TBA, myORB_TBA)) * sqrt(dot_product(myORB_DFT(:,ie),myORB_DFT(:,ie))))
! enddo
! return
! endfunction
endmodule