fleur.F90 30.6 KB
Newer Older
1
2
3
4
5
!--------------------------------------------------------------------------------
! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
! This file is part of FLEUR and available as free software under the conditions
! of the MIT license as expressed in the LICENSE file in more detail.
!--------------------------------------------------------------------------------
6
MODULE m_fleur
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
  IMPLICIT NONE
CONTAINS
  SUBROUTINE fleur_execute(mpi_comm)

    !     ***************************************************************
    !
    !     based on flapw7 (c.l.fu, m.weinert, e.wimmer):
    !     full potential linearized augmented plane wave method for thin
    !     films and superlattices (version 7 ---- general symmetry)
    !     symmetry part       ---  e.wimmer
    !     potential generator ---  c.l.fu,r.podloucky
    !     matrix elements     ---  m.weinert
    !     charge density      ---  c.l.fu
    !                                c.l.fu        1987
    !     2nd variation diagon.  --- r.-q. wu      1992
    !     forces a la Yu et al   --- r.podloucky   1995
    !     full relativistic core --- a.shick       1996
    !     broyden mixing         --- r.pentcheva   1996
    !     gga (pw91, pbe)        --- t.asada       1997
    !     local orbitals         --- p.kurz        1997
    !     automatic symmetry     --- w.hofer       1997
    !     core tails & start     --- r.abt         1998
    !     spin orbit coupling    --- a.shick,x.nie 1998
    !     non-colinear magnet.   --- p.kurz        1999
    !     one-dimensional        --- y.mokrousov   2002
    !     exchange parameters    --- m.lezaic      2004
    !
    !                       g.bihlmayer, s.bluegel 1999
    !     ***************************************************************
    !----------------------------------------
    ! this routine is the main PROGRAM

    USE m_types
    USE m_constants
    USE m_fleur_init
    USE m_pldngen
    USE m_optional
44
45
    USE m_cdn_io
    USE m_qfix
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
    USE m_vgen
    USE m_rhodirgen
    USE m_writexcstuff
    USE m_vmatgen
    USE m_eigen
    USE m_eigenso
    USE m_fermie
    USE m_force0
    USE m_cdngen
    USE m_totale
    USE m_potdis
    USE m_mix
    USE m_xmlOutput
    USE m_juDFT_time
    USE m_calc_hybrid
    !          USE m_jcoff
    !          USE m_jcoff2
    !          USE m_ssomat
    USE m_wann_optional
    USE m_wannier
    USE m_bs_comfort

    USE m_gen_map
    USE m_dwigner
    !          USE m_generate_pntgpt
    !          USE m_rotate_eig
    USE m_ylm
73
#ifdef CPP_MPI
74
    USE m_mpi_bc_all,  ONLY : mpi_bc_all
75
    USE m_mpi_bc_potden
76
#endif
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
    USE m_eig66_io,   ONLY : open_eig, close_eig
    IMPLICIT NONE

    INTEGER,INTENT(IN) :: mpi_comm

    !     Types, these variables contain a lot of data!
    TYPE(t_input)    :: input
    TYPE(t_dimension):: DIMENSION
    TYPE(t_atoms)    :: atoms
    TYPE(t_sphhar)   :: sphhar
    TYPE(t_cell)     :: cell
    TYPE(t_stars)    :: stars
    TYPE(t_sym)      :: sym
    TYPE(t_noco)     :: noco
    TYPE(t_vacuum)   :: vacuum
    TYPE(t_sliceplot):: sliceplot
    TYPE(t_banddos)  :: banddos
    TYPE(t_obsolete) :: obsolete
    TYPE(t_enpara)   :: enpara
    TYPE(t_xcpot)    :: xcpot
    TYPE(t_results)  :: results
    TYPE(t_jij)      :: jij
    TYPE(t_kpts)     :: kpts
    TYPE(t_hybrid)   :: hybrid
    TYPE(t_oneD)     :: oneD
    TYPE(t_mpi)      :: mpi
103
    TYPE(t_coreSpecInput) :: coreSpecInput
104
105
    TYPE(t_wann)     :: wann
    TYPE(t_potden)   :: v,vx
106
    TYPE(t_potden)   :: inDen, outDen, mixDen
107
108

    !     .. Local Scalars ..
109
    INTEGER:: eig_id, archiveType
110
    INTEGER:: n,it,ithf,pc
111
    LOGICAL:: stop80,reap,l_endit,l_opti,l_cont,l_qfix, l_error
112
    REAL   :: fermiEnergyTemp, fix
113
114
115
116
117
118
119
    !--- J<
    INTEGER             :: phn
    REAL, PARAMETER     :: tol = 1.e-8
    INTEGER             :: qcount ,imt,i_J,j_J
    !--- J>
    !     HF/hybrid-functionals/EXX
    LOGICAL               ::  l_restart
