TABLE OF CONTENTS


ABINIT/m_m1geo [ Modules ]

[ Top ] [ Modules ]

NAME

 m_m1geo

FUNCTION

 This module contains definition the m1geo type (="move unique geometry")
 and its related init and destroy routines

COPYRIGHT

 Copyright (C) 2018-2024 ABINIT group (XG)
 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_m1geo
24 
25  use defs_basis
26  use m_profiling_abi
27  use m_abimover
28  use m_abihist
29  use m_dtset
30  use m_dtfil
31 
32  implicit none
33 
34  private
35 
36  public :: m1geo_init
37  public :: m1geo_destroy

m_gstateimg/m1geo_destroy [ Functions ]

[ Top ] [ m_gstateimg ] [ Functions ]

NAME

 m1geo_destroy

FUNCTION

 Destroy the m1geo structure

INPUTS

OUTPUT

SOURCE

217  subroutine m1geo_destroy(m1geo_param)
218 
219 !Arguments ------------------------------------
220 !scalars
221  type(m1geo_type),intent(inout) :: m1geo_param
222 
223 !************************************************************************
224 
225  ABI_SFREE(m1geo_param%mixesimgf)
226  ABI_SFREE(m1geo_param%ab_xfh_1geo%xfhist)
227 
228  call abimover_destroy(m1geo_param%ab_mover)
229  call delocint_fin(m1geo_param%deloc)
230  call abihist_free(m1geo_param%hist_1geo)
231  call mttk_fin(m1geo_param%mttk_vars)
232 
233  end subroutine m1geo_destroy
234 
235 !----------------------------------------------------------------------
236 
237 end module m_m1geo

m_gstateimg/m1geo_init [ Functions ]

[ Top ] [ m_gstateimg ] [ Functions ]

NAME

 m1geo_init

FUNCTION

 Initializes the m1geo structure

INPUTS

OUTPUT

SOURCE

102  subroutine m1geo_init(dtfil,dtset,m1geo_param)
103 
104 !Arguments ------------------------------------
105 !scalars
106  type(datafiles_type),target,intent(in) :: dtfil
107  type(dataset_type),target,intent(in) :: dtset
108  type(m1geo_type),intent(inout) :: m1geo_param
109 
110 !Local variables
111 !integer
112  integer :: ionmov,natom,ncycle,nhisttot,nimage,ntimimage,ntypat
113  real(dp), allocatable :: amu_curr(:)
114 
115 !************************************************************************
116 
117 
118  ionmov=dtset%ionmov
119  natom=dtset%natom
120  nimage=dtset%nimage
121  ntimimage=dtset%ntimimage
122  ntypat=dtset%ntypat
123 
124 !Straight initialization from dtset
125  m1geo_param%dt_chkdilatmx  =dtset%chkdilatmx
126  m1geo_param%dilatmx        =dtset%dilatmx
127  m1geo_param%hmctt          =dtset%hmctt
128  m1geo_param%nctime         =dtset%nctime
129  m1geo_param%nimage         =dtset%nimage
130  m1geo_param%npsp           =dtset%npsp
131  m1geo_param%usewvl         =dtset%usewvl
132 
133  ABI_MALLOC(m1geo_param%mixesimgf,(nimage))
134  m1geo_param%mixesimgf(:)=dtset%mixesimgf(1:nimage)
135 
136 
137 !From dtset, with change or adjustment of meaning
138  m1geo_param%ntime          =dtset%ntimimage  ! Beware ntime vs ntimimage
139  m1geo_param%rprimd_orig    =dtset%rprimd_orig(:,:,1)
140 
141  ABI_MALLOC(amu_curr,(ntypat))
142  amu_curr(:)=dtset%amu_orig(:,1)
143 
144 !From dtfil
145  m1geo_param%filnam_ds4     =dtfil%filnam_ds(4)
146 
147 !###########################################################
148 !Init sub-datastructures
149 
150 !1) ab_mover
151  call abimover_ini(m1geo_param%ab_mover,amu_curr,dtfil,dtset,m1geo_param%specs)
152  ABI_FREE(amu_curr)
153 
154 !2) deloc
155  if(ionmov==10 .or. ionmov==11)then
156    call delocint_ini(m1geo_param%deloc)
157  end if
158 
159 !3) hist
160 !In a simple approach, hist_1geo is initialized here.
161 !Likely one should also make an interface with the hist(nimage) from the calling routine, and average it over the images...
162 
163  ncycle=m1geo_param%specs%ncycle
164  if(ionmov==25.and.dtset%hmctt>=0) then
165    ncycle=dtset%hmctt
166    if(dtset%hmcsst>0.and.dtset%optcell/=0) then
167      ncycle=ncycle+dtset%hmcsst
168    endif
169  endif
170  nhisttot=ncycle*ntimimage  ! WARNING, here ntime is used instead of ntimimage
171  if (dtset%nctime>0) nhisttot=nhisttot+1
172 
173 !We just store the needed history step not all of them.
174  if(m1geo_param%specs%nhist/=-1) nhisttot = m1geo_param%specs%nhist ! We don't need to store all the history
175 
176  call abihist_init(m1geo_param%hist_1geo,natom,ntimimage,m1geo_param%specs%isVused,m1geo_param%specs%isARused)
177 
178 !4) ab_xfh
179  m1geo_param%ab_xfh_1geo%mxfh=ntimimage
180  m1geo_param%ab_xfh_1geo%nxfh=0
181  m1geo_param%ab_xfh_1geo%nxfhr=0
182  ABI_MALLOC(m1geo_param%ab_xfh_1geo%xfhist,(3,natom+4,2,ntimimage))
183 
184 !5) mttk
185  if (ionmov==13)then
186    call mttk_ini(m1geo_param%mttk_vars,m1geo_param%ab_mover%nnos)
187  end if
188 
189 !###########################################################
190 !Init other values
191 
192  m1geo_param%ncycle=m1geo_param%specs%ncycle
193  m1geo_param%icycle=0
194  m1geo_param%iexit=0
195  m1geo_param%itime=0
196  m1geo_param%nerr_dilatmx=0
197  m1geo_param%skipcycle=.FALSE.
198 
199  end subroutine m1geo_init

