atomes 1.1.16
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
prepdata.F90
Go to the documentation of this file.
1! This file is part of the 'atomes' software.
2!
3! 'atomes' is free software: you can redistribute it and/or modify it under the terms
4! of the GNU Affero General Public License as published by the Free Software Foundation,
5! either version 3 of the License, or (at your option) any later version.
6!
7! 'atomes' is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
8! without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
9! See the GNU General Public License for more details.
10!
11! You should have received a copy of the GNU Affero General Public License along with 'atomes'.
12! If not, see <https://www.gnu.org/licenses/>
13!
14! Copyright (C) 2022-2024 by CNRS and University of Strasbourg
15!
20
21#if defined (HAVE_CONFIG_H)
22# include <config.h>
23#endif
24
25INTEGER (KIND=c_int) FUNCTION alloc_data (N1, N2, N3) bind (C,NAME='alloc_data_')
26
27USE parameters
28
29IMPLICIT NONE
30
31INTEGER (KIND=c_int), INTENT(IN) :: n1, n2, n3
32
33INTERFACE
34 INTEGER FUNCTION allochem()
35 END FUNCTION
36END INTERFACE
37
38na=n1
39nsp=n2
40ns=n3
41
42if (allocated(fullpos)) deallocate(fullpos)
43allocate (fullpos(na,3,ns), stat=err)
44if (err .ne. 0) then
45 call show_error ("Impossible to allocate memory"//char(0), &
46 "Function: alloc_data"//char(0), "Table: FULLPOS"//char(0))
47 alloc_data = 0
48 goto 001
49endif
50if (allocated(lot)) deallocate(lot)
51allocate(lot(na), stat=err)
52if (err .ne. 0) then
53 call show_error ("Impossible to allocate memory"//char(0), &
54 "Function: alloc_data"//char(0), "Table: LOT"//char(0))
55 alloc_data = 0
56 goto 001
57endif
59
60001 continue
61
62END FUNCTION
63
64SUBROUTINE prep_spec (idatoms, nsps, open_apf) bind (C,NAME='prep_spec_')
65
66USE parameters
67USE mendeleiev
68
69IMPLICIT NONE
70
71INTEGER (KIND=c_int), DIMENSION(NSP), INTENT(IN) :: nsps
72real(kind=c_double), DIMENSION(NSP), INTENT(IN) :: idatoms
73INTEGER (KIND=c_int), INTENT(IN) :: open_apf
74
75CHARACTER (LEN=14) :: ELEM
76! Now we are calling the GTK+ routines
77
78do i=1, nsp
79 nbspbs(i) = nsps(i)
80 j = int(idatoms(i))
81 atomid(i) = j
82 if (open_apf .eq. 1) then
83 if (j .gt. 0) then
84 tl(i) = atsym(j)
85 elem = element(j)
86 else
87 tl(i) = "X "
88 elem = "Unknown"
89 endif
90 ! In C all string must be terminated by a CHAR(0)
91 ! To spec_data_
92 call spec_data (0, i-1, atomid(i), nbspbs(i), &
93 tl(i)//char(0), elem//char(0), &
94 0.0d0, 0.0d0, 0.0d0, 0.0d0)
95 endif
96enddo
97nbspbs(nsp+1) = na
98
99END SUBROUTINE
100
101SUBROUTINE read_chem (PMASS, PRAD, PNSCATT, PXSCATT) bind (C,NAME='read_chem_')
102
103USE parameters
104
105IMPLICIT NONE
106
107real(kind=c_double), DIMENSION(NSP) :: pmass, prad, pnscatt, pxscatt
108
109do i=1, nsp
110 mass(i) = pmass(i)
111 rvdw(i) = prad(i)
112 nscattl(i) = pnscatt(i)
113 xscattl(i) = pxscatt(i)
114enddo
115
116END SUBROUTINE
117
118SUBROUTINE read_data (PLOT, PBS) bind (C,NAME='read_data_')
119
120USE parameters
121
122IMPLICIT NONE
123
124INTEGER (KIND=c_int), DIMENSION(NA), INTENT(IN) :: PLOT
125INTEGER (KIND=c_int), DIMENSION(NSP), INTENT(IN) :: PBS
126
127do i=1, na
128 lot(i)=plot(i)+1
129enddo
130do i=1, nsp
131 nbspbs(i) = pbs(i)
132 xi(i) = dble(nbspbs(i)) / dble(na)
133enddo
134nbspbs(nsp+1) = na
135
136END SUBROUTINE
137
138SUBROUTINE read_pos (PCX, PCY, PCZ) bind (C,NAME='read_pos_')
139
140USE parameters
141
142real(kind=c_double), DIMENSION(NA*NS), INTENT(IN) :: pcx, pcy, pcz
143
144k=0
145do i=1, ns
146 do j=1, na
147 k=k+1
148 fullpos(j,1,i)=pcx(k)
149 fullpos(j,2,i)=pcy(k)
150 fullpos(j,3,i)=pcz(k)
151 enddo
152enddo
153
154END SUBROUTINE
155
156INTEGER FUNCTION send_pos (NPA, NPS, NLOT, POSTAB)
157
158INTEGER :: i, j, k, err
159INTEGER, INTENT(IN) :: npa, nps
160INTEGER, DIMENSION(:), INTENT(IN) :: nlot
161DOUBLE PRECISION, DIMENSION(:,:,:), INTENT(IN) :: postab
162DOUBLE PRECISION, DIMENSION(:), ALLOCATABLE :: xpos, ypos, zpos
163
164if (allocated(xpos)) deallocate(xpos)
165allocate(xpos(npa*nps), stat=err)
166if (err .ne. 0) then
167 call show_error ("Impossible to allocate memory"//char(0), &
168 "Function: SEND_POS"//char(0), "Table: XPOS"//char(0))
169 send_pos = 0
170 goto 001
171endif
172if (allocated(ypos)) deallocate(ypos)
173allocate(ypos(npa*nps), stat=err)
174if (err .ne. 0) then
175 call show_error ("Impossible to allocate memory"//char(0), &
176 "Function: SEND_POS"//char(0), "Table: YPOS"//char(0))
177 send_pos = 0
178 goto 001
179endif
180if (allocated(zpos)) deallocate(zpos)
181allocate(zpos(npa*nps), stat=err)
182if (err .ne. 0) then
183 call show_error ("Impossible to allocate memory"//char(0), &
184 "Function: SEND_POS"//char(0), "Table: ZPOS"//char(0))
185 send_pos = 0
186 goto 001
187endif
188
189k=0
190do i=1, nps
191do j=1, npa
192 k=k+1
193 xpos(k)=postab(j,1,i)
194 ypos(k)=postab(j,2,i)
195 zpos(k)=postab(j,3,i)
196enddo
197enddo
198
199! To 'save_pos_'
200call save_pos (npa, nlot, k, xpos, ypos, zpos)
201
202send_pos = 1
203
204001 continue
205
206if (allocated(xpos)) deallocate (xpos)
207if (allocated(ypos)) deallocate (ypos)
208if (allocated(zpos)) deallocate (zpos)
209
210END FUNCTION
211
212INTEGER (KIND=c_int) FUNCTION prep_data () bind (C,NAME='prep_data_')
213
214!
215! Data initialization
216!
217
218USE parameters
219USE mendeleiev
220
221IMPLICIT NONE
222
223CHARACTER (LEN=2), DIMENSION(MAXE) :: ttype
224
225INTERFACE
226 INTEGER FUNCTION send_pos(NPA, NPS, NLOT, POSTAB)
227 INTEGER, INTENT(IN) :: npa, nps
228 INTEGER, DIMENSION(:), INTENT(IN) :: nlot
229 DOUBLE PRECISION, DIMENSION(:,:,:), INTENT(IN) :: postab
230 END FUNCTION
231 INTEGER FUNCTION allochem()
232 END FUNCTION
233 INTEGER FUNCTION findid(NAMEAT)
234 CHARACTER (LEN=2), INTENT(IN) :: nameat
235 END FUNCTION
236 INTEGER FUNCTION chemistry ()
237 END FUNCTION
238END INTERFACE
239
240prep_data = 0
241
242if (allocated(lot)) deallocate(lot)
243allocate(lot(na), stat=err)
244if (err .ne. 0) then
245 call show_error ("Impossible to allocate memory"//char(0), &
246 "Function: prep_data"//char(0), "Table: LOT"//char(0))
247 goto 001
248endif
249lot(:)=0 ! Initialisation du tableau
250
251noc=1
252do noa=1, na
253 if (noa .eq. 1) then
254 lot(noa)=noc
255 ttype(noc)=tab_of_type(noa)
256 else if (lot(noa) .eq. 0) then
257 do nob=1, noc
258 if (tab_of_type(noa) .eq. ttype(nob)) then
259 lot(noa)= nob
260 endif
261 enddo
262 if (lot(noa) .eq. 0) then
263 noc=noc+1
264 lot(noa)=noc
265 ttype(noc)=tab_of_type(noa)
266 endif
267 endif
268enddo
269nsp=noc
270
271if (allochem() == 0) then
272 prep_data = 0
273 goto 001
274endif
275
276if (allocated(elemid)) deallocate(elemid)
277allocate(elemid(nsp), stat=err)
278if (err .ne. 0) then
279 call show_error ("Impossible to allocate memory"//char(0), &
280 "Function: prep_data"//char(0), "Table: ELEMID"//char(0))
281 goto 001
282endif
283
284nbspbs(:) = 0
285do noa=1, na
286 nob=lot(noa)
287 nbspbs(nob)=nbspbs(nob)+1
288enddo
289
290nbspbs(nsp+1)=na
291
292!
293! Link between index and label
294! Correspondance entre indice de type et label
295!
296
297do nob=1, nsp
298 tl(nob)=ttype(nob)
299enddo
300
301do noa=1, nsp
302 atomid(noa) = findid(tl(noa))
303 select case (atomid(noa))
304 case (-1)
305 prep_data = 0
306 goto 001
307 case (0)
308 elemid(noa) = "Dummy "//tl(noa)
309 mass(noa) = 1.0
310 rvdw(noa) = 0.5
311 nscattl(noa) = 0.0
312 xscattl(noa) = 0.0
313 case default
315 mass(noa) = amass(atomid(noa))
317 if(atomid(noa) < 105) then
319 ! Covalent radius are the defaut values
320 rvdw(noa) = arcov(atomid(noa))
321 ! To use ionic radius uncomment the following line
322 ! RVDW(NOA) = ARION(ATOMID(NOA))
323 ! To use Wander Waals radius uncomment the following line
324 ! RVDW(NOA) = ARVDW(ATOMID(NOA))
325 else
326 nscattl(noa) = 0.0d0
327 rvdw(noa) = 0.0d0
328 endif
329 end select
330 if (nscattl(noa) .eq. 0.00) then
331 call show_warning ("Element "//tl(noa)//" does not have neutron scattering length "//char(0), &
332 "If this is a bug please report it to"//char(0), package_bugreport//char(0))
333 endif
334enddo
335
336! Now we are calling the GTK+ routines
337! To init_data_
338call init_data (na, nsp, ns, 1)
339
340do i=1, nsp
341! In C all string must be terminated by a CHAR(0)
342! To spec_data_
343 call spec_data (1, i-1, atomid(i), nbspbs(i), &
344 tl(i)//char(0), elemid(i)//char(0), &
345 mass(i), rvdw(i), nscattl(i), xscattl(i))
346enddo
347
348if (chemistry() .eq. 1) then
350endif
351
352001 continue
353
354if (allocated(elemid)) deallocate(elemid)
355
356END FUNCTION
357
358INTEGER FUNCTION findid(NAMEAT)
359
360USE parameters
361USE mendeleiev
362
363INTEGER (KIND=c_int), EXTERNAL :: dummy_ask
364INTEGER :: elemt
365CHARACTER (LEN=2), INTENT(IN) :: nameat
366
367findid=0
368do elemt=1, maxe
369
370 if (atsym(elemt) .eq. nameat) then
371 findid=elemt
372 exit
373 endif
374
375enddo
376
377if (findid .eq. 0) then
378 call show_warning ("Problem with the atomic coordinates"//char(0), &
379 "Element "//nameat//" does not exist "//char(0), " "//char(0))
380 findid = dummy_ask("Do you want to use dummy atom(s) for unknown species "//nameat//" ?"//char(0));
381endif
382
383END FUNCTION
384
385REAL (KIND=c_double) FUNCTION set_mass (SW) bind (C,NAME='set_mass_')
386
387USE mendeleiev
388USE parameters
389
390INTEGER (KIND=c_int), INTENT(IN) :: sw
391
392set_mass = amass(sw)
393
394END FUNCTION
395
396REAL (KIND=c_double) FUNCTION set_radius (SW, RD) bind (C,NAME='set_radius_')
397
398USE mendeleiev
399USE parameters
400
401INTEGER (KIND=c_int), INTENT(IN) :: sw, rd
402DOUBLE PRECISION :: rad
403
404if (sw .lt. 108) then
405 if (rd .eq. 0) then
406 rad = arcov(sw)
407 else if (rd .eq. 1) then
408 rad = arion(sw)
409 else if (rd .eq. 2) then
410 rad = arvdw(sw)
411 else if (rd .eq. 3) then
412 rad = arcry(sw)
413 endif
414else
415 rad = 0.0
416endif
417set_radius = rad
418
419END FUNCTION
420
421REAL (KIND=c_double) FUNCTION set_neutron (SW) bind (C,NAME='set_neutron_')
422
423USE mendeleiev
424USE parameters
425
426INTEGER (KIND=c_int), INTENT(IN) :: sw
427
428if (sw .lt. 105) then
429 set_neutron = coheb(sw)
430else
431 set_neutron = 0.0
432endif
433
434END FUNCTION
435
436!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
437! Reconstruction des trajectoires réelles
438!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
439SUBROUTINE prep_pos (PINFO, FINFO) bind (C,NAME='prep_pos_')
440
441USE parameters
442
443IMPLICIT NONE
444
445INTEGER (KIND=c_int), INTENT(IN) :: PINFO, FINFO
446
447INTERFACE
448 INTEGER FUNCTION send_pos(NPA, NPS, NLOT, POSTAB)
449 INTEGER, INTENT(IN) :: NPA, NPS
450 INTEGER, DIMENSION(:), INTENT(IN) :: NLOT
451 DOUBLE PRECISION, DIMENSION(:,:,:), INTENT(IN) :: POSTAB
452 END FUNCTION
453END INTERFACE
454
455pbc=.false.
456if (pinfo .eq. 1) then
457 pbc=.true.
458endif
459
460frac=.false.
461if (finfo > 0) frac=.true.
462
463if (frac) then
464 do noc=1, ns
465 do noa=1, na
466 if (ncells .gt. 1) then
467 fullpos(noa,:,noc) = matmul(fullpos(noa,:,noc),the_box(noc)%fractocart)
468 else
469 fullpos(noa,:,noc) = matmul(fullpos(noa,:,noc),the_box(1)%fractocart)
470 endif
471 enddo
472 enddo
473endif
474
476
477END SUBROUTINE
integer function allochem()
Definition allochem.F90:22
integer(kind=c_int) function chemistry()
Definition chemistry.F90:22
void show_warning(char *warning, GtkWidget *win)
show warning
Definition interface.c:260
void show_error(char *error, int val, GtkWidget *win)
show error message
Definition interface.c:293
double precision, dimension(118), parameter amass
double precision, dimension(105), parameter arvdw
double precision, dimension(105), parameter coheb
double precision, dimension(105), parameter arcov
double precision, dimension(105), parameter arion
character(len=2), dimension(118), parameter atsym
double precision, dimension(105), parameter arcry
integer maxe
character(len=14), dimension(118), parameter element
double precision, dimension(:,:,:), allocatable fullpos
integer ncells
double precision, dimension(:), allocatable nscattl
double precision, dimension(:), allocatable mass
character(len=2), dimension(:), allocatable tab_of_type
character(len=15), dimension(:), allocatable elemid
double precision, dimension(:), allocatable xi
character(len=2), dimension(:), allocatable tl
integer nob
integer, dimension(:), allocatable nbspbs
integer noc
integer noa
integer err
type(lattice), dimension(:), allocatable, target the_box
logical frac
double precision, dimension(:), allocatable xscattl
integer, dimension(:), allocatable atomid
integer, dimension(:), allocatable lot
double precision, dimension(:), allocatable rvdw
logical pbc
integer nsp
integer(kind=c_int) function alloc_data(n1, n2, n3)
Definition prepdata.F90:26
subroutine read_data(plot, pbs)
Definition prepdata.F90:119
subroutine prep_pos(pinfo, finfo)
Definition prepdata.F90:440
integer(kind=c_int) function prep_data()
Definition prepdata.F90:213
subroutine read_chem(pmass, prad, pnscatt, pxscatt)
Definition prepdata.F90:102
integer function send_pos(npa, nps, nlot, postab)
Definition prepdata.F90:157
subroutine read_pos(pcx, pcy, pcz)
Definition prepdata.F90:139
real(kind=c_double) function set_mass(sw)
Definition prepdata.F90:386
integer function findid(nameat)
Definition prepdata.F90:359
real(kind=c_double) function set_radius(sw, rd)
Definition prepdata.F90:397
real(kind=c_double) function set_neutron(sw)
Definition prepdata.F90:422
subroutine prep_spec(idatoms, nsps, open_apf)
Definition prepdata.F90:65