TABLE OF CONTENTS


ABINIT/m_precpred_1geo [ Modules ]

[ Top ] [ Modules ]

NAME

  m_precpred_1geo

FUNCTION

 Single geometry: apply force and stress preconditioner followed by geometry predictor.
 Choose among the whole set of geometry predictors defined by iommov.

COPYRIGHT

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

SOURCE

17 #if defined HAVE_CONFIG_H
18 #include "config.h"
19 #endif
20 
21 #include "abi_common.h"
22 
23 module m_precpred_1geo
24 
25  use defs_basis
26  use m_errors
27  use m_abicore
28  use m_abimover
29  use m_abihist
30  use m_xmpi
31  use m_nctk
32 #ifdef HAVE_NETCDF
33  use netcdf
34 #endif
35 #if defined HAVE_LOTF
36  use lotfpath
37  use m_pred_lotf
38 #endif
39 
40  use m_fstrings,           only : strcat, sjoin, itoa
41  use m_geometry,           only : chkdilatmx
42  use m_crystal,            only : crystal_init, crystal_t
43  use m_pred_bfgs,          only : pred_bfgs, pred_lbfgs
44  use m_pred_delocint,      only : pred_delocint
45  use m_pred_fire,          only : pred_fire
46  use m_pred_isokinetic,    only : pred_isokinetic
47  use m_pred_diisrelax,     only : pred_diisrelax
48  use m_pred_nose,          only : pred_nose
49  use m_pred_srkhna14,      only : pred_srkna14
50  use m_pred_isothermal,    only : pred_isothermal
51  use m_pred_verlet,        only : pred_verlet
52  use m_pred_velverlet,     only : pred_velverlet
53  use m_pred_moldyn,        only : pred_moldyn
54  use m_pred_langevin,      only : pred_langevin
55  use m_pred_steepdesc,     only : pred_steepdesc
56  use m_pred_simple,        only : pred_simple, prec_simple
57  use m_pred_hmc,           only : pred_hmc
58  use m_ipi,                only : ipi_pred
59 !use m_generate_training_set, only : generate_training_set
60 
61  implicit none
62 
63  private

ABINIT/precpred_1geo [ Functions ]

[ Top ] [ Functions ]

NAME

 mover

FUNCTION

 Single geometry: apply force and stress preconditioner followed by geometry predictor.
 Choose among the whole set of geometry predictors defined by ionmov.

INPUTS

  (TO BE DESCRIBED)

OUTPUT

  (TO BE DESCRIBED)

SIDE EFFECTS

  hist is the main quantity updated, contains xred, rprimd and acell.
  (TO BE DESCRIBED)

NOTES

SOURCE

 94 subroutine precpred_1geo(ab_mover,ab_xfh,amu_curr,deloc,dt_chkdilatmx,comm_cell,dilatmx,filnam_ds4,hist,hmctt,&
 95 & icycle,iexit,itime,mttk_vars,nctime,ncycle,nerr_dilatmx,npsp,ntime,rprimd_orig,skipcycle,usewvl)
 96 
 97 !Arguments ------------------------------------
 98 !scalars
 99  integer, intent(in) :: comm_cell,dt_chkdilatmx,hmctt,icycle,itime,nctime,npsp,ntime
