TABLE OF CONTENTS


ABINIT/m_scf_history [ Modules ]

[ Top ] [ Modules ]

NAME

  m_scf_history

FUNCTION

  This module provides the definition of the scf_history_type used to store
  various arrays obtained from previous SCF cycles (density, positions, wavefunctions ...),
  as needed by the specific SCF algorithm.

COPYRIGHT

 Copyright (C) 2011-2018 ABINIT group (MT)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

OUTPUT

PARENTS

CHILDREN

SOURCE

27 #if defined HAVE_CONFIG_H
28 #include "config.h"
29 #endif
30 
31 #include "abi_common.h"
32 
33 MODULE m_scf_history
34 
35  use defs_basis
36  use defs_abitypes
37  use m_profiling_abi
38  use m_errors
39 
40  use m_pawcprj,  only : pawcprj_type, pawcprj_free
41  use m_pawrhoij, only : pawrhoij_type, pawrhoij_nullify, pawrhoij_free
42 
43  implicit none
44 
45  private
46 
47 ! public procedures.
48  public :: scf_history_init
49  public :: scf_history_free
50  public :: scf_history_nullify

m_scf_history/scf_history_free [ Functions ]

[ Top ] [ m_scf_history ] [ Functions ]

NAME

  scf_history_free

FUNCTION

  Clean and destroy a scf_history datastructure

SIDE EFFECTS

  scf_history(:)=<type(scf_history_type)>=scf_history datastructure

PARENTS

      gstateimg,scfcv

CHILDREN

SOURCE

344 subroutine scf_history_free(scf_history)
345 
346 
347 !This section has been created automatically by the script Abilint (TD).
348 !Do not modify the following lines by hand.
349 #undef ABI_FUNC
350 #define ABI_FUNC 'scf_history_free'
351 !End of the abilint section
352 
353  implicit none
354 
355 !Arguments ------------------------------------
356 !arrays
357  type(scf_history_type),intent(inout) :: scf_history
358 
359 !Local variables-------------------------------
360 !scalars
361  integer :: jj
362 
363 !************************************************************************
364 
365  !@scf_history_type
366 
367  if (allocated(scf_history%pawrhoij_last)) then
368    call pawrhoij_free(scf_history%pawrhoij_last)
369    ABI_DATATYPE_DEALLOCATE(scf_history%pawrhoij_last)
370  end if
371  if (allocated(scf_history%pawrhoij)) then
372    do jj=1,size(scf_history%pawrhoij,2)
373      call pawrhoij_free(scf_history%pawrhoij(:,jj))
374    end do
375    ABI_DATATYPE_DEALLOCATE(scf_history%pawrhoij)
376  end if
377  if (allocated(scf_history%cprj)) then
378    do jj=1,size(scf_history%cprj,3)
379      call pawcprj_free(scf_history%cprj(:,:,jj))
380    end do
381    ABI_DATATYPE_DEALLOCATE(scf_history%cprj)
382  end if
383 
384  if (allocated(scf_history%hindex))       then
385    ABI_DEALLOCATE(scf_history%hindex)
386  end if
387  if (allocated(scf_history%deltarhor))    then
388    ABI_DEALLOCATE(scf_history%deltarhor)
389  end if
390  if (allocated(scf_history%xreddiff))     then
391    ABI_DEALLOCATE(scf_history%xreddiff)
392  end if
393  if (allocated(scf_history%atmrho_last))  then
394    ABI_DEALLOCATE(scf_history%atmrho_last)
395  end if
396  if (allocated(scf_history%xred_last))    then
397    ABI_DEALLOCATE(scf_history%xred_last)
398  end if
399  if (allocated(scf_history%rhor_last))    then
400    ABI_DEALLOCATE(scf_history%rhor_last)
401  end if
402  if (allocated(scf_history%taur_last))    then
403    ABI_DEALLOCATE(scf_history%taur_last)
404  end if
405  if (allocated(scf_history%cg))           then
406    ABI_DEALLOCATE(scf_history%cg)
407  end if
408  if (allocated(scf_history%dotprod_sumdiag_cgcprj_ij))           then
409    ABI_DEALLOCATE(scf_history%dotprod_sumdiag_cgcprj_ij)
410  end if
411 
412  scf_history%history_size=-1
413  scf_history%usecg=0
414  scf_history%icall=0
415  scf_history%mcprj=0
416  scf_history%mcg=0
417 
418 end subroutine scf_history_free

