Statistiques
| Révision :

root / src / ReadGeom_vasp.f90 @ 12

Historique | Voir | Annoter | Télécharger (12,22 ko)

1 5 pfleura2
SUBROUTINE ReadGeom_vasp
2 5 pfleura2
3 12 pfleura2
!----------------------------------------------------------------------
4 12 pfleura2
!  Copyright 2003-2014 Ecole Normale Supérieure de Lyon,
5 12 pfleura2
!  Centre National de la Recherche Scientifique,
6 12 pfleura2
!  Université Claude Bernard Lyon 1. All rights reserved.
7 12 pfleura2
!
8 12 pfleura2
!  This work is registered with the Agency for the Protection of Programs
9 12 pfleura2
!  as IDDN.FR.001.100009.000.S.P.2014.000.30625
10 12 pfleura2
!
11 12 pfleura2
!  Authors: P. Fleurat-Lessard, P. Dayal
12 12 pfleura2
!  Contact: optnpath@gmail.com
13 12 pfleura2
!
14 12 pfleura2
! This file is part of "Opt'n Path".
15 12 pfleura2
!
16 12 pfleura2
!  "Opt'n Path" is free software: you can redistribute it and/or modify
17 12 pfleura2
!  it under the terms of the GNU Affero General Public License as
18 12 pfleura2
!  published by the Free Software Foundation, either version 3 of the License,
19 12 pfleura2
!  or (at your option) any later version.
20 12 pfleura2
!
21 12 pfleura2
!  "Opt'n Path" is distributed in the hope that it will be useful,
22 12 pfleura2
!  but WITHOUT ANY WARRANTY; without even the implied warranty of
23 12 pfleura2
!
24 12 pfleura2
!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
25 12 pfleura2
!  GNU Affero General Public License for more details.
26 12 pfleura2
!
27 12 pfleura2
!  You should have received a copy of the GNU Affero General Public License
28 12 pfleura2
!  along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>.
29 12 pfleura2
!
30 12 pfleura2
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr,
31 12 pfleura2
! for commercial licensing opportunities.
32 12 pfleura2
!----------------------------------------------------------------------
33 12 pfleura2
34 5 pfleura2
  Use Path_module
35 5 pfleura2
  Use Io_module
36 5 pfleura2
37 5 pfleura2
  IMPLICIT NONE
38 5 pfleura2
39 5 pfleura2
40 12 pfleura2
  CHARACTER(LCHARS) :: LineTmp,Line
41 5 pfleura2
42 5 pfleura2
  INTEGER(KINT), ALLOCATABLE :: NbAtType(:) !na
43 5 pfleura2
  INTEGER(KINT) :: NbType, NbTypeUser
44 5 pfleura2
  INTEGER(KINT) :: I, J, Iat, NAtP, Idx
45 5 pfleura2
  INTEGER(KINT) :: IGeom
46 5 pfleura2
  REAL(KREAL) :: Xt,Yt,Zt
47 5 pfleura2
  LOGICAL :: Debug, TChk
48 5 pfleura2
  CHARACTER(4), ALLOCATABLE :: FFFt(:,:) ! 3,Na
49 5 pfleura2
  LOGICAL, ALLOCATABLE :: FFFl(:,:) ! 3,Na
50 5 pfleura2
51 5 pfleura2
!!! for VASP inputs
52 5 pfleura2
! latx_loc used to check that the lattice parameters are the same for all geometries
53 5 pfleura2
  REAL(KREAL) :: lata_loc(3), latb_loc(3), latc_loc(3)
54 5 pfleura2
! VaspPar_loc is the local value of the vasp parametre
55 5 pfleura2
  REAL(KREAL) :: VaspPar_loc
56 5 pfleura2
! Vdir_loc is used to convert the read geometry to Cartesian.
57 5 pfleura2
! V_direct is set by the first geometry.
58 5 pfleura2
  CHARACTER(LCHARS) :: Vdir_loc
59 5 pfleura2
60 5 pfleura2
61 5 pfleura2
  INTERFACE
62 5 pfleura2
     function valid(string) result (isValid)
63 5 pfleura2
       CHARACTER(*), intent(in) :: string
64 5 pfleura2
       logical                  :: isValid
65 5 pfleura2
     END function VALID
66 5 pfleura2
  END INTERFACE
67 5 pfleura2
68 5 pfleura2
  debug=valid('Read_geom').or.valid('ReadGeom_vasp')
