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
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 if (parse_is_defined(namespace, 'XYZ'//trim(what))) then ! read an xyz file
188 call check_duplicated(done)
190 gf%source = read_coords_xyz
191 ! no default, since we do not do this unless the input tag is present
192 call parse_variable(namespace, 'XYZ'//trim(what), '', str)
193
194 message(1) = "Reading " // trim(what) // " from " // trim(str)
195 call messages_info(1, namespace=namespace)
196
197 iunit = io_open(str, global_namespace, action='read', status='old')
198 read(iunit, *) gf%n
199
200 if (gf%n <= 0) then
201 write(message(1),'(a,i6)') "Invalid number of atoms ", gf%n
202 call messages_fatal(1, namespace=namespace)
203 end if
204
205 read(iunit, *) ! skip comment line
206
207 safe_allocate(gf%atom(1:gf%n))
209 do ia = 1, gf%n
210 safe_allocate(gf%atom(ia)%x(1:space%dim))
211 read(iunit,*) gf%atom(ia)%label, gf%atom(ia)%x
212 end do
213
214 call io_close(iunit)
215 end if
216
217 !%Variable XSFCoordinates
218 !%Type string
219 !%Section System::Coordinates
220 !%Description
221 !% Another option besides PDB and XYZ coordinates formats is XSF, as <a href=http://www.xcrysden.org/doc/XSF.html>defined</a>
222 !% by the XCrySDen visualization program. Specify the filename with this variable.
223 !% The lattice vectors will also be read from this file and the value of
224 !% <tt>PeriodicDimensions</tt> needs to be compatible with the first line
225 !% (<tt>CRYSTAL</tt>, <tt>SLAB</tt>, <tt>POLYMER</tt>, or <tt>MOLECULE</tt>).
226 !% The file should not contain <tt>ATOMS</tt>, <tt>CONVVEC</tt>, or <tt>PRIMCOORD</tt>.
227 !% NOTE: The coordinates are treated in the units specified by <tt>Units</tt> and/or <tt>UnitsInput</tt>.
228 !%End
229
230 if (parse_is_defined(namespace, 'XSF'//trim(what))) then ! read an xsf file
231 call check_duplicated(done)
232
233 gf%source = read_coords_xsf
234 ! no default, since we do not do this unless the input tag is present
235 call parse_variable(namespace, 'XSF'//trim(what), '', str)
236
237 message(1) = "Reading " // trim(what) // " from " // trim(str)
238 call messages_info(1, namespace=namespace)
239
240 iunit = io_open(str, namespace, action='read', status='old')
241
242 read(iunit, '(a256)') str
243 if (str(1:9) == 'ANIMSTEPS') then
244 read(str(10:), *) nsteps
245 read(iunit, *) str
246
247 !%Variable XSFCoordinatesAnimStep
248 !%Type integer
249 !%Default 1
250 !%Section System::Coordinates
251 !%Description
252 !% If an animated file is given with <tt>XSFCoordinates</tt>, this variable selects which animation step
253 !% will be used. The <tt>PRIMVEC</tt> block must be written for each step.
254 !%End
255 call parse_variable(namespace, 'XSFCoordinatesAnimStep', 1, step_to_use)
256 if (step_to_use < 1) then
257 message(1) = "XSFCoordinatesAnimStep must be > 0."
258 call messages_fatal(1, namespace=namespace)
259 else if (step_to_use > nsteps) then
260 write(message(1),'(a,i9)') "XSFCoordinatesAnimStep must be <= available number of steps ", nsteps
261 call messages_fatal(1, namespace=namespace)
262 end if
263 write(message(1),'(a,i9,a,i9)') 'Using animation step ', step_to_use, ' out of ', nsteps
264 call messages_info(1, namespace=namespace)
265 else
266 nsteps = 0
267 step_to_use = 1
268 end if
269
270 ! periodicity = 'CRYSTAL', 'SLAB', 'POLYMER', 'MOLECULE'
271 select case (trim(str))
272 case ('CRYSTAL')
273 periodic_dim = 3
274 case ('SLAB')
275 periodic_dim = 2
276 case ('POLYMER')
277 periodic_dim = 1
278 case ('MOLECULE')
279 periodic_dim = 0
280 case ('ATOMS')
281 call messages_not_implemented("Input from XSF file beginning with ATOMS", namespace=namespace)
282 case default
283 write(message(1),'(3a)') 'Line in file was "', trim(str), '" instead of CRYSTAL/SLAB/POLYMER/MOLECULE.'
284 call messages_fatal(1, namespace=namespace)
285 end select
286 if (periodic_dim /= space%periodic_dim) then
287 message(1) = "Periodicity in XSF input is incompatible with the value of PeriodicDimensions."
288 call messages_fatal(1, namespace=namespace)
289 end if
290
291 do istep = 1, step_to_use - 1
292 read(iunit, *) ! PRIMVEC
293 do jdir = 1, space%dim
294 read(iunit, *) ! lattice vectors
295 end do
296 read(iunit, *) ! PRIMCOORD istep
297 read(iunit, *) gf%n, int_one
298 do ia = 1, gf%n
299 read(iunit, *) ! atoms
300 end do
301 end do
302
303 read(iunit, '(a256)') str
304 if (str(1:7) /= 'PRIMVEC') then
305 write(message(1),'(3a)') 'Line in file was "', trim(str), '" instead of "PRIMVEC".'
306 call messages_warning(1, namespace=namespace)
307 end if
308 if (nsteps > 0) then
309 read(str(8:), *) istep
310 if (istep /= step_to_use) then
311 write(message(1), '(a,i9,a,i9)') 'Found PRIMVEC index ', istep, ' instead of ', step_to_use
312 call messages_fatal(1, namespace=namespace)
313 end if
314 end if
315
316 safe_allocate(gf%latvec(1:space%dim, 1:space%dim))
317 gf%latvec = m_zero
318 do jdir = 1, space%dim
319 read(iunit, *) gf%latvec(1:space%dim, jdir)
320 end do
321 gf%latvec = units_to_atomic(units_inp%length, gf%latvec)
322
323 read(iunit, '(a256)') str
324 if (str(1:9) /= 'PRIMCOORD') then
325 write(message(1),'(3a)') 'Line in file was "', trim(str), '" instead of "PRIMCOORD".'
326 call messages_warning(1, namespace=namespace)
327 end if
328 if (nsteps > 0) then
329 read(str(10:), *) istep
330 if (istep /= step_to_use) then
331 write(message(1), '(a,i9,a,i9)') 'Found PRIMCOORD index ', istep, ' instead of ', step_to_use
332 call messages_fatal(1, namespace=namespace)
333 end if
334 end if
335
336 read(iunit, *) gf%n, int_one
337 if (gf%n <= 0) then
338 write(message(1),'(a,i6)') "Invalid number of atoms ", gf%n
339 call messages_fatal(1, namespace=namespace)
340 end if
341 if (int_one /= 1) then
342 write(message(1),'(a,i6,a)') 'Number in file was ', int_one, ' instead of 1.'
343 call messages_warning(1, namespace=namespace)
344 end if
345 safe_allocate(gf%atom(1:gf%n))
346
347 ! TODO: add support for velocities as vectors here?
348 do ia = 1, gf%n
349 safe_allocate(gf%atom(ia)%x(1:space%dim))
350 read(iunit,*) gf%atom(ia)%label, gf%atom(ia)%x(1:space%dim)
351 end do
352
353 call io_close(iunit)
354 end if
355
356 !%Variable Coordinates
357 !%Type block
358 !%Section System::Coordinates
359 !%Description
360 !% If <tt>XYZCoordinates</tt>, <tt>PDBCoordinates</tt>, and <tt>XSFCoordinates</tt> were not found,
361 !% <tt>Octopus</tt> tries to read the coordinates for the atoms from the block <tt>Coordinates</tt>. The
362 !% format is quite straightforward:
363 !%
364 !% <tt>%Coordinates
365 !% <br>&nbsp;&nbsp;'C' | -0.56415 | 0.0 | 0.0 | no
366 !% <br>&nbsp;&nbsp;'O' | &nbsp;0.56415 | 0.0 | 0.0 | no
367 !% <br>%</tt>
368 !%
369 !% The first line defines a carbon atom at coordinates (-0.56415, 0.0, 0.0),
370 !% that is <b>not</b> allowed to move during dynamical simulations. The second line has
371 !% a similar meaning. This block obviously defines a carbon monoxide molecule, if the
372 !% input units are <tt>eV_Angstrom</tt>. The number of coordinates for each species
373 !% must be equal to the dimension of your space (generally 3).
374 !% Note that in this way it is possible to fix some of the atoms (this
375 !% is not possible when specifying the coordinates through a <tt>PDBCoordinates</tt> or
376 !% <tt>XYZCoordinates</tt> file). The last column is optional, and the default is yes.
377 !% It is always possible to fix <b>all</b> atoms using the <tt>MoveIons</tt> directive.
378 !%End
379
380 if (parse_block(namespace, trim(what), blk) == 0) then
381 call check_duplicated(done)
382
383 gf%n = parse_block_n(blk)
384
385 gf%source = read_coords_inp
386 gf%flags = ior(gf%flags, xyz_flags_move)
387
388 message(1) = "Reading " // trim(what) // " from " // trim(what) // " block"
389 call messages_info(1, namespace=namespace)
390
391 safe_allocate(gf%atom(1:gf%n))
392
393 do ia = 1, gf%n
394 ncol = parse_block_cols(blk, ia - 1)
395 if ((ncol < space%dim + 1) .or. (ncol > space%dim + 2)) then
396 write(message(1), '(3a,i2)') 'Error in block ', what, ' line #', ia
397 call messages_fatal(1, namespace=namespace)
398 end if
399 call parse_block_string (blk, ia - 1, 0, gf%atom(ia)%label)
400
401 safe_allocate(gf%atom(ia)%x(1:space%dim))
402 do jdir = 1, space%dim
403 call parse_block_float (blk, ia - 1, jdir, gf%atom(ia)%x(jdir))
404 end do
405 if (ncol == space%dim + 2) then
406 call parse_block_logical(blk, ia - 1, space%dim + 1, gf%atom(ia)%move)
407 else
408 gf%atom(ia)%move = .true.
409 end if
410 end do
411
412 call parse_block_end(blk)
413 end if
414
415 !%Variable ReducedCoordinates
416 !%Type block
417 !%Section System::Coordinates
418 !%Description
419 !% This block gives the atomic coordinates relative to the real
420 !% space unit cell. The format is the same as the
421 !% <tt>Coordinates</tt> block.
422 !%
423 !% Note that in Octopus the origin of coordinates is in the center
424 !% of the cell, so the coordinates inside the cell are in the
425 !% range [-0.5, 0.5).
426 !%
427 !% This block cannot be used with the <tt>minimum</tt> box shape.
428 !%End
429
430 ! This is valid only for Coordinates.
431 if (trim(what) == 'Coordinates') then
432 if (parse_block(namespace, 'Reduced'//trim(what), blk) == 0) then
433 call check_duplicated(done)
434
435 gf%n = parse_block_n(blk)
436
437 gf%source = read_coords_reduced
438 gf%flags = ior(gf%flags, xyz_flags_move)
439
440 message(1) = "Reading " // trim(what) // " from Reduced" // trim(what) // " block"
441 call messages_info(1, namespace=namespace)
442
443 safe_allocate(gf%atom(1:gf%n))
444
445 do ia = 1, gf%n
446 ncol = parse_block_cols(blk, ia - 1)
447 if ((ncol < space%dim + 1) .or. (ncol > space%dim + 2)) then
448 write(message(1), '(3a,i2)') 'Error in block ', what, ' line #', ia
449 call messages_fatal(1, namespace=namespace)
450 end if
451 call parse_block_string (blk, ia - 1, 0, gf%atom(ia)%label)
452
453 safe_allocate(gf%atom(ia)%x(1:space%dim))
454 do jdir = 1, space%dim
455 call parse_block_float (blk, ia - 1, jdir, gf%atom(ia)%x(jdir))
456 end do
457 if (ncol == space%dim + 2) then
458 call parse_block_logical(blk, ia - 1, space%dim + 1, gf%atom(ia)%move)
459 else
460 gf%atom(ia)%move = .true.
461 end if
462 end do
463
464 call parse_block_end(blk)
465 end if
466 end if
467
468 if (gf%source /= read_coords_reduced) then
469 ! adjust units
470 do ia = 1, gf%n
471 if (gf%source == read_coords_xyz) then
472 gf%atom(ia)%x = units_to_atomic(units_inp%length_xyz_file, gf%atom(ia)%x)
473 else
474 gf%atom(ia)%x = units_to_atomic(units_inp%length, gf%atom(ia)%x)
475 end if
476
477 end do
478 end if
479
480 pop_sub(read_coords_read)
481
482 contains
483
484 subroutine check_duplicated(done)
485 logical, intent(inout) :: done
486
488
489 if (.not. done) then
490 done = .true.
491 else
492 message(1) = 'Multiple definitions of '//trim(what)//' in the input file.'
493 call messages_fatal(1, namespace=namespace)
494 end if
495
497 end subroutine check_duplicated
498
499 end subroutine read_coords_read
500
501
502 ! ---------------------------------------------------------
503 subroutine read_coords_read_pdb(what, iunit, gf)
504 character(len=*), intent(in) :: what
505 integer, intent(in) :: iunit
506 type(read_coords_info), intent(inout) :: gf
507
508 character(len=80) :: record
509 character(len=6) :: record_name
510 integer :: na
511
512 push_sub(read_coords_read_pdb)
513
514 ! TODO: read the specification of the crystal structure.
515
516 ! First count number of atoms
517 rewind(iunit)
518 do
519 read(iunit, '(a80)', err=990, end=990) record
520 read(record, '(a6)') record_name
521 if (trim(record_name) == 'ATOM' .or. trim(record_name) == 'HETATM') then
522 gf%n = gf%n + 1
523 end if
524 end do
525990 continue
526
527 safe_allocate(gf%atom(1:gf%n))
528
529 ! read in the data
530 rewind(iunit)
531 na = 1
532 do
533 read(iunit, '(a80)', err=991, end=991) record
534 read(record, '(a6)') record_name
535 if (trim(record_name) == 'ATOM' .or. trim(record_name) == 'HETATM') then
536 read(record, '(12x,a4,1x,a3)') gf%atom(na)%label, gf%atom(na)%residue
537 call str_trim(gf%atom(na)%label)
538 gf%atom(na)%label = gf%atom(na)%label(1:1)
539 call str_trim(gf%atom(na)%residue)
540
541 read(record, '(30x,3f8.3)') gf%atom(na)%x(1:3)
542 ! PDB files are always in angstrom
543 ! See the format definition
544 ! https://www.wwpdb.org/documentation/file-format-content/format33/sect9.html#HETATM
545 gf%atom(na)%x = units_to_atomic(unit_angstrom, gf%atom(na)%x)
546
547 na = na + 1
548 end if
549 end do
550991 continue
551
552 pop_sub(read_coords_read_pdb)
553 end subroutine read_coords_read_pdb
554
555end module read_coords_oct_m
556
557!! Local Variables:
558!! mode: f90
559!! coding: utf-8
560!! End:
real(real64), parameter, public m_zero
Definition: global.F90:187
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
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)