120
#ifdef CPP_MPI
121
122
    INCLUDE 'mpif.h'
    INTEGER:: ierr(2)
123
#endif
124
125
    INTEGER, ALLOCATABLE :: eig_idList(:)
    mpi%mpi_comm = mpi_comm
Daniel Wortmann's avatar
Daniel Wortmann committed
126

127
128
    CALL timestart("Initialization")
    CALL fleur_init(mpi,input,DIMENSION,atoms,sphhar,cell,stars,sym,noco,vacuum,&
129
130
                    sliceplot,banddos,obsolete,enpara,xcpot,results,jij,kpts,hybrid,&
                    oneD,coreSpecInput,wann,l_opti)
131
    CALL timestop("Initialization")
Daniel Wortmann's avatar
Daniel Wortmann committed
132
133


134

135
136
137
138
    IF (l_opti) THEN
       IF (sliceplot%iplot .AND. (mpi%irank==0) ) THEN
          IF (noco%l_noco) THEN
             CALL pldngen(sym,stars,atoms,sphhar,vacuum,&
139
                          cell,input,noco,oneD,sliceplot)
140
          ENDIF
141
142
       ENDIF
       CALL OPTIONAL(mpi,atoms,sphhar,vacuum,DIMENSION,&
143
                     stars,input,sym,cell,sliceplot,obsolete,xcpot,noco,oneD)
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
    ENDIF
    !
    IF (sliceplot%iplot)      CALL juDFT_end("density plot o.k.",mpi%irank)
    IF (input%strho)          CALL juDFT_end("starting density generated",mpi%irank)
    IF (input%swsp)           CALL juDFT_end("spin polarised density generated",mpi%irank)
    IF (input%lflip)          CALL juDFT_end("magnetic moments flipped",mpi%irank)
    IF (input%l_bmt)          CALL juDFT_end('"cdnbmt" written',mpi%irank)


    !+Wannier
    input%l_wann = .FALSE.
    INQUIRE (file='wann_inp',exist=input%l_wann)
    IF (input%l_wann.AND.(mpi%irank==0).AND.(.NOT.wann%l_bs_comf)) THEN
       IF(mpi%isize.NE.1) CALL juDFT_error('No Wannier+MPI at the moment',calledby = 'fleur')
       CALL wann_optional(input,atoms,sym,cell,oneD,noco,wann)
    END IF
    IF (wann%l_gwf) input%itmax = 1
161
    
162
163
164
165
166
167
168
169
170
171
172
173
    !-Wannier

    ALLOCATE (eig_idList(wann%nparampts))

    l_restart = .TRUE.

    it     = 0
    ithf   = 0
    l_cont = ( it < input%itmax )
    results%last_distance = -1.0
    IF (mpi%irank.EQ.0) CALL openXMLElementNoAttributes('scfLoop')

174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
    ! Initialize and load inDen density without density matrix(start)
    CALL inDen%init(stars,atoms,sphhar,vacuum,oneD,DIMENSION%jspd,.FALSE.)
    IF (noco%l_noco) THEN
       ALLOCATE (inDen%cdom(stars%ng3),inDen%cdomvz(vacuum%nmzd,2))
       ALLOCATE (inDen%cdomvxy(vacuum%nmzxyd,oneD%odi%n2d-1,2))
       archiveType = CDN_ARCHIVE_TYPE_NOCO_const
    ELSE
       ALLOCATE (inDen%cdom(1),inDen%cdomvz(1,1),inDen%cdomvxy(1,1,1))
       archiveType = CDN_ARCHIVE_TYPE_CDN1_const
    END IF
    IF(mpi%irank.EQ.0) THEN
       CALL readDensity(stars,vacuum,atoms,cell,sphhar,input,sym,oneD,archiveType,CDN_INPUT_DEN_const,&
                        0,fermiEnergyTemp,l_qfix,inDen%iter,inDen%mt,inDen%pw,inDen%vacz,inDen%vacxy,&
                        inDen%cdom,inDen%cdomvz,inDen%cdomvxy)
       CALL timestart("Qfix")
       CALL qfix(stars,atoms,sym,vacuum, sphhar,input,cell,oneD,inDen%pw,inDen%vacxy,inDen%mt,inDen%vacz,&
                 .FALSE.,.false.,fix)
       CALL timestop("Qfix")
       CALL writeDensity(stars,vacuum,atoms,cell,sphhar,input,sym,oneD,archiveType,CDN_INPUT_DEN_const,&
                         0,-1.0,0.0,.FALSE.,inDen%iter,inDen%mt,inDen%pw,inDen%vacz,inDen%vacxy,inDen%cdom,&
                         inDen%cdomvz,inDen%cdomvxy)
    END IF
    ! Initialize and load inDen density without density matrix(end)