69 5 pfleura2
70 5 pfleura2
 if (debug) Call Header("Entering ReadGeom_vasp")
71 5 pfleura2
72 5 pfleura2
     ALLOCATE(FFF(3,nat),FFFl(3,Nat),FFFt(3,Nat))
73 5 pfleura2
     ! First geometry is a bit special for VASP as we have to set
74 5 pfleura2
     ! many things
75 5 pfleura2
     IGeom=1
76 5 pfleura2
     IF (DEBUG) WRITE(*,*) "Reading Geom :",IGeom
77 5 pfleura2
     READ(IOIN,'(A)') Vasp_Title
78 5 pfleura2
     READ(IOIN,*) Vasp_param
79 5 pfleura2
80 5 pfleura2
     READ(IOIN,*) Lat_a
81 5 pfleura2
     READ(IOIN,*) Lat_b
82 5 pfleura2
     READ(IOIN,*) Lat_c
83 5 pfleura2
84 5 pfleura2
     Lat_a=Lat_a*Vasp_param
85 5 pfleura2
     Lat_b=Lat_b*Vasp_param
86 5 pfleura2
     Lat_c=Lat_c*Vasp_param
87 5 pfleura2
88 5 pfleura2
     ALLOCATE(NbAtType(nat))
89 12 pfleura2
! PFL 2013.12.07 ->
90 12 pfleura2
! From VASP 5.3, POSCAR can have an additional line containing
91 12 pfleura2
! the names of the atoms
92 12 pfleura2
     VASP5=.FALSE.
93 5 pfleura2
     READ(IOIN,'(A)') Vasp_types
94 12 pfleura2
     LineTMp=AdjustL(Trim(Vasp_types))
95 12 pfleura2
     if (ICHAR(LineTmp(1:1))<48) THEN
96 12 pfleura2
        Call Warning("ReadGeom_Vasp","Illegal caracter found in reading Vasp types")
97 12 pfleura2
        STOP
98 12 pfleura2
     END IF
99 12 pfleura2
     If (ICHAR(LineTmp(1:1))>58) THEN
100 12 pfleura2
! we have a letter, the user has surely provided the names.
101 12 pfleura2
! We follow VASP and do not use them.
102 12 pfleura2
        VASP5=.TRUE.
103 12 pfleura2
        Vasp_types_user=Vasp_types
104 12 pfleura2
! we read the usual "number of atom by type" line
105 12 pfleura2
        READ(IOIN,'(A)') Vasp_types
106 12 pfleura2
     END IF
107 12 pfleura2
! <- PFL 2013.12.07
108 12 pfleura2
109 5 pfleura2
     ! First, how many different types ?
110 5 pfleura2
     NbAtType=0
111 5 pfleura2
     READ(Vasp_types,*,END=999) NbAtType
112 5 pfleura2
999  CONTINUE
113 5 pfleura2
     NbType=0
114 5 pfleura2
     NatP=0
115 5 pfleura2
     DO WHILE (NbAtType(NbType+1).NE.0)
116 5 pfleura2
        NbType=NbType+1
117 5 pfleura2
        Natp=Natp+NbAtType(NbType)
118 5 pfleura2
     END DO
119 5 pfleura2
120 5 pfleura2
121 5 pfleura2
     if (NAtp.NE.Nat) THEN
122 5 pfleura2
        WRITE(IOOUT,*) 'WARNING Number of atoms not consistent between NAMELIST &path and First geom'
123 5 pfleura2
        WRITE(IOOUT,*) "Using:",Natp
124 5 pfleura2
        DEALLOCATE(XyzGeomI, AtName)
125 5 pfleura2
        Nat=Natp
126 5 pfleura2
        ALLOCATE(XyzGeomI(NGeomI,3,Nat), AtName(NAt))
127 5 pfleura2
     END if
128 5 pfleura2
129 5 pfleura2
     ! Do we know the atom types ?
130 5 pfleura2
     IF (AtTypes(1).EQ.'  ') THEN
131 5 pfleura2
        ! user has not provided atom types... we have to find them by ourselves
132 5 pfleura2
        ! by looking into the POTCAR file...
133 5 pfleura2
        INQUIRE(File="POTCAR", EXIST=Tchk)
134 5 pfleura2
        IF (.NOT.Tchk) THEN