m_scf_history/scf_history_init [ Functions ]

[ Top ] [ m_scf_history ] [ Functions ]

NAME

  scf_history_init

FUNCTION

  Init all scalars and pointers in a scf_history datastructure
  according to scf_history%history_size value which has to be
  defined before calling this routine

INPUTS

  dtset <type(dataset_type)>=all input variables in this dataset
  mpi_enreg=MPI-parallelisation information
  usecg= if ==0 => no handling of wfs, if==1 => handling of density/potential AND wfs, if==2 => ONLY handling of wfs

SIDE EFFECTS

  scf_history=<type(scf_history_type)>=scf_history datastructure
    hindex is always allocated
    The density/potential arrays that are possibly allocated are : atmrho_last, deltarhor, 
      pawrhoij, pawrhoij_last, rhor_last, taur_last, xreddiff, xred_last.
    The wfs arrays that are possibly allocated are : cg and cprj

PARENTS

      gstate,scfcv

CHILDREN

SOURCE

216 subroutine scf_history_init(dtset,mpi_enreg,usecg,scf_history)
217 
218 
219 !This section has been created automatically by the script Abilint (TD).
220 !Do not modify the following lines by hand.
221 #undef ABI_FUNC
222 #define ABI_FUNC 'scf_history_init'
223 !End of the abilint section
224 
225  implicit none
226 
227 !Arguments ------------------------------------
228 !scalars
229  integer, intent(in) :: usecg
230  type(dataset_type),intent(in) :: dtset
231  type(MPI_type),intent(in) :: mpi_enreg
232 !arrays
233  type(scf_history_type),intent(inout) :: scf_history
234 
235 !Local variables-------------------------------
236 !scalars
237  integer :: jj,mband_cprj,my_natom,my_nspinor,nfft
238 !arrays
239 
240 !************************************************************************
241 
242  !@scf_history_type
243 
244  if (scf_history%history_size<0) then
245    call scf_history_nullify(scf_history)
246  else
247 
248    scf_history%usecg=usecg
249    scf_history%wfmixalg=dtset%fockoptmix/100
250 
251    nfft=dtset%nfft
252    if (dtset%usepaw==1.and.(dtset%pawecutdg>=1.0000001_dp*dtset%ecut)) nfft=dtset%nfftdg
253    my_natom=mpi_enreg%my_natom
254 
255    if (scf_history%history_size>=0 .and. usecg<2) then
256      ABI_ALLOCATE(scf_history%rhor_last,(nfft,dtset%nspden))
257      ABI_ALLOCATE(scf_history%taur_last,(nfft,dtset%nspden*dtset%usekden))
258      ABI_ALLOCATE(scf_history%xred_last,(3,dtset%natom))
259      ABI_DATATYPE_ALLOCATE(scf_history%pawrhoij_last,(my_natom*dtset%usepaw))
260      if (dtset%usepaw==1) then
261        call pawrhoij_nullify(scf_history%pawrhoij_last)
262      end if
263    end if
264 
265    if (scf_history%history_size>0) then
266 
267      scf_history%natom=dtset%natom
268      scf_history%nfft=nfft
269      scf_history%nspden=dtset%nspden
270      scf_history%beta=zero
271      scf_history%icall=0
272 
273      scf_history%mcg=0
274      scf_history%mcprj=0
275      if (usecg>0) then
276        my_nspinor=max(1,dtset%nspinor/mpi_enreg%nproc_spinor)
277        scf_history%mcg=dtset%mpw*my_nspinor*dtset%nbandhf*dtset%mkmem*dtset%nsppol ! This is for scf_history_wf
278        if(usecg==1)scf_history%mcg=dtset%mpw*my_nspinor*dtset%mband*dtset%mkmem*dtset%nsppol ! This is for scf_history (when extrapwf==1)
279        if (dtset%usepaw==1) then
280          mband_cprj=dtset%nbandhf
281          if(usecg==1)mband_cprj=dtset%mband
282          if (dtset%paral_kgb/=0) mband_cprj=mband_cprj/mpi_enreg%nproc_band
283          scf_history%mcprj=my_nspinor*mband_cprj*dtset%mkmem*dtset%nsppol
284        end if
285      end if
286 
287      if (usecg<2) then
288        ABI_ALLOCATE(scf_history%hindex,(scf_history%history_size))
289        scf_history%alpha=zero
290      else
291        ABI_ALLOCATE(scf_history%hindex,(2*scf_history%history_size+1))
292        scf_history%alpha=dtset%wfmix
293      endif
294      scf_history%hindex(:)=0
295 
296      if (usecg<2) then 
297        ABI_ALLOCATE(scf_history%deltarhor,(nfft,dtset%nspden,scf_history%history_size))
298        ABI_ALLOCATE(scf_history%xreddiff,(3,dtset%natom,scf_history%history_size))
299        ABI_ALLOCATE(scf_history%atmrho_last,(nfft))
300        if (dtset%usepaw==1) then
301          ABI_DATATYPE_ALLOCATE(scf_history%pawrhoij,(my_natom,scf_history%history_size))
302          do jj=1,scf_history%history_size
303            call pawrhoij_nullify(scf_history%pawrhoij(:,jj))
304          end do
305        endif
306      end if
307 
308      if (scf_history%usecg>0) then
309        ABI_ALLOCATE(scf_history%cg,(2,scf_history%mcg,scf_history%history_size))
310 !      Note that the allocation is made even when usepaw==0. Still, scf_history%mcprj=0 ...
311        ABI_DATATYPE_ALLOCATE(scf_history%cprj,(dtset%natom,scf_history%mcprj,scf_history%history_size))
312      end if
313 
314      if (scf_history%usecg==2)then
315 !      This relatively small matrix is always allocated when usecg==1, even if not used
316        ABI_ALLOCATE(scf_history%dotprod_sumdiag_cgcprj_ij,(2,scf_history%history_size,scf_history%history_size))
317      endif
318 
319    end if
320  end if
321 
322 end subroutine scf_history_init