198
199
200
201
202
203
204
205
206
207
    DO WHILE (l_cont)

       it = it + 1
       !+t3e
       IF (input%alpha.LT.10.0) THEN

          IF (it.GT.1) THEN
             obsolete%pot8 = .FALSE.
             input%alpha = input%alpha - NINT(input%alpha)
          END IF
208

209
210
211
212
213
214
215
216
          CALL resetIterationDependentTimers()
          CALL timestart("Iteration")
          IF (mpi%irank.EQ.0) THEN
             WRITE(*,"(a,i3)",advance="no") "Iteration:",it
             !-t3e
             WRITE (6,FMT=8100) it
             WRITE (16,FMT=8100) it
8100         FORMAT (/,10x,'   it=    ',i5)
217

218
219
220
             IF (.NOT.obsolete%pot8) THEN

                !      ----> potential generator
221
                !
222
223
224
225
226
227
228
229
                !---> pk non-collinear
                !--->        reload the density matrix from file rhomat_in
                !--->        calculate spin-up and -down density for USE in the
                !--->        potential generator and store the direction of
                !--->        magnetization on file dirofmag
                IF (noco%l_noco) THEN
                   CALL timestart("gen. spin-up and -down density")
                   CALL rhodirgen(DIMENSION,sym,stars,atoms,sphhar,&
230
                                  vacuum,22,cell,input,oneD)
231
232
233
                   CALL timestop("gen. spin-up and -down density")
                ENDIF
                !---> pk non-collinear
234

235
236
237
238
                reap=.NOT.obsolete%disp
                input%total = .TRUE.
             ENDIF!(obsolete%pot8)
          ENDIF !mpi%irank.eq.0
239
#ifdef CPP_MPI
240
          CALL MPI_BCAST(input%total,1,MPI_LOGICAL,0,mpi%mpi_comm,ierr)
241
#endif
242

243
244
245
246
247
248
249
250
251
252
253
          !--- J<
          IF(jij%l_jenerg) GOTO 234

          jij%alph1(:)=noco%alph(:)
          stop80= .FALSE.
          IF ((noco%l_soc.AND.noco%l_ss)) THEN
             IF ((.NOT.wann%l_gwf).AND.&
                  ((jij%l_J).OR.(jij%nqpt/=1).OR.(jij%nmagn/=1).OR.(jij%phnd/=1))) THEN
                CALL juDFT_error("fleur: J-loop with ss+soc", calledby ="fleur")
             END IF
          END IF
254

255
256
257
          ! Initialize and load inDen density matrix and broadcast inDen(start)
          WRITE(*,*) 'test-1'
          IF (ALLOCATED(inDen%mmpMat)) DEALLOCATE (inDen%mmpMat)
258
          IF(mpi%irank.EQ.0) THEN
259
260
261
262
             IF (isDensityMatrixPresent().AND.atoms%n_u>0) THEN
                ALLOCATE (inDen%mmpMat(-lmaxU_const:lmaxU_const,-lmaxU_const:lmaxU_const,atoms%n_u,input%jspins))
                CALL readDensityMatrix(input,atoms,inDen%mmpMat,l_error)
                IF(l_error) CALL juDFT_error('Error in reading density matrix!',calledby='fleur')
263
264
265
             ELSE
                ALLOCATE (inDen%mmpMat(-lmaxU_const:-lmaxU_const,-lmaxU_const:-lmaxU_const,1,input%jspins))
                inDen%mmpMat = CMPLX(0.0,0.0)
266
             END IF
267
          END IF
268
          WRITE(*,*) 'test-2'
269
270
271
#ifdef CPP_MPI
          CALL mpi_bc_potden(mpi,stars,sphhar,atoms,input,vacuum,oneD,noco,inDen)
#endif
272
          ! Initialize and load inDen density matrix and broadcast inDen(end)
273

274
275
276
277
278
279
280
281
          DO qcount=1,jij%nqpt
             IF (jij%l_J) THEN
                noco%qss(:)=jij%qj(:,qcount)
                jij%qn = ( noco%qss(1)**2 + noco%qss(2)**2 + noco%qss(3)**2 )
             END IF
             IF ((input%l_wann.OR.jij%l_J).AND.(mpi%irank.EQ.0)) THEN
                WRITE(6,*) 'qss=(',noco%qss(1),',',noco%qss(2),',',noco%qss(3),')'
                CALL timestart("Q-point for J_ij(total)")
282

283
             ENDIF