135 5 pfleura2
           WRITE(*,*) "ERROR! No AtTypes provided, and POTCAR file not found"
136 5 pfleura2
           STOP
137 5 pfleura2
        END IF
138 5 pfleura2
        OPEN(IOTMP,File="POTCAR")
139 5 pfleura2
        DO I=1,NbType
140 5 pfleura2
           Line='Empty'
141 5 pfleura2
           DO WHILE (Line(1:5).NE.'TITEL')
142 5 pfleura2
              READ(IOTMP,'(A)') Line
143 5 pfleura2
              Line=AdjustL(Line)
144 5 pfleura2
           END DO
145 5 pfleura2
           Line=adjustl(Line(9:))
146 5 pfleura2
           Idx=index(Line," ")
147 5 pfleura2
           Line=adjustl(Line(Idx+1:))
148 5 pfleura2
           AtTypes(I)=Line(1:2)
149 5 pfleura2
        END DO
150 5 pfleura2
        if (debug) WRITE(*,'(A,100(1X,A2))') "ReadG:VASP AtTypes",AtTypes(1:NbType)
151 5 pfleura2
        CLOSE(IOTMP)
152 5 pfleura2
153 5 pfleura2
     ELSE  !AtTypes(1).EQ.'  '
154 5 pfleura2
        ! user has provided atom types
155 5 pfleura2
        NbTypeUser=0
156 5 pfleura2
        DO WHILE (AtTypes(NbTypeUser+1).NE.'  ')
157 5 pfleura2
           NbTypeUser=NbTypeUser+1
158 5 pfleura2
        END DO
159 5 pfleura2
        IF (NbType.NE.NbTypeUser) THEN
160 5 pfleura2
           WRITE(*,*) "ERROR Read_Geom : NbType in POSCAR do not match AtTypes"
161 5 pfleura2
           STOP
162 5 pfleura2
        END IF
163 5 pfleura2
     END IF
164 5 pfleura2
165 5 pfleura2
     IAt=1
166 5 pfleura2
     DO I=1,NbType
167 5 pfleura2
        DO J=1,NbAtType(I)
168 5 pfleura2
           AtName(Iat)=AtTypes(I)
169 5 pfleura2
           Iat=Iat+1
170 5 pfleura2
        END DO
171 5 pfleura2
     END DO
172 5 pfleura2
     DEALLOCATE(NbAtType)
173 5 pfleura2
174 5 pfleura2
     NbTypes=NbType
175 12 pfleura2
! PFL 2013.12.07 ->
176 12 pfleura2
! A comment: in fact, this line might just not be 'Selective Dynamics'
177 5 pfleura2
     READ(IOIN,'(A)') Vasp_comment
178 12 pfleura2
     LineTMp=AdjustL(Trim(Vasp_comment))
179 12 pfleura2
     Call upcase(LineTMP)
180 12 pfleura2
     if (LineTMP(1:1)=="S") THEN
181 12 pfleura2
        Vasp_SelectD=.TRUE.
182 12 pfleura2
        READ(IOIN,'(A)') Vasp_direct
183 12 pfleura2
     ELSE
184 12 pfleura2
        Vasp_SelectD=.FALSE.
185 12 pfleura2
        Vasp_direct=Vasp_comment
186 12 pfleura2
     END IF
187 12 pfleura2
! <- PFL 2013.12.07
188 5 pfleura2
189 5 pfleura2
     V_direct=adjustl(Vasp_direct)
190 5 pfleura2
     Call upcase(V_direct)
191 5 pfleura2
192 5 pfleura2
     if (debug) WRITE(*,*) "V_direct=",Trim(V_Direct)
193 5 pfleura2
194 5 pfleura2
     DO I=1,Nat
195 12 pfleura2
        IF (Vasp_SelectD) THEN
196 12 pfleura2
           READ(IOIN,*) XyzGeomI(IGeom,1:3,I),FFF(1:3,I)
197 5 pfleura2
           DO J=1,3
198 5 pfleura2
              FFF(J,I)=AdjustL(FFF(J,I))
199 5 pfleura2
              FFFl(j,i)=(FFF(j,i)(1:1).EQ.'T').OR.(FFF(j,i)(1:1).EQ.'t')
200 5 pfleura2
           END DO
201 12 pfleura2
        ELSE
202 12 pfleura2
           READ(IOIN,*) XyzGeomI(IGeom,1:3,I)