m_scf_history/scf_history_nullify [ Functions ]

[ Top ] [ m_scf_history ] [ Functions ]

NAME

  scf_history_nullify

FUNCTION

  Nullify (set to null) an scf_history datastructure

SIDE EFFECTS

  scf_history(:)=<type(scf_history_type)>=scf_history datastructure

PARENTS

      gstateimg,m_scf_history

CHILDREN

SOURCE

440 subroutine scf_history_nullify(scf_history)
441 
442 
443 !This section has been created automatically by the script Abilint (TD).
444 !Do not modify the following lines by hand.
445 #undef ABI_FUNC
446 #define ABI_FUNC 'scf_history_nullify'
447 !End of the abilint section
448 
449  implicit none
450 
451 !Arguments ------------------------------------
452 !arrays
453  type(scf_history_type),intent(inout) :: scf_history
454 !Local variables-------------------------------
455 !scalars
456 
457 !************************************************************************
458 
459  !@scf_history_type
460  scf_history%history_size=-1
461  scf_history%icall=0
462  scf_history%mcprj=0
463  scf_history%mcg=0
464 
465 end subroutine scf_history_nullify

m_scf_history/scf_history_type [ Types ]

[ Top ] [ m_scf_history ] [ Types ]

NAME

 scf_history_type

FUNCTION

 This structured datatype contains various arrays obtained from
 previous SCF cycles (density, positions...)

SOURCE

 63  type, public :: scf_history_type
 64 
 65 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
 66 ! declared in another part of ABINIT, that might need to take into account your modification.
 67 
 68 ! Integer scalar
 69 
 70   integer :: history_size
 71    ! Number of previous SCF cycles stored in history
 72    ! If history_size<0, scf_history is not used
 73    ! If history_size=0, scf_history only contains
 74    !    current values of data (rhor, taur, pawrhoih, xred)
 75    ! If history_size>0, scf_history contains
 76    !    current values of data and also
 77    !    history_size previous values of these data
 78 
 79   integer :: icall
 80    ! Number of call for the routine extraprho or wf_mixing
 81 
 82   integer :: mcg
 83    ! Size of cg array
 84 
 85   integer :: mcprj
 86    ! Size of cprj datsatructure array
 87 
 88   integer :: natom
 89    ! Number of atoms in cell
 90 
 91   integer :: nfft
 92    ! Size of FFT grid (for density)
 93 
 94   integer :: nspden
 95    ! Number of independant spin components for density
 96 
 97   integer :: usecg
 98    ! usecg=0 if the extrapolation/mixing of the density/potential is active but not the one of the wavefunction
 99    ! usecg=1 if the extrapolation/mixing of the density/potential and wavefunctions is active