Daniel Wortmann's avatar
Daniel Wortmann committed
284

285

286
287
288
289
290
291
             IF ( noco%l_soc ) THEN
                dimension%neigd2 = dimension%neigd*2
             ELSE
                dimension%neigd2 = dimension%neigd
             END IF

292

293
294
             !HF
             IF (hybrid%l_hybrid) CALL  calc_hybrid(hybrid,kpts,atoms,input,DIMENSION,mpi,noco,&
295
                                                    cell,vacuum,oneD,banddos,results,sym,xcpot,v,it)
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
             !#endif

             DO pc = 1, wann%nparampts
                !---> gwf
                IF (wann%l_sgwf.OR.wann%l_ms) THEN
                   noco%qss(:) = wann%param_vec(:,pc)
                   jij%qn = (noco%qss(1)**2 + noco%qss(2)**2 + noco%qss(3)**2)
                   noco%alph(:) = wann%param_alpha(:,pc)
                ELSE IF (wann%l_socgwf) THEN
                   IF(wann%l_dim(2)) noco%phi   = tpi_const * wann%param_vec(2,pc)
                   IF(wann%l_dim(3)) noco%theta = tpi_const * wann%param_vec(3,pc)
                END IF
                !---< gwf

                IF (.NOT.obsolete%pot8) THEN
                   CALL timestart("generation of potential")
                   IF (mpi%irank==0) WRITE(*,"(a)",advance="no") " * Potential generation "
                   CALL vgen(hybrid,reap,input,xcpot,DIMENSION, atoms,sphhar,stars,vacuum,&
314
                        sym,obsolete,cell, oneD,sliceplot,mpi ,results,noco,inDen,v,vx)
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
                   CALL timestop("generation of potential")

                   IF (mpi%irank.EQ.0) THEN
                      !---> pk non-collinear
                      !--->          generate the four component matrix potential from spin up
                      !--->          and down potentials and direction of the magnetic field
                      IF (noco%l_noco) THEN
                         CALL timestart("generation of potential-matrix")
                         CALL vmatgen(stars, atoms,sphhar,vacuum,sym,input,oneD,8,22,26)
                         CALL timestop("generation of potential-matrix")
                      ENDIF
                      !---> end pk non-collinear
                      !---> do some output for the tddft calculations:
                      IF (input%gw /= 0) THEN
                         CALL write_xcstuff(sphhar,atoms,DIMENSION,sym, stars,vacuum,input)
                      ENDIF
331
                      !
332
333
334
335
336
                   ENDIF ! mpi%irank.eq.0

                   !
                   !+t3e
                ENDIF ! .not.obsolete%pot8
Daniel Wortmann's avatar
Daniel Wortmann committed
337

Daniel Wortmann's avatar
Daniel Wortmann committed
338
339


340
#ifdef CPP_MPI
341
                CALL MPI_BARRIER(mpi%mpi_comm,ierr)
342
343
#endif

344
345
346
347
348
349
350
351
352
353
354
355
356
                !
                !          ----> eigenvalues and eigenfunctions
                !
                !--- J<
                IF(jij%l_disp)THEN
                   jij%mtypes=1
                   jij%nmagn=1
                   jij%phnd=1
                ENDIF
                i_J=1
                DO imt=1,jij%mtypes
                   DO j_J=i_J,jij%nmagn
                      DO phn=1,jij%phnd
357

Daniel Wortmann's avatar
Daniel Wortmann committed
358

359
360
361
                         input%eigvar(1)= .TRUE.
                         input%eigvar(2)= .TRUE.
                         input%eigvar(3)= .TRUE.
Daniel Wortmann's avatar
Daniel Wortmann committed
362

363
364
                         input%eigvar(2)= input%eigvar(2) .AND. ( noco%l_soc .AND. (.NOT.noco%l_noco) )
                         ! eigvar(1/2)= 1st/2nd var. ; eigvar(3)= calc density,etc
365

366
367
368
369
370
371
372
                         IF ( noco%l_soc .AND. (.NOT.noco%l_noco) ) THEN
                            input%evonly(1)= .FALSE.
                            input%evonly(2)= input%eonly
                         ELSE
                            input%evonly(1)= input%eonly
                            input%evonly(2)= .FALSE.
                         ENDIF
373

374
375
376
                         IF ( input%eigvar(1).OR.input%eigvar(2) ) THEN
                            IF (jij%l_J) THEN
                               input%tkb=0.
377
#ifdef CPP_NEVER
378
                               CALL jcoff(i_J,j_J,phn,mpi,atoms,atoms, noco,jij)