203 12 pfleura2
           FFFl(1:3,i)=.TRUE.
204 12 pfleura2
           FFF(1:3,i)="T"
205 12 pfleura2
        END IF
206 12 pfleura2
     END DO
207 12 pfleura2
! PFL 2013.12.07 ->
208 12 pfleura2
! Reading the VASP manual more carefully, it is stated:
209 12 pfleura2
! ! The seventh line (or eighth line if 'selective dynamics' is switched on)
210 12 pfleura2
! ! specifies whether the atomic positions are provided in cartesian coordinates
211 12 pfleura2
! ! or in direct coordinates (respectively fractional coordinates). As in the
212 12 pfleura2
! ! file KPOINTS only the first character on the line is significant and the
213 12 pfleura2
! ! only key characters recognized by VASP are 'C', 'c', 'K' or 'k' for
214 12 pfleura2
! ! switching to the cartesian mode.
215 12 pfleura2
!
216 12 pfleura2
! We thus change V_direct to Direct if its first letter is not 'C' or 'K'
217 12 pfleura2
! (we have it all in upper case)
218 12 pfleura2
     IF ((V_direct(1:1).NE.'C').AND.(V_direct(1:1).NE.'K')) THEN
219 12 pfleura2
        V_direct='DIRECT'
220 12 pfleura2
     END IF
221 12 pfleura2
! <- PFL 2013.12.07
222 5 pfleura2
223 5 pfleura2
     IF (V_direct(1:6).EQ.'DIRECT') THEN
224 5 pfleura2
        DO I=1,Nat
225 5 pfleura2
           Xt=XyzGeomI(IGeom,1,I)
226 5 pfleura2
           Yt=XyzGeomI(IGeom,2,I)
227 5 pfleura2
           Zt=XyzGeomI(IGeom,3,I)
228 5 pfleura2
           XyzGeomI(IGeom,1,I)=Xt*lat_a(1)+Yt*lat_b(1)+Zt*lat_c(1)
229 5 pfleura2
           XyzGeomI(IGeom,2,I)=Xt*lat_a(2)+Yt*lat_b(2)+Zt*lat_c(2)
230 5 pfleura2
           XyzGeomI(IGeom,3,I)=Xt*lat_a(3)+Yt*lat_b(3)+Zt*lat_c(3)
231 5 pfleura2
        END DO
232 5 pfleura2
233 5 pfleura2
     END IF
234 5 pfleura2
235 5 pfleura2
     ALLOCATE(X0_vasp(Nat),Y0_vasp(Nat),Z0_vasp(Nat))
236 5 pfleura2
     X0_vasp=XyzGeomI(IGeom,1,:)
237 5 pfleura2
     Y0_vasp=XyzGeomI(IGeom,2,:)
238 5 pfleura2
     Z0_vasp=XyzGeomI(IGeom,3,:)
239 5 pfleura2
240 5 pfleura2
     if (debug)  THEN
241 5 pfleura2
        WRITE(*,*) Nat
242 5 pfleura2
        WRITE(*,*) "Geom ",IGeom,"/",NGeomI
243 5 pfleura2
        DO I=1,Nat
244 5 pfleura2
           WRITE(*,'(1X,A2,3(1X,F15.6))') AtName(I),XyzGeomI(IGeom,1:3,I)
245 5 pfleura2
        END DO
246 5 pfleura2
     END IF
247 5 pfleura2
248 5 pfleura2
249 5 pfleura2
     DO IGeom=2,NGeomI
250 5 pfleura2
!!! PFL 16th March 2008 ->
251 5 pfleura2
! Users might want to use different conditions for different geometries
252 5 pfleura2
! i.e. to mix Cartesian and Direct :-p
253 5 pfleura2
! thus we check for each geom which condition is used.
254 5 pfleura2
!
255 12 pfleura2
!! PFL 2013 Dec 07 ->
256 12 pfleura2
!!!!! For now we assume that all "POSCAR" used follow the same
257 12 pfleura2
!!!! pattern as the first one
258 12 pfleura2
!!! this has to be changed: reading the VASP header should be in
259 12 pfleura2
!!! a subroutine by itself and combined afterwards.
260 12 pfleura2
!!!! <- PFL 2013 Dec 07
261 5 pfleura2
262 5 pfleura2
! Title
263 5 pfleura2
        READ(IOIN,*) LineTmp