100  integer, intent(inout) :: iexit,ncycle,nerr_dilatmx
101  integer, intent(in) :: usewvl
102  real(dp), intent(in) :: dilatmx
103  logical, intent(inout) :: skipcycle
104  character(len=fnlen), intent(in) :: filnam_ds4
105  type(ab_xfh_type),intent(inout) :: ab_xfh
106  type(abihist), intent(inout) :: hist
107  type(abimover), intent(in) :: ab_mover
108  type(delocint), intent(inout) :: deloc
109  type(mttk_type), intent(inout) :: mttk_vars
110 !arrays
111  real(dp), intent(in) :: amu_curr(ab_mover%ntypat)
112  real(dp), intent(in) :: rprimd_orig(3,3)
113 
114 !Local variables-------------------------------
115 !scalars
116  integer,parameter :: master=0
117  integer :: ii,me,nloop
118 !integer :: iatom
119  logical,parameter :: DEBUG=.FALSE.
120 !character(len=500) :: message
121  character(len=500) :: dilatmx_errmsg
122  character(len=fnlen) :: filename
123  type(abiforstr) :: preconforstr ! Preconditioned forces and stress
124  type(crystal_t) :: crystal
125 !arrays
126  real(dp) :: acell(3),rprimd(3,3)
127  real(dp), allocatable :: xred(:,:)
128 
129 ! ***************************************************************
130 
131 !DEBUG
132 !write(std_out,'(a,i4)')' m_precpred_1geo, enter : ab_mover%ionmov=',ab_mover%ionmov
133 !ABI_MALLOC(xred,(3,ab_mover%natom))
134 !call hist2var(acell,hist,ab_mover%natom,rprimd,xred,DEBUG)
135 !do iatom=1,ab_mover%natom
136 !  write(std_out,'(i4,3es14.6)') iatom,xred(1:3,iatom)
137 !enddo
138 !ABI_FREE(xred)
139 !ENDDEBUG
140 
141  me=xmpi_comm_rank(comm_cell)
142 
143 !Precondition forces, stress and energy
144  call abiforstr_ini(preconforstr,ab_mover%natom)
145 
146  if (ab_mover%goprecon>0 .or. iexit==1)then
147    call prec_simple(ab_mover,preconforstr,hist,icycle,itime,iexit)
148  end if
149 
150 !Call to each predictor
151 !MT->GAF: dirty trick to predict vel(t)
152 !do a double loop: 1- compute vel, 2- exit
153  nloop=1
154 
155  if (nctime>0.and.iexit==1) then
156    iexit=0;nloop=2
157  end if
158 
159  do ii=1,nloop
160    if (ii==2) iexit=1
161 
162    select case (ab_mover%ionmov)
163    case (1)
164      call pred_moldyn(ab_mover,hist,icycle,itime,ncycle,ntime,DEBUG,iexit)
165    case (2,3)
166      call pred_bfgs(ab_mover,ab_xfh,preconforstr,hist,ab_mover%ionmov,itime,DEBUG,iexit)
167    case (4,5)
168      call pred_simple(ab_mover,hist,iexit)
169    case (6,7)
170      call pred_verlet(ab_mover,hist,ab_mover%ionmov,itime,ntime,DEBUG,iexit)
171    case (8)
172      call pred_nose(ab_mover,hist,itime,ntime,DEBUG,iexit)
173    case (9)
174      call pred_langevin(ab_mover,hist,icycle,itime,ncycle,ntime,DEBUG,iexit,skipcycle)
175    case (10,11)
176      call pred_delocint(ab_mover,ab_xfh,deloc,preconforstr,hist,ab_mover%ionmov,itime,DEBUG,iexit)
177    case (12)
178      call pred_isokinetic(ab_mover,hist,itime,ntime,DEBUG,iexit)
179    case (13)
180      call pred_isothermal(ab_mover,hist,itime,mttk_vars,ntime,DEBUG,iexit)
181    case (14)
182      call pred_srkna14(ab_mover,hist,icycle,DEBUG,iexit,skipcycle)
183    case (15)
184      call pred_fire(ab_mover, ab_xfh,preconforstr,hist,ab_mover%ionmov,itime,DEBUG,iexit)
185    case (20)
186      call pred_diisrelax(ab_mover,hist,itime,ntime,DEBUG,iexit)
187    case (21)
188      call pred_steepdesc(ab_mover,preconforstr,hist,itime,DEBUG,iexit)
189    case (22)
190      call pred_lbfgs(ab_mover,ab_xfh,preconforstr,hist,ab_mover%ionmov,itime,DEBUG,iexit)
191 #if defined HAVE_LOTF
192    case (23)
193      call pred_lotf(ab_mover,hist,itime,icycle,DEBUG,iexit)
194 #endif
195    case (24)
196      call pred_velverlet(ab_mover,hist,itime,ntime,DEBUG,iexit)
197    case (25)
198      call pred_hmc(ab_mover,hist,itime,icycle,ntime,hmctt,mttk_vars,DEBUG,iexit)
199    case (27)
200      !In case of ionmov 27, all the atomic configurations have been computed at the
201      !beginning of the routine in generate_training_set, thus we just need to increase the indexes
202      !in the hist
203      hist%ihist = abihist_findIndex(hist,+1)
204    case (28)
205      call ipi_pred(ab_mover, hist, itime, ntime, DEBUG, iexit, comm_cell)
206    case default
207      ABI_ERROR(sjoin("Wrong value of ionmov:", itoa(ab_mover%ionmov)))
208    end select
209 
210  end do
211 
212  ABI_MALLOC(xred,(3,ab_mover%natom))
213  call hist2var(acell,hist,ab_mover%natom,rprimd,xred,DEBUG)
214 
215  ! check dilatmx here and correct if necessary
216  if (usewvl == 0) then
217    call chkdilatmx(dt_chkdilatmx,dilatmx,rprimd,rprimd_orig,dilatmx_errmsg)
218    if (LEN_TRIM(dilatmx_errmsg) /= 0) then
219      ABI_WARNING(dilatmx_errmsg)
220      nerr_dilatmx = nerr_dilatmx+1
221      if (nerr_dilatmx > 3) then
222        ! Write last structure before aborting, so that we can restart from it.
223        ! zion is not available, but it's not useful here.
224        if (me == master) then
225          ! Init crystal
226          hist%ihist = abihist_findIndex(hist,-1)
227          call hist2var(acell,hist,ab_mover%natom,rprimd,xred,DEBUG)
228          call crystal_init(amu_curr,crystal,0,ab_mover%natom,&
229            npsp,ab_mover%ntypat,ab_mover%nsym,rprimd,ab_mover%typat,xred,&
230            [(-one, ii=1,ab_mover%ntypat)],ab_mover%znucl,2,.False.,.False.,"dilatmx_structure",&
231            symrel=ab_mover%symrel,tnons=ab_mover%tnons,symafm=ab_mover%symafm)
232 
233 #ifdef HAVE_NETCDF
234          ! Write netcdf file
235          filename = strcat(filnam_ds4, "_DILATMX_STRUCT.nc")
236          NCF_CHECK(crystal%ncwrite_path(filename))
237 #endif
238          call crystal%free()
239        end if
240        call xmpi_barrier(comm_cell)
241        write (dilatmx_errmsg, '(a,i0,9a)') &
242         'Dilatmx has been exceeded too many times (', nerr_dilatmx, ')',ch10, &
243         'See the description of dilatmx and chkdilatmx input variables.',ch10, &
244         'Action: either first do a calculation with chkdilatmx=0, or ',ch10,&
245         'restart your calculation with a larger dilatmx, or larger lattice vectors.',ch10,&
246         'Warning: With chkdilatmx = 0 the final computation of lattice parameters might be inaccurate.'
247        ABI_ERROR_CLASS(dilatmx_errmsg, "DilatmxError")
248      end if
249    else
250      nerr_dilatmx=0
251    end if
252  end if
253 
254 !DEBUG
255 !write(std_out,'(a,i4)')' m_precpred_1geo, exit '
256 !do iatom=1,ab_mover%natom
257 !  write(std_out,'(i4,3es14.6)') iatom,xred(1:3,iatom)
258 !enddo
259 !ENDDEBUG
260 
261  call abiforstr_fin(preconforstr)
262  ABI_FREE(xred)
263 
264 end subroutine precpred_1geo