379
#endif
380
381
382
383
384
385
386
387
388
                            ENDIF
                            IF (input%eigvar(1)) THEN
                               CALL timestart("generation of hamiltonian and diagonalization (total)")

                               ! WRITE(6,fmt='(A)') 'Starting 1st variation ...'
                               CALL timestart("eigen")
                               IF (mpi%irank==0) WRITE(*,"(a)",advance="no") "* Eigenvalue problem "
                               CALL eigen(mpi,stars,sphhar,atoms,obsolete,xcpot,&
                                    sym,kpts,DIMENSION,vacuum,input,cell,enpara,banddos,noco,jij,oneD,hybrid,&
389
                                    it,eig_id,inDen,results,v,vx)
390
391
392
393
394
                               eig_idList(pc) = eig_id
                               CALL timestop("eigen")
                               !
                               !                   add all contributions to total energy
                               !
395
#ifdef CPP_MPI
396
397
398
399
400
401
402
403
404
405
406
407
408
                               ! send all result of local total energies to the r
                               IF (mpi%irank==0) THEN
                                  CALL MPI_Reduce(MPI_IN_PLACE,results%te_hfex%valence,&
                                       1,MPI_REAL8,MPI_SUM,0,mpi%mpi_comm,ierr(1))
                                  CALL MPI_Reduce(MPI_IN_PLACE,results%te_hfex%core,&
                                       1,MPI_REAL8,MPI_SUM,0,mpi%mpi_comm,ierr(1))
                               ELSE
                                  CALL MPI_Reduce(results%te_hfex%valence,MPI_IN_PLACE,&
                                       1,MPI_REAL8,MPI_SUM,0, mpi%mpi_comm,ierr(1))
                                  CALL MPI_Reduce(results%te_hfex%core,MPI_IN_PLACE,&
                                       1,MPI_REAL8,MPI_SUM,0, mpi%mpi_comm,ierr(1))
                               ENDIF
                               !                                  END IF
409
410
#endif

Daniel Wortmann's avatar
Daniel Wortmann committed
411

412
413
414
415
416
417
418
419
420

                            ENDIF
                            IF (input%eigvar(2)) THEN
                               ! RS: open unit for SOC vectors for GW
                               IF(noco%l_soc.AND.input%gw.EQ.2) THEN
                                  WRITE(6,'(A)') 'RS: open SOCVEC unit 4649'
                                  OPEN(4649,file='SOCVEC',form='unformatted')
                               ENDIF
                               ! WRITE(6,fmt='(A)') 'Starting 2nd variation ...'
421
                               CALL eigenso(eig_id,mpi,DIMENSION,stars,vacuum,atoms,sphhar,&
422
423
424
425
426
427
428
                                    obsolete,sym,cell,noco,input,kpts, oneD)
                               IF(noco%l_soc.AND.input%gw.EQ.2) THEN
                                  CLOSE(4649)
                                  INQUIRE(1014,opened=l_endit)
                                  IF(l_endit) CLOSE(1014)
                                  INQUIRE(667,opened=l_endit)
                                  IF(l_endit) CLOSE(667)
429
                                  CALL juDFT_end("GW+SOC finished",mpi%irank)
430
431
432
433
434
435
436
437
438
439
440
441
                               ENDIF
                            ENDIF
                            CALL timestop("generation of hamiltonian and diagonalization (total)")

#ifdef CPP_MPI
                            CALL MPI_BARRIER(mpi%mpi_comm,ierr)
#endif

                         ENDIF ! ( input%eigvar(1) .OR. input%eigvar(2) )

                         IF ( input%eigvar(3) .AND. noco%l_soc .AND. noco%l_ss ) THEN
#ifdef CPP_NEVER
442
                            CALL ssomat(eig_id, mpi,DIMENSION,stars,vacuum,atoms,sphhar,&
443
444
445
446
447
448
449
450
                                 sym,cell,noco,input,obsolete,kpts,oneD,MPI_DOUBLE_PRECISION )
#endif
                            stop80= .TRUE.
                         ENDIF
                         !-t3e
                         !
                         !              ----> fermi level and occupancies
                         !
451
                         IF (mpi%irank==0) WRITE(*,"(a)",advance="no") "* Fermi-energy "
452
453
454
                         IF ( input%eigvar(3) .AND. ( .NOT.(noco%l_soc .AND. noco%l_ss) ) ) THEN
                            IF (jij%l_J) THEN
                               CALL timestart("determination of fermi energy")
Daniel Wortmann's avatar
Daniel Wortmann committed
455
                               ALLOCATE ( results%w_iks(dimension%neigd,kpts%nkpt,dimension%jspd) )
456
                               CALL fermie(eig_id, mpi,kpts,obsolete,input,&
Daniel Wortmann's avatar
Daniel Wortmann committed
457
                                    noco,enpara%epara_min,jij,cell,results)