264 5 pfleura2
        !        if (debug) WRITE(*,*) "1 LineTmp:",Trim(LineTmp)
265 5 pfleura2
! Vasp_Param
266 5 pfleura2
        READ(IOIN,*) VaspPar_loc
267 5 pfleura2
        !        if (debug) WRITE(*,*) "2 LineTmp:",Trim(LineTmp)
268 5 pfleura2
! Lattice parameters
269 5 pfleura2
        READ(IOIN,*) Lata_loc
270 5 pfleura2
        !        if (debug) WRITE(*,*) "3 LineTmp:",Trim(LineTmp)
271 5 pfleura2
        READ(IOIN,*) latb_loc
272 5 pfleura2
        !        if (debug) WRITE(*,*) "4 LineTmp:",Trim(LineTmp)
273 5 pfleura2
        READ(IOIN,*) Latc_loc
274 5 pfleura2
        !        if (debug) WRITE(*,*) "5 LineTmp:",Trim(LineTmp)
275 12 pfleura2
! PFL 2013 Dec 07 ->
276 12 pfleura2
! From VASP 5.3, POSCAR can have an additional line containing
277 12 pfleura2
! the names of the atoms
278 12 pfleura2
! Even thought VASP5 flag has been set on the first geometry,
279 12 pfleura2
! we are quite flexible and let the user mix old and new POSCAR:
280 12 pfleura2
        READ(IOIN,'(A)') LineTmp
281 12 pfleura2
        LineTMp=AdjustL(Trim(Vasp_types))
282 12 pfleura2
        if (ICHAR(LineTmp(1:1))<48) THEN
283 12 pfleura2
           Call Warning("ReadGeom_Vasp","Illegal caracter found in reading Vasp types")
284 12 pfleura2
           STOP
285 12 pfleura2
        END IF
286 12 pfleura2
        If (ICHAR(LineTmp(1:1))>58) THEN
287 12 pfleura2
! we have a letter, the user has surely provided the names.
288 12 pfleura2
! we read the usual "number of atom by type" line
289 12 pfleura2
           READ(IOIN,*) LineTmp
290 5 pfleura2
        !        if (debug) WRITE(*,*) "6 LineTmp:",Trim(LineTmp)
291 12 pfleura2
        END IF
292 12 pfleura2
! Do we have Selective Dynamics ?
293 12 pfleura2
        READ(IOIN,'(A)') Line
294 12 pfleura2
        LineTMp=AdjustL(Trim(Line))
295 12 pfleura2
        Call upcase(LineTMP)
296 12 pfleura2
        if (LineTMP(1:1)=="S") THEN
297 12 pfleura2
           if (.NOT.Vasp_SelectD) Vasp_comment=AdjustL(Line)
298 12 pfleura2
           Vasp_SelectD=.TRUE.
299 12 pfleura2
           READ(IOIN,'(A)') Vdir_loc
300 12 pfleura2
        ELSE
301 12 pfleura2
           Vasp_SelectD=.FALSE.
302 12 pfleura2
           Vdir_loc=Line
303 12 pfleura2
     END IF
304 12 pfleura2
! <- PFL 2013.12.07
305 12 pfleura2
306 12 pfleura2
307 5 pfleura2
! V_direct
308 5 pfleura2
        Vdir_loc=adjustl(Vdir_loc)
309 5 pfleura2
        Call upcase(Vdir_loc)
310 12 pfleura2
! PFL 2013.12.07 ->
311 12 pfleura2
     IF ((Vdir_loc(1:1).NE.'C').AND.(Vdir_loc(1:1).NE.'K')) THEN
312 12 pfleura2
        Vdir_loc='DIRECT'
313 12 pfleura2
     END IF
314 12 pfleura2
! <- PFL 2013.12.07
315 12 pfleura2
316 12 pfleura2
317 5 pfleura2
        !        if (debug) WRITE(*,*) "8 LineTmp:",Trim(LineTmp)
318 5 pfleura2
319 5 pfleura2
        Lata_loc=Lata_loc*VaspPar_loc-Lat_a
320 5 pfleura2
        if (dot_product(lata_loc,Lata_loc).GE.1e-6) THEN
321 5 pfleura2
           WRITE(*,*) "Lattice parameter a is different between geometries 1 and ",IGeom
322 5 pfleura2
           WRITE(IOOUT,*) "Lattice parameter a is different between geometries 1 and ",IGeom
