TABLE OF CONTENTS


ABINIT/bond_lotf [ Modules ]

[ Top ] [ Modules ]

NAME

 bond_lotf

FUNCTION

  Define BOND variables and the procedure to
  set them.

COPYRIGHT

 Copyright (C) 2005-2024 ABINIT group (MMancini)
 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 bond_lotf
24 
25  use defs_basis
26  use m_errors
27  use m_abicore
28 
29  implicit none
30 
31  public
32 
33  !--Fittedatoms variables
34  integer :: nfit !--dimension needed for fits
35  integer :: nfitmax  !--dimension of the fit arrays
36  integer,allocatable :: ifit(:)
37  logical,allocatable :: tafit(:)
38 
39  !--Fittedbonds variables
40  integer :: nbondex
41  integer :: ibn_tot,ibn_tot2,ibn_tots
42  integer,allocatable :: imat(:)
43  integer,dimension(:,:),allocatable :: ibnd_mat,ibmat,ibmat_large
44 
45  public ::             &
46    bond_atom_init,     &
47    bond_tafit_init,    &
48    bond_matrix_alloc,  &
49    bond_matrix_set,    &
50    bond_compute,       &
51    bond_fit_set,       &
52    bond_dealloc
53 
54 
55 contains

bond_lotf/bond_atom_init [ Functions ]

[ Top ] [ bond_lotf ] [ Functions ]

NAME

 bond_atom_init

FUNCTION

INPUTS

SOURCE

 95  subroutine bond_atom_init(nneigx,nneig,neighl)
 96 
 97   implicit none
 98 
 99   !Arguments ------------------------
100   integer,intent(in) :: nneigx
101   integer,intent(in) :: nneig(:)
102   integer,intent(in) :: neighl(:,:)
103 
104   ! nneigx : max number of neighbours
105   ! tau0(3,natom) : atomic positions
106   ! neighl(nneigx,natom) : list of neighbours
107   ! nneig(natom) : number of neighbours
108   ! niter  : iteration number (itime)
109 
110   !Local ---------------------------
111   integer :: i,j,iat,jat, ibn
112   integer,allocatable :: ibnd_dum(:,:)
113   character(len=500) :: msg
114 
115 ! *************************************************************************
116 
117   !--Initialize nbondex
118    nbondex = ((nfitmax/2 * nneigx)+1)/2
119 
120   !--Now find initial numbers of active/border bonds :
121   !  (0)        CLEARS THE BOND(atom) MATRIX :
122    ABI_MALLOC(ibnd_mat,(2,nbondex))
123    ABI_MALLOC(ibnd_dum,(2,nfitmax*6))
124    ibnd_mat = 0
125    ibnd_dum = 0
126 
127    ibn_tot  = 0    !-- bonds between fitted atoms
128    ibn_tot2 = 0    !--existing but non optimized bonds with border atoms
129    ibn_tots = 0    !--total number of bonds in the fit + border zone
130 
131    do i =1,nfit
132      iat = ifit(i)
133      do j = 1,nneig(iat)
134        jat = neighl(j,iat)
135        if(tafit(jat)) then
136          if(jat > iat) then !--jat is a fitted atom
137            ibn_tot = ibn_tot + 1
138            ibnd_mat(:,ibn_tot) = (/iat,jat/)
139          end if
140        else    !--jat is a border atom
141          ibn_tot2 = ibn_tot2 + 1
142          ibnd_dum(:,ibn_tot2) = (/iat,jat/)
143        end if
144      end do
145    end do
146 
147    if(ibn_tot2 > (6*nfitmax)) then
148      write(msg,'(3a,i8,2a)')&
149 &     'ERROR: BOND_ATOM_INIT',ch10,&
150 &     'IBN_TOT2 =  ',ibn_tot2,ch10,&
151 &     ' ibnd_dum out of bounds, ibn_tot2 too large '
152      ABI_ERROR(msg)
153 
154    end if
155 
156   !--Reorder to keep 'variational' bonds first :
157    ibn_tots = ibn_tot + ibn_tot2
158    do ibn=ibn_tot+1,ibn_tots
159      ibnd_mat(:,ibn) = ibnd_dum(:,ibn-ibn_tot)
160    end do
161 
162 
163    ABI_FREE(ibnd_dum)
164  end subroutine bond_atom_init

