-
Notifications
You must be signed in to change notification settings - Fork 25
/
io_module.f90
3682 lines (3376 loc) · 125 KB
/
io_module.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
! -*- mode: F90; mode: font-lock -*-
! ------------------------------------------------------------------------------
! $Id$
! ------------------------------------------------------------------------------
! Module io_module
! ------------------------------------------------------------------------------
! Code area 9: general
! ------------------------------------------------------------------------------
!!****h* Conquest/io_module
!! NAME
!! io_module -- governs input and output of atom positions, types,
!! hopping and matrix writing
!! PURPOSE
!! Locate all input and output (bar comment writing) in one place so
!! that MPI-type calls can be isolated from main body of code
!! USES
!! GenComms, basic_types, common, datatypes, global_module,
!! group_module, matrix_module, species_module
!! AUTHOR
!! D.R.Bowler
!! CREATION DATE
!! 25/01/00 by D.R.Bowler
!! MODIFICATION HISTORY
!! 02/03/00 by DRB
!! Changed the unit chosen to write out C-matrix
!! 12/04/00 by DRB
!! write_c -> write_matrix
!! 18->19/04/00 by DRB
!! Implementing the new generic group_set and primary_set types
!! 1/05/00 by DRB
!! Incorporation into Conquest - removed the parameters read in
!! from stdin, and left it with reading make_prt.dat
!! 04/05/01 by DRB
!! Used the ParaDens tidied output to improve readability
!! 21/06/2001 dave
!! Added GenComms and RCS Id and Log tags and removed mpi
!! 18/03/2002 dave
!! Added static tag for object file id
!! 24/06/2002 dave
!! Changed read_atoms to use fdf io_assign routines (not raw
!! open/close)
!! 14:00, 04/02/2003 drb
!! Added dump/grab routines for charge, matrices, blips, local and
!! non-local pseudos for force checking
!! 16:43, 2003/04/15 dave
!! Added read_atomic_coordinates and changed read_mult
!! 16:55, 2003/06/09 dave
!! Further changes based on TM's work: now use atom_coord instead
!! of various x_atom etc variables Changes also to mapping and
!! numbering
!! 13:30, 22/09/2003 drb
!! Added constrained atoms
!! 10:09, 13/02/2006 drb
!! Removed all explicit references to data_ variables and rewrote
!! in terms of new matrix routines
!! 11:40, 10/11/2006 Veronika
!! Added option for reading pdb files
!! 20/11/2006 Veronika
!! Added option for writing out pdb files
!! 28/01/2008 Veronika
!! Fixed writing out pdb files - the pdb-style atom names are now
!! also written out
!! 2008/02/06 08:07 dave
!! Changed for output to file not stdout
!! 2008/07/31 ast
!! Added timers
!! 2008/09/01 08:20 dave
!! Added new Conquest input routines
!! 2009/07/24 ast
!! Filenames based on node/iteration number are now not limited to
!! 3 figures New routine get_file_name takes care of creating a
!! sufficiently long string
!! 2015/06/25 17:18 dave
!! Changed all file name lengths to 50 characters
!! 2016/03/15 14:02 dave
!! Added tests for existence of files to all grab_ routines and abort if
!! they don't exist
!! 2019/11/04 11:36 dave
!! Removed redundant code (old SFC routines)
module io_module
use datatypes, only: double
use numbers, only: zero
use global_module, only: io_lun
use GenComms, only: cq_abort, gcopy
use timer_module, only: start_timer, stop_timer, cq_timer
use timer_module, only: start_backtrace, stop_backtrace
use timer_stdclocks_module, only: tmr_std_matrices
use input_module, only: leqi, io_assign, io_close
implicit none
logical :: pdb_format, append_coords, pdb_output
character(len=1) :: pdb_altloc
character(len=80) :: pdb_template
!Maximum of wallclock time (in seconds): See subroutine 'check_stop' 2018.Jan.17 TM
real(double) :: time_max =zero
integer :: atom_output_threshold
!Name and Format of MatrixFile used in store_matrix
logical :: flag_MatrixFile_RankFromZero ! Starting from 0 in ##### (*matrix2.i**.p#####)
logical :: flag_MatrixFile_BinaryFormat ! Binary or Ascii
logical :: flag_MatrixFile_BinaryFormat_Grab ! Binary or Ascii for reading
logical :: flag_MatrixFile_BinaryFormat_Dump ! Binary or Ascii for writing
! flag_MatrixFile_BinaryFormat_Dump is set as flag_MatrixFile_BinaryFormat_Grab in the beginning,
! then, if we get a signal to finalise CONQUEST by calling "check_stop", it is set ast
! flag_MatrixFile_BinaryFormat_Dump_END.
! (during the job, we need to keep the format of the file, prepared in the last job.)
logical :: flag_MatrixFile_BinaryFormat_Dump_END ! set from Conquest_input
! System signature
! Moved here from read_and_write so that it can be used for extended XYZ output
! Moved here from initial_read_module to slove the dependence problem
character(len=80), save :: titles
logical :: flag_coords_xyz
!!***
contains
! --------------------------------------------------------------------
! Subroutine read_atomic_positions
! --------------------------------------------------------------------
!!****f* io_module/read_atomic_positions *
!!
!! NAME
!! read_atomic_positions
!! USAGE
!!
!! PURPOSE
!! Reads positions, species and constraint type for atoms
!! INPUTS
!!
!!
!! USES
!! datatypes, dimens, GenComms, global_module, species_module
!! AUTHOR
!! D.R.Bowler
!! CREATION DATE
!! 16:26, 2003/04/15
!! MODIFICATION HISTORY
!! 16:46, 2003/06/09 tm
!! id_glob_inv, atom_coord, species_glob is added.
!! and labelling in make_prt.dat is changed
!! 13:31, 22/09/2003 drb
!! Added reading of flags to constrain atom movement in three directions
!! 2006/10/09 08:21 dave
!! Tidying up output
!! 20/11/2006 Veronika
!! Corrected reading of constraints from pdb file
!! 15:05, 27/04/2007 drb
!! Check to ensure right number of species added
!! 2007/06/28, 21:37 mt + drb
!! Added coordinate wrapping to non-pdb coordinates.
!! 2007/10/11 Veronika
!! Added coordinate wrapping to pdb coordinates
!! Added coordinate wrapping for coordinates outside
!! the [-1*cell parameter; 2*cell parameter] range
!! 2013/07/01 M.Arita & T.Miyazaki
!! Added shift_in_bohr when wrapping atoms and allocation of atom_coord_diff
!! 2013/08/20 M.Arita
!! Bug fix & correct the if-statement
!! 2015/06/08 lat
!! Added experimental backtrace
!! 2018/01/22 tsuyoshi (with dave)
!! Allocate atom_coord_diff for all calculations
!! 2019/04/04 14:17 dave
!! Correct bug in wrapping with non-fractional coordinates and Angstroms
!! 2020/07/27 tsuyoshi
!! Added atom_vels
!! 2020/10/07 tsuyoshi
!! Removed allocation of atom_vels (moved to "control")
!! 2022/12/09 16:52 dave
!! Added simple check for Cartesian coordinates when fractional flag set and vice versa
!! SOURCE
!!
subroutine read_atomic_positions(filename)
use datatypes
use dimens, only: r_super_x, r_super_y, r_super_z
use global_module, only: x_atom_cell, y_atom_cell, z_atom_cell, &
ni_in_cell, numprocs, &
flag_fractional_atomic_coords, rcellx, &
rcelly, rcellz, id_glob, iprint_init, &
id_glob_inv, atom_coord, species_glob, &
flag_move_atom, area_init, shift_in_bohr, &
runtype,atom_coord_diff,id_glob_old,id_glob_inv_old
use species_module, only: species, species_label, n_species
use GenComms, only: inode, ionode, cq_abort, cq_warn
use memory_module, only: reg_alloc_mem, type_dbl, type_int
use units, only: AngToBohr
use units, only: dist_units, ang
use numbers, only: RD_ERR, zero, one, two
! Passed variables
character(len=*) :: filename
! Local variables
type(cq_timer) :: backtrace_timer
integer :: lun, i, spec, stat, n_wrong_coords
logical :: movex, movey, movez
real(double) :: x, y, z
real(double), dimension(3) :: cell
! Local variables for reading pdb files
integer :: ios, j
integer :: num_move_atom
real(double) :: num_move_atom_real
real(double), dimension(3) :: angle
character(len=80) :: pdb_line
character(len=2) :: atom_name
!****lat<$
call start_backtrace(t=backtrace_timer,who='read_atomic_positions',where=9,level=2)
!****lat>$
if(inode==ionode) then
! read a pdb file if General.pdb == .true.)
pdb: if (pdb_format) then
if(flag_fractional_atomic_coords) &
call cq_abort('Pdb file format does not support &
&fractional coordinates.')
if (iprint_init>2) &
write(io_lun,fmt='(4x,a)') 'Entering read_atomic_positions, pdb file'
call io_assign(lun)
! Go through the file and count the atoms
open( unit=lun, file=filename, status='old', iostat=ios)
if ( ios > 0 ) call cq_abort('Reading pdb file: file error')
if (iprint_init > 2) &
write (io_lun,'(5x,a)') &
'Counting atoms, checking for alternate locations'
ni_in_cell = 0
first: do
read (lun, '(a)', iostat=ios) pdb_line
! End of file?
if (ios < 0) exit first
select case (pdb_line(1:6))
case ('ATOM ','HETATM')
if (lgt( pdb_line(17:17),' ')) then
! File with alternate locations and no
! General.PdbAltLoc specified?
!
! This check fails if General.PdbAltLoc is in the
! input file but has no value. In this case, ALL
! entries with altloc will be ignored. User
! beware. (Of course we can always sell it as a
! feature.)
if (leqi( pdb_altloc,' ') ) then
call io_close(lun)
call cq_abort('PDB file with alternate &
&locations, specify General.PdbAltLoc')
endif
! Count the alternate location entry?
if (leqi( pdb_line(17:17), pdb_altloc)) &
ni_in_cell = ni_in_cell + 1
else
ni_in_cell = ni_in_cell + 1
end if
end select
end do first
call io_close(lun)
if (iprint_init>0) &
write(io_lun,'(5x,a,i5)') 'Number of atoms: ', ni_in_cell
! Now read the file again and extracts the coordinates
open(unit=lun, file=filename, status='old', iostat=ios)
if (ios > 0) call cq_abort('Reading pdb file: file error')
if (iprint_init > 2) write (*,*) 'Reading atomic positions'
! Allocate array for coordinates
allocate(flag_move_atom(3,ni_in_cell),&
atom_coord(3,ni_in_cell),&
species_glob(ni_in_cell),STAT=stat)
if(stat/=0) &
call cq_abort("Failure to allocate coordinates: ",ni_in_cell)
! Are the following two lines necessary?
call reg_alloc_mem(area_init, 3*ni_in_cell,type_dbl)
call reg_alloc_mem(area_init, 4*ni_in_cell,type_int)
i = 0
second: do
read (lun, '(a)', iostat=ios) pdb_line
! End of file?
if (ios < 0) exit second
select case (pdb_line(1:6))
case ('ATOM ','HETATM')
if (lgt( pdb_line(17:17),' ')) then
! Accept the alternate location?
if ( leqi( pdb_line(17:17), pdb_altloc)) then
i = i + 1
read (pdb_line(31:38),'(f8.3)') atom_coord(1,i)
read (pdb_line(39:46),'(f8.3)') atom_coord(2,i)
read (pdb_line(47:54),'(f8.3)') atom_coord(3,i)
atom_coord(1,i) = AngToBohr * atom_coord(1,i)
atom_coord(2,i) = AngToBohr * atom_coord(2,i)
atom_coord(3,i) = AngToBohr * atom_coord(3,i)
read (pdb_line(13:14),'(a2)') atom_name
atom_name = adjustl (atom_name)
do j = 1, n_species
if (leqi (atom_name,species_label(j) ) ) then
species_glob(i) = j
if(species_glob(i)>n_species) then
write(io_lun,&
fmt='(4x,"** WARNING ! ** &
&Species incompatibility between &
&coordinates and input")')
call cq_abort("Species specified &
&greater than number in input &
&file: ",species_glob(i),&
n_species)
end if
exit
end if
end do
flag_move_atom(:,i) = .false.
read (pdb_line(61:66),'(f6.0)') num_move_atom_real
num_move_atom = floor(num_move_atom_real)
if (num_move_atom_real - (real(num_move_atom)) > &
RD_ERR ) &
call cq_abort('The constraints on atom(s) &
&have not been set. Check the pdb file.')
if (num_move_atom > 7) &
call cq_abort('Out of range value for &
&flag_move_atom. Check the pdb file.')
if (num_move_atom >= 4) then
flag_move_atom(3,i) = .true.
num_move_atom = num_move_atom - 4
end if
if (num_move_atom >= 2) then
flag_move_atom(2,i) = .true.
num_move_atom = num_move_atom - 2
end if
if (num_move_atom == 1) flag_move_atom(1,i) = .true.
!if (iprint_init > 0) &
! write (io_lun,fmt='(3x, i7, 3f15.8, i3, 3L2)') &
! i, atom_coord(1:3,i), species_glob(i), &
! flag_move_atom(1:3,i)
end if
else
i = i + 1
read (pdb_line(31:38),'(f8.3)') atom_coord(1,i)
read (pdb_line(39:46),'(f8.3)') atom_coord(2,i)
read (pdb_line(47:54),'(f8.3)') atom_coord(3,i)
atom_coord(1,i) = AngToBohr * atom_coord(1,i)
atom_coord(2,i) = AngToBohr * atom_coord(2,i)
atom_coord(3,i) = AngToBohr * atom_coord(3,i)
read (pdb_line(13:14),'(a2)') atom_name
atom_name = adjustl (atom_name)
do j = 1, n_species
if (leqi (atom_name,species_label(j) ) ) then
species_glob(i) = j
if(species_glob(i)>n_species) then
write (io_lun,&
fmt='(2x,"** WARNING ! ** Species &
&incompatibility between &
&coordinates and input")')
call cq_abort("Species specified greater &
&than number in input file: ",&
species_glob(i),n_species)
end if
exit
end if
end do
flag_move_atom(:,i) = .false.
read (pdb_line(61:66),'(f6.0)') num_move_atom_real
num_move_atom = floor(num_move_atom_real)
if (num_move_atom_real - (real(num_move_atom)) > RD_ERR ) &
call cq_abort('The constraints on atom(s) &
&have not been set. Check the pdb file.')
if (num_move_atom > 7) &
call cq_abort('Out of range value for &
&flag_move_atom. Check the pdb file.')
if (num_move_atom >= 4) then
flag_move_atom(3,i) = .true.
num_move_atom = num_move_atom - 4
end if
if (num_move_atom >= 2) then
flag_move_atom(2,i) = .true.
num_move_atom = num_move_atom - 2
end if
if (num_move_atom == 1) flag_move_atom(1,i) = .true.
!if (iprint_init > 0) &
! write (io_lun,fmt='(3x, i7, 3f15.8, i3, 3L2)')&
! i, atom_coord(1:3,i), species_glob(i), &
! flag_move_atom(1:3,i)
end if
case ('CRYST1')
read (pdb_line,'(6x,3f9.3,3f7.2,15x)') &
r_super_x, r_super_y, r_super_z, angle
r_super_x = AngToBohr * r_super_x
r_super_y = AngToBohr * r_super_y
r_super_z = AngToBohr * r_super_z
! The following write statement is also in
! initial_read_module, but I think it's better to have
! it (also) here, because then the cell size will be
! printed together with the coordinates
!write(io_lun,4) r_super_x, r_super_y, r_super_z
!4 format(/10x,'The simulation box has the following dimensions',/, &
! 10x,'a = ',f9.5,' b = ',f9.5,' c = ',f9.5,' a.u.')
end select
end do second
! Wrap coordinates
cell(1) = r_super_x
cell(2) = r_super_y
cell(3) = r_super_z
do i = 1, ni_in_cell
do j = 1, 3
! Introduce shift_in_bohr: (small shift for fractional coordinates may
! cause a problem for very large systems.)
!
! Originally, 'if'-statement in the next line was active, but since we need
! a common and strict rule to treat the atoms on the boundary of the unit cell,
! we need to do this wrapping with shift_in_bohr for all cases.
!
! Thus, I have commented out the next line. 23/Jul/2014 TM
! (coordinate of the atoms on the boundary of the unit cell must be 0.000*** not 0.999***.)
!
!if ((atom_coord(j,i) < zero) .or. (atom_coord(j,i) > cell(j))) &
atom_coord(j,i) = &
atom_coord(j,i) - floor((atom_coord(j,i)+shift_in_bohr)/cell(j)) * cell(j)
end do
end do
call io_close(lun)
else ! Read normal coordinate file
if(iprint_init>2) write(io_lun,'(10x,a40,a20)') 'Entering read_atomic_positions; reading ', filename
call io_assign(lun)
open(unit=lun,file=filename,status='old')
! Read supercell vector - for now it must be orthorhombic so
! we use x and y as dummy variables
read(lun,*) r_super_x, x, y
if(abs(x)>RD_ERR.OR.abs(y)>RD_ERR) call cq_warn('read_atomic_positions', &
'Non-orthorhombic simulation cells are not supported by CONQUEST')
read(lun,*) x,r_super_y, y
if(abs(x)>RD_ERR.OR.abs(y)>RD_ERR) call cq_warn('read_atomic_positions', &
'Non-orthorhombic simulation cells are not supported by CONQUEST')
read(lun,*) x,y,r_super_z
if(abs(x)>RD_ERR.OR.abs(y)>RD_ERR) call cq_warn('read_atomic_positions', &
'Non-orthorhombic simulation cells are not supported by CONQUEST')
read(lun,*) ni_in_cell
!2010.06.25 TM (Angstrom Units in coords file, but not pdb)
if(dist_units == ang) then
r_super_x = r_super_x*AngToBohr
r_super_y = r_super_y*AngToBohr
r_super_z = r_super_z*AngToBohr
endif
!2010.06.25 TM (Angstrom Units in coords file, but not pdb)
allocate(flag_move_atom(3,ni_in_cell),atom_coord(3,ni_in_cell),&
species_glob(ni_in_cell),STAT=stat)
if(stat/=0) &
call cq_abort("Failure to allocate coordinates: ",ni_in_cell)
call reg_alloc_mem(area_init, 6*ni_in_cell,type_dbl)
call reg_alloc_mem(area_init, 4*ni_in_cell,type_int)
n_wrong_coords = 0
do i=1,ni_in_cell
read(lun,*) x,y,z,species_glob(i),movex,movey,movez
if(species_glob(i)>n_species) then
write(io_lun,fmt='(6x,"** WARNING ! ** Species &
&incompatibility between coordinates and &
&input")')
call cq_abort("Species specified greater than number &
&in input file: ", species_glob(i), n_species)
end if
if(flag_fractional_atomic_coords) then
atom_coord(1,i) = x*r_super_x
atom_coord(2,i) = y*r_super_y
atom_coord(3,i) = z*r_super_z
! Simple check for Cartesian coordinates
if(abs(x)>two .and. abs(y)>two .and. abs(z)>two) n_wrong_coords = n_wrong_coords + 1
else
atom_coord(1,i) = x
atom_coord(2,i) = y
atom_coord(3,i) = z
!2010.06.25 TM (Angstrom Units in coords file, but not pdb)
! 2019/04/04 14:16 dave
! Moved here so that distances are corrected *before* wrapping below
if(dist_units == ang) atom_coord(:,i)=atom_coord(:,i)*AngToBohr
! Simple check for fractional coordinates
if(abs(x)<one .and. abs(y)<one .and. abs(z)<one) n_wrong_coords = n_wrong_coords + 1
end if
! Wrap coordinates
cell(1) = r_super_x
cell(2) = r_super_y
cell(3) = r_super_z
do j = 1, 3
! Introduce shift_in_bohr: (small shift for fractional coordinates may
! cause a problem for very large systems.)
!
! Originally, 'if'-statement in the next line was active, but since we need
! a common and strict rule to treat the atoms on the boundary of the unit cell,
! we need to do this wrapping with shift_in_bohr for all cases.
!
! Thus, I have commented out the next line. 23/Jul/2014 TM
! (coordinate of the atoms on the boundary of the unit cell must be 0.000*** not 0.999***.)
!
!if ((atom_coord(j,i) < zero) .or. (atom_coord(j,i) > cell(j))) &
atom_coord(j,i) = &
atom_coord(j,i) - floor((atom_coord(j,i)+shift_in_bohr)/cell(j)) * cell(j)
end do
flag_move_atom(1,i) = movex
flag_move_atom(2,i) = movey
flag_move_atom(3,i) = movez
end do
! Warn user if it looks like Cartesian coordinates with fractional flag set or vice versa
if(flag_fractional_atomic_coords .and. n_wrong_coords>2) then
call cq_warn("read_atomic_positions", &
"Expected fractional coordinates but many are greater than one: ",n_wrong_coords)
else if(n_wrong_coords>2) then
call cq_warn("read_atomic_positions", &
"Expected Cartesian coordinates but many are less than one: ",n_wrong_coords)
end if
call io_close(lun)
end if pdb
! Check for sensible number of processes
if(ni_in_cell<numprocs) call cq_abort("We must have at least one atom per process: ",ni_in_cell,numprocs)
end if
if((iprint_init>0) .or. (iprint_init==0.AND.ni_in_cell<atom_output_threshold)) &
call print_atomic_positions
call gcopy(ni_in_cell)
if(inode/=ionode) &
allocate(flag_move_atom(3,ni_in_cell), atom_coord(3,ni_in_cell),&
species_glob(ni_in_cell), STAT=stat)
if(stat/=0) call cq_abort("Failure to allocate coordinates: ",ni_in_cell)
allocate(id_glob(ni_in_cell), id_glob_inv(ni_in_cell),&
x_atom_cell(ni_in_cell), y_atom_cell(ni_in_cell),&
z_atom_cell(ni_in_cell),species(ni_in_cell),STAT=stat)
if(stat/=0) call cq_abort("Failure to allocate coordinates(2): ",ni_in_cell)
call gcopy(r_super_x)
call gcopy(r_super_y)
call gcopy(r_super_z)
call gcopy(species_glob,ni_in_cell)
call gcopy(atom_coord, 3, ni_in_cell)
call gcopy(flag_move_atom, 3, ni_in_cell)
rcellx = r_super_x
rcelly = r_super_y
rcellz = r_super_z
allocate(atom_coord_diff(3,ni_in_cell), STAT=stat)
if (stat.NE.0) call cq_abort('Error allocating atom_coord_diff: ', 3, ni_in_cell)
allocate(id_glob_old(ni_in_cell),id_glob_inv_old(ni_in_cell), STAT=stat)
if (stat.NE.0) call cq_abort('Error allocating id_glob_old/id_glob_inv_old: ', &
ni_in_cell)
atom_coord_diff=zero
!call gcopy(atom_coord_diff,3,ni_in_cell)
id_glob_old=0
id_glob_inv_old=0
!****lat<$
call stop_backtrace(t=backtrace_timer,who='read_atomic_positions')
!****lat>$
return
end subroutine read_atomic_positions
!!***
! --------------------------------------------------------------------
! Subroutine write_atomic_positions
! --------------------------------------------------------------------
!!****f* io_module/write_atomic_positions *
!!
!! NAME
!! write_atomic_positions
!! USAGE
!!
!! PURPOSE
!! writes positions, species and constraint type for atoms
!! INPUTS
!!
!!
!! USES
!! datatypes, dimens, GenComms, global_module, species_module
!! AUTHOR
!! D.R.Bowler
!! CREATION DATE
!! 2006/11/07 08:02 dave
!! MODIFICATION HISTORY
!! 20/11/2006 Veronika
!! Added option for writing out pdb files
!! 2008/07/18
!! Added timers
!! 2013/07/01 M.Arita & T.Miyazaki
!! Changed formats from f18.10 to f25.17
!! SOURCE
!!
subroutine write_atomic_positions(filename, pdb_temp)
use datatypes
use numbers, only: zero
use dimens, only: r_super_x, r_super_y, r_super_z
use global_module, only: x_atom_cell, y_atom_cell, z_atom_cell, &
ni_in_cell, &
flag_fractional_atomic_coords, rcellx, &
rcelly, rcellz, iprint_init, atom_coord, &
species_glob, flag_move_atom, area_init, &
IPRINT_TIME_THRES3, min_layer
use species_module, only: species, species_label
use GenComms, only: inode, ionode, cq_abort
use units, only: BohrToAng
use timer_module
! Passed variables
character(len=*) :: filename, pdb_temp
! Local variables
integer :: lun, i, spec, stat, template, ios, j
real(double) :: x, y, z, beta
logical :: movex, movey, movez
character(len=80) :: pdb_line
character(len=2) :: atom_name
real(double), dimension(3) :: coords
type(cq_timer) :: tmr_l_tmp1
if(inode==ionode) then
call start_timer(tmr_l_tmp1,WITH_LEVEL)
if (pdb_output) then
if (iprint_init + min_layer > 2) write(io_lun,fmt='(6x,a)') 'Writing atomic positions'
call io_assign(lun)
! No appending of coords for a pdb file
open (unit = lun, file = 'output.pdb', status = 'replace', iostat = ios)
if ( ios > 0 ) call cq_abort('Cannot create file.')
! Open the template file
call io_assign(template)
open (unit = template, file = pdb_temp, status = 'old', iostat = ios)
if ( ios > 0 ) then
write(io_lun,*) 'Filename was: |',pdb_temp,'|'
call cq_abort('Reading template file: file error ',ios)
end if
! Read the template file, reprint all lines that are not ATOM, HETATM, or CRYST1
! to the output file
! If the line contains coords/cell parameters, replace them by the calculated ones.
! NB We do not update the connectivity information (the CONECT entries).
i = 0
do
read (template, '(a)', iostat=ios) pdb_line
! End of file?
if (ios < 0) exit
select case (pdb_line(1:6))
case ('ATOM ','HETATM')
if ( leqi( pdb_line(17:17), pdb_altloc) .or. leqi( pdb_line(17:17), ' ')) then
i = i + 1
coords(1) = BohrToAng * atom_coord(1,i)
coords(2) = BohrToAng * atom_coord(2,i)
coords(3) = BohrToAng * atom_coord(3,i)
atom_name = adjustr (species_label(species_glob(i))(1:2))
beta = 0.0
if (flag_move_atom(1,i)) beta = beta + 1
if (flag_move_atom(2,i)) beta = beta + 2
if (flag_move_atom(3,i)) beta = beta + 4
write (lun,'(a12,a2,a16,3f8.3,a6,f6.2,a14)') &
pdb_line(1:12), atom_name, pdb_line(15:30), &
coords(:), pdb_line(55:60), beta, pdb_line(67:80)
else
! If the entry is ATOM or HETATM but the alternate location
! does not match, just copy the line from the template
write (lun,'(a)') pdb_line
endif
case ('CRYST1')
coords(1) = BohrToAng * r_super_x
coords(2) = BohrToAng * r_super_y
coords(3) = BohrToAng * r_super_z
! Once we get to non-orthorhombic cells the line below has to be updated
write (lun,'(a6,3f9.3,3f7.2,a)') pdb_line(1:6), coords(:), &
90.0, 90.0, 90.0, pdb_line(56:80)
case default
write (lun,'(a)') pdb_line
end select
end do
call io_close(lun)
call io_close(template)
else
if (iprint_init + min_layer > 2) write(io_lun,fmt='(6x,a)') 'Writing atomic positions'
call io_assign(lun)
if(append_coords) then
open(unit=lun,file=filename,position='append')
write(lun,*)
else
open(unit=lun,file=filename)
end if
! Read supercell vector - for now it must be orthorhombic so we use x and y as dummy variables
write(lun,fmt='(3f25.17)') r_super_x, zero, zero
write(lun,fmt='(3f25.17)') zero, r_super_y, zero
write(lun,fmt='(3f25.17)') zero, zero, r_super_z
write(lun,fmt='(i8)') ni_in_cell
do i=1,ni_in_cell
if(flag_fractional_atomic_coords) then
write(lun,fmt='(3f25.17,i6,3L3)') atom_coord(1,i)/r_super_x,atom_coord(2,i)/r_super_y,&
atom_coord(3,i)/r_super_z, species_glob(i),flag_move_atom(1,i),flag_move_atom(2,i), &
flag_move_atom(3,i)
else
write(lun,fmt='(3f25.17,i6,3L3)') atom_coord(1,i),atom_coord(2,i),atom_coord(3,i), &
species_glob(i),flag_move_atom(1,i),flag_move_atom(2,i), &
flag_move_atom(3,i)
end if
end do
call io_close(lun)
end if ! pdb_output
call stop_print_timer(tmr_l_tmp1,"writing atomic positions",IPRINT_TIME_THRES3)
end if
return
end subroutine write_atomic_positions
!!***
! --------------------------------------------------------------------
! Subroutine read_mult
! --------------------------------------------------------------------
!!****f* io_module/read_mult *
!!
!! NAME
!! read_mult -- reads in data for multiplications
!! USAGE
!! read_mult(Processor id, partition data structure)
!! read_mult(myid,parts)
!! PURPOSE
!! Reads in data for partitions and distributes to all processors
!! INPUTS
!! integer :: myid - processor id
!! type(group_set) :: parts - contains info on partitions
!! USES
!! datatypes, global_module, group_module, basic_types
!! AUTHOR
!! D.R.Bowler
!! CREATION DATE
!!
!! MODIFICATION HISTORY
!! 21/06/2001 dave
!! Added cq_abort and gcopy
!! 2006/10/19 16:41 dave
!! Added filename for partition reading
!! 2013/04/03 L.Tong
!! Swapped the original create_sfc_partitions with the new
!! sfc_partitions_to_processors from sfc_partitions_module
!! 2013/08/20 M.Arita
!! Correct the if-statement
!! 2017/08/29 jack baker & dave
!! Removed r_super_x references (redundant)
!! SOURCE
!!
subroutine read_mult(myid,parts,part_file)
use datatypes
use global_module, only: numprocs, iprint_init, id_glob, &
ni_in_cell, x_atom_cell, y_atom_cell, &
z_atom_cell, atom_coord, id_glob_inv, &
species_glob, runtype, id_glob_old, id_glob_inv_old
use group_module, only: make_cc2, part_method, PYTHON, HILBERT
use basic_types
use species_module, only: species
use maxima_module, only: maxpartsproc, maxatomspart, &
maxatomsproc, maxpartscell
use construct_module, only: init_group
use sfc_partitions_module, only: sfc_partitions_to_processors
implicit none
! Passed variables
integer :: myid
type(group_set) :: parts
character(len=80) :: part_file
integer :: np_in_cell, mx_tmp_edge
integer :: irc,ierr,nnd,nnd1,np,np1,ind_part
integer :: n_cont,n_beg,ni,ni1,i
integer :: id_global
integer :: ind_file, stat
integer :: isendbuf(6)
real(double) :: rsendbuf(6)
type(cq_timer) :: backtrace_timer
!****lat<$
call start_backtrace(t=backtrace_timer,who='read_mult',where=9,level=2)
!****lat>$
!parts%ng_on_node = 0
if(part_method == PYTHON) then
if(myid==0) then
if(iprint_init>1) then
write(io_lun,'(10x,a29)') ' Reading partition data '
write(io_lun,'(10x,a29/)') '-----------------------------'
write(io_lun,'(10x,a29)') '*** Partition Data Starts ***'
end if
call read_partitions(parts,part_file)
if(iprint_init>1) write(io_lun,'(10x,a29)') '*** Partition Data Ends ***'
end if
isendbuf = 0
if(myid==0) then
isendbuf(1) = parts%ngcellx
isendbuf(2) = parts%ngcelly
isendbuf(3) = parts%ngcellz
isendbuf(4) = maxpartsproc
isendbuf(5) = maxatomspart
isendbuf(6) = maxatomsproc
np_in_cell = parts%ngcellx*parts%ngcelly*parts%ngcellz
endif
call gcopy(isendbuf,6)
if(myid/=0) then
parts%ngcellx = isendbuf(1)
parts%ngcelly = isendbuf(2)
parts%ngcellz = isendbuf(3)
maxpartsproc = isendbuf(4)
maxatomspart = isendbuf(5)
maxatomsproc = isendbuf(6)
np_in_cell = parts%ngcellx*parts%ngcelly*parts%ngcellz
maxpartscell = np_in_cell
mx_tmp_edge = max(parts%ngcellx,parts%ngcelly,parts%ngcellz)
call init_group(parts, maxpartsproc, mx_tmp_edge, &
np_in_cell, maxatomspart, numprocs)
endif
call gcopy(parts%ng_on_node,numprocs)
call gcopy(parts%inode_beg,numprocs)
call gcopy(parts%ngnode,parts%mx_gcell)
call gcopy(parts%nm_group,parts%mx_gcell)
call gcopy(parts%icell_beg,parts%mx_gcell)
call gcopy(id_glob, ni_in_cell)
call gcopy(id_glob_inv, ni_in_cell)
!if (leqi(runtype,'md') .OR. leqi(runtype,'cg')) then
if (.NOT. leqi(runtype,'static')) then
call gcopy(id_glob_old,ni_in_cell)
call gcopy(id_glob_inv_old,ni_in_cell)
endif
else if (part_method == HILBERT) then
if (iprint_init > 1.AND.myid==0) then
write(io_lun,'(10x,a33)' ) 'Partitioning using Hilbert curves'
write(io_lun,'(10x,a33/)') '---------------------------------'
end if
!call create_sfc_partitions(myid, parts)
call sfc_partitions_to_processors(parts)
np_in_cell = parts%ngcellx*parts%ngcelly*parts%ngcellz
if (iprint_init > 2.AND.myid==0) write(io_lun,fmt='(10x,a)') 'Finished partitioning'
end if
! inverse table to npnode
do np=1,np_in_cell
parts%inv_ngnode(parts%ngnode(np))=np
end do
call make_cc2(parts,numprocs)
do ni = 1, ni_in_cell
id_global= id_glob(ni)
x_atom_cell(ni) = atom_coord(1,id_global)
y_atom_cell(ni) = atom_coord(2,id_global)
z_atom_cell(ni) = atom_coord(3,id_global)
species(ni) = species_glob(id_global)
end do
!****lat<$
call stop_backtrace(t=backtrace_timer,who='read_mult')
!****lat>$
return
end subroutine read_mult
!!***
! --------------------------------------------------------------------
! Subroutine read_partitions
! --------------------------------------------------------------------
!!****f* io_module/read_partitions *
!!
!! NAME
!! read_atoms -- reads in partition data
!! USAGE
!! read_partitions(Partition data structure, number of processors)
!! PURPOSE
!! Reads in all data associated with partitions
!! (numbers within partition) on master processor and distributes to
!! all other processors
!! INPUTS
!! type(group_set) :: parts - partition data structure
!! integer :: nnode - number of processors
!! USES
!! datatypes, global_module, basic_types, common, ham
!! AUTHOR
!! D.R.Bowler
!! CREATION DATE
!! 16:50, 2003/04/15 dave
!! MODIFICATION HISTORY
!! 2006/10/19 16:42 dave
!! Added filename for partitions (defaults to make_prt.dat)
!! 2007/06/28 21:16 mt + drb
!! Added check for partition file access error.
!! 2011/11/17 10:19 dave
!! Bug fix to format for 1141
!! 2013/08/2013
!! Added storage of id_glob_old & id_glob_inv_old
!! SOURCE
!!
subroutine read_partitions(parts, part_file)
! Module usage
use datatypes
use global_module, only: numprocs, iprint_init, ni_in_cell, &
id_glob, id_glob_inv, id_glob_old, &
id_glob_inv_old, runtype
use maxima_module, only: maxpartsproc, maxatomspart, &
maxatomsproc, maxpartscell
use basic_types
use construct_module, only: init_group
use GenComms, only: cq_abort, myid
implicit none
! Passed variables
type(group_set) :: parts
character(len=80) :: part_file
! Local variables
type(cq_timer) :: backtrace_timer
integer :: nnode
integer :: nnd,nnd1,np,np1,ind_part,n_cont,n_beg,ni,ni1
integer :: irc,ierr,np_in_cell, lun, glob_count, map, ios
integer :: ind_global, ntmpx, ntmpy, ntmpz, ntmp1, ntmp2, ntmp3, mx_tmp_edge
!****lat<$
call start_backtrace(t=backtrace_timer,who='read_partitions',where=1,level=4)
!****lat>$
if(iprint_init>2.AND.myid==0) write(io_lun,fmt='(4x,a)') 'Entering read_partitions'
call io_assign(lun)
open(unit=lun, file=part_file, status='old', iostat=ios)
if ( ios > 0 ) call cq_abort('Reading partition file: file error')
!open(unit=lun,file=part_file)
! Read and check partitions along cell sides
read(lun,*) ntmpx,ntmpy,ntmpz
mx_tmp_edge = max(ntmpx,ntmpy,ntmpz)
! Find and check total number of partitions
np_in_cell=ntmpx*ntmpy*ntmpz
if(np_in_cell<=0) &
call cq_abort('read_partitions: no partitions ! ',np_in_cell)
! Read and check number of processors
read(lun,*) nnode
if((nnode<1).or.(nnode>numprocs)) &
call cq_abort('read_partitions: nodes wrong ',nnode,numprocs)
! Loop over nodes, reading bundle information
maxpartsproc = 0
maxatomspart = 0
maxatomsproc = 0
do nnd=1,nnode
ntmp3 = 0 ! count atoms on processor
! First node number, partitions and starting pointer
read(lun,*) nnd1,ntmp1, ntmp2
if(nnd1/=nnd) call cq_abort('read_partitions: node label wrong ',nnd1,nnd)
! Check for correct numbers of partitions
if(ntmp1<0) call cq_abort('read_partitions: parts on proc ', ntmp1)
if(ntmp1>maxpartsproc) maxpartsproc = ntmp1
! If there are partitions, read them in
if(ntmp1>0) then
! Check for correct maxima
if(ntmp2 + ntmp1 -1>np_in_cell) &
call cq_abort('read_partitions: too many parts in cell ',&
ntmp1 + ntmp2-1, np_in_cell)
do np=1,ntmp1
read(lun,*) np1,ind_part,n_cont,n_beg
if(np1/=np) &
call cq_abort('read_partitions: part label error ',np1,np)
if(n_cont>maxatomspart) maxatomspart = n_cont
ntmp3 = ntmp3 + n_cont
if(n_cont>0) then
do ni=1,n_cont
read(lun,*) ni1,map
if(ni1/=ni) &
call cq_abort('read_partitions: atom label wrong ',ni1,ni)
enddo
endif
enddo ! Loop over np = parts%ng_on_node
endif ! if(parts%ng_on_node>0)
if(ntmp3>maxatomsproc) maxatomsproc = ntmp3
enddo ! nnd = 1,nnode
if(iprint_init>3.AND.myid==0) &
write(io_lun,fmt='(6x,a,2i5)') 'Atoms proc max: ',maxatomsproc, maxpartsproc
call init_group(parts, maxpartsproc, mx_tmp_edge, np_in_cell, &
maxatomspart, numprocs)
maxpartscell = np_in_cell
close(unit=lun)
open(unit=lun, file=part_file, status='old', iostat=ios)
if ( ios > 0 ) call cq_abort('Reading partition file: file error')
!open(unit=lun,file=part_file)
! Read and check partitions along cell sides
read(lun,*) parts%ngcellx,parts%ngcelly,parts%ngcellz
if(iprint_init>0.AND.myid==0) &
write(io_lun,102) parts%ngcellx,parts%ngcelly,parts%ngcellz
if(max(parts%ngcellx,parts%ngcelly,parts%ngcellz)>parts%mx_gedge) then
call cq_abort('read_partitions: too many parts edge ', &
max(parts%ngcellx,parts%ngcelly,parts%ngcellz),parts%mx_gedge)
endif
! Read and check number of processors
read(lun,*) nnode
if(iprint_init>0.AND.myid==0) write(io_lun,110) nnode
if((nnode<1).or.(nnode>numprocs)) then
call cq_abort('read_partitions: nodes wrong ',nnode,numprocs)
endif
! Output header for atoms
! if(iprint_init==0) write(io_lun,*) ' Partn x y z Species Proc'
! Loop over nodes, reading bundle information
glob_count = 0
do nnd=1,nnode
! First node number, partitions and starting pointer