323 5 pfleura2
           STOP
324 5 pfleura2
        END IF
325 5 pfleura2
326 5 pfleura2
        Latb_loc=Latb_loc*VaspPar_loc-Lat_b
327 5 pfleura2
        if (dot_product(latb_loc,Latb_loc).GE.1e-6) THEN
328 5 pfleura2
           WRITE(*,*) "Lattice parameter b is different between geometries 1 and ",IGeom
329 5 pfleura2
           WRITE(IOOUT,*) "Lattice parameter b is different between geometries 1 and ",IGeom
330 5 pfleura2
           STOP
331 5 pfleura2
        END IF
332 5 pfleura2
333 5 pfleura2
        Latc_loc=Latc_loc*VaspPar_loc-Lat_c
334 5 pfleura2
        if (dot_product(latc_loc,Latc_loc).GE.1e-6) THEN
335 5 pfleura2
           WRITE(*,*) "Lattice parameter c is different between geometries 1 and ",IGeom
336 5 pfleura2
           WRITE(IOOUT,*) "Lattice parameter c is different between geometries 1 and ",IGeom
337 5 pfleura2
           STOP
338 5 pfleura2
        END IF
339 5 pfleura2
340 5 pfleura2
        DO I=1,Nat
341 12 pfleura2
           IF (Vasp_SelectD) THEN
342 12 pfleura2
              READ(IOIN,*) XyzGeomI(IGeom,1:3,I),FFFt(1:3,I)
343 12 pfleura2
              DO J=1,3
344 5 pfleura2
! A frozen atom in any geom is considered frozen in all geom
345 12 pfleura2
                 FFFt(j,I)=AdjustL(FFFt(j,i))
346 12 pfleura2
                 FFFl(j,i)=FFFl(j,i).AND.(FFFt(j,i).EQ.'T')
347 12 pfleura2
              END DO
348 12 pfleura2
           ELSE
349 12 pfleura2
              READ(IOIN,*) XyzGeomI(IGeom,1:3,I)
350 12 pfleura2
           END IF
351 12 pfleura2
352 5 pfleura2
        END DO
353 5 pfleura2
354 5 pfleura2
        IF (Vdir_loc(1:6).EQ.'DIRECT') THEN
355 5 pfleura2
           DO I=1,Nat
356 5 pfleura2
              Xt=XyzGeomI(IGeom,1,I)
357 5 pfleura2
              Yt=XyzGeomI(IGeom,2,I)
358 5 pfleura2
              Zt=XyzGeomI(IGeom,3,I)
359 5 pfleura2
              XyzGeomI(IGeom,1,I)=Xt*lat_a(1)+Yt*lat_b(1)+Zt*lat_c(1)
360 5 pfleura2
              XyzGeomI(IGeom,2,I)=Xt*lat_a(2)+Yt*lat_b(2)+Zt*lat_c(2)
361 5 pfleura2
              XyzGeomI(IGeom,3,I)=Xt*lat_a(3)+Yt*lat_b(3)+Zt*lat_c(3)
362 5 pfleura2
           END DO
363 5 pfleura2
        END IF
364 5 pfleura2
365 5 pfleura2
        if (debug)  THEN
366 5 pfleura2
           WRITE(*,*) Nat
367 5 pfleura2
           WRITE(*,*) "Geom ",IGeom,"/",NGeomI
368 5 pfleura2
           DO I=1,Nat
369 5 pfleura2
              WRITE(*,'(1X,A2,3(1X,F15.6))') AtName(I),XyzGeomI(IGeom,1:3,I)
370 5 pfleura2
           END DO
371 5 pfleura2
        END IF
372 5 pfleura2
373 5 pfleura2
     END DO
374 5 pfleura2
375 5 pfleura2
     DO I=1,Nat
376 5 pfleura2
        DO J=1,3
377 5 pfleura2
           if (FFFl(j,i)) FFF(j,i)="T"
378 5 pfleura2
        END DO
379 5 pfleura2
     END DO
380 5 pfleura2
381 5 pfleura2
 Deallocate(FFFl,FFFt)
382 5 pfleura2
383 5 pfleura2
384 5 pfleura2
  if (debug) Call Header("Exiting ReadGeom_vasp")
385 5 pfleura2
386 5 pfleura2
END SUBROUTINE ReadGeom_vasp