458
459
460
461
462
                               DEALLOCATE ( results%w_iks )
                               CALL timestop("determination of fermi energy")
                            ENDIF
                            IF ( noco%l_soc .AND. (.NOT. noco%l_noco) ) dimension%neigd = 2*dimension%neigd
                            IF( .NOT. ALLOCATED(results%w_iks) )&
Daniel Wortmann's avatar
Daniel Wortmann committed
463
                                 ALLOCATE ( results%w_iks(dimension%neigd,kpts%nkpt,dimension%jspd) )
464
465
466
467
                            IF ( (mpi%irank.EQ.0).AND.(.NOT.jij%l_J) ) THEN
                               CALL timestart("determination of fermi energy")

                               IF ( noco%l_soc .AND. (.NOT. noco%l_noco) ) THEN
468
                                  input%zelec = input%zelec*2
469
                                  CALL fermie(eig_id,mpi,kpts,obsolete,&
470
                                       input,noco,enpara%epara_min,jij,cell,results)
471
472
                                  results%seigscv = results%seigscv/2
                                  results%ts = results%ts/2
473
                                  input%zelec = input%zelec/2
474
475
                               ELSE
                                  CALL fermie(eig_id,mpi,kpts,obsolete,&
476
                                       input,noco,enpara%epara_min,jij,cell,results)
477
478
                               ENDIF
                               CALL timestop("determination of fermi energy")
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
                               !+Wannier
                            IF(wann%l_bs_comf)THEN
                               IF(pc.EQ.1) THEN
                                  OPEN(777,file='out_eig.1')
                                  OPEN(778,file='out_eig.2')
                                  OPEN(779,file='out_eig.1_diag')
                                  OPEN(780,file='out_eig.2_diag')
                               END IF

                               CALL bs_comfort(eig_id,DIMENSION,input,noco,kpts%nkpt,pc)

                               IF(pc.EQ.wann%nparampts)THEN
                                  CLOSE(777)
                                  CLOSE(778)
                                  CLOSE(779)
                                  CLOSE(780)
                               END IF
                            END IF
                            !-Wannier
498
499
500
501
502
503
504
505
506
507
508
                            ENDIF

                            IF (input%eonly) THEN
                               CALL close_eig(eig_id)

                               IF (.NOT. jij%l_J) THEN
                                  DEALLOCATE( results%w_iks )
#ifdef CPP_MPI
                                  CALL MPI_BARRIER(mpi%mpi_comm,ierr)
#endif
                                  IF (mpi%irank==0) THEN
Daniel Wortmann's avatar
Daniel Wortmann committed
509
                                     WRITE (*,fmt='(A)') 'eigenvalues written, program stops'
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
                                  ENDIF
                                  stop80= .TRUE.
                               ENDIF
                            ENDIF ! input%eonly

                         ENDIF ! ( input%eigvar(3) .and. .not.(noco%l_soc .and. noco%l_ss) )
                         !--- J<
                         IF(jij%l_J) THEN
                            IF (.NOT. input%eonly) THEN
                               DEALLOCATE ( results%w_iks )
                            ENDIF
                            IF (((i_J.EQ.j_J)).OR.(sym%invs.AND.(jij%qn.GT.tol))) GOTO 33
                         ENDIF
                      ENDDO !phn
33                    CONTINUE
                   ENDDO !j_J
                   i_J=i_J+jij%nmagtype(imt)
                ENDDO !imt
528
                IF ((mpi%irank.EQ.0).AND.(input%l_wann.OR.jij%l_J)) THEN
529
530
                   CALL timestop("Q-point for J_ij(total)")
                ENDIF
531
532
533
534
535
536
537
538
             END DO !pc
          ENDDO !qcount
          IF (stop80) THEN
             IF ((mpi%irank.EQ.0).AND.(isCurrentXMLElement("iteration"))) THEN
                CALL closeXMLElement('iteration')
             END IF
             EXIT ! it
          ENDIF
539

540
234       CONTINUE
541

542
543
544
545
546
          IF (mpi%irank.EQ.0) THEN
             IF(jij%l_J) THEN
                IF(.NOT.jij%l_disp)THEN
                   REWIND(113)
                   REWIND(114)
547
#ifdef CPP_NEVER
548
                   CALL jcoff2(atoms,sym,cell,jij,input)
549
550
#endif
                ENDIF
551
552
                CLOSE(113)
                CLOSE(114)
553
             ENDIF
554
          ENDIF
555

556
557
          IF (.NOT.jij%l_J) THEN
             !--- J>
