forked from m3g/packmol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
getinp.f90
1108 lines (1022 loc) · 36.5 KB
/
getinp.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
!
! Written by Leandro Martínez, 2009-2011.
! Copyright (c) 2009-2018, Leandro Martínez, Jose Mario Martinez,
! Ernesto G. Birgin.
!
! Subroutine getinp: subroutine that reads the input file
!
subroutine getinp()
use sizes
use compute_data, only : ntype, natoms, idfirst, nmols, ityperest, coor, restpars
use input
use usegencan
implicit none
integer :: i, k, ii, iarg, iline, idatom, iatom, in, lixo, irest, itype, itest,&
imark, ioerr, nloop0, iread, idfirstatom
double precision :: clen
character(len=strl) :: record, blank
logical :: inside_structure
! Clearing the blank character arrays
do i = 1, strl
blank(i:i) = ' '
end do
! Getting random seed and optional optimization parameters if set
seed = 1234567
randini = .false.
check = .false.
chkgrad = .false.
iprint1 = 2
iprint2 = 2
discale = 1.1d0
writeout = 10
maxit = 20
nloop = 0
nloop0 = 0
movefrac = 0.05
movebadrandom = .false.
precision = 1.d-2
writebad = .false.
add_amber_ter = .false.
add_box_sides = .false.
add_sides_fix = 0.d0
sidemax = 1000.d0
ioerr = 0
avoidoverlap = .true.
packall = .false.
use_short_tol = .false.
crd = .false.
inside_structure = .false.
do i = 1, nlines
if ( keyword(i,1).eq.'structure') inside_structure = .true.
if ( keyword(i,1).eq.'end' .and. &
keyword(i,2).eq.'structure') inside_structure = .false.
if(keyword(i,1).eq.'seed') then
read(keyword(i,2),*,iostat=ioerr) seed
if ( ioerr /= 0 ) exit
if ( seed == -1 ) call seed_from_time(seed)
else if(keyword(i,1).eq.'randominitialpoint') then
randini = .true.
else if(keyword(i,1).eq.'check') then
check = .true.
else if(keyword(i,1).eq.'writebad') then
writebad = .true.
else if(keyword(i,1).eq.'precision') then
read(keyword(i,2),*,iostat=ioerr) precision
if ( ioerr /= 0 ) exit
write(*,*) ' Optional precision set: ', precision
else if(keyword(i,1).eq.'movefrac') then
read(keyword(i,2),*,iostat=ioerr) movefrac
if ( ioerr /= 0 ) exit
write(*,*) ' Optional movefrac set: ', movefrac
else if(keyword(i,1).eq.'movebadrandom') then
movebadrandom = .true.
write(*,*) ' Will move randomly bad molecues (movebadrandom) '
else if(keyword(i,1).eq.'chkgrad') then
chkgrad = .true.
else if(keyword(i,1).eq.'writeout') then
read(keyword(i,2),*,iostat=ioerr) writeout
if ( ioerr /= 0 ) exit
write(*,*) ' Output frequency: ', writeout
else if(keyword(i,1).eq.'maxit') then
read(keyword(i,2),*,iostat=ioerr) maxit
if ( ioerr /= 0 ) exit
write(*,*) ' User defined GENCAN number of iterations: ', maxit
else if(keyword(i,1).eq.'nloop') then
if( .not. inside_structure ) then
read(keyword(i,2),*,iostat=ioerr) nloop
if ( ioerr /= 0 ) exit
end if
else if(keyword(i,1).eq.'nloop0') then
if( .not. inside_structure ) then
read(keyword(i,2),*,iostat=ioerr) nloop0
if ( ioerr /= 0 ) exit
end if
else if(keyword(i,1).eq.'discale') then
read(keyword(i,2),*,iostat=ioerr) discale
if ( ioerr /= 0 ) exit
write(*,*) ' Optional initial tolerance scale: ', discale
else if(keyword(i,1).eq.'sidemax') then
read(keyword(i,2),*,iostat=ioerr) sidemax
if ( ioerr /= 0 ) exit
write(*,*) ' User set maximum system dimensions: ', sidemax
else if(keyword(i,1).eq.'fbins') then
read(keyword(i,2),*,iostat=ioerr) fbins
if ( ioerr /= 0 ) exit
write(*,*) ' User set linked-cell bin parameter: ', fbins
else if(keyword(i,1).eq.'add_amber_ter') then
add_amber_ter = .true.
write(*,*) ' Will add the TER flag between molecules. '
else if(keyword(i,1).eq.'avoid_overlap') then
if ( keyword(i,2).eq.'yes') then
avoidoverlap = .true.
write(*,*) ' Will avoid overlap to fixed molecules at initial point. '
else
avoidoverlap = .false.
write(*,*) ' Will NOT avoid overlap to fixed molecules at initial point. '
end if
else if(keyword(i,1).eq.'packall') then
packall = .true.
write(*,*) ' Will pack all molecule types from the beginning. '
else if(keyword(i,1).eq.'use_short_tol') then
use_short_tol = .true.
write(*,*) ' Will use a short distance penalty for all atoms. '
else if(keyword(i,1).eq.'writecrd') then
crd = .true.
write(*,*) ' Will write output also in CRD format '
read(keyword(i,2),*,iostat=ioerr) crdfile
else if(keyword(i,1).eq.'add_box_sides') then
add_box_sides = .true.
write(*,*) ' Will print BOX SIDE informations. '
read(keyword(i,2),*,iostat=ioerr) add_sides_fix
if ( ioerr /= 0 ) then
ioerr = 0
cycle
end if
write(*,*) ' Will sum ', add_sides_fix,' to each side length on print'
else if(keyword(i,1).eq.'iprint1') then
read(keyword(i,2),*,iostat=ioerr) iprint1
if ( ioerr /= 0 ) exit
write(*,*) ' Optional printvalue 1 set: ', iprint1
else if(keyword(i,1).eq.'iprint2') then
read(keyword(i,2),*,iostat=ioerr) iprint2
if ( ioerr /= 0 ) exit
write(*,*) ' Optional printvalue 2 set: ', iprint2
else if( keyword(i,1) /= 'tolerance' .and. &
keyword(i,1) /= 'short_tol_dist' .and. &
keyword(i,1) /= 'short_tol_scale' .and. &
keyword(i,1) /= 'structure' .and. &
keyword(i,1) /= 'end' .and. &
keyword(i,1) /= 'atoms' .and. &
keyword(i,1) /= 'output' .and. &
keyword(i,1) /= 'filetype' .and. &
keyword(i,1) /= 'number' .and. &
keyword(i,1) /= 'inside' .and. &
keyword(i,1) /= 'outside' .and. &
keyword(i,1) /= 'fixed' .and. &
keyword(i,1) /= 'center' .and. &
keyword(i,1) /= 'centerofmass' .and. &
keyword(i,1) /= 'over' .and. &
keyword(i,1) /= 'above' .and. &
keyword(i,1) /= 'below' .and. &
keyword(i,1) /= 'constrain_rotation' .and. &
keyword(i,1) /= 'radius' .and. &
keyword(i,1) /= 'fscale' .and. &
keyword(i,1) /= 'short_radius' .and. &
keyword(i,1) /= 'short_radius_scale' .and. &
keyword(i,1) /= 'resnumbers' .and. &
keyword(i,1) /= 'connect' .and. &
keyword(i,1) /= 'changechains' .and. &
keyword(i,1) /= 'chain' .and. &
keyword(i,1) /= 'discale' .and. &
keyword(i,1) /= 'maxit' .and. &
keyword(i,1) /= 'movebadrandom' .and. &
keyword(i,1) /= 'maxmove' .and. &
keyword(i,1) /= 'add_amber_ter' .and. &
keyword(i,1) /= 'sidemax' .and. &
keyword(i,1) /= 'seed' .and. &
keyword(i,1) /= 'randominitialpoint' .and. &
keyword(i,1) /= 'restart_from' .and. &
keyword(i,1) /= 'restart_to' .and. &
keyword(i,1) /= 'nloop' .and. &
keyword(i,1) /= 'nloop0' .and. &
keyword(i,1) /= 'writeout' .and. &
keyword(i,1) /= 'writebad' .and. &
keyword(i,1) /= 'check' .and. &
keyword(i,1) /= 'iprint1' .and. &
keyword(i,1) /= 'iprint2' .and. &
keyword(i,1) /= 'writecrd' .and. &
keyword(i,1) /= 'segid' .and. &
keyword(i,1) /= 'chkgrad' ) then
write(*,*) ' ERROR: Keyword not recognized: ', trim(keyword(i,1))
stop
end if
end do
if ( ioerr /= 0 ) then
write(*,*) ' ERROR: Some optional keyword was not used correctly: ', trim(keyword(i,1))
stop
end if
write(*,*) ' Seed for random number generator: ', seed
call init_random_number(seed)
! Checking for the name of the output file to be created
xyzout = '####'
do iline = 1, nlines
if(keyword(iline,1).eq.'output') then
xyzout = keyword(iline,2)
xyzout = trim(adjustl(xyzout))
end if
end do
if(xyzout(1:4) == '####') then
write(*,*)' ERROR: Output file not (correctly?) specified. '
stop
end if
write(*,*)' Output file: ', trim(adjustl(xyzout))
! Reading structure files
itype = 0
do iline = 1, nlines
if(keyword(iline,1).eq.'structure') then
itype = itype + 1
record = keyword(iline,2)
write(*,*) ' Reading coordinate file: ', trim(adjustl(record))
! Reading pdb input files
if(pdb) then
name(itype) = trim(adjustl(record))
record = keyword(iline,2)
pdbfile(itype) = trim(record)
idfirst(itype) = 1
idfirstatom = 0
do ii = itype - 1, 1, -1
idfirst(itype) = idfirst(itype) + natoms(ii)
end do
open(10,file=keyword(iline,2),status='old',iostat=ioerr)
if ( ioerr /= 0 ) call failopen(keyword(iline,2))
! Read coordinates
record(1:6) = '######'
do while(record(1:4).ne.'ATOM'.and.record(1:6).ne.'HETATM')
read(10,str_format) record
end do
idatom = idfirst(itype) - 1
do while(idatom.lt.natoms(itype)+idfirst(itype)-1)
if(record(1:4).eq.'ATOM'.or.record(1:6).eq.'HETATM') then
idatom = idatom + 1
amass(idatom) = 1.d0
maxcon(idatom) = 0
! Read the index of the first atom, to adjust connectivities, if any
if(idfirstatom == 0) read(record(7:11),*,iostat=ioerr) idfirstatom
read(record,"( t31,f8.3,t39,f8.3,t47,f8.3 )",iostat=ioerr) &
(coor(idatom,k),k=1,3)
if( ioerr /= 0 ) then
record = keyword(iline,2)
write(*,*) ' ERROR: Failed to read coordinates from', &
' file: ', trim(adjustl(record))
write(*,*) ' Probably the coordinates are not in', &
' standard PDB file format. '
write(*,*) ' Standard PDB format specifications', &
' can be found at: '
write(*,*) ' www.rcsb.org/pdb '
stop
end if
! This only tests if residue numbers can be read, they are used
! only for output
read(record(23:26),*,iostat=ioerr) itest
if( ioerr /= 0 ) then
record = pdbfile(itype)
write(*,*) ' ERROR: Failed reading residue number',&
' from PDB file: ', trim(adjustl(record))
write(*,*) ' Residue numbers are integers that',&
' must be within columns 23 and 26. '
write(*,*) ' Other characters within these columns',&
' will cause input/output errors. '
write(*,*) ' Standard PDB format specifications',&
' can be found at: '
write(*,*) ' www.rcsb.org/pdb '
stop
end if
end if
read(10,str_format,iostat=ioerr) record
end do
!
! Read connectivity, if there is any specified
!
do while(.true.)
if ( ioerr /= 0 ) exit
if(record(1:6).eq.'CONECT') then
iread = 7
read(record(iread:iread+4),*,iostat=ioerr) iatom
iatom = iatom - idfirstatom + 1
idatom = idfirst(itype) - 1 + iatom
if(ioerr /= 0) then
write(*,*) " ERROR: Could not read atom index from CONECT line: "
write(*,*) trim(adjustl(record))
stop
end if
iread = iread + 5
read(record(iread:iread+4),*,iostat=ioerr) nconnect(idatom,1)
if(ioerr /= 0) then
write(*,*) " ERROR: Could not read any connection index from CONECT line: "
write(*,*) trim(adjustl(record))
stop
end if
nconnect(idatom,1) = nconnect(idatom,1) - idfirstatom + 1
maxcon(idatom) = 1
do while(.true.)
iread = iread + 5
read(record(iread:iread+4),*,iostat=ioerr) nconnect(idatom,maxcon(idatom)+1)
if(ioerr == 0) then
maxcon(idatom) = maxcon(idatom) + 1
nconnect(idatom,maxcon(idatom)) = nconnect(idatom,maxcon(idatom)) - idfirstatom + 1
else
exit
end if
end do
end if
read(10,str_format,iostat=ioerr) record
end do
close(10)
end if
! Reading tinker input files
if(tinker) then
open(10,file=keyword(iline,2),status='old',iostat=ioerr)
if ( ioerr /= 0 ) call failopen(keyword(iline,2))
idfirst(itype) = 1
do ii = itype - 1, 1, -1
idfirst(itype) = idfirst(itype) + natoms(ii)
end do
record = keyword(iline,2)
call setcon(record(1:64),idfirst(itype))
open(10,file = keyword(iline,2), status = 'old')
record = blank
do while(record.le.blank)
read(10,str_format) record
end do
i = 1
do while(record(i:i).le.' ')
i = i + 1
if ( i > strl ) exit
end do
iarg = i
if ( i < strl ) then
do while(record(i:i).gt.' ')
i = i + 1
if ( i > strl ) exit
end do
end if
read(record(iarg:i-1),*) natoms(itype)
if ( i < strl ) then
do while(record(i:i).le.' ')
i = i + 1
if ( i > strl ) exit
end do
end if
iarg = i
if ( i < strl ) then
do while(record(i:i).gt.' ')
i = i + 1
if ( i > strl ) exit
end do
end if
read(record(iarg:i-1),str_format) name(itype)
record = name(itype)
name(itype) = trim(adjustl(record))
if(name(itype).lt.' ') name(itype) = 'Without_title'
idatom = idfirst(itype) - 1
do iatom = 1, natoms(itype)
idatom = idatom + 1
record = blank
do while(record.le.blank)
read(10,str_format) record
end do
i = 1
do while(record(i:i).le.' ')
i = i + 1
if ( i > strl ) exit
end do
iarg = i
if ( i < strl ) then
do while(record(i:i).gt.' ')
i = i + 1
if ( i > strl ) exit
end do
end if
read(record(iarg:i-1),*) in
if ( i < strl ) then
do while(record(i:i).le.' ')
i = i + 1
if ( i > strl ) exit
end do
end if
iarg = i
if ( i < strl ) then
do while(record(i:i).gt.' ')
i = i + 1
if ( i > strl ) exit
end do
end if
read(record(iarg:i-1),*) ele(idatom)
read(record(i:strl),*) (coor(idatom,k), k = 1, 3),&
(nconnect(idatom, k), k = 1, maxcon(idatom))
amass(idatom) = 1.d0
end do
close(10)
end if
! Reading xyz input files
if(xyz) then
open(10,file=keyword(iline,2),status='old',iostat=ioerr)
if ( ioerr /= 0 ) call failopen(keyword(iline,2))
read(10,*) natoms(itype)
read(10,str_format) name(itype)
if(name(itype).lt.' ') name(itype) = 'Without_title'
idfirst(itype) = 1
do ii = itype - 1, 1, -1
idfirst(itype) = idfirst(itype) + natoms(ii)
end do
idatom = idfirst(itype) - 1
do iatom = 1, natoms(itype)
idatom = idatom + 1
record = blank
read(10,str_format) record
read(record,*) ele(idatom), (coor(idatom,k),k=1,3)
amass(idatom) = 1.d0
end do
close(10)
end if
! Reading moldy input files
if(moldy) then
open(10,file=keyword(iline,2), status ='old',iostat=ioerr)
if ( ioerr /= 0 ) call failopen(keyword(iline,2))
read(10,*) name(itype), nmols(itype)
natoms(itype) = 0
do while(.true.)
read(10,str_format,iostat=ioerr) record
if ( ioerr /= 0 ) exit
if(record.gt.' '.and.record(1:3).ne.'end') &
natoms(itype) = natoms(itype) + 1
end do
close(10)
idfirst(itype) = 1
do ii = itype - 1, 1, -1
idfirst(itype) = idfirst(itype) + natoms(ii)
end do
open(10,file=keyword(iline,2),status='old')
read(10,str_format) record
idatom = idfirst(itype) - 1
do iatom = 1, natoms(itype)
idatom = idatom + 1
read(10,str_format) record
read(record,*) lixo, (coor(idatom,k), k = 1, 3),&
amass(idatom), charge(idatom), ele(idatom)
end do
close(10)
end if
end if
end do
ntype = itype
write(*,*) ' Number of independent structures: ', ntype
write(*,*) ' The structures are: '
do itype = 1, ntype
record = name(itype)
write(*,*) ' Structure ', itype, ':', trim(adjustl(record)),&
'(',natoms(itype),' atoms)'
end do
! Setting the vectors for the number of GENCAN loops
if(nloop.eq.0) then
nloop_all = 200*ntype
nloop = nloop_all
else
nloop_all = nloop
end if
write(*,*) ' Maximum number of GENCAN loops for all molecule packing: ', nloop_all
do itype = 1, ntype
if ( nloop_type(itype) == 0 ) then
nloop_type(itype) = nloop_all
else
write(*,*) ' Maximum number of GENCAN loops for type: ', itype, ': ', nloop_type(itype)
end if
end do
! nloop0 are the number of loops for the initial phase packing
if(nloop0.eq.0) then
nloop0 = 20*ntype
else
write(*,*) ' Maximum number of GENCAN loops-0 for all molecule packing: ', nloop0
end if
do itype = 1, ntype
if ( nloop0_type(itype) == 0 ) then
nloop0_type(itype) = nloop0
else
write(*,*) ' Maximum number of GENCAN loops-0 for type: ', itype, ': ', nloop0_type(itype)
end if
end do
! Reading the restrictions that were set
irest = 0
ioerr = 0
do iline = 1, nlines
if(keyword(iline,1).eq.'fixed') then
irest = irest + 1
irestline(irest) = iline
ityperest(irest) = 1
read(keyword(iline,2),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,4)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,5)
read(keyword(iline,7),*,iostat=ioerr) restpars(irest,6)
end if
if(keyword(iline,1).eq.'inside') then
irest = irest + 1
irestline(irest) = iline
if(keyword(iline,2).eq.'cube') then
ityperest(irest) = 2
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
else if(keyword(iline,2).eq.'box') then
ityperest(irest) = 3
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
read(keyword(iline,7),*,iostat=ioerr) restpars(irest,5)
read(keyword(iline,8),*,iostat=ioerr) restpars(irest,6)
else if(keyword(iline,2).eq.'sphere') then
ityperest(irest) = 4
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
else if(keyword(iline,2).eq.'ellipsoid') then
ityperest(irest) = 5
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
read(keyword(iline,7),*,iostat=ioerr) restpars(irest,5)
read(keyword(iline,8),*,iostat=ioerr) restpars(irest,6)
read(keyword(iline,9),*,iostat=ioerr) restpars(irest,7)
else if(keyword(iline,2).eq.'cylinder') then
ityperest(irest) = 12
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
read(keyword(iline,7),*,iostat=ioerr) restpars(irest,5)
read(keyword(iline,8),*,iostat=ioerr) restpars(irest,6)
read(keyword(iline,9),*,iostat=ioerr) restpars(irest,7)
read(keyword(iline,10),*,iostat=ioerr) restpars(irest,9)
restpars(irest,8) = restpars(irest,4)**2 + &
restpars(irest,5)**2 + &
restpars(irest,6)**2
if(restpars(irest,8).lt.1.d-10) then
write(*,*) ' ERROR: The norm of the director vector', &
' of the cylinder constraint cannot be zero.'
ioerr = 1
else
clen = dsqrt(restpars(irest,8))
restpars(irest,4) = restpars(irest,4) / clen
restpars(irest,5) = restpars(irest,5) / clen
restpars(irest,6) = restpars(irest,6) / clen
end if
else
ioerr = 1
end if
end if
if(keyword(iline,1).eq.'outside') then
irest = irest + 1
irestline(irest) = iline
if(keyword(iline,2).eq.'cube') then
ityperest(irest) = 6
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
else if(keyword(iline,2).eq.'box') then
ityperest(irest) = 7
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
read(keyword(iline,7),*,iostat=ioerr) restpars(irest,5)
read(keyword(iline,8),*,iostat=ioerr) restpars(irest,6)
else if(keyword(iline,2).eq.'sphere') then
ityperest(irest) = 8
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
else if(keyword(iline,2).eq.'ellipsoid') then
ityperest(irest) = 9
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
read(keyword(iline,7),*,iostat=ioerr) restpars(irest,5)
read(keyword(iline,8),*,iostat=ioerr) restpars(irest,6)
read(keyword(iline,9),*,iostat=ioerr) restpars(irest,7)
else if(keyword(iline,2).eq.'cylinder') then
ityperest(irest) = 13
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
read(keyword(iline,7),*,iostat=ioerr) restpars(irest,5)
read(keyword(iline,8),*,iostat=ioerr) restpars(irest,6)
read(keyword(iline,9),*,iostat=ioerr) restpars(irest,7)
read(keyword(iline,10),*,iostat=ioerr) restpars(irest,9)
restpars(irest,8) = restpars(irest,4)**2 + &
restpars(irest,5)**2 + &
restpars(irest,6)**2
if(restpars(irest,8).lt.1.d-10) then
write(*,*) ' ERROR: The norm of the director vector',&
' of the cylinder constraint cannot be zero.'
ioerr = 1
else
clen = dsqrt(restpars(irest,8))
restpars(irest,4) = restpars(irest,4) / clen
restpars(irest,5) = restpars(irest,5) / clen
restpars(irest,6) = restpars(irest,6) / clen
end if
else
ioerr = 1
end if
end if
if(keyword(iline,1).eq.'over' .or. keyword(iline,1).eq.'above') then
irest = irest + 1
irestline(irest) = iline
ityperest(irest) = 10
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
if(keyword(iline,2).ne.'plane') ioerr = 1
end if
if(keyword(iline,1).eq.'below') then
irest = irest + 1
irestline(irest) = iline
ityperest(irest) = 11
read(keyword(iline,3),*,iostat=ioerr) restpars(irest,1)
read(keyword(iline,4),*,iostat=ioerr) restpars(irest,2)
read(keyword(iline,5),*,iostat=ioerr) restpars(irest,3)
read(keyword(iline,6),*,iostat=ioerr) restpars(irest,4)
if(keyword(iline,2).ne.'plane') ioerr = 1
end if
if ( ioerr /= 0 ) then
write(*,*) ' ERROR: Some restriction is not set correctly. '
stop
end if
end do
nrest = irest
write(*,*) ' Total number of restrictions: ', nrest
! Getting the tolerance
ioerr = 1
dism = -1.d0
do iline = 1, nlines
if(keyword(iline,1).eq.'tolerance') then
read(keyword(iline,2),*,iostat=ioerr) dism
if ( ioerr /= 0 ) then
write(*,*) ' ERROR: Failed reading tolerance. '
stop
end if
exit
end if
end do
if ( ioerr /= 0 ) then
write(*,*) ' ERROR: Overall tolerance not set. Use, for example: tolerance 2.0 '
stop
end if
write(*,*) ' Distance tolerance: ', dism
! Reading, if defined, the short distance penalty parameters
ioerr = 1
short_tol_dist = dism/2.d0
! Reading short_tol_dist
do iline = 1, nlines
if(keyword(iline,1).eq.'short_tol_dist') then
read(keyword(iline,2),*,iostat=ioerr) short_tol_dist
if ( ioerr /= 0 ) then
write(*,*) ' ERROR: Failed reading short_tol_dist. '
stop
end if
if ( short_tol_dist > dism ) then
write(*,*) ' ERROR: The short_tol_dist parameter must be smaller than the tolerance. '
stop
end if
write(*,*) ' User defined short tolerance distance: ', short_tol_dist
short_tol_dist = short_tol_dist**2
exit
end if
end do
! Reading short_tol_scale
short_tol_scale = 3.d0
do iline = 1, nlines
if(keyword(iline,1).eq.'short_tol_scale') then
read(keyword(iline,2),*,iostat=ioerr) short_tol_scale
if ( ioerr /= 0 ) then
write(*,*) ' ERROR: Failed reading short_tol_scale. '
stop
end if
if ( short_tol_dist <= 0.d0 ) then
write(*,*) ' ERROR: The short_tol_scale parameter must be positive. '
stop
end if
write(*,*) ' User defined short tolerance scale: ', short_tol_scale
exit
end if
end do
! Assigning the input lines that correspond to each structure
itype = 0
iline = 0
do while(iline < nlines)
iline = iline + 1
if(keyword(iline,1).eq.'structure') then
itype = itype + 1
linestrut(itype,1) = iline
iline = iline + 1
do while(keyword(iline,1).ne.'end'.or.&
keyword(iline,2).ne.'structure')
if(keyword(iline,1) == 'structure'.or.&
iline == nlines) then
write(*,*) ' ERROR: Structure specification not ending with "end structure"'
stop
end if
iline = iline + 1
end do
linestrut(itype,2) = iline
end if
end do
! If pdb files, get the type of residue numbering output for each
! molecule
if(pdb) then
do itype = 1, ntype
connect(itype) = .true.
resnumbers(itype) = -1
changechains(itype) = .false.
chain(itype) = "#"
segid(itype) = ""
maxmove(itype) = nmols(itype)
do iline = 1, nlines
if(iline.gt.linestrut(itype,1).and.&
iline.lt.linestrut(itype,2)) then
if(keyword(iline,1).eq.'changechains') then
changechains(itype) = .true.
end if
if(keyword(iline,1).eq.'maxmove') then
read(keyword(iline,2),*) maxmove(itype)
end if
if(keyword(iline,1).eq.'resnumbers') then
read(keyword(iline,2),*) resnumbers(itype)
end if
if(keyword(iline,1).eq.'connect') then
if(keyword(iline,2) == "no") then
connect(itype) = .false.
end if
end if
if(keyword(iline,1).eq.'chain') then
read(keyword(iline,2),*) chain(itype)
end if
if(keyword(iline,1).eq.'segid') then
read(keyword(iline,2),*) segid(itype)
end if
end if
end do
if (crd) then
if (itype.gt.1 .and. segid(itype)=="") then
if (segid(itype-1) /= "") then
write(*,*) ' Warning: Type of segid not defined for ', itype,'. Keeping it same as previous'
endif
segid(itype) = segid(itype-1)
endif
endif
if ( resnumbers(itype) == -1 ) then
write(*,*) ' Warning: Type of residue numbering not',&
' set for structure ',itype
call setrnum(pdbfile(itype),imark)
if(imark.eq.1) resnumbers(itype) = 0
if(imark.gt.1) resnumbers(itype) = 1
end if
write(*,*) ' Residue numbering set for structure ',itype,':',&
resnumbers(itype)
write(*,*) ' Swap chains of molecules of structure ',&
itype,':', changechains(itype)
if ( chain(itype) /= "#" ) then
write(*,*) ' Specific chain identifier set for structure ',itype,':',chain(itype)
end if
if ( chain(itype) /= "#" .and. changechains(itype) ) then
write(*,*) " ERROR: 'changechains' and 'chain' input parameters are not compatible "
write(*,*) " for a single structure. "
stop
end if
end do
end if
! Write the number of molecules of each type
do itype = 1, ntype
write(*,*) ' Number of molecules of type ', itype, ': ', nmols(itype)
if(pdb.and.nmols(itype).gt.9999) then
write(*,*) ' Warning: There will be more than 9999 molecules of type ',itype
if (.not. crd) write(*,*) ' Residue numbering is reset after 9999. '
if (crd) write(*,*) ' Residue numbering is reset after 9999 in pdb but not in crd. '
if ( chain(itype) == "#" ) then
write(*,*) ' Each set be will be assigned a different chain in the PDB output file. '
end if
end if
if(crd.and.nmols(itype).gt.99999999) then
write(*,*) ' Warning: There will be more than 99999999 molecules of type ',itype
write(*,*) ' Residue numbering is reset after 99999999 in crd. '
endif
end do
! Checking if restart files will be used for each structure or for the whole system
restart_from(0) = "none"
restart_to(0) = "none"
do itype = 1, ntype
restart_from(itype) = "none"
restart_to(itype) = "none"
end do
lines: do iline = 1, nlines
if ( keyword(iline,1) == 'restart_from' ) then
do itype = 1, ntype
if(iline.gt.linestrut(itype,1).and.&
iline.lt.linestrut(itype,2)) then
restart_from(itype) = keyword(iline,2)
cycle lines
end if
end do
if( restart_from(0) == 'none' ) then
restart_from(0) = keyword(iline,2)
else
write(*,*) ' ERROR: More than one definition of restart_from file for all system. '
stop
end if
end if
if ( keyword(iline,1) == 'restart_to' ) then
do itype = 1, ntype
if(iline.gt.linestrut(itype,1).and.&
iline.lt.linestrut(itype,2)) then
restart_to(itype) = keyword(iline,2)
cycle lines
end if
end do
if( restart_to(0) == 'none' ) then
restart_to(0) = keyword(iline,2)
else
write(*,*) ' ERROR: More than one definition of restart_to file for all system. '
stop
end if
end if
end do lines
return
end subroutine getinp
!
! Subroutine that stops if failed to open file
!
subroutine failopen(record)
use sizes
character(len=strl) :: record
write(*,*)
write(*,*) ' ERROR: Could not open file. '
write(*,*) ' Could not find file: ',trim(record)
write(*,*) ' Please check if all the input and structure '
write(*,*) ' files are in the current directory or if the'
write(*,*) ' correct paths are provided.'
write(*,*)
stop
end subroutine failopen
!
! Subroutine that checks if a pdb structure has one or more than
! one residue
!
subroutine setrnum(file,nres)
use sizes
implicit none
integer :: iread, ires, ireslast, nres, ioerr
character(len=strl) :: file
character(len=strl) :: record
open(10,file=file,status='old')
iread = 0
nres = 1
do while(nres.eq.1)
read(10,str_format,iostat=ioerr) record
if ( ioerr /= 0 ) exit
if(record(1:4).eq.'ATOM'.or.record(1:6).eq.'HETATM') then
read(record(23:26),*,iostat=ioerr) ires
if ( ioerr /= 0 ) cycle
iread = iread + 1
if(iread.gt.1) then
if(ires.ne.ireslast) then
nres = 2
close(10)
return
end if
end if
ireslast = ires
end if
end do
close(10)
return
end subroutine setrnum
!
! Subroutine that computes de number of connections of each atom
! for tinker xyz files
!
subroutine setcon(xyzfile,idfirst)
use sizes
use input, only : maxcon
implicit none
integer :: idfirst
integer :: natoms, idatom, iatom, ic, i
character(len=64) :: xyzfile
character(len=120) :: record
open(10, file = xyzfile, status='old')
read(10,*) natoms
idatom = idfirst - 1
do iatom = 1, natoms
idatom = idatom + 1
read(10,"( a120 )") record
ic = 0
do i = 1, 119
if(record(i:i).gt.' '.and.record(i+1:i+1).le.' ') ic = ic + 1
end do
maxcon(idatom) = ic - 5
end do
close(10)
return
end subroutine setcon
!
! Subroutine getkeywords: gets keywords and values from the input
! file in a more robust way
!
subroutine getkeywords()
use sizes
use input, only : keyword, nlines, inputfile, forbidden_char
implicit none
character(len=strl) :: record
integer :: iline, i, j, k, ilast, ival, ioerr