bond_lotf/bond_compute [ Functions ]

[ Top ] [ bond_lotf ] [ Functions ]

NAME

 bond_compute

FUNCTION

  Updates bond matrix, associates the bond to the atom neighlists

INPUTS

  nneig
  neighl

SOURCE

266  subroutine bond_compute(nneig,neighl)
267 
268   implicit none
269 
270   !Arguments ------------------------
271   integer,intent(in) :: nneig(:)
272   integer,intent(in) :: neighl(:,:)
273   !Local ---------------------------
274   integer :: ii,jj,iat,jat
275   integer, allocatable, dimension(:,:) :: ibnd_dum
276   character(len=500) :: msg
277 
278 ! *************************************************************************
279 
280    ABI_MALLOC(ibnd_dum,(2,nfitmax*6))
281    ibnd_dum(:,:) = 0
282    ibn_tot  = 0    ! bonds between the fitted atoms
283    ibn_tot2 = 0    ! existing but non optimized bonds with border atoms
284    ibn_tots = 0    ! total number of bonds
285 
286 
287    if(nfit <= 0) then
288      write(msg,'(a,i8)')' UPDLIS WARNING : nfit <= 0 = ', nfit
289      ABI_WARNING(msg)
290    end if
291 
292    do ii = 1,nfit
293      iat = ifit(ii)
294      do jj = 1,nneig(iat)
295        jat = neighl(jj,iat)
296        if(tafit(jat))then
297          if(jat > iat) then  ! jat is a fitted atom
298            ibn_tot = ibn_tot + 1
299            ibnd_mat(:,ibn_tot) = (/ iat, jat /)
300          end if
301        else  ! jat is a border atom
302          ibn_tot2 = ibn_tot2 + 1
303          ibnd_dum(:,ibn_tot2) = (/ iat, jat /)
304        end if
305      end do
306    end do
307 
308    if(ibn_tot2 > (6*nfitmax)) then
309      write(msg,'(3a,i8,2a)')&
310 &     'ERROR: BOND_ATOM_INIT',ch10,&
311 &     'IBN_TOT2 =  ',ibn_tot2,ch10,&
312 &     ' ibnd_dum out of bounds, ibn_tot2 too large '
313      ABI_ERROR(msg)
314    end if
315 
316   !--reorder to keep 'variational' bonds first :
317    ibn_tots = ibn_tot + ibn_tot2
318    do ii = ibn_tot+1,ibn_tots
319      ibnd_mat(:,ii) = ibnd_dum(:,ii-ibn_tot)
320    end do
321 
322    ABI_FREE (ibnd_dum)
323  end subroutine bond_compute

bond_lotf/bond_dealloc [ Functions ]

[ Top ] [ bond_lotf ] [ Functions ]

NAME

 bond_dealloc

FUNCTION

  deallocate variables

INPUTS

SOURCE

380  subroutine  bond_dealloc()
381 
382 ! *************************************************************************
383    ABI_FREE(ibnd_mat)
384    ABI_FREE(tafit)
385    ABI_FREE(ifit)
386    ABI_FREE(ibmat)
387    ABI_FREE(ibmat_large)
388    ABI_FREE(imat)
389  end subroutine bond_dealloc
390 
391 end module bond_lotf

bond_lotf/bond_fit_set [ Functions ]

[ Top ] [ bond_lotf ] [ Functions ]

NAME

 bond_fit_set

FUNCTION

  set nfitmax (or control it), ifit,nfit

INPUTS

SOURCE