m_m1geo/m1geo [ Types ]

[ Top ] [ m_m1geo ] [ Types ]

NAME

 m1geo type

FUNCTION

 This datatype has the purpose to store all the data taken
 usually from dtset (but not only) needed to move a single geometry
 opbtained by the weighting of different images,
 using the algorithms defined by ionmov.

SOURCE

55 type, public :: m1geo_type
56 ! Scalars
57   integer  :: dt_chkdilatmx  ! chkdilatmx from dtset
58   integer  :: hmctt          ! hmctt from dtset
59   integer  :: icycle         ! cycle index
60   integer  :: iexit          ! if nonzero, will have to exit
61   integer  :: itime          ! time index for the move1geo algorithms
62   integer  :: nctime         !
63   integer  :: ncycle         ! foreseen number of cycles
64   integer  :: nerr_dilatmx   ! number of dilatmx trespass errors, if chkdilatmx
65   integer  :: nimage         ! nimage from dtset
66   integer  :: npsp           ! npsp from dtset
67   integer  :: ntime          ! number of time steps for the move1geo algorithms
68   integer  :: usewvl         ! usewvl from dtset
69   real(dp) :: dilatmx        ! dilatmx from dtset
70   logical  :: skipcycle      ! .TRUE. when the remaining of the cycle has to be skipped. .FALSE. otherwise.
71 ! Arrays
72   real(dp) :: rprimd_orig(3,3)   ! rprimd from dtset
73   real(dp),allocatable :: mixesimgf(:)   ! mixesimgf from dtset
74 ! Character string
75   character(len=fnlen) :: filnam_ds4
76 ! Datatypes
77   type(abimover) :: ab_mover
78   type(abimover_specs) :: specs
79   type(ab_xfh_type) :: ab_xfh_1geo
80   type(delocint) :: deloc
81   type(abihist) :: hist_1geo
82   type(mttk_type) :: mttk_vars
83 
84  end type m1geo_type
85 
86 contains  !=============================================================