-
Notifications
You must be signed in to change notification settings - Fork 1
/
topology_generate_util.F
2057 lines (1947 loc) · 92.4 KB
/
topology_generate_util.F
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
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
! **************************************************************************************************
!> \brief Collection of subroutine needed for topology related things
!> \par History
!> Teodor Laino 09.2006 - Major rewriting with linear scaling routines
! **************************************************************************************************
MODULE topology_generate_util
USE atomic_kind_types, ONLY: atomic_kind_type,&
deallocate_atomic_kind_set,&
set_atomic_kind
USE cell_types, ONLY: pbc
USE cp_log_handling, ONLY: cp_get_default_logger,&
cp_logger_get_default_io_unit,&
cp_logger_type,&
cp_to_string
USE cp_output_handling, ONLY: cp_print_key_finished_output,&
cp_print_key_unit_nr
USE cp_units, ONLY: cp_unit_to_cp2k
USE fist_neighbor_list_types, ONLY: fist_neighbor_deallocate,&
fist_neighbor_type
USE fist_neighbor_lists, ONLY: build_fist_neighbor_lists
USE input_constants, ONLY: do_add,&
do_bondparm_covalent,&
do_bondparm_vdw,&
do_conn_off,&
do_conn_user,&
do_remove
USE input_section_types, ONLY: section_vals_get,&
section_vals_get_subs_vals,&
section_vals_type,&
section_vals_val_get
USE kinds, ONLY: default_string_length,&
dp
USE memory_utilities, ONLY: reallocate
USE message_passing, ONLY: mp_para_env_type
USE particle_types, ONLY: allocate_particle_set,&
deallocate_particle_set,&
particle_type
USE periodic_table, ONLY: get_ptable_info
USE qmmm_types_low, ONLY: qmmm_env_mm_type
USE string_table, ONLY: id2str,&
s2s,&
str2id
USE string_utilities, ONLY: integer_to_string,&
uppercase
USE topology_types, ONLY: atom_info_type,&
connectivity_info_type,&
topology_parameters_type
USE topology_util, ONLY: array1_list_type,&
array2_list_type,&
find_molecule,&
give_back_molecule,&
reorder_list_array,&
reorder_structure
USE util, ONLY: find_boundary,&
sort
#include "./base/base_uses.f90"
IMPLICIT NONE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'topology_generate_util'
PRIVATE
LOGICAL, PARAMETER :: debug_this_module = .FALSE.
PUBLIC :: topology_generate_bend, &
topology_generate_bond, &
topology_generate_dihe, &
topology_generate_impr, &
topology_generate_onfo, &
topology_generate_ub, &
topology_generate_molecule, &
topology_generate_molname
CONTAINS
! **************************************************************************************************
!> \brief Generates molnames: useful when the connectivity on file does not
!> provide them
!> \param conn_info ...
!> \param natom ...
!> \param natom_prev ...
!> \param nbond_prev ...
!> \param id_molname ...
!> \author Teodoro Laino [tlaino] - University of Zurich 10.2008
! **************************************************************************************************
SUBROUTINE topology_generate_molname(conn_info, natom, natom_prev, nbond_prev, &
id_molname)
TYPE(connectivity_info_type), POINTER :: conn_info
INTEGER, INTENT(IN) :: natom, natom_prev, nbond_prev
INTEGER, DIMENSION(:), INTENT(INOUT) :: id_molname
CHARACTER(LEN=default_string_length), PARAMETER :: basename = "MOL"
CHARACTER(LEN=default_string_length) :: molname
INTEGER :: i, N, nmol
LOGICAL :: check
TYPE(array1_list_type), ALLOCATABLE, DIMENSION(:) :: atom_bond_list
! convert a simple list of bonds to a list of bonds per atom
! (each bond is present in the forward and backward direction
ALLOCATE (atom_bond_list(natom))
DO I = 1, natom
ALLOCATE (atom_bond_list(I)%array1(0))
END DO
N = 0
IF (ASSOCIATED(conn_info%bond_a)) N = SIZE(conn_info%bond_a) - nbond_prev
CALL reorder_structure(atom_bond_list, conn_info%bond_a(nbond_prev + 1:) - natom_prev, &
conn_info%bond_b(nbond_prev + 1:) - natom_prev, N)
nmol = 0
check = ALL(id_molname == str2id(s2s("__UNDEF__"))) .OR. ALL(id_molname /= str2id(s2s("__UNDEF__")))
CPASSERT(check)
DO i = 1, natom
IF (id2str(id_molname(i)) == "__UNDEF__") THEN
molname = TRIM(basename)//ADJUSTL(cp_to_string(nmol))
CALL generate_molname_low(i, atom_bond_list, molname, id_molname)
nmol = nmol + 1
END IF
END DO
DO I = 1, natom
DEALLOCATE (atom_bond_list(I)%array1)
END DO
DEALLOCATE (atom_bond_list)
END SUBROUTINE topology_generate_molname
! **************************************************************************************************
!> \brief Generates molnames: useful when the connectivity on file does not
!> provide them
!> \param i ...
!> \param atom_bond_list ...
!> \param molname ...
!> \param id_molname ...
!> \author Teodoro Laino [tlaino] - University of Zurich 10.2008
! **************************************************************************************************
RECURSIVE SUBROUTINE generate_molname_low(i, atom_bond_list, molname, id_molname)
INTEGER, INTENT(IN) :: i
TYPE(array1_list_type), DIMENSION(:) :: atom_bond_list
CHARACTER(LEN=default_string_length), INTENT(IN) :: molname
INTEGER, DIMENSION(:), INTENT(INOUT) :: id_molname
INTEGER :: j, k
IF (debug_this_module) THEN
WRITE (*, *) "Entered with :", i
WRITE (*, *) TRIM(molname)//": entering with i:", i, " full series to test:: ", atom_bond_list(i)%array1
IF ((TRIM(id2str(id_molname(i))) /= "__UNDEF__") .AND. &
(TRIM(id2str(id_molname(i))) /= TRIM(molname))) THEN
WRITE (*, *) "Atom (", i, ") has already a molecular name assigned ! ("//TRIM(id2str(id_molname(i)))//")."
WRITE (*, *) "New molecular name would be: ("//TRIM(molname)//")."
WRITE (*, *) "Detecting something wrong in the molecular setup!"
CPABORT("")
END IF
END IF
id_molname(i) = str2id(molname)
DO j = 1, SIZE(atom_bond_list(i)%array1)
k = atom_bond_list(i)%array1(j)
IF (debug_this_module) WRITE (*, *) "entering with i:", i, "testing :", k
IF (k == -1) CYCLE
atom_bond_list(i)%array1(j) = -1
WHERE (atom_bond_list(k)%array1 == i) atom_bond_list(k)%array1 = -1
CALL generate_molname_low(k, atom_bond_list, molname, id_molname)
END DO
END SUBROUTINE generate_molname_low
! **************************************************************************************************
!> \brief Use information from bond list to generate molecule. (ie clustering)
!> \param topology ...
!> \param qmmm ...
!> \param qmmm_env ...
!> \param subsys_section ...
! **************************************************************************************************
SUBROUTINE topology_generate_molecule(topology, qmmm, qmmm_env, subsys_section)
TYPE(topology_parameters_type), INTENT(INOUT) :: topology
LOGICAL, INTENT(in), OPTIONAL :: qmmm
TYPE(qmmm_env_mm_type), OPTIONAL, POINTER :: qmmm_env
TYPE(section_vals_type), POINTER :: subsys_section
CHARACTER(len=*), PARAMETER :: routineN = 'topology_generate_molecule'
INTEGER, PARAMETER :: nblock = 100
INTEGER :: atom_in_kind, atom_in_mol, first, handle, handle2, i, iatm, iatom, iend, ifirst, &
ilast, inum, istart, itype, iw, j, jump1, jump2, last, max_mol_num, mol_num, mol_res, &
mol_typ, myind, N, natom, nlocl, ntype, resid
INTEGER, DIMENSION(:), POINTER :: qm_atom_index, wrk1, wrk2
LOGICAL :: do_again, found, my_qmmm
TYPE(array1_list_type), ALLOCATABLE, DIMENSION(:) :: atom_bond_list
TYPE(atom_info_type), POINTER :: atom_info
TYPE(connectivity_info_type), POINTER :: conn_info
TYPE(cp_logger_type), POINTER :: logger
NULLIFY (logger)
logger => cp_get_default_logger()
iw = cp_print_key_unit_nr(logger, subsys_section, "PRINT%TOPOLOGY_INFO/UTIL_INFO", &
extension=".subsysLog")
CALL timeset(routineN, handle)
NULLIFY (qm_atom_index)
NULLIFY (wrk1)
NULLIFY (wrk2)
atom_info => topology%atom_info
conn_info => topology%conn_info
!
! QM/MM coordinate_control
!
my_qmmm = .FALSE.
IF (PRESENT(qmmm) .AND. PRESENT(qmmm_env)) my_qmmm = qmmm
natom = topology%natoms
IF (ASSOCIATED(atom_info%map_mol_typ)) DEALLOCATE (atom_info%map_mol_typ)
ALLOCATE (atom_info%map_mol_typ(natom))
IF (ASSOCIATED(atom_info%map_mol_num)) DEALLOCATE (atom_info%map_mol_num)
ALLOCATE (atom_info%map_mol_num(natom))
IF (ASSOCIATED(atom_info%map_mol_res)) DEALLOCATE (atom_info%map_mol_res)
ALLOCATE (atom_info%map_mol_res(natom))
! Initialisation
atom_info%map_mol_typ(:) = 0
atom_info%map_mol_num(:) = -1
atom_info%map_mol_res(:) = 1
! Parse the atom list to find the different molecule types and residues
ntype = 1
atom_info%map_mol_typ(1) = 1
resid = 1
CALL reallocate(wrk1, 1, nblock)
wrk1(1) = atom_info%id_molname(1)
DO iatom = 2, natom
IF (topology%conn_type == do_conn_off) THEN
! No connectivity: each atom becomes a molecule of its own molecule kind
ntype = ntype + 1
atom_info%map_mol_typ(iatom) = ntype
ELSE IF (topology%conn_type == do_conn_user) THEN
! User-defined connectivity: 5th column of COORD section or molecule
! or residue name in the case of PDB files
IF (atom_info%id_molname(iatom) == atom_info%id_molname(iatom - 1)) THEN
atom_info%map_mol_typ(iatom) = atom_info%map_mol_typ(iatom - 1)
IF (atom_info%id_resname(iatom) == atom_info%id_resname(iatom - 1)) THEN
atom_info%map_mol_res(iatom) = atom_info%map_mol_res(iatom - 1)
ELSE
resid = resid + 1
atom_info%map_mol_res(iatom) = resid
END IF
ELSE
! Check if the type is already known
found = .FALSE.
DO itype = 1, ntype
IF (atom_info%id_molname(iatom) == wrk1(itype)) THEN
atom_info%map_mol_typ(iatom) = itype
found = .TRUE.
EXIT
END IF
END DO
IF (.NOT. found) THEN
ntype = ntype + 1
atom_info%map_mol_typ(iatom) = ntype
IF (ntype > SIZE(wrk1)) CALL reallocate(wrk1, 1, 2*SIZE(wrk1))
wrk1(ntype) = atom_info%id_molname(iatom)
END IF
resid = resid + 1
atom_info%map_mol_res(iatom) = resid
END IF
ELSE
IF (atom_info%id_molname(iatom - 1) == atom_info%id_molname(iatom)) THEN
atom_info%map_mol_typ(iatom) = ntype
ELSE
ntype = ntype + 1
atom_info%map_mol_typ(iatom) = ntype
END IF
END IF
END DO
DEALLOCATE (wrk1)
IF (iw > 0) WRITE (iw, '(/,T2,A)') "Start of molecule generation"
! convert a simple list of bonds to a list of bonds per atom
! (each bond is present in the forward and backward direction
ALLOCATE (atom_bond_list(natom))
DO I = 1, natom
ALLOCATE (atom_bond_list(I)%array1(0))
END DO
N = 0
IF (ASSOCIATED(conn_info%bond_a)) N = SIZE(conn_info%bond_a)
CALL reorder_structure(atom_bond_list, conn_info%bond_a, conn_info%bond_b, N)
CALL find_molecule(atom_bond_list, atom_info%map_mol_num, atom_info%id_molname)
DO I = 1, natom
DEALLOCATE (atom_bond_list(I)%array1)
END DO
DEALLOCATE (atom_bond_list)
IF (iw > 0) WRITE (iw, '(/,T2,A)') "End of molecule generation"
! Modify according map_mol_typ the array map_mol_num
IF (iw > 0) WRITE (iw, '(/,T2,A)') "Checking for non-continuous generated molecules"
! Check molecule number
ALLOCATE (wrk1(natom))
ALLOCATE (wrk2(natom))
wrk1 = atom_info%map_mol_num
IF (debug_this_module) THEN
DO i = 1, natom
WRITE (*, '(2I10)') i, atom_info%map_mol_num(i)
END DO
END IF
CALL sort(wrk1, natom, wrk2)
istart = 1
mol_typ = wrk1(istart)
DO i = 2, natom
IF (mol_typ /= wrk1(i)) THEN
iend = i - 1
first = MINVAL(wrk2(istart:iend))
last = MAXVAL(wrk2(istart:iend))
nlocl = last - first + 1
IF (iend - istart + 1 /= nlocl) THEN
IF (debug_this_module) WRITE (*, *) iend, istart, iend - istart + 1, first, last, nlocl
CALL cp_abort(__LOCATION__, &
"CP2K requires molecules to be contiguous and we have detected a non contiguous one!! "// &
"In particular a molecule defined from index ("//cp_to_string(first)//") to ("// &
cp_to_string(last)//") contains other molecules, not connected! "// &
"Too late at this stage everything should be already ordered! "// &
"If you have not yet employed the REORDERING keyword, please do so. "// &
"It may help to fix this issue.")
END IF
istart = i
mol_typ = wrk1(istart)
END IF
END DO
iend = i - 1
first = MINVAL(wrk2(istart:iend))
last = MAXVAL(wrk2(istart:iend))
nlocl = last - first + 1
IF (iend - istart + 1 /= nlocl) THEN
IF (debug_this_module) WRITE (*, *) iend, istart, iend - istart + 1, first, last, nlocl
CALL cp_abort(__LOCATION__, &
"CP2K requires molecules to be contiguous and we have detected a non contiguous one!! "// &
"In particular a molecule defined from index ("//cp_to_string(first)//") to ("// &
cp_to_string(last)//") contains other molecules, not connected! "// &
"Too late at this stage everything should be already ordered! "// &
"If you have not yet employed the REORDERING keyword, please do so. "// &
"It may help to fix this issue.")
END IF
DEALLOCATE (wrk1)
DEALLOCATE (wrk2)
IF (iw > 0) WRITE (iw, '(/,T2,A)') "End of check"
IF (iw > 0) WRITE (UNIT=iw, FMT="(/,T2,A)") "Start of renumbering molecules"
IF (topology%conn_type == do_conn_user) THEN
mol_num = 1
atom_info%map_mol_num(1) = 1
DO iatom = 2, natom
IF (atom_info%id_molname(iatom) /= atom_info%id_molname(iatom - 1)) THEN
mol_num = 1
ELSE IF (atom_info%map_mol_res(iatom) /= atom_info%map_mol_res(iatom - 1)) THEN
mol_num = mol_num + 1
END IF
atom_info%map_mol_num(iatom) = mol_num
END DO
ELSE
mol_typ = atom_info%map_mol_typ(1)
mol_num = atom_info%map_mol_num(1)
DO i = 2, natom
IF (atom_info%map_mol_typ(i) /= mol_typ) THEN
myind = atom_info%map_mol_num(i) - mol_num + 1
CPASSERT(myind /= atom_info%map_mol_num(i - 1))
mol_typ = atom_info%map_mol_typ(i)
mol_num = atom_info%map_mol_num(i)
END IF
atom_info%map_mol_num(i) = atom_info%map_mol_num(i) - mol_num + 1
END DO
END IF
IF (iw > 0) WRITE (UNIT=iw, FMT="(/,T2,A)") "End of renumbering molecules"
! Optionally, use the residues as molecules
CALL timeset(routineN//"_PARA_RES", handle2)
IF (iw > 0) WRITE (UNIT=iw, FMT="(/,T2,A,L2)") "Starting PARA_RES: ", topology%para_res
IF (topology%para_res) THEN
IF (topology%conn_type == do_conn_user) THEN
atom_info%id_molname(:) = atom_info%id_resname(:)
ntype = 1
atom_info%map_mol_typ(1) = 1
mol_num = 1
atom_info%map_mol_num(1) = 1
DO iatom = 2, natom
IF (atom_info%id_molname(iatom) /= atom_info%id_molname(iatom - 1)) THEN
ntype = ntype + 1
mol_num = 1
ELSE IF (atom_info%map_mol_res(iatom) /= atom_info%map_mol_res(iatom - 1)) THEN
mol_num = mol_num + 1
END IF
atom_info%map_mol_typ(iatom) = ntype
atom_info%map_mol_num(iatom) = mol_num
END DO
ELSE
mol_res = 1
mol_typ = atom_info%map_mol_typ(1)
mol_num = atom_info%map_mol_num(1)
atom_info%map_mol_res(1) = mol_res
DO i = 2, natom
IF ((atom_info%resid(i - 1) /= atom_info%resid(i)) .OR. &
(atom_info%id_resname(i - 1) /= atom_info%id_resname(i))) THEN
mol_res = mol_res + 1
END IF
IF ((atom_info%map_mol_typ(i) /= mol_typ) .OR. &
(atom_info%map_mol_num(i) /= mol_num)) THEN
mol_typ = atom_info%map_mol_typ(i)
mol_num = atom_info%map_mol_num(i)
mol_res = 1
END IF
atom_info%map_mol_res(i) = mol_res
END DO
END IF
END IF
IF (iw > 0) WRITE (UNIT=iw, FMT="(/,T2,A)") "End of PARA_RES"
CALL timestop(handle2)
IF (iw > 0) THEN
DO iatom = 1, natom
WRITE (iw, '(4(1X,A,":",I0),2(1X,A,1X,A))') "iatom", iatom, &
"map_mol_typ", atom_info%map_mol_typ(iatom), &
"map_mol_num", atom_info%map_mol_num(iatom), &
"map_mol_res", atom_info%map_mol_res(iatom), &
"mol_name:", TRIM(id2str(atom_info%id_molname(iatom))), &
"res_name:", TRIM(id2str(atom_info%id_resname(iatom)))
END DO
END IF
IF (my_qmmm) THEN
do_again = .FALSE.
IF (iw > 0) WRITE (iw, *) "MAP_MOL_NUM ", atom_info%map_mol_num
IF (iw > 0) WRITE (iw, *) "MAP_MOL_TYP ", atom_info%map_mol_typ
IF (iw > 0) WRITE (iw, *) "MAP_MOL_RES ", atom_info%map_mol_res
ALLOCATE (qm_atom_index(SIZE(qmmm_env%qm_atom_index)))
qm_atom_index = qmmm_env%qm_atom_index
CPASSERT(ALL(qm_atom_index /= 0))
DO myind = 1, SIZE(qm_atom_index)
IF (qm_atom_index(myind) == 0) CYCLE
CALL find_boundary(atom_info%map_mol_typ, natom, ifirst, ilast, &
atom_info%map_mol_typ(qm_atom_index(myind)))
CALL find_boundary(atom_info%map_mol_typ, atom_info%map_mol_num, natom, ifirst, ilast, &
atom_info%map_mol_typ(qm_atom_index(myind)), atom_info%map_mol_num(qm_atom_index(myind)))
IF (iw > 0) WRITE (iw, *) "qm fragment:: ifirst, ilast", ifirst, ilast
CPASSERT(((ifirst /= 0) .OR. (ilast /= natom)))
DO iatm = ifirst, ilast
atom_info%id_molname(iatm) = str2id(s2s("_QM_"// &
TRIM(id2str(atom_info%id_molname(iatm)))))
IF (iw > 0) WRITE (iw, *) "QM Molecule name :: ", id2str(atom_info%id_molname(iatm))
WHERE (qm_atom_index == iatm) qm_atom_index = 0
END DO
DO iatm = 1, ifirst - 1
IF (ANY(qm_atom_index == iatm)) do_again = .TRUE.
END DO
DO iatm = ilast + 1, natom
IF (ANY(qm_atom_index == iatm)) do_again = .TRUE.
END DO
IF (iw > 0) WRITE (iw, *) " Another QM fragment? :: ", do_again
IF (ifirst /= 1) THEN
jump1 = atom_info%map_mol_typ(ifirst) - atom_info%map_mol_typ(ifirst - 1)
CPASSERT(jump1 <= 1 .AND. jump1 >= 0)
jump1 = ABS(jump1 - 1)
ELSE
jump1 = 0
END IF
IF (ilast /= natom) THEN
jump2 = atom_info%map_mol_typ(ilast + 1) - atom_info%map_mol_typ(ilast)
CPASSERT(jump2 <= 1 .AND. jump2 >= 0)
jump2 = ABS(jump2 - 1)
ELSE
jump2 = 0
END IF
! Changing mol_type consistently
DO iatm = ifirst, natom
atom_info%map_mol_typ(iatm) = atom_info%map_mol_typ(iatm) + jump1
END DO
DO iatm = ilast + 1, natom
atom_info%map_mol_typ(iatm) = atom_info%map_mol_typ(iatm) + jump2
END DO
IF (jump1 == 1) THEN
DO iatm = ifirst, ilast
atom_info%map_mol_num(iatm) = 1
END DO
END IF
IF (jump2 == 1) THEN
CALL find_boundary(atom_info%map_mol_typ, natom, first, last, atom_info%map_mol_typ(ilast + 1))
CALL find_boundary(atom_info%map_mol_typ, atom_info%map_mol_num, natom, ifirst, ilast, &
atom_info%map_mol_typ(ilast + 1), atom_info%map_mol_num(ilast + 1))
atom_in_mol = ilast - ifirst + 1
inum = 1
DO iatm = first, last, atom_in_mol
atom_info%map_mol_num(iatm:iatm + atom_in_mol - 1) = inum
inum = inum + 1
END DO
END IF
IF (.NOT. do_again) EXIT
END DO
DEALLOCATE (qm_atom_index)
IF (iw > 0) THEN
WRITE (iw, *) "After the QM/MM Setup:"
DO iatom = 1, natom
WRITE (iw, *) " iatom,map_mol_typ,map_mol_num ", iatom, &
atom_info%map_mol_typ(iatom), atom_info%map_mol_num(iatom)
END DO
END IF
END IF
!
! Further check : see if the number of atoms belonging to same molecule kinds
! are equal
IF (iw > 0) THEN
WRITE (iw, *) "SUMMARY:: Number of molecule kinds found:", ntype
ntype = MAXVAL(atom_info%map_mol_typ)
DO i = 1, ntype
atom_in_kind = COUNT(atom_info%map_mol_typ == i)
WRITE (iw, *) "Molecule kind:", i, " contains", atom_in_kind, " atoms"
IF (atom_in_kind <= 1) CYCLE
CALL find_boundary(atom_info%map_mol_typ, natom, first, last, i)
WRITE (iw, *) "Boundary atoms:", first, last
CPASSERT(last - first + 1 == atom_in_kind)
max_mol_num = MAXVAL(atom_info%map_mol_num(first:last))
WRITE (iw, *) "Number of molecules of kind", i, "is ::", max_mol_num
atom_in_mol = atom_in_kind/max_mol_num
WRITE (iw, *) "Number of atoms per each molecule:", atom_in_mol
WRITE (iw, *) "MAP_MOL_TYP::", atom_info%map_mol_typ(first:last)
WRITE (iw, *) "MAP_MOL_NUM::", atom_info%map_mol_num(first:last)
WRITE (iw, *) "MAP_MOL_RES::", atom_info%map_mol_res(first:last)
!
DO j = 1, max_mol_num
IF (COUNT(atom_info%map_mol_num(first:last) == j) /= atom_in_mol) THEN
WRITE (iw, *) "molecule type:", i, "molecule num:", j, " has ", &
COUNT(atom_info%map_mol_num(first:last) == j), &
" atoms instead of ", atom_in_mol, " ."
CALL cp_abort(__LOCATION__, &
"Two molecules of the same kind have "// &
"been created with different numbers of atoms!")
END IF
END DO
END DO
END IF
CALL cp_print_key_finished_output(iw, logger, subsys_section, &
"PRINT%TOPOLOGY_INFO/UTIL_INFO")
CALL timestop(handle)
END SUBROUTINE topology_generate_molecule
! **************************************************************************************************
!> \brief Use info from periodic table and assumptions to generate bonds
!> \param topology ...
!> \param para_env ...
!> \param subsys_section ...
!> \author Teodoro Laino 09.2006
! **************************************************************************************************
SUBROUTINE topology_generate_bond(topology, para_env, subsys_section)
TYPE(topology_parameters_type), INTENT(INOUT) :: topology
TYPE(mp_para_env_type), POINTER :: para_env
TYPE(section_vals_type), POINTER :: subsys_section
CHARACTER(len=*), PARAMETER :: routineN = 'topology_generate_bond'
CHARACTER(LEN=2) :: upper_sym_1
INTEGER :: cbond, handle, handle2, i, iatm1, iatm2, iatom, ibond, idim, iw, j, jatom, k, &
n_bonds, n_heavy_bonds, n_hydr_bonds, n_rep, natom, npairs, output_unit
INTEGER, ALLOCATABLE, DIMENSION(:) :: bond_a, bond_b, list, map_nb
INTEGER, DIMENSION(:), POINTER :: isolated_atoms, tmp_v
LOGICAL :: connectivity_ok, explicit, print_info
LOGICAL, ALLOCATABLE, DIMENSION(:) :: h_list
REAL(KIND=dp) :: bondparm_factor, cell_v(3), dr(3), &
ksign, my_maxrad, r2, r2_min, rbond, &
rbond2, tmp
REAL(KIND=dp), DIMENSION(1, 1) :: r_max, r_minsq
REAL(KIND=dp), DIMENSION(:), POINTER :: radius
REAL(KIND=dp), DIMENSION(:, :), POINTER :: pbc_coord
TYPE(array2_list_type), DIMENSION(:), POINTER :: bond_list
TYPE(atom_info_type), POINTER :: atom_info
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(atomic_kind_type), POINTER :: atomic_kind
TYPE(connectivity_info_type), POINTER :: conn_info
TYPE(cp_logger_type), POINTER :: logger
TYPE(fist_neighbor_type), POINTER :: nonbonded
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(section_vals_type), POINTER :: bond_section, generate_section, &
isolated_section
NULLIFY (logger, particle_set, atomic_kind_set, nonbonded, bond_section, generate_section)
NULLIFY (isolated_atoms, tmp_v)
CALL timeset(routineN, handle)
logger => cp_get_default_logger()
output_unit = cp_logger_get_default_io_unit(logger)
iw = cp_print_key_unit_nr(logger, subsys_section, "PRINT%TOPOLOGY_INFO/GENERATE_INFO", &
extension=".subsysLog")
! Get atoms that one considers isolated (like ions in solution)
ALLOCATE (isolated_atoms(0))
generate_section => section_vals_get_subs_vals(subsys_section, "TOPOLOGY%GENERATE")
isolated_section => section_vals_get_subs_vals(generate_section, "ISOLATED_ATOMS")
CALL section_vals_get(isolated_section, explicit=explicit)
IF (explicit) THEN
CALL section_vals_val_get(isolated_section, "LIST", n_rep_val=n_rep)
DO i = 1, n_rep
CALL section_vals_val_get(isolated_section, "LIST", i_vals=tmp_v, i_rep_val=i)
CALL reallocate(isolated_atoms, 1, SIZE(isolated_atoms) + SIZE(tmp_v))
isolated_atoms(SIZE(isolated_atoms) - SIZE(tmp_v) + 1:SIZE(isolated_atoms)) = tmp_v
END DO
END IF
atom_info => topology%atom_info
conn_info => topology%conn_info
bondparm_factor = topology%bondparm_factor
cbond = 0
natom = topology%natoms
NULLIFY (radius)
! Allocate temporary arrays
ALLOCATE (radius(natom))
ALLOCATE (list(natom))
ALLOCATE (h_list(natom))
ALLOCATE (pbc_coord(3, natom))
h_list = .FALSE.
CALL timeset(TRIM(routineN)//"_1", handle2)
DO iatom = 1, natom
list(iatom) = iatom
upper_sym_1 = id2str(atom_info%id_element(iatom))
IF (topology%bondparm_type == do_bondparm_covalent) THEN
CALL get_ptable_info(symbol=upper_sym_1, covalent_radius=radius(iatom))
ELSE IF (topology%bondparm_type == do_bondparm_vdw) THEN
CALL get_ptable_info(symbol=upper_sym_1, vdw_radius=radius(iatom))
ELSE
CPABORT("Illegal bondparm_type")
END IF
IF (upper_sym_1 == "H ") h_list(iatom) = .TRUE.
! isolated atoms? put the radius to 0.0_dp
IF (ANY(isolated_atoms == iatom)) radius(iatom) = 0.0_dp
radius(iatom) = cp_unit_to_cp2k(radius(iatom), "angstrom")
IF (iw > 0) WRITE (iw, '(T2,"GENERATE|",5X,A,T50,A5,T60,A,T69,F12.6)') &
"In topology_generate_bond :: iatom = ", upper_sym_1, &
"radius:", radius(iatom)
END DO
CALL timestop(handle2)
CALL timeset(TRIM(routineN)//"_2", handle2)
! Initialize fake particle_set and atomic_kinds to generate the bond list
! using the neighboring list routine
ALLOCATE (atomic_kind_set(1))
CALL allocate_particle_set(particle_set, natom)
!
my_maxrad = MAXVAL(radius)*2.0_dp
atomic_kind => atomic_kind_set(1)
CALL set_atomic_kind(atomic_kind=atomic_kind, kind_number=1, &
name="XXX", element_symbol="XXX", mass=0.0_dp, atom_list=list)
CALL section_vals_val_get(subsys_section, "TOPOLOGY%GENERATE%BONDLENGTH_MAX", r_val=tmp)
r_max = tmp
IF (my_maxrad*bondparm_factor > r_max(1, 1) .AND. (.NOT. topology%molname_generated)) THEN
IF (output_unit > 0) THEN
WRITE (output_unit, '(T2,"GENERATE|",A)') &
" ERROR in connectivity generation!", &
" The THRESHOLD to select possible bonds is larger than the max. bondlength", &
" used to build the neighbors lists. Increase the BONDLENGTH_MAX parameter"
WRITE (output_unit, '(T2,"GENERATE|",2(A,F11.6),A)') &
" Present THRESHOLD (", my_maxrad*bondparm_factor, " )."// &
" Present BONDLENGTH_MAX (", r_max(1, 1), " )"
END IF
CPABORT("")
END IF
DO i = 1, natom
particle_set(i)%atomic_kind => atomic_kind_set(1)
particle_set(i)%r(1) = atom_info%r(1, i)
particle_set(i)%r(2) = atom_info%r(2, i)
particle_set(i)%r(3) = atom_info%r(3, i)
pbc_coord(:, i) = pbc(atom_info%r(:, i), topology%cell)
END DO
CALL section_vals_val_get(subsys_section, "TOPOLOGY%GENERATE%BONDLENGTH_MIN", r_val=tmp)
r_minsq = tmp*tmp
CALL timestop(handle2)
CALL timeset(TRIM(routineN)//"_3", handle2)
CALL build_fist_neighbor_lists(atomic_kind_set, particle_set, &
cell=topology%cell, r_max=r_max, r_minsq=r_minsq, &
ei_scale14=1.0_dp, vdw_scale14=1.0_dp, nonbonded=nonbonded, &
para_env=para_env, build_from_scratch=.TRUE., geo_check=.TRUE., &
mm_section=generate_section)
IF (iw > 0) THEN
WRITE (iw, '(T2,"GENERATE| Number of prescreened bonds (neighbors):",T71,I10)') &
nonbonded%neighbor_kind_pairs(1)%npairs
END IF
npairs = 0
DO i = 1, SIZE(nonbonded%neighbor_kind_pairs)
npairs = npairs + nonbonded%neighbor_kind_pairs(i)%npairs
END DO
ALLOCATE (bond_a(npairs))
ALLOCATE (bond_b(npairs))
ALLOCATE (map_nb(npairs))
idim = 0
DO j = 1, SIZE(nonbonded%neighbor_kind_pairs)
DO i = 1, nonbonded%neighbor_kind_pairs(j)%npairs
idim = idim + 1
bond_a(idim) = nonbonded%neighbor_kind_pairs(j)%list(1, i)
bond_b(idim) = nonbonded%neighbor_kind_pairs(j)%list(2, i)
map_nb(idim) = j
END DO
END DO
CALL timestop(handle2)
CALL timeset(TRIM(routineN)//"_4", handle2)
! We have a list of neighbors let's order the list w.r.t. the particle number
ALLOCATE (bond_list(natom))
DO I = 1, natom
ALLOCATE (bond_list(I)%array1(0))
ALLOCATE (bond_list(I)%array2(0))
END DO
CALL reorder_structure(bond_list, bond_a, bond_b, map_nb, SIZE(bond_a))
DEALLOCATE (bond_a)
DEALLOCATE (bond_b)
DEALLOCATE (map_nb)
! Find the Real bonds in the system
! Let's start with heavy atoms.. hydrogens will be treated only later on...
! Heavy atoms loop
CALL reallocate(conn_info%bond_a, 1, 1)
CALL reallocate(conn_info%bond_b, 1, 1)
connectivity_ok = .FALSE.
! No need to check consistency between provided molecule name and
! generated connectivity since we overrided the molecule definition.
IF (topology%create_molecules) THEN
atom_info%id_molname = str2id(s2s("TO_DEFINE_LATER"))
! A real name assignment will then be performed in the reorder module..
END IF
! It may happen that the connectivity created is fault for the missing
! of one bond.. this external loop ensures that everything was created
! fits exactly with the definition of molecules..
DO WHILE (.NOT. connectivity_ok)
n_heavy_bonds = 0
n_bonds = 0
DO iatm1 = 1, natom
IF (h_list(iatm1)) CYCLE
DO j = 1, SIZE(bond_list(iatm1)%array1)
iatm2 = bond_list(iatm1)%array1(j)
IF (atom_info%id_molname(iatm1) /= atom_info%id_molname(iatm2)) CYCLE
IF (h_list(iatm2) .OR. (iatm2 <= iatm1)) CYCLE
k = bond_list(iatm1)%array2(j)
ksign = SIGN(1.0_dp, REAL(k, KIND=dp))
k = ABS(k)
cell_v = MATMUL(topology%cell%hmat, &
REAL(nonbonded%neighbor_kind_pairs(k)%cell_vector, KIND=dp))
dr = pbc_coord(:, iatm1) - pbc_coord(:, iatm2) - ksign*cell_v
r2 = DOT_PRODUCT(dr, dr)
IF (r2 <= r_minsq(1, 1)) THEN
CALL cp_abort(__LOCATION__, &
"bond distance between atoms less then the smallest distance provided "// &
"in input "//cp_to_string(tmp)//" [bohr]")
END IF
! Screen isolated atoms
IF ((ANY(isolated_atoms == iatm1)) .OR. (ANY(isolated_atoms == iatm2))) CYCLE
! Screen neighbors
IF (topology%bondparm_type == do_bondparm_covalent) THEN
rbond = radius(iatm1) + radius(iatm2)
ELSE IF (topology%bondparm_type == do_bondparm_vdw) THEN
rbond = MAX(radius(iatm1), radius(iatm2))
END IF
rbond2 = rbond*rbond
rbond2 = rbond2*(bondparm_factor)**2
!Test the distance to the sum of the covalent radius
IF (r2 <= rbond2) THEN
n_heavy_bonds = n_heavy_bonds + 1
CALL add_bonds_list(conn_info, iatm1, iatm2, n_heavy_bonds)
END IF
END DO
END DO
n_hydr_bonds = 0
n_bonds = n_heavy_bonds
! Now check bonds formed by hydrogens...
! The hydrogen valence is 1 so we can choose the closest atom..
IF (output_unit > 0) WRITE (output_unit, *)
DO iatm1 = 1, natom
IF (.NOT. h_list(iatm1)) CYCLE
r2_min = HUGE(0.0_dp)
ibond = -1
print_info = .TRUE.
DO j = 1, SIZE(bond_list(iatm1)%array1)
iatm2 = bond_list(iatm1)%array1(j)
print_info = .FALSE.
IF (atom_info%id_molname(iatm1) /= atom_info%id_molname(iatm2)) CYCLE
IF (h_list(iatm2) .AND. (iatm2 <= iatm1)) CYCLE
! Screen isolated atoms
IF ((ANY(isolated_atoms == iatm1)) .OR. (ANY(isolated_atoms == iatm2))) CYCLE
k = bond_list(iatm1)%array2(j)
ksign = SIGN(1.0_dp, REAL(k, KIND=dp))
k = ABS(k)
cell_v = MATMUL(topology%cell%hmat, &
REAL(nonbonded%neighbor_kind_pairs(k)%cell_vector, KIND=dp))
dr = pbc_coord(:, iatm1) - pbc_coord(:, iatm2) - ksign*cell_v
r2 = DOT_PRODUCT(dr, dr)
IF (r2 <= r_minsq(1, 1)) THEN
CALL cp_abort(__LOCATION__, &
"bond distance between atoms less then the smallest distance provided "// &
"in input "//cp_to_string(tmp)//" [bohr]")
END IF
IF (r2 <= r2_min) THEN
r2_min = r2
ibond = iatm2
END IF
END DO
IF (ibond == -1) THEN
IF (output_unit > 0 .AND. print_info) THEN
WRITE (output_unit, '(T2,"GENERATE|",1X,A,I10,A)') &
"WARNING:: No connections detected for Hydrogen - Atom Nr:", iatm1, " !"
END IF
ELSE
n_hydr_bonds = n_hydr_bonds + 1
n_bonds = n_bonds + 1
CALL add_bonds_list(conn_info, MIN(iatm1, ibond), MAX(iatm1, ibond), n_bonds)
END IF
END DO
IF (output_unit > 0) THEN
WRITE (output_unit, '(T2,"GENERATE|",1X,A,T71,I10)') &
" Preliminary Number of Bonds generated:", n_bonds
END IF
! External defined bonds (useful for complex connectivity)
bond_section => section_vals_get_subs_vals(generate_section, "BOND")
CALL connectivity_external_control(section=bond_section, &
Iarray1=conn_info%bond_a, &
Iarray2=conn_info%bond_b, &
nvar=n_bonds, &
topology=topology, &
output_unit=output_unit)
! Resize arrays to their proper size..
CALL reallocate(conn_info%bond_a, 1, n_bonds)
CALL reallocate(conn_info%bond_b, 1, n_bonds)
IF (topology%create_molecules) THEN
! Since we create molecule names we're not sure that all atoms are contiguous
! so we need to reorder them on the basis of the generated name
IF (.NOT. topology%reorder_atom) THEN
topology%reorder_atom = .TRUE.
IF (output_unit > 0) WRITE (output_unit, '(T2,"GENERATE|",A)') &
" Molecules names have been generated. Now reordering particle set in order to have ", &
" atoms belonging to the same molecule in a sequential order."
END IF
connectivity_ok = .TRUE.
ELSE
! Check created connectivity and possibly give the OK to proceed
connectivity_ok = check_generate_mol(conn_info%bond_a, conn_info%bond_b, &
atom_info, bondparm_factor, output_unit)
END IF
IF (my_maxrad*bondparm_factor > r_max(1, 1) .AND. (.NOT. topology%molname_generated)) THEN
IF (output_unit > 0) THEN
WRITE (output_unit, '(T2,"GENERATE|",A)') &
" ERROR in connectivity generation!", &
" The THRESHOLD to select possible bonds is bigger than the MAX bondlength", &
" used to build the neighbors lists. Increase the BONDLENGTH_MAX patameter"
WRITE (output_unit, '(T2,"GENERATE|",2(A,F11.6),A)') &
" Present THRESHOLD (", my_maxrad*bondparm_factor, " )."// &
" Present BONDLENGTH_MAX (", r_max(1, 1), " )"
END IF
CPABORT("")
END IF
END DO
IF (connectivity_ok .AND. (output_unit > 0)) THEN
WRITE (output_unit, '(T2,"GENERATE|",A)') &
" Achieved consistency in connectivity generation."
END IF
CALL fist_neighbor_deallocate(nonbonded)
CALL timestop(handle2)
CALL timeset(TRIM(routineN)//"_6", handle2)
! Deallocate temporary working arrays
DO I = 1, natom
DEALLOCATE (bond_list(I)%array1)
DEALLOCATE (bond_list(I)%array2)
END DO
DEALLOCATE (bond_list)
DEALLOCATE (pbc_coord)
DEALLOCATE (radius)
DEALLOCATE (list)
CALL deallocate_particle_set(particle_set)
CALL deallocate_atomic_kind_set(atomic_kind_set)
!
CALL timestop(handle2)
IF (output_unit > 0 .AND. n_bonds > 0) THEN
WRITE (output_unit, '(T2,"GENERATE|",1X,A,T71,I10)') " Number of Bonds generated:", &
n_bonds
END IF
CALL timeset(TRIM(routineN)//"_7", handle2)
! If PARA_RES then activate RESIDUES
CALL reallocate(conn_info%c_bond_a, 1, 0)
CALL reallocate(conn_info%c_bond_b, 1, 0)
IF (topology%para_res) THEN
DO ibond = 1, SIZE(conn_info%bond_a)
iatom = conn_info%bond_a(ibond)
jatom = conn_info%bond_b(ibond)
IF ((atom_info%id_molname(iatom) /= atom_info%id_molname(jatom)) .OR. &
(atom_info%resid(iatom) /= atom_info%resid(jatom)) .OR. &
(atom_info%id_resname(iatom) /= atom_info%id_resname(jatom))) THEN
IF (iw > 0) WRITE (iw, *) " PARA_RES, bond between molecules atom ", &
iatom, jatom
cbond = cbond + 1
CALL reallocate(conn_info%c_bond_a, 1, cbond)
CALL reallocate(conn_info%c_bond_b, 1, cbond)
conn_info%c_bond_a(cbond) = iatom
conn_info%c_bond_b(cbond) = jatom
ELSE
IF (atom_info%id_molname(iatom) /= atom_info%id_molname(jatom)) THEN
CPABORT("Bonds between different molecule types?")
END IF
END IF
END DO
END IF
CALL timestop(handle2)
DEALLOCATE (isolated_atoms)
CALL timestop(handle)
CALL cp_print_key_finished_output(iw, logger, subsys_section, &
"PRINT%TOPOLOGY_INFO/GENERATE_INFO")
END SUBROUTINE topology_generate_bond
! **************************************************************************************************
!> \brief Performs a check on the generated connectivity
!> \param bond_a ...
!> \param bond_b ...
!> \param atom_info ...
!> \param bondparm_factor ...
!> \param output_unit ...
!> \return ...
!> \author Teodoro Laino 09.2006
! **************************************************************************************************
FUNCTION check_generate_mol(bond_a, bond_b, atom_info, bondparm_factor, output_unit) &
RESULT(conn_ok)
INTEGER, DIMENSION(:), POINTER :: bond_a, bond_b
TYPE(atom_info_type), POINTER :: atom_info
REAL(KIND=dp), INTENT(INOUT) :: bondparm_factor
INTEGER, INTENT(IN) :: output_unit
LOGICAL :: conn_ok
CHARACTER(len=*), PARAMETER :: routineN = 'check_generate_mol'
CHARACTER(LEN=10) :: ctmp1, ctmp2, ctmp3
INTEGER :: handle, i, idim, itype, j, mol_natom, &
natom, nsize
INTEGER, ALLOCATABLE, DIMENSION(:, :) :: mol_info_tmp
INTEGER, DIMENSION(:), POINTER :: mol_map, mol_map_o, wrk
INTEGER, DIMENSION(:, :), POINTER :: mol_info
LOGICAL, DIMENSION(:), POINTER :: icheck
TYPE(array1_list_type), DIMENSION(:), POINTER :: bond_list
CALL timeset(routineN, handle)
conn_ok = .TRUE.
natom = SIZE(atom_info%id_atmname)
ALLOCATE (bond_list(natom))
DO I = 1, natom
ALLOCATE (bond_list(I)%array1(0))
END DO
CALL reorder_structure(bond_list, bond_a, bond_b, SIZE(bond_a))
ALLOCATE (mol_map(natom))
ALLOCATE (mol_map_o(natom))
ALLOCATE (wrk(natom))
DO i = 1, natom
mol_map(i) = atom_info%id_molname(i)
END DO
mol_map_o = mol_map
CALL sort(mol_map, natom, wrk)
!
! mol(i,1) : stores id of the molecule
! mol(i,2) : stores the total number of atoms forming that kind of molecule
! mol(i,3) : contains the number of molecules generated for that kind
! mol(i,4) : contains the number of atoms forming one molecule of that kind
! Connectivity will be considered correct only if for each i:
!
! mol(i,2) = mol(i,3)*mol(i,4)
!
! If not, very probably, a bond is missing increase bondparm by 10% and let's
! check if the newest connectivity is bug free..
!
ALLOCATE (mol_info_tmp(natom, 2))
itype = mol_map(1)
nsize = 1
idim = 1
mol_info_tmp(1, 1) = itype
DO i = 2, natom
IF (mol_map(i) /= itype) THEN
nsize = nsize + 1
itype = mol_map(i)
mol_info_tmp(nsize, 1) = itype
mol_info_tmp(nsize - 1, 2) = idim
idim = 1
ELSE
idim = idim + 1
END IF
END DO
mol_info_tmp(nsize, 2) = idim
ALLOCATE (mol_info(nsize, 4))
mol_info(1:nsize, 1:2) = mol_info_tmp(1:nsize, 1:2)
DEALLOCATE (mol_info_tmp)
DO i = 1, nsize