Octopus
read_coords.F90
Go to the documentation of this file.
1
2!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23 use debug_oct_m
24 use global_oct_m
25 use io_oct_m
28 use parser_oct_m
30 use space_oct_m
31 use string_oct_m
32 use unit_oct_m
34
35 implicit none
36
37 private
38 public :: &
43
45 integer, public, parameter :: &
46 READ_COORDS_ERR = 0, &
47 read_coords_pdb = 1, &
48 read_coords_xyz = 2, &
49 read_coords_inp = 3, &
52
54 integer, public, parameter :: &
55 XYZ_FLAGS_RESIDUE = 1, &
56 xyz_flags_charge = 2, &
58
60 character(len=15) :: label
61 real(real64), allocatable :: x(:)
62
63 real(real64) :: charge = m_huge
64 character(len=3) :: residue
65 logical :: move
66 end type read_coords_atom
67
69 integer :: source
70 integer :: flags
71
72 integer :: n
73 type(read_coords_atom), allocatable :: atom(:)
74
75 real(real64), allocatable :: latvec(:,:)
76 end type read_coords_info
77
78
79contains
80
81 ! ---------------------------------------------------------
82 subroutine read_coords_init(gf)
83 type(read_coords_info), intent(out) :: gf
84
85 push_sub(read_coords_init)
86
87 gf%source = read_coords_err
88 gf%flags = 0
89 gf%n = 0
90
91 pop_sub(read_coords_init)
92 end subroutine read_coords_init
93
94
95 ! ---------------------------------------------------------
96 subroutine read_coords_end(gf)
97 type(read_coords_info), intent(inout) :: gf
98
99 integer :: ia
100
101 push_sub(read_coords_end)
102
103 do ia = 1, gf%n
104 safe_deallocate_a(gf%atom(ia)%x)
105 end do
106 safe_deallocate_a(gf%atom)
107 call read_coords_init(gf)
108 safe_deallocate_a(gf%latvec)
109
110 pop_sub(read_coords_end)
111 end subroutine read_coords_end
112
113
114 ! ---------------------------------------------------------
115 subroutine read_coords_read(what, gf, space, namespace)
116 character(len=*), intent(in) :: what
117 type(read_coords_info), intent(inout) :: gf
118 class(space_t), intent(in) :: space
119 type(namespace_t), intent(in) :: namespace
120
121 integer :: ia, ncol, iunit, jdir, int_one, nsteps, istep, step_to_use, periodic_dim
122 type(block_t) :: blk
123 character(len=256) :: str
124 logical :: done
125
126 push_sub(read_coords_read)
127
128 done = .false.
129
130 !%Variable PDBCoordinates
131 !%Type string
132 !%Section System::Coordinates
133 !%Description
134 !% If this variable is present, the program tries to read the atomic coordinates
135 !% from the file specified by its value. The PDB (<a href=http://www.rcsb.org/pdb>Protein Data Bank</a>)
136 !% format is quite complicated, and it goes
137 !% well beyond the scope of this manual. You can find a comprehensive
138 !% description <a href=http://www.wwpdb.org/docs.html>here</a>.
139 !% From the plethora of instructions defined in the PDB standard, <tt>Octopus</tt>
140 !% only reads two, <tt>ATOM</tt> and <tt>HETATOM</tt>. From these fields, it reads:
141 !% <ul>
142 !% <li> columns 13-16: The species; in fact <tt>Octopus</tt> only cares about the
143 !% first letter - <tt>CA</tt> and <tt>CB</tt> will both refer to carbon - so elements whose
144 !% chemical symbol has more than one letter cannot be represented in this way.
145 !% So, if you want to run mercury (Hg), please use one of the other methods
146 !% to input the coordinates.
147 !% <li> columns 18-21: The residue. Ignored.
148 !% <li> columns 31-54: The Cartesian coordinates. The Fortran format is <tt>(3f8.3)</tt>.</li>
149 !% <li> columns 61-65: Classical charge of the atom. Required if reading classical atoms, ignored otherwise.
150 !% The Fortran format is <tt>(f6.2)</tt>.</li>
151 !% </ul>
152 !% NOTE: The coordinates are treated in the units specified by <tt>Units</tt> and/or <tt>UnitsInput</tt>.
153 !%End
155 if (parse_is_defined(namespace, 'PDB'//trim(what))) then
158 gf%source = read_coords_pdb
159 gf%flags = ior(gf%flags, xyz_flags_residue)
160 gf%flags = ior(gf%flags, xyz_flags_charge)
162 ! no default, since we do not do this unless the input tag is present
163 call parse_variable(namespace, 'PDB'//trim(what), '', str)
164
165 message(1) = "Reading " // trim(what) // " from " // trim(str)
166 call messages_info(1, namespace=namespace)
167
168 iunit = io_open(str, namespace, action='read')
169 call read_coords_read_pdb(what, iunit, gf)
170 call io_close(iunit)
171 end if
172
173 !%Variable XYZCoordinates
174 !%Type string
175 !%Section System::Coordinates
176 !%Description
177 !% If <tt>PDBCoordinates</tt> is not present, the program reads the atomic coordinates from
178 !% the XYZ file specified by the variable <tt>XYZCoordinates</tt> -- in case this variable
179 !% is present. The XYZ format is very simple: The first line of the file has an integer
180 !% indicating the number of atoms. The second can contain comments that are simply ignored by
181 !% <tt>Octopus</tt>. Then there follows one line per atom, containing the chemical species and
182 !% the Cartesian coordinates of the atom.
183 !%
184 !% If you want to specify the unit of the XYZ file, you can use the variable <tt>UnitsXYZFiles</tt>.
185 !%End
186
187 !Here we use the global namespace to look for geometry.
188 if (parse_is_defined(namespace, 'XYZ'//trim(what))) then ! read an xyz file
190
191 gf%source = read_coords_xyz
192 ! no default, since we do not do this unless the input tag is present
193 call parse_variable(namespace, 'XYZ'//trim(what), '', str)
194
195 message(1) = "Reading " // trim(what) // " from " // trim(str)
196 call messages_info(1, namespace=namespace)
197
198 iunit = io_open(str, global_namespace, action='read', status='old')
199 read(iunit, *) gf%n
200
201 if (gf%n <= 0) then
202 write(message(1),'(a,i6)') "Invalid number of atoms ", gf%n
203 call messages_fatal(1, namespace=namespace)
204 end if
205
206 read(iunit, *) ! skip comment line
207
208 safe_allocate(gf%atom(1:gf%n))
209
210 do ia = 1, gf%n
211 safe_allocate(gf%atom(ia)%x(1:space%dim))
212 read(iunit,*) gf%atom(ia)%label, gf%atom(ia)%x
213 end do
214
215 call io_close(iunit)
216 end if
217
218 !%Variable XSFCoordinates
219 !%Type string
220 !%Section System::Coordinates
221 !%Description
222 !% Another option besides PDB and XYZ coordinates formats is XSF, as <a href=http://www.xcrysden.org/doc/XSF.html>defined</a>
223 !% by the XCrySDen visualization program. Specify the filename with this variable.
224 !% The lattice vectors will also be read from this file and the value of
225 !% <tt>PeriodicDimensions</tt> needs to be compatible with the first line
226 !% (<tt>CRYSTAL</tt>, <tt>SLAB</tt>, <tt>POLYMER</tt>, or <tt>MOLECULE</tt>).
227 !% The file should not contain <tt>ATOMS</tt>, <tt>CONVVEC</tt>, or <tt>PRIMCOORD</tt>.
228 !% NOTE: The coordinates are treated in the units specified by <tt>Units</tt> and/or <tt>UnitsInput</tt>.
229 !%End
230
231 if (parse_is_defined(namespace, 'XSF'//trim(what))) then ! read an xsf file
232 call check_duplicated(done)
233
234 gf%source = read_coords_xsf
235 ! no default, since we do not do this unless the input tag is present
236 call parse_variable(namespace, 'XSF'//trim(what), '', str)
237
238 message(1) = "Reading " // trim(what) // " from " // trim(str)
239 call messages_info(1, namespace=namespace)
240
241 iunit = io_open(str, namespace, action='read', status='old')
242
243 read(iunit, '(a256)') str
244 if (str(1:9) == 'ANIMSTEPS') then
245 read(str(10:), *) nsteps
246 read(iunit, *) str
247
248 !%Variable XSFCoordinatesAnimStep
249 !%Type integer
250 !%Default 1
251 !%Section System::Coordinates
252 !%Description
253 !% If an animated file is given with <tt>XSFCoordinates</tt>, this variable selects which animation step
254 !% will be used. The <tt>PRIMVEC</tt> block must be written for each step.
255 !%End
256 call parse_variable(namespace, 'XSFCoordinatesAnimStep', 1, step_to_use)
257 if (step_to_use < 1) then
258 message(1) = "XSFCoordinatesAnimStep must be > 0."
259 call messages_fatal(1, namespace=namespace)
260 else if (step_to_use > nsteps) then
261 write(message(1),'(a,i9)') "XSFCoordinatesAnimStep must be <= available number of steps ", nsteps
262 call messages_fatal(1, namespace=namespace)
263 end if
264 write(message(1),'(a,i9,a,i9)') 'Using animation step ', step_to_use, ' out of ', nsteps
265 call messages_info(1, namespace=namespace)
266 else
267 nsteps = 0
268 step_to_use = 1
269 end if
270
271 ! periodicity = 'CRYSTAL', 'SLAB', 'POLYMER', 'MOLECULE'
272 select case (trim(str))
273 case ('CRYSTAL')
274 periodic_dim = 3
275 case ('SLAB')
276 periodic_dim = 2
277 case ('POLYMER')
278 periodic_dim = 1
279 case ('MOLECULE')
280 periodic_dim = 0
281 case ('ATOMS')
282 call messages_not_implemented("Input from XSF file beginning with ATOMS", namespace=namespace)
283 case default
284 write(message(1),'(3a)') 'Line in file was "', trim(str), '" instead of CRYSTAL/SLAB/POLYMER/MOLECULE.'
285 call messages_fatal(1, namespace=namespace)
286 end select
287 if (periodic_dim /= space%periodic_dim) then
288 message(1) = "Periodicity in XSF input is incompatible with the value of PeriodicDimensions."
289 call messages_fatal(1, namespace=namespace)
290 end if
291
292 do istep = 1, step_to_use - 1
293 read(iunit, *) ! PRIMVEC
294 do jdir = 1, space%dim
295 read(iunit, *) ! lattice vectors
296 end do
297 read(iunit, *) ! PRIMCOORD istep
298 read(iunit, *) gf%n, int_one
299 do ia = 1, gf%n
300 read(iunit, *) ! atoms
301 end do
302 end do
303
304 read(iunit, '(a256)') str
305 if (str(1:7) /= 'PRIMVEC') then
306 write(message(1),'(3a)') 'Line in file was "', trim(str), '" instead of "PRIMVEC".'
307 call messages_warning(1, namespace=namespace)
308 end if
309 if (nsteps > 0) then
310 read(str(8:), *) istep
311 if (istep /= step_to_use) then
312 write(message(1), '(a,i9,a,i9)') 'Found PRIMVEC index ', istep, ' instead of ', step_to_use
313 call messages_fatal(1, namespace=namespace)
314 end if
315 end if
316
317 safe_allocate(gf%latvec(1:space%dim, 1:space%dim))
318 gf%latvec = m_zero
319 do jdir = 1, space%dim
320 read(iunit, *) gf%latvec(1:space%dim, jdir)
321 end do
322 gf%latvec = units_to_atomic(units_inp%length, gf%latvec)
323
324 read(iunit, '(a256)') str
325 if (str(1:9) /= 'PRIMCOORD') then
326 write(message(1),'(3a)') 'Line in file was "', trim(str), '" instead of "PRIMCOORD".'
327 call messages_warning(1, namespace=namespace)
328 end if
329 if (nsteps > 0) then
330 read(str(10:), *) istep
331 if (istep /= step_to_use) then
332 write(message(1), '(a,i9,a,i9)') 'Found PRIMCOORD index ', istep, ' instead of ', step_to_use
333 call messages_fatal(1, namespace=namespace)
334 end if
335 end if
336
337 read(iunit, *) gf%n, int_one
338 if (gf%n <= 0) then
339 write(message(1),'(a,i6)') "Invalid number of atoms ", gf%n
340 call messages_fatal(1, namespace=namespace)
341 end if
342 if (int_one /= 1) then
343 write(message(1),'(a,i6,a)') 'Number in file was ', int_one, ' instead of 1.'
344 call messages_warning(1, namespace=namespace)
345 end if
346 safe_allocate(gf%atom(1:gf%n))
347
348 ! TODO: add support for velocities as vectors here?
349 do ia = 1, gf%n
350 safe_allocate(gf%atom(ia)%x(1:space%dim))
351 read(iunit,*) gf%atom(ia)%label, gf%atom(ia)%x(1:space%dim)
352 end do
353
354 call io_close(iunit)
355 end if
356
357 !%Variable Coordinates
358 !%Type block
359 !%Section System::Coordinates
360 !%Description
361 !% If <tt>XYZCoordinates</tt>, <tt>PDBCoordinates</tt>, and <tt>XSFCoordinates</tt> were not found,
362 !% <tt>Octopus</tt> tries to read the coordinates for the atoms from the block <tt>Coordinates</tt>. The
363 !% format is quite straightforward:
364 !%
365 !% <tt>%Coordinates
366 !% <br>&nbsp;&nbsp;'C' | -0.56415 | 0.0 | 0.0 | no
367 !% <br>&nbsp;&nbsp;'O' | &nbsp;0.56415 | 0.0 | 0.0 | no
368 !% <br>%</tt>
369 !%
370 !% The first line defines a carbon atom at coordinates (-0.56415, 0.0, 0.0),
371 !% that is <b>not</b> allowed to move during dynamical simulations. The second line has
372 !% a similar meaning. This block obviously defines a carbon monoxide molecule, if the
373 !% input units are <tt>eV_Angstrom</tt>. The number of coordinates for each species
374 !% must be equal to the dimension of your space (generally 3).
375 !% Note that in this way it is possible to fix some of the atoms (this
376 !% is not possible when specifying the coordinates through a <tt>PDBCoordinates</tt> or
377 !% <tt>XYZCoordinates</tt> file). The last column is optional, and the default is yes.
378 !% It is always possible to fix <b>all</b> atoms using the <tt>MoveIons</tt> directive.
379 !%End
380
381 if (parse_block(namespace, trim(what), blk) == 0) then
382 call check_duplicated(done)
383
384 gf%n = parse_block_n(blk)
385
386 gf%source = read_coords_inp
387 gf%flags = ior(gf%flags, xyz_flags_move)
388
389 message(1) = "Reading " // trim(what) // " from " // trim(what) // " block"
390 call messages_info(1, namespace=namespace)
391
392 safe_allocate(gf%atom(1:gf%n))
393
394 do ia = 1, gf%n
395 ncol = parse_block_cols(blk, ia - 1)
396 if ((ncol < space%dim + 1) .or. (ncol > space%dim + 2)) then
397 write(message(1), '(3a,i2)') 'Error in block ', what, ' line #', ia
398 call messages_fatal(1, namespace=namespace)
399 end if
400 call parse_block_string (blk, ia - 1, 0, gf%atom(ia)%label)
401
402 safe_allocate(gf%atom(ia)%x(1:space%dim))
403 do jdir = 1, space%dim
404 call parse_block_float (blk, ia - 1, jdir, gf%atom(ia)%x(jdir))
405 end do
406 if (ncol == space%dim + 2) then
407 call parse_block_logical(blk, ia - 1, space%dim + 1, gf%atom(ia)%move)
408 else
409 gf%atom(ia)%move = .true.
410 end if
411 end do
412
413 call parse_block_end(blk)
414 end if
415
416 !%Variable ReducedCoordinates
417 !%Type block
418 !%Section System::Coordinates
419 !%Description
420 !% This block gives the atomic coordinates relative to the real
421 !% space unit cell. The format is the same as the
422 !% <tt>Coordinates</tt> block.
423 !%
424 !% Note that in Octopus the origin of coordinates is in the center
425 !% of the cell, so the coordinates inside the cell are in the
426 !% range [-0.5, 0.5).
427 !%
428 !% This block cannot be used with the <tt>minimum</tt> box shape.
429 !%End
430
431 ! This is valid only for Coordinates.
432 if (trim(what) == 'Coordinates') then
433 if (parse_block(namespace, 'Reduced'//trim(what), blk) == 0) then
434 call check_duplicated(done)
435
436 gf%n = parse_block_n(blk)
437
438 gf%source = read_coords_reduced
439 gf%flags = ior(gf%flags, xyz_flags_move)
440
441 message(1) = "Reading " // trim(what) // " from Reduced" // trim(what) // " block"
442 call messages_info(1, namespace=namespace)
443
444 safe_allocate(gf%atom(1:gf%n))
445
446 do ia = 1, gf%n
447 ncol = parse_block_cols(blk, ia - 1)
448 if ((ncol < space%dim + 1) .or. (ncol > space%dim + 2)) then
449 write(message(1), '(3a,i2)') 'Error in block ', what, ' line #', ia
450 call messages_fatal(1, namespace=namespace)
451 end if
452 call parse_block_string (blk, ia - 1, 0, gf%atom(ia)%label)
453
454 safe_allocate(gf%atom(ia)%x(1:space%dim))
455 do jdir = 1, space%dim
456 call parse_block_float (blk, ia - 1, jdir, gf%atom(ia)%x(jdir))
457 end do
458 if (ncol == space%dim + 2) then
459 call parse_block_logical(blk, ia - 1, space%dim + 1, gf%atom(ia)%move)
460 else
461 gf%atom(ia)%move = .true.
462 end if
463 end do
464
465 call parse_block_end(blk)
466 end if
467 end if
468
469 if (gf%source /= read_coords_reduced) then
470 ! adjust units
471 do ia = 1, gf%n
472 if (gf%source == read_coords_xyz) then
473 gf%atom(ia)%x = units_to_atomic(units_inp%length_xyz_file, gf%atom(ia)%x)
474 else
475 gf%atom(ia)%x = units_to_atomic(units_inp%length, gf%atom(ia)%x)
476 end if
477
478 end do
479 end if
480
481 pop_sub(read_coords_read)
482
483 contains
484
485 subroutine check_duplicated(done)
486 logical, intent(inout) :: done
487
489
490 if (.not. done) then
491 done = .true.
492 else
493 message(1) = 'Multiple definitions of '//trim(what)//' in the input file.'
494 call messages_fatal(1, namespace=namespace)
495 end if
496
498 end subroutine check_duplicated
499
500 end subroutine read_coords_read
501
502
503 ! ---------------------------------------------------------
504 subroutine read_coords_read_pdb(what, iunit, gf)
505 character(len=*), intent(in) :: what
506 integer, intent(in) :: iunit
507 type(read_coords_info), intent(inout) :: gf
508
509 character(len=80) :: record
510 character(len=6) :: record_name
511 integer :: na, ia
512
513 push_sub(read_coords_read_pdb)
514
515 ! TODO: read the specification of the crystal structure.
516
517 ! First count number of atoms
518 rewind(iunit)
519 gf%n = 0
520 do
521 read(iunit, '(a80)', err=990, end=990) record
522 read(record, '(a6)') record_name
523 if (trim(record_name) == 'ATOM' .or. trim(record_name) == 'HETATM') then
524 gf%n = gf%n + 1
525 end if
526 end do
527990 continue
528
529 safe_allocate(gf%atom(1:gf%n))
530 do ia = 1, gf%n
531 safe_allocate(gf%atom(ia)%x(1:3))
532 end do
533
534 ! read in the data
535 rewind(iunit)
536 na = 1
537 do
538 read(iunit, '(a80)', err=991, end=991) record
539 read(record, '(a6)') record_name
540 if (trim(record_name) == 'ATOM' .or. trim(record_name) == 'HETATM') then
541 read(record, '(12x,a4,1x,a3)') gf%atom(na)%label, gf%atom(na)%residue
542 gf%atom(na)%label = trim(adjustl(gf%atom(na)%label))
543 gf%atom(na)%label = gf%atom(na)%label(1:1)
544 gf%atom(na)%residue = trim(adjustl(gf%atom(na)%residue))
545
546 read(record, '(30x,3f8.3)') gf%atom(na)%x(1:3)
547 ! PDB files are always in angstrom
548 ! See the format definition
549 ! https://www.wwpdb.org/documentation/file-format-content/format33/sect9.html#HETATM
550 gf%atom(na)%x = units_to_atomic(unit_angstrom, gf%atom(na)%x)
551
552 na = na + 1
553 end if
554 end do
555991 continue
556
557 pop_sub(read_coords_read_pdb)
558 end subroutine read_coords_read_pdb
559
560end module read_coords_oct_m
561
562!! Local Variables:
563!! mode: f90
564!! coding: utf-8
565!! End:
real(real64), parameter, public m_huge
Definition: global.F90:206
real(real64), parameter, public m_zero
Definition: global.F90:188
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
type(namespace_t), public global_namespace
Definition: namespace.F90:132
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public parse_block_logical(blk, l, c, res)
Definition: parser.F90:638
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
integer, parameter, public read_coords_xsf
integer, parameter, public xyz_flags_charge
integer, parameter, public read_coords_reduced
subroutine, public read_coords_init(gf)
integer, parameter, public xyz_flags_move
subroutine, public read_coords_end(gf)
integer, parameter, public read_coords_xyz
subroutine read_coords_read_pdb(what, iunit, gf)
subroutine, public read_coords_read(what, gf, space, namespace)
integer, parameter, public read_coords_pdb
integer, parameter, public read_coords_inp
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
subroutine check_duplicated(done)
int true(void)