100    ! usecg=2 if the extrapolation/mixing of the wavefunctions is active but not the one of the density/potential
101 
102   integer :: wfmixalg
103    ! algorithm used to mix the wavefunctions (in case usecg=2)
104 
105   real(dp) :: alpha
106    ! alpha mixing coefficient for the prediction of density and wavefunctions
107    ! In the case of wavefunction simple mixing, contain wfmix factor
108 
109   real(dp) :: beta
110    ! beta mixing coefficient for the prediction of density and wavefunctions
111 
112 ! Integer arrays
113 
114   integer,allocatable :: hindex(:)
115    ! Indexes of SCF cycles in the history
116    !
117    ! For the density-based schemes (with or without wavefunctions) : 
118    ! hindex(history_size)
119    ! hindex(1) is the newest SCF cycle
120    ! hindex(history_size) is the oldest SCF cycle
121    !
122    ! For wavefunction-based schemes (outer loop of a double loop SCF):
123    ! hindex(2*history_size+1)
124    ! The odd indices refer to the out wavefunction,
125    ! the even indices refer to the in wavefunction (not all such wavefunctions being stored, though). 
126    ! hindex(1:2) is the newest SCF cycle, hindex(3:4) is the SCF cycle before the newest one ... In case of an 
127    ! algorithm based on a biorthogonal ensemble of wavefunctions, the reference is stored in hindex(2*history_size+1)
128    ! When the index points to a location beyond history_size, the corresponding wavefunction set must be reconstructed
129    ! from the existing wavefunctions sets (to be implemented)
130 
131 ! Real (real(dp)) arrays
132 
133    real(dp),allocatable :: cg(:,:,:)
134     ! cg(2,mcg,history_size)
135     ! wavefunction coefficients needed for each SCF cycle of history
136     ! Might also contain the wf residuals
137 
138    real(dp),allocatable :: deltarhor(:,:,:)
139     ! deltarhor(nfft,nspden,history_size)
140     ! Difference between electronic density (in real space)
141     ! and sum of atomic densities at the end of each SCF cycle of history
142 
143    real(dp),allocatable :: atmrho_last(:)
144     ! atmrho_last(nfft)
145     ! Sum of atomic densities at the end of the LAST SCF cycle
146 
147    real(dp),allocatable :: rhor_last(:,:)
148     ! rhor_last(nfft,nspden)
149     ! Last computed electronic density (in real space)
150 
151    real(dp),allocatable :: taur_last(:,:)
152     ! taur_last(nfft,nspden*usekden)
153     ! Last computed kinetic energy density (in real space)
154 
155    real(dp),allocatable :: xreddiff(:,:,:)
156     ! xreddiff(3,natom,history_size)
157     ! Difference of reduced coordinates of atoms between a
158     ! SCF cycle and the previous
159 
160    real(dp),allocatable :: xred_last(:,:)
161     ! xred_last(3,natom)
162     ! Last computed atomic positions (reduced coordinates)
163 
164    real(dp),allocatable :: dotprod_sumdiag_cgcprj_ij(:,:,:)
165     ! dotprod_sumdiag_cgcprj_mn(2,history_size,history_size)
166     ! Container for the scalar products between aligned sets of wavefunctions or their residuals
167     ! S_ij=Sum_nk <wf_nk(set i)|wf_nk(set j)> possibly with some weighting factor that might depend on nk.
168 
169 ! Structured datatypes arrays
170 
171   type(pawrhoij_type), allocatable :: pawrhoij(:,:)
172     ! pawrhoij(natom,history_size)
173     ! PAW only: occupancies matrix at the end of each SCF cycle of history
174 
175   type(pawrhoij_type), allocatable :: pawrhoij_last(:)
176     ! pawrhoij_last(natom)
177     ! PAW only: last computed occupancies matrix
178 
179   type(pawcprj_type),allocatable :: cprj(:,:,:)
180     !cprj(natom,nspinor*mband*mkmem*nsppol,history_size)
181 
182  end type scf_history_type