336  subroutine bond_fit_set(nax,nfitdum)
337 
338   implicit none
339 
340   !Arguments ------------------------
341   integer,intent(in) :: nax
342   integer,intent(in) :: nfitdum
343   !Local -----------------------------
344   integer  :: i
345   character(len=500) :: msg
346 
347 ! *************************************************************************
348 
349    if(.not. allocated(ifit)) then
350      nfitmax = nfitdum + 1000
351      ABI_MALLOC(ifit,(nfitmax))
352    elseif(nfitdum > nfitmax) then
353      write(msg,'(a)')' BOND_FIT_SET : PROBLEM OF dimensionS !! '
354      ABI_ERROR(msg)
355    end if
356 
357    ifit = 0
358    nfit = 0
359    do i = 1, nax
360      if(tafit(i)) then
361        nfit = nfit + 1
362        ifit(nfit) = i
363      end if
364    end do
365 
366  end subroutine bond_fit_set

bond_lotf/bond_matrix_alloc [ Functions ]

[ Top ] [ bond_lotf ] [ Functions ]

NAME

 bond_matrix_alloc

FUNCTION

  allocate imat,ibmat,ibmat_large

INPUTS

SOURCE

179  subroutine bond_matrix_alloc(nax,nneigx)
180 
181   implicit none
182 
183   !Arguments ------------------------
184   integer,intent(in) :: nax
185   integer,intent(in) :: nneigx
186 
187 ! *************************************************************************
188 
189    ABI_MALLOC(ibmat,(nneigx,0:nfitmax))
190    ABI_MALLOC(ibmat_large,(nfitmax,nfitmax))
191    ABI_MALLOC(imat,(nax))
192  end subroutine bond_matrix_alloc

bond_lotf/bond_matrix_set [ Functions ]

[ Top ] [ bond_lotf ] [ Functions ]

NAME

 bond_matrix_set

FUNCTION

  Set or update bond matrix imat,ibmat,ibmat_large
  associates the bond to the atom neighlists

INPUTS

SOURCE

206  subroutine bond_matrix_set(nneig,neighl)
207 
208   implicit none
209 
210   !Arguments ------------------------
211   integer,intent(in) :: nneig(:)
212   integer,intent(in) :: neighl(:,:)
213   !Local ---------------------------
214   integer :: ibn,ibn2,iat,jat,ii,jj,j2
215 
216 ! *************************************************************************
217 
218    ibmat(:,:) = 0
219    ibmat_large(:,:) = 0
220    imat(:) = 0
221 
222    ibn = 0
223    do ii =1,nfit
224      iat = ifit(ii)
225      do jj = 1,nneig(iat)
226        jat = neighl(jj,iat)
227        if(jat > iat.AND.tafit(jat)) then
228          ibn  = ibn  + 1
229          ibmat(jj,ii) = ibn
230        end if
231      end do
232      imat(iat) = ii
233    end do
234 
235   !--ibmat_large: useful for the glue potential
236    ibn2 = 0
237    do ii =1,nfit
238      iat = ifit(ii)
239      do jj = 1,nneig(iat)
240        jat = neighl(jj,iat)
241        if(jat > iat.AND.tafit(jat)) then
242          ibn2  = ibn2  + 1
243          j2 = imat(jat)
244          ibmat_large(j2,ii) = ibn2
245          ibmat_large(ii,j2) = ibn2
246        end if
247      end do
248    end do
249  end subroutine bond_matrix_set

bond_lotf/bond_tafit_init [ Functions ]

[ Top ] [ bond_lotf ] [ Functions ]

NAME

 bond_tafit_init

FUNCTION

INPUTS

SOURCE

68  subroutine bond_tafit_init(nax)
69 
70   implicit none
71 
72   !Arguments ------------------------
73   integer,intent(in) :: nax
74 
75 ! *************************************************************************
76 
77    ABI_MALLOC(tafit,(nax))
78 
79   !--MMANCINI strange!!!!!
80   !  tafit(:nax) = tquant(:nax)
81    tafit(:nax) = .true.
82 
83  end subroutine bond_tafit_init