558
#ifdef CPP_MPI
559
560
             CALL MPI_BCAST(results%ef,1,MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
             CALL MPI_BCAST(results%w_iks,SIZE(results%w_iks),MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
561
#endif
562
563
564
565
566
567
568
569
570
571
572
573
574
575
             !-t3e
             !
             !        ----> initialise force_old
             !
             CALL force_0(results)
             !
             !        ----> charge density
             !
             !+Wannier functions
             input%l_wann = .FALSE.
             INQUIRE (file='wann_inp',exist=input%l_wann)
             IF ((input%l_wann).AND.(.NOT.wann%l_bs_comf)) THEN
                CALL wannier(DIMENSION,mpi,input,sym,atoms,stars,vacuum,sphhar,oneD,&
                     wann,noco,cell,enpara,banddos,sliceplot,results,&
576
                     eig_idList,(sym%invs).AND.(.NOT.noco%l_noco),kpts%nkpt)
577
578
579
580
581
582
583
             END IF
             IF (wann%l_gwf) CALL juDFT_error("provide wann_inp if l_gwf=T", calledby = "fleur")
             !-Wannier

             CALL timestart("generation of new charge density (total)")
             IF (mpi%irank==0) WRITE(*,"(a)",advance="no") "* New Charge "
             CALL cdngen(eig_id,mpi,input,banddos,sliceplot,vacuum,&
584
                         DIMENSION,kpts,atoms,sphhar,stars,sym,obsolete,&
585
                         enpara,cell,noco,jij,results,oneD,coreSpecInput,&
586
                         inDen%iter,outDen)
587
588
589

             IF ( noco%l_soc .AND. (.NOT. noco%l_noco) ) dimension%neigd=dimension%neigd/2
             !+t3e
590
#ifdef CPP_MPI
591
592
593
594
595
596
597
598
599
600
601
602
603
             CALL MPI_BCAST(enpara%evac0,SIZE(enpara%evac0),MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
             CALL MPI_BCAST(enpara%el0,SIZE(enpara%el0),MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
             CALL MPI_BCAST(enpara%ello0,SIZE(enpara%ello0),MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)

             IF (noco%l_noco) THEN
                DO n= 1,atoms%ntype
                   IF (noco%l_relax(n)) THEN
                      CALL MPI_BCAST(noco%alph(n),1,MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
                      CALL MPI_BCAST(noco%beta(n),1,MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
                   ENDIF
                ENDDO
                IF (noco%l_constr) THEN
                   CALL MPI_BCAST(noco%b_con,SIZE(noco%b_con),MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
604
                ENDIF
605
             ENDIF
606
#endif
607
608
609
             CALL timestop("generation of new charge density (total)")
             IF (mpi%irank.EQ.0) THEN
                !-t3e
610

611
612
613
                IF (banddos%ndir.GT.0) THEN
                   CALL juDFT_end("NDIR",mpi%irank)
                END IF
614

615
                !----> output potential and potential difference
616
617
618
619
620
                IF (obsolete%disp) THEN
                   reap = .FALSE.
                   input%total = .FALSE.
                   CALL timestart("generation of potential (total)")
                   CALL vgen(hybrid,reap,input,xcpot,DIMENSION, atoms,sphhar,stars,vacuum,sym,&
621
                        obsolete,cell,oneD,sliceplot,mpi, results,noco,outDen,v,vx)
622
                   CALL timestop("generation of potential (total)")
Daniel Wortmann's avatar
Daniel Wortmann committed
623

624
625
                   CALL potdis(stars,vacuum,atoms,sphhar, input,cell,sym)
                END IF
626

627
                !----> total energy
628
629
630
                CALL timestart('determination of total energy')
                CALL totale(atoms,sphhar,stars,vacuum,DIMENSION,&
                     sym,input,noco,cell,oneD,xcpot,hybrid,it,results)
631

632
                CALL timestop('determination of total energy')
633

634
635
636
637
638
639
640
641
                ! in case of parallel processing, the total energy calculation is done
                ! only for irank.eq.0, since no parallelization is required here. once
                ! a force calculation is applied, however, the irank.eq.0 process is
                ! led into a MPI_FINALIZE after convergence, while the other processes
                ! are not, resulting in fleur not terminating despite having finished
                ! the calculation. the next 7 lines correct that issue.
                ! (other files subject to this correction: geo.F, force_w.F)
                ! Schlipf/Klueppelberg Jun 2012
642
#ifdef CPP_MPI
643
644
645
646
647
648
649
             ELSEIF (input%l_f) THEN ! forces, but mpi%irank.ne.0
                !This does not work, you can not call MPI_BCAST within a
                !else part of irank==0 as PE=0 will not call this!
                !CALL MPI_BCAST(lconv,1,MPI_LOGICAL,0,mpi_comm,ierr)
                !IF (lconv) THEN
                !  CALL MPI_FINALIZE(ierr)
                !END IF
650
651
#endif

Daniel Wortmann's avatar
Daniel Wortmann committed
652
             ENDIF ! mpi%irank.EQ.0
653
654
             !Close file if not a hybrid calculation
             IF ( hybrid%l_hybrid ) CALL close_eig(eig_id)
655

Daniel Wortmann's avatar
Daniel Wortmann committed
656
          ENDIF !(if not jij%l_J)
657
658
659
660
661
662
663
664
665
666
667
668
       ELSE
          input%alpha = input%alpha - 10.
       END IF !(if input%alpha <10.)
       IF (.NOT.jij%l_J) THEN

          IF (mpi%irank.EQ.0) THEN
             !-t3e
             !
             !          ----> mix input and output densities
             !
             CALL timestart("mixing")
             IF (mpi%irank==0) WRITE(*,"(a)",advance="no") "* Mixing distance: "
669
670
             CALL mix(stars,atoms,sphhar,vacuum,input,sym,cell,noco,oneD,hybrid,inDen,outDen,results,mixDen)
             inDen = mixDen
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
             !
             CALL timestop("mixing")
             WRITE (6,FMT=8130) it
             WRITE (16,FMT=8130) it
             IF (mpi%irank==0) THEN
                WRITE(*,"(f11.7)",advance='no') results%last_distance
                WRITE(*,*)
             ENDIF
8130         FORMAT (/,5x,'******* it=',i3,'  is completed********',/,/)
             WRITE(*,*) "Iteration:",it," Distance:",results%last_distance
             CALL timestop("Iteration")
             !+t3e
          ENDIF ! mpi%irank.EQ.0

          !--- J<
       ELSE
       ENDIF !(if not jij%l_J)
       !--- J>
689
690

#ifdef CPP_MPI
691
692
       CALL MPI_BCAST(results%last_distance,1,MPI_DOUBLE_PRECISION,0,mpi%mpi_comm,ierr)
       CALL MPI_BARRIER(mpi%mpi_comm,ierr)
693
#endif
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
       !-t3e

       !+fo
       CALL priv_geo_end(mpi)

       !-fo
       IF ( hybrid%l_calhf ) ithf = ithf + 1
       IF ( hybrid%l_subvxc ) THEN
          l_cont = ( ithf < input%itmax )
          results%te_hfex%core    = 0
          results%te_hfex%valence = 0
       ELSE
          l_cont = ( it < input%itmax )
       END IF
       CALL writeTimesXML()
       CALL check_time_for_next_iteration(it,l_cont)

       l_cont = l_cont.AND.((input%mindistance<=results%last_distance).OR.input%l_f)

       IF ((mpi%irank.EQ.0).AND.(isCurrentXMLElement("iteration"))) THEN
          CALL closeXMLElement('iteration')
       END IF

    END DO ! DO WHILE (l_cont)

    IF (mpi%irank.EQ.0) CALL closeXMLElement('scfLoop')
    CALL juDFT_end("all done",mpi%irank)

  CONTAINS
    SUBROUTINE priv_geo_end(mpi)
      TYPE(t_mpi),INTENT(IN)::mpi
      LOGICAL :: l_exist
      !Check if a new input was generated
      INQUIRE (file='inp_new',exist=l_exist)
      IF (l_exist) THEN
         CALL juDFT_end(" GEO new inp created ! ",mpi%irank)
      END IF
      !check for inp.xml
      INQUIRE (file='inp_new.xml',exist=l_exist)
      IF (.NOT.l_exist) RETURN
      IF (mpi%irank==0) THEN
         CALL system('mv inp.xml inp_old.xml')
         CALL system('mv inp_new.xml inp.xml')
         INQUIRE (file='qfix',exist=l_exist)
         IF (l_exist) THEN
            OPEN(2,file='qfix')
            WRITE(2,*)"F"
            CLOSE(2)
            PRINT *,"qfix set to F"
         ENDIF
         INQUIRE(file='broyd',exist=l_exist)
         IF (l_exist) THEN
            CALL system('rm broyd')
            PRINT *,"broyd file deleted"
         ENDIF
         INQUIRE(file='broyd.7',exist=l_exist)
750
         IF (l_exist) THEN
751
752
753
754
755
756
757
758
759
            CALL system('rm broyd.7')
            PRINT *,"broyd.7 file deleted"
         ENDIF
      ENDIF
      CALL juDFT_end(" GEO new inp.xml created ! ",mpi%irank)
    END SUBROUTINE priv_geo_end

  END SUBROUTINE fleur_execute
END MODULE m_fleur