-
Notifications
You must be signed in to change notification settings - Fork 1
/
topology_util.F
1257 lines (1160 loc) · 55.1 KB
/
topology_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
!> jgh (23-05-2004) Last atom of molecule information added
! **************************************************************************************************
MODULE topology_util
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 graphcon, ONLY: graph_type,&
hash_molecule,&
reorder_graph,&
vertex
USE input_section_types, ONLY: section_vals_get,&
section_vals_get_subs_vals,&
section_vals_type,&
section_vals_val_get,&
section_vals_val_set
USE kinds, ONLY: default_string_length,&
dp
USE mm_mapping_library, ONLY: amber_map,&
charmm_map,&
gromos_map
USE periodic_table, ONLY: get_ptable_info
USE qmmm_types_low, ONLY: qmmm_env_mm_type
USE string_table, ONLY: id2str,&
str2id
USE string_utilities, ONLY: uppercase
USE topology_types, ONLY: atom_info_type,&
connectivity_info_type,&
topology_parameters_type
USE util, ONLY: sort
#include "./base/base_uses.f90"
IMPLICIT NONE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'topology_util'
! **************************************************************************************************
TYPE array1_list_type
INTEGER, POINTER, DIMENSION(:) :: array1 => NULL()
END TYPE array1_list_type
! **************************************************************************************************
TYPE array2_list_type
INTEGER, POINTER, DIMENSION(:) :: array1 => NULL()
INTEGER, POINTER, DIMENSION(:) :: array2 => NULL()
END TYPE array2_list_type
PRIVATE
PUBLIC :: topology_set_atm_mass, &
topology_reorder_atoms, &
topology_molecules_check, &
check_subsys_element, &
reorder_structure, &
find_molecule, &
array1_list_type, &
array2_list_type, &
give_back_molecule, &
reorder_list_array, &
tag_molecule
INTERFACE reorder_structure
MODULE PROCEDURE reorder_structure1d, reorder_structure2d
END INTERFACE
CONTAINS
! **************************************************************************************************
!> \brief ...
!> \param topology ...
!> \param qmmm ...
!> \param qmmm_env_mm ...
!> \param subsys_section ...
!> \param force_env_section ...
!> \par History
!> Teodoro Laino 09.2006 - Rewritten with a graph matching algorithm
! **************************************************************************************************
SUBROUTINE topology_reorder_atoms(topology, qmmm, qmmm_env_mm, subsys_section, force_env_section)
TYPE(topology_parameters_type), INTENT(INOUT) :: topology
LOGICAL, INTENT(in), OPTIONAL :: qmmm
TYPE(qmmm_env_mm_type), OPTIONAL, POINTER :: qmmm_env_mm
TYPE(section_vals_type), POINTER :: subsys_section, force_env_section
CHARACTER(len=*), PARAMETER :: routineN = 'topology_reorder_atoms'
CHARACTER(LEN=default_string_length) :: mol_id
CHARACTER(LEN=default_string_length), POINTER :: molname(:), telement(:), &
tlabel_atmname(:), tlabel_molname(:), &
tlabel_resname(:)
INTEGER :: handle, i, iatm, iindex, ikind, imol, imol_ref, iref, iund, iw, j, k, location, &
max_mol_num, mm_index, n, n_rep, n_var, natom, natom_loc, nkind, nlinks, old_hash, &
old_mol, output_unit, qm_index, unique_mol
INTEGER, DIMENSION(:), POINTER :: mm_link_atoms, qm_atom_index
INTEGER, POINTER :: atm_map1(:), atm_map2(:), atm_map3(:), map_atom_type(:), &
map_mol_hash(:), mm_indexes_n(:), mm_indexes_v(:), mol_bnd(:, :), mol_hash(:), &
mol_num(:), new_position(:), order(:), tmp_n(:), tmp_v(:), wrk(:)
LOGICAL :: explicit, matches, my_qmmm
REAL(KIND=dp), DIMENSION(:, :), POINTER :: tr
REAL(KIND=dp), POINTER :: tatm_charge(:), tatm_mass(:)
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
TYPE(graph_type), DIMENSION(:), POINTER :: reference_set
TYPE(section_vals_type), POINTER :: generate_section, isolated_section, &
qm_kinds, qmmm_link_section, &
qmmm_section
TYPE(vertex), DIMENSION(:), POINTER :: reference, unordered
NULLIFY (logger, generate_section, isolated_section, tmp_v, tmp_n)
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)
output_unit = cp_logger_get_default_io_unit(logger)
IF (output_unit > 0) WRITE (output_unit, '(T2,"REORDER | ")')
my_qmmm = .FALSE.
IF (PRESENT(qmmm) .AND. PRESENT(qmmm_env_mm)) my_qmmm = qmmm
atom_info => topology%atom_info
conn_info => topology%conn_info
natom = topology%natoms
NULLIFY (new_position, reference_set)
NULLIFY (tlabel_atmname, telement, mol_num, tlabel_molname, tlabel_resname)
NULLIFY (tr, tatm_charge, tatm_mass, atm_map1, atm_map2, atm_map3)
! This routine can be called only at a very high level where these structures are still
! not even taken into account...
CPASSERT(.NOT. ASSOCIATED(atom_info%map_mol_num))
CPASSERT(.NOT. ASSOCIATED(atom_info%map_mol_typ))
CPASSERT(.NOT. ASSOCIATED(atom_info%map_mol_res))
!ALLOCATE all the temporary arrays needed for this routine
ALLOCATE (new_position(natom))
ALLOCATE (mol_num(natom))
ALLOCATE (molname(natom))
ALLOCATE (tlabel_atmname(natom))
ALLOCATE (tlabel_molname(natom))
ALLOCATE (tlabel_resname(natom))
ALLOCATE (tr(3, natom))
ALLOCATE (tatm_charge(natom))
ALLOCATE (tatm_mass(natom))
ALLOCATE (telement(natom))
ALLOCATE (atm_map1(natom))
ALLOCATE (atm_map2(natom))
! The only information we have at this level is the connectivity of the system.
! 0. Build a list of mapping atom types
ALLOCATE (map_atom_type(natom))
! 1. Build a list of bonds
ALLOCATE (atom_bond_list(natom))
DO I = 1, natom
map_atom_type(I) = atom_info%id_atmname(i)
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)
ALLOCATE (atom_info%map_mol_num(natom))
atom_info%map_mol_num = -1
CALL find_molecule(atom_bond_list, atom_info%map_mol_num, atom_info%id_molname)
max_mol_num = MAXVAL(atom_info%map_mol_num)
! In atom_info%map_mol_num have already been mapped molecules
ALLOCATE (mol_bnd(2, max_mol_num))
ALLOCATE (mol_hash(max_mol_num))
ALLOCATE (map_mol_hash(max_mol_num))
! 2. Sort the map_mol_num array.. atm_map1 contains now the mapped index
! of the reordered array
CALL sort(atom_info%map_mol_num, natom, atm_map1)
old_mol = 0
iindex = 0
imol = 0
DO i = 1, natom
IF (old_mol .NE. atom_info%map_mol_num(I)) THEN
old_mol = atom_info%map_mol_num(I)
iindex = 0
IF (imol > 0) THEN
mol_bnd(2, imol) = i - 1
END IF
imol = imol + 1
mol_bnd(1, imol) = i
END IF
iindex = iindex + 1
atm_map2(atm_map1(i)) = iindex
END DO
mol_bnd(2, imol) = natom
! Indexes of the two molecules to check
iref = 1
iund = max_mol_num/2 + 1
! Allocate reference and unordered
NULLIFY (reference, unordered)
! This is the real matching of graphs
DO j = 1, max_mol_num
CALL setup_graph(j, reference, map_atom_type, &
atom_bond_list, mol_bnd, atm_map1, atm_map2)
ALLOCATE (order(SIZE(reference)))
CALL hash_molecule(reference, order, mol_hash(j))
DEALLOCATE (order)
DO I = 1, SIZE(reference)
DEALLOCATE (reference(I)%bonds)
END DO
DEALLOCATE (reference)
END DO
! Reorder molecules hashes
CALL sort(mol_hash, max_mol_num, map_mol_hash)
! Now find unique molecules and reorder atoms too (if molecules match)
old_hash = -1
unique_mol = 0
natom_loc = 0
IF (output_unit > 0) THEN
WRITE (output_unit, '(T2,"REORDER | ",A)') &
"Reordering Molecules. The Reordering of molecules will override all", &
"information regarding molecule names and residue names.", &
"New ones will be provided on the basis of the connectivity!"
END IF
DO j = 1, max_mol_num
IF (mol_hash(j) .NE. old_hash) THEN
unique_mol = unique_mol + 1
old_hash = mol_hash(j)
CALL setup_graph_set(reference_set, unique_mol, map_mol_hash(j), &
map_atom_type, atom_bond_list, mol_bnd, atm_map1, atm_map2, &
atm_map3)
! Reorder Last added reference
mol_id = TRIM(ADJUSTL(cp_to_string(unique_mol)))
DO i = 1, SIZE(atm_map3)
natom_loc = natom_loc + 1
new_position(natom_loc) = atm_map3(i)
molname(natom_loc) = mol_id
mol_num(natom_loc) = unique_mol
END DO
DEALLOCATE (atm_map3)
ELSE
CALL setup_graph(map_mol_hash(j), unordered, map_atom_type, &
atom_bond_list, mol_bnd, atm_map1, atm_map2, atm_map3)
DO imol_ref = 1, unique_mol
!
reference => reference_set(imol_ref)%graph
ALLOCATE (order(SIZE(reference)))
CALL reorder_graph(reference, unordered, order, matches)
IF (matches) EXIT
DEALLOCATE (order)
END DO
IF (matches) THEN
! Reorder according the correct index
ALLOCATE (wrk(SIZE(order)))
CALL sort(order, SIZE(order), wrk)
DO i = 1, SIZE(order)
natom_loc = natom_loc + 1
new_position(natom_loc) = atm_map3(wrk(i))
molname(natom_loc) = mol_id
mol_num(natom_loc) = unique_mol
END DO
!
DEALLOCATE (order)
DEALLOCATE (wrk)
ELSE
unique_mol = unique_mol + 1
CALL setup_graph_set(reference_set, unique_mol, map_mol_hash(j), &
map_atom_type, atom_bond_list, mol_bnd, atm_map1, atm_map2, &
atm_map3)
! Reorder Last added reference
mol_id = TRIM(ADJUSTL(cp_to_string(unique_mol)))
DO i = 1, SIZE(atm_map3)
natom_loc = natom_loc + 1
new_position(natom_loc) = atm_map3(i)
molname(natom_loc) = mol_id
mol_num(natom_loc) = unique_mol
END DO
DEALLOCATE (atm_map3)
END IF
DO I = 1, SIZE(unordered)
DEALLOCATE (unordered(I)%bonds)
END DO
DEALLOCATE (unordered)
DEALLOCATE (atm_map3)
END IF
END DO
IF (output_unit > 0) THEN
WRITE (output_unit, '(T2,"REORDER | ",A,I7,A)') "Number of unique molecules found:", unique_mol, "."
END IF
CPASSERT(natom_loc == natom)
DEALLOCATE (map_atom_type)
DEALLOCATE (atm_map1)
DEALLOCATE (atm_map2)
DEALLOCATE (mol_bnd)
DEALLOCATE (mol_hash)
DEALLOCATE (map_mol_hash)
! Deallocate working arrays
DO I = 1, natom
DEALLOCATE (atom_bond_list(I)%array1)
END DO
DEALLOCATE (atom_bond_list)
DEALLOCATE (atom_info%map_mol_num)
! Deallocate reference_set
DO i = 1, SIZE(reference_set)
DO j = 1, SIZE(reference_set(i)%graph)
DEALLOCATE (reference_set(i)%graph(j)%bonds)
END DO
DEALLOCATE (reference_set(i)%graph)
END DO
DEALLOCATE (reference_set)
!Lets swap the atoms now
DO iatm = 1, natom
location = new_position(iatm)
tlabel_molname(iatm) = id2str(atom_info%id_molname(location))
tlabel_resname(iatm) = id2str(atom_info%id_resname(location))
tlabel_atmname(iatm) = id2str(atom_info%id_atmname(location))
telement(iatm) = id2str(atom_info%id_element(location))
tr(1, iatm) = atom_info%r(1, location)
tr(2, iatm) = atom_info%r(2, location)
tr(3, iatm) = atom_info%r(3, location)
tatm_charge(iatm) = atom_info%atm_charge(location)
tatm_mass(iatm) = atom_info%atm_mass(location)
END DO
IF (topology%create_molecules) THEN
DO iatm = 1, natom
tlabel_molname(iatm) = "MOL"//TRIM(molname(iatm))
tlabel_resname(iatm) = "R"//TRIM(molname(iatm))
END DO
topology%create_molecules = .FALSE.
END IF
DO iatm = 1, natom
atom_info%id_molname(iatm) = str2id(tlabel_molname(iatm))
atom_info%id_resname(iatm) = str2id(tlabel_resname(iatm))
atom_info%id_atmname(iatm) = str2id(tlabel_atmname(iatm))
atom_info%id_element(iatm) = str2id(telement(iatm))
atom_info%resid(iatm) = mol_num(iatm)
atom_info%r(1, iatm) = tr(1, iatm)
atom_info%r(2, iatm) = tr(2, iatm)
atom_info%r(3, iatm) = tr(3, iatm)
atom_info%atm_charge(iatm) = tatm_charge(iatm)
atom_info%atm_mass(iatm) = tatm_mass(iatm)
END DO
! Let's reorder all the list provided in the input..
ALLOCATE (wrk(SIZE(new_position)))
CALL sort(new_position, SIZE(new_position), wrk)
! NOTE: In the future the whole list of possible integers should be updated at this level..
! Let's update all the integer lists in the qmmm_env_mm section and in the input sections
! from where qmmm_env_qm will be read.
IF (my_qmmm) THEN
! Update the qmmm_env_mm
CPASSERT(SIZE(qmmm_env_mm%qm_atom_index) /= 0)
CPASSERT(ALL(qmmm_env_mm%qm_atom_index /= 0))
ALLOCATE (qm_atom_index(SIZE(qmmm_env_mm%qm_atom_index)))
DO iatm = 1, SIZE(qmmm_env_mm%qm_atom_index)
qm_atom_index(iatm) = wrk(qmmm_env_mm%qm_atom_index(iatm))
END DO
qmmm_env_mm%qm_atom_index = qm_atom_index
CPASSERT(ALL(qmmm_env_mm%qm_atom_index /= 0))
DEALLOCATE (qm_atom_index)
! Update the qmmm_section: MM_INDEX of the QM_KIND
qmmm_section => section_vals_get_subs_vals(force_env_section, "QMMM")
qm_kinds => section_vals_get_subs_vals(qmmm_section, "QM_KIND")
CALL section_vals_get(qm_kinds, n_repetition=nkind)
DO ikind = 1, nkind
CALL section_vals_val_get(qm_kinds, "MM_INDEX", i_rep_section=ikind, n_rep_val=n_var)
DO k = 1, n_var
CALL section_vals_val_get(qm_kinds, "MM_INDEX", i_rep_section=ikind, i_rep_val=k, &
i_vals=mm_indexes_v)
ALLOCATE (mm_indexes_n(SIZE(mm_indexes_v)))
DO j = 1, SIZE(mm_indexes_v)
mm_indexes_n(j) = wrk(mm_indexes_v(j))
END DO
CALL section_vals_val_set(qm_kinds, "MM_INDEX", i_rep_section=ikind, &
i_vals_ptr=mm_indexes_n, i_rep_val=k)
END DO
END DO
! Handle the link atoms
IF (qmmm_env_mm%qmmm_link) THEN
! Update the qmmm_env_mm
CPASSERT(SIZE(qmmm_env_mm%mm_link_atoms) > 0)
ALLOCATE (mm_link_atoms(SIZE(qmmm_env_mm%mm_link_atoms)))
DO iatm = 1, SIZE(qmmm_env_mm%mm_link_atoms)
mm_link_atoms(iatm) = wrk(qmmm_env_mm%mm_link_atoms(iatm))
END DO
qmmm_env_mm%mm_link_atoms = mm_link_atoms
CPASSERT(ALL(qmmm_env_mm%mm_link_atoms /= 0))
DEALLOCATE (mm_link_atoms)
! Update the qmmm_section: MM_INDEX,QM_INDEX of the LINK atom list
qmmm_link_section => section_vals_get_subs_vals(qmmm_section, "LINK")
CALL section_vals_get(qmmm_link_section, n_repetition=nlinks)
CPASSERT(nlinks /= 0)
DO ikind = 1, nlinks
CALL section_vals_val_get(qmmm_link_section, "QM_INDEX", i_rep_section=ikind, i_val=qm_index)
CALL section_vals_val_get(qmmm_link_section, "MM_INDEX", i_rep_section=ikind, i_val=mm_index)
mm_index = wrk(mm_index)
qm_index = wrk(qm_index)
CALL section_vals_val_set(qmmm_link_section, "MM_INDEX", i_rep_section=ikind, i_val=mm_index)
CALL section_vals_val_set(qmmm_link_section, "QM_INDEX", i_rep_section=ikind, i_val=qm_index)
END DO
END IF
END IF
!
!LIST of ISOLATED atoms
!
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)
ALLOCATE (tmp_n(SIZE(tmp_v)))
DO j = 1, SIZE(tmp_v)
tmp_n(j) = wrk(tmp_v(j))
END DO
CALL section_vals_val_set(isolated_section, "LIST", i_vals_ptr=tmp_n, i_rep_val=i)
END DO
END IF
DEALLOCATE (wrk)
!DEALLOCATE all the temporary arrays needed for this routine
DEALLOCATE (new_position)
DEALLOCATE (tlabel_atmname)
DEALLOCATE (tlabel_molname)
DEALLOCATE (tlabel_resname)
DEALLOCATE (telement)
DEALLOCATE (tr)
DEALLOCATE (tatm_charge)
DEALLOCATE (tatm_mass)
DEALLOCATE (molname)
DEALLOCATE (mol_num)
! DEALLOCATE the bond structures in the connectivity
DEALLOCATE (conn_info%bond_a)
DEALLOCATE (conn_info%bond_b)
IF (output_unit > 0) WRITE (output_unit, '(T2,"REORDER | ")')
CALL timestop(handle)
CALL cp_print_key_finished_output(iw, logger, subsys_section, &
"PRINT%TOPOLOGY_INFO/UTIL_INFO")
END SUBROUTINE topology_reorder_atoms
! **************************************************************************************************
!> \brief Set up a SET of graph kind
!> \param graph_set ...
!> \param idim ...
!> \param ind ...
!> \param array2 ...
!> \param atom_bond_list ...
!> \param map_mol ...
!> \param atm_map1 ...
!> \param atm_map2 ...
!> \param atm_map3 ...
!> \author Teodoro Laino 10.2006
! **************************************************************************************************
SUBROUTINE setup_graph_set(graph_set, idim, ind, array2, atom_bond_list, map_mol, &
atm_map1, atm_map2, atm_map3)
TYPE(graph_type), DIMENSION(:), POINTER :: graph_set
INTEGER, INTENT(IN) :: idim, ind
INTEGER, DIMENSION(:), POINTER :: array2
TYPE(array1_list_type), DIMENSION(:), INTENT(IN) :: atom_bond_list
INTEGER, DIMENSION(:, :), POINTER :: map_mol
INTEGER, DIMENSION(:), POINTER :: atm_map1, atm_map2, atm_map3
INTEGER :: ldim
TYPE(graph_type), DIMENSION(:), POINTER :: tmp_graph_set
ldim = 0
NULLIFY (tmp_graph_set)
IF (ASSOCIATED(graph_set)) THEN
ldim = SIZE(graph_set)
CPASSERT(ldim + 1 == idim)
MARK_USED(idim)
NULLIFY (tmp_graph_set)
CALL allocate_graph_set(graph_set, tmp_graph_set)
END IF
CALL allocate_graph_set(tmp_graph_set, graph_set, ldim, ldim + 1)
CALL setup_graph(ind, graph_set(ldim + 1)%graph, array2, &
atom_bond_list, map_mol, atm_map1, atm_map2, atm_map3)
END SUBROUTINE setup_graph_set
! **************************************************************************************************
!> \brief Allocate a new graph_set deallocating an old one..
!> \param graph_set_in ...
!> \param graph_set_out ...
!> \param ldim_in ...
!> \param ldim_out ...
!> \author Teodoro Laino 10.2006
! **************************************************************************************************
SUBROUTINE allocate_graph_set(graph_set_in, graph_set_out, ldim_in, ldim_out)
TYPE(graph_type), DIMENSION(:), POINTER :: graph_set_in, graph_set_out
INTEGER, INTENT(IN), OPTIONAL :: ldim_in, ldim_out
INTEGER :: b_dim, i, j, mydim_in, mydim_out, v_dim
CPASSERT(.NOT. ASSOCIATED(graph_set_out))
mydim_in = 0
mydim_out = 0
IF (ASSOCIATED(graph_set_in)) THEN
mydim_in = SIZE(graph_set_in)
mydim_out = SIZE(graph_set_in)
END IF
IF (PRESENT(ldim_in)) mydim_in = ldim_in
IF (PRESENT(ldim_out)) mydim_out = ldim_out
ALLOCATE (graph_set_out(mydim_out))
DO i = 1, mydim_out
NULLIFY (graph_set_out(i)%graph)
END DO
! Copy graph structure into the temporary array
DO i = 1, mydim_in
v_dim = SIZE(graph_set_in(i)%graph)
ALLOCATE (graph_set_out(i)%graph(v_dim))
DO j = 1, v_dim
graph_set_out(i)%graph(j)%kind = graph_set_in(i)%graph(j)%kind
b_dim = SIZE(graph_set_in(i)%graph(j)%bonds)
ALLOCATE (graph_set_out(i)%graph(j)%bonds(b_dim))
graph_set_out(i)%graph(j)%bonds = graph_set_in(i)%graph(j)%bonds
DEALLOCATE (graph_set_in(i)%graph(j)%bonds)
END DO
DEALLOCATE (graph_set_in(i)%graph)
END DO
IF (ASSOCIATED(graph_set_in)) THEN
DEALLOCATE (graph_set_in)
END IF
END SUBROUTINE allocate_graph_set
! **************************************************************************************************
!> \brief Set up a graph kind
!> \param ind ...
!> \param graph ...
!> \param array2 ...
!> \param atom_bond_list ...
!> \param map_mol ...
!> \param atm_map1 ...
!> \param atm_map2 ...
!> \param atm_map3 ...
!> \author Teodoro Laino 09.2006
! **************************************************************************************************
SUBROUTINE setup_graph(ind, graph, array2, atom_bond_list, map_mol, &
atm_map1, atm_map2, atm_map3)
INTEGER, INTENT(IN) :: ind
TYPE(vertex), DIMENSION(:), POINTER :: graph
INTEGER, DIMENSION(:), POINTER :: array2
TYPE(array1_list_type), DIMENSION(:), INTENT(IN) :: atom_bond_list
INTEGER, DIMENSION(:, :), POINTER :: map_mol
INTEGER, DIMENSION(:), POINTER :: atm_map1, atm_map2
INTEGER, DIMENSION(:), OPTIONAL, POINTER :: atm_map3
INTEGER :: i, idim, ifirst, ilast, j, nbonds, &
nelement
IF (PRESENT(atm_map3)) THEN
CPASSERT(.NOT. ASSOCIATED(atm_map3))
END IF
CPASSERT(.NOT. ASSOCIATED(graph))
! Setup reference graph
idim = 0
ifirst = map_mol(1, ind)
ilast = map_mol(2, ind)
nelement = ilast - ifirst + 1
ALLOCATE (graph(nelement))
IF (PRESENT(atm_map3)) THEN
ALLOCATE (atm_map3(nelement))
END IF
DO i = ifirst, ilast
idim = idim + 1
graph(idim)%kind = array2(atm_map1(i))
nbonds = SIZE(atom_bond_list(atm_map1(i))%array1)
ALLOCATE (graph(idim)%bonds(nbonds))
DO j = 1, nbonds
graph(idim)%bonds(j) = atm_map2(atom_bond_list(atm_map1(i))%array1(j))
END DO
IF (PRESENT(atm_map3)) THEN
atm_map3(idim) = atm_map1(i)
END IF
END DO
END SUBROUTINE setup_graph
! **************************************************************************************************
!> \brief Order arrays of lists
!> \param Ilist1 ...
!> \param Ilist2 ...
!> \param Ilist3 ...
!> \param Ilist4 ...
!> \param nsize ...
!> \param ndim ...
!> \author Teodoro Laino 09.2006
! **************************************************************************************************
RECURSIVE SUBROUTINE reorder_list_array(Ilist1, Ilist2, Ilist3, Ilist4, nsize, ndim)
INTEGER, DIMENSION(:), POINTER :: Ilist1
INTEGER, DIMENSION(:), OPTIONAL, POINTER :: Ilist2, Ilist3, Ilist4
INTEGER, INTENT(IN) :: nsize, ndim
INTEGER :: i, iend, istart, ldim
INTEGER, DIMENSION(:), POINTER :: tmp, wrk
CPASSERT(nsize > 0)
ALLOCATE (wrk(Ndim))
CALL sort(Ilist1, Ndim, wrk)
IF (nsize /= 1) THEN
ALLOCATE (tmp(Ndim))
tmp = Ilist2(1:Ndim)
DO i = 1, Ndim
Ilist2(i) = tmp(wrk(i))
END DO
SELECT CASE (nsize)
CASE (3)
tmp = Ilist3(1:Ndim)
DO i = 1, Ndim
Ilist3(i) = tmp(wrk(i))
END DO
CASE (4)
tmp = Ilist3(1:Ndim)
DO i = 1, Ndim
Ilist3(i) = tmp(wrk(i))
END DO
tmp = Ilist4(1:Ndim)
DO i = 1, Ndim
Ilist4(i) = tmp(wrk(i))
END DO
END SELECT
DEALLOCATE (tmp)
istart = 1
DO i = 1, Ndim
IF (Ilist1(i) /= Ilist1(istart)) THEN
iend = i - 1
ldim = iend - istart + 1
CALL reorder_list_array_low(Ilist2, Ilist3, Ilist4, nsize, &
ldim, istart, iend)
istart = i
END IF
END DO
! Last term to sort
iend = Ndim
ldim = iend - istart + 1
CALL reorder_list_array_low(Ilist2, Ilist3, Ilist4, nsize, &
ldim, istart, iend)
END IF
DEALLOCATE (wrk)
END SUBROUTINE reorder_list_array
! **************************************************************************************************
!> \brief Low level routine for ordering arrays of lists
!> \param Ilist2 ...
!> \param Ilist3 ...
!> \param Ilist4 ...
!> \param nsize ...
!> \param ldim ...
!> \param istart ...
!> \param iend ...
!> \author Teodoro Laino 09.2006
! **************************************************************************************************
RECURSIVE SUBROUTINE reorder_list_array_low(Ilist2, Ilist3, Ilist4, nsize, &
ldim, istart, iend)
INTEGER, DIMENSION(:), OPTIONAL, POINTER :: Ilist2, Ilist3, Ilist4
INTEGER, INTENT(IN) :: nsize, ldim, istart, iend
INTEGER, DIMENSION(:), POINTER :: tmp_2, tmp_3, tmp_4
SELECT CASE (nsize)
CASE (2)
ALLOCATE (tmp_2(ldim))
tmp_2(:) = Ilist2(istart:iend)
CALL reorder_list_array(tmp_2, nsize=nsize - 1, ndim=ldim)
Ilist2(istart:iend) = tmp_2(:)
DEALLOCATE (tmp_2)
CASE (3)
ALLOCATE (tmp_2(ldim))
ALLOCATE (tmp_3(ldim))
tmp_2(:) = Ilist2(istart:iend)
tmp_3(:) = Ilist3(istart:iend)
CALL reorder_list_array(tmp_2, tmp_3, nsize=nsize - 1, ndim=ldim)
Ilist2(istart:iend) = tmp_2(:)
Ilist3(istart:iend) = tmp_3(:)
DEALLOCATE (tmp_2)
DEALLOCATE (tmp_3)
CASE (4)
ALLOCATE (tmp_2(ldim))
ALLOCATE (tmp_3(ldim))
ALLOCATE (tmp_4(ldim))
tmp_2(:) = Ilist2(istart:iend)
tmp_3(:) = Ilist3(istart:iend)
tmp_4(:) = Ilist4(istart:iend)
CALL reorder_list_array(tmp_2, tmp_3, tmp_4, nsize=nsize - 1, ndim=ldim)
Ilist2(istart:iend) = tmp_2(:)
Ilist3(istart:iend) = tmp_3(:)
Ilist4(istart:iend) = tmp_4(:)
DEALLOCATE (tmp_2)
DEALLOCATE (tmp_3)
DEALLOCATE (tmp_4)
END SELECT
END SUBROUTINE reorder_list_array_low
! **************************************************************************************************
!> \brief ...
!> \param icheck ...
!> \param bond_list ...
!> \param i ...
!> \param mol_natom ...
!> \param mol_map ...
!> \param my_mol ...
!> \author Teodoro Laino 09.2006
! **************************************************************************************************
RECURSIVE SUBROUTINE give_back_molecule(icheck, bond_list, i, mol_natom, mol_map, my_mol)
LOGICAL, DIMENSION(:), POINTER :: icheck
TYPE(array1_list_type), DIMENSION(:), POINTER :: bond_list
INTEGER, INTENT(IN) :: i
INTEGER, INTENT(INOUT) :: mol_natom
INTEGER, DIMENSION(:), POINTER :: mol_map
INTEGER, INTENT(IN) :: my_mol
INTEGER :: j, k
IF (mol_map(i) == my_mol) THEN
icheck(i) = .TRUE.
DO j = 1, SIZE(bond_list(i)%array1)
k = bond_list(i)%array1(j)
IF (icheck(k)) CYCLE
mol_natom = mol_natom + 1
CALL give_back_molecule(icheck, bond_list, k, mol_natom, mol_map, my_mol)
END DO
ELSE
! Do nothing means only that bonds were found between molecules
! as we defined them.. This could easily be a bond detected but not
! physically present.. so just skip this part and go on..
END IF
END SUBROUTINE give_back_molecule
! **************************************************************************************************
!> \brief gives back a mapping of molecules.. icheck needs to be initialized with -1
!> \param icheck ...
!> \param bond_list ...
!> \param i ...
!> \param my_mol ...
!> \author Teodoro Laino 04.2007 - Zurich University
! **************************************************************************************************
RECURSIVE SUBROUTINE tag_molecule(icheck, bond_list, i, my_mol)
INTEGER, DIMENSION(:), POINTER :: icheck
TYPE(array1_list_type), DIMENSION(:), POINTER :: bond_list
INTEGER, INTENT(IN) :: i, my_mol
INTEGER :: j, k
icheck(i) = my_mol
DO j = 1, SIZE(bond_list(i)%array1)
k = bond_list(i)%array1(j)
IF (k <= i) CYCLE
CALL tag_molecule(icheck, bond_list, k, my_mol)
END DO
END SUBROUTINE tag_molecule
! **************************************************************************************************
!> \brief Given two lists of corresponding element a complex type is built in
!> which each element of the type has a list of mapping elements
!> \param work ...
!> \param list1 ...
!> \param list2 ...
!> \param N ...
!> \author Teodoro Laino 08.2006
! **************************************************************************************************
SUBROUTINE reorder_structure1d(work, list1, list2, N)
TYPE(array1_list_type), DIMENSION(:), &
INTENT(INOUT) :: work
INTEGER, DIMENSION(:), INTENT(IN) :: list1, list2
INTEGER, INTENT(IN) :: N
INTEGER :: I, index1, index2, Nsize
INTEGER, DIMENSION(:), POINTER :: wrk_tmp
DO I = 1, N
index1 = list1(I)
index2 = list2(I)
wrk_tmp => work(index1)%array1
Nsize = SIZE(wrk_tmp)
ALLOCATE (work(index1)%array1(Nsize + 1))
work(index1)%array1(1:Nsize) = wrk_tmp
work(index1)%array1(Nsize + 1) = index2
DEALLOCATE (wrk_tmp)
wrk_tmp => work(index2)%array1
Nsize = SIZE(wrk_tmp)
ALLOCATE (work(index2)%array1(Nsize + 1))
work(index2)%array1(1:Nsize) = wrk_tmp
work(index2)%array1(Nsize + 1) = index1
DEALLOCATE (wrk_tmp)
END DO
END SUBROUTINE reorder_structure1d
! **************************************************************************************************
!> \brief Given two lists of corresponding element a complex type is built in
!> which each element of the type has a list of mapping elements
!> \param work ...
!> \param list1 ...
!> \param list2 ...
!> \param list3 ...
!> \param N ...
!> \author Teodoro Laino 09.2006
! **************************************************************************************************
SUBROUTINE reorder_structure2d(work, list1, list2, list3, N)
TYPE(array2_list_type), DIMENSION(:), &
INTENT(INOUT) :: work
INTEGER, DIMENSION(:), INTENT(IN) :: list1, list2, list3
INTEGER, INTENT(IN) :: N
INTEGER :: I, index1, index2, index3, Nsize
INTEGER, DIMENSION(:), POINTER :: wrk_tmp
DO I = 1, N
index1 = list1(I)
index2 = list2(I)
index3 = list3(I)
wrk_tmp => work(index1)%array1
Nsize = SIZE(wrk_tmp)
ALLOCATE (work(index1)%array1(Nsize + 1))
work(index1)%array1(1:Nsize) = wrk_tmp
work(index1)%array1(Nsize + 1) = index2
DEALLOCATE (wrk_tmp)
wrk_tmp => work(index2)%array1
Nsize = SIZE(wrk_tmp)
ALLOCATE (work(index2)%array1(Nsize + 1))
work(index2)%array1(1:Nsize) = wrk_tmp
work(index2)%array1(Nsize + 1) = index1
DEALLOCATE (wrk_tmp)
wrk_tmp => work(index1)%array2
Nsize = SIZE(wrk_tmp)
ALLOCATE (work(index1)%array2(Nsize + 1))
work(index1)%array2(1:Nsize) = wrk_tmp
work(index1)%array2(Nsize + 1) = index3
DEALLOCATE (wrk_tmp)
wrk_tmp => work(index2)%array2
Nsize = SIZE(wrk_tmp)
ALLOCATE (work(index2)%array2(Nsize + 1))
work(index2)%array2(1:Nsize) = wrk_tmp
work(index2)%array2(Nsize + 1) = -index3
DEALLOCATE (wrk_tmp)
END DO
END SUBROUTINE reorder_structure2d
! **************************************************************************************************
!> \brief each atom will be assigned a molecule number based on bonded fragments
!> The array mol_info should be initialized with -1 before calling the
!> find_molecule routine
!> \param atom_bond_list ...
!> \param mol_info ...
!> \param mol_name ...
!> \author Joost 05.2006
! **************************************************************************************************
SUBROUTINE find_molecule(atom_bond_list, mol_info, mol_name)
TYPE(array1_list_type), DIMENSION(:), INTENT(IN) :: atom_bond_list
INTEGER, DIMENSION(:), POINTER :: mol_info, mol_name
INTEGER :: I, my_mol_name, N, nmol
N = SIZE(atom_bond_list)
nmol = 0
DO I = 1, N
IF (mol_info(I) == -1) THEN
nmol = nmol + 1
my_mol_name = mol_name(I)
CALL spread_mol(atom_bond_list, mol_info, i, nmol, my_mol_name, &
mol_name)
END IF
END DO
END SUBROUTINE find_molecule
! **************************************************************************************************
!> \brief spreads the molnumber over the bonded list
!> \param atom_bond_list ...
!> \param mol_info ...
!> \param iatom ...
!> \param imol ...
!> \param my_mol_name ...
!> \param mol_name ...
!> \author Joost 05.2006
! **************************************************************************************************
RECURSIVE SUBROUTINE spread_mol(atom_bond_list, mol_info, iatom, imol, &
my_mol_name, mol_name)
TYPE(array1_list_type), DIMENSION(:), INTENT(IN) :: atom_bond_list
INTEGER, DIMENSION(:), POINTER :: mol_info
INTEGER, INTENT(IN) :: iatom, imol, my_mol_name
INTEGER, DIMENSION(:), POINTER :: mol_name
INTEGER :: atom_b, i
mol_info(iatom) = imol
DO I = 1, SIZE(atom_bond_list(iatom)%array1)
atom_b = atom_bond_list(iatom)%array1(I)
! In this way we're really sure that all atoms belong to the same
! molecule. This should correct possible errors in the generation of
! the bond list..
IF (mol_info(atom_b) == -1 .AND. my_mol_name == mol_name(atom_b)) THEN
CALL spread_mol(atom_bond_list, mol_info, atom_b, imol, my_mol_name, mol_name)
IF (mol_info(atom_b) /= imol) CPABORT("internal error")
END IF
END DO
END SUBROUTINE spread_mol
! **************************************************************************************************
!> \brief Use info from periodic table and set atm_mass
!> \param topology ...
!> \param subsys_section ...
! **************************************************************************************************
SUBROUTINE topology_set_atm_mass(topology, subsys_section)
TYPE(topology_parameters_type), INTENT(INOUT) :: topology
TYPE(section_vals_type), POINTER :: subsys_section
CHARACTER(len=*), PARAMETER :: routineN = 'topology_set_atm_mass'
CHARACTER(LEN=2) :: upper_sym_1
CHARACTER(LEN=default_string_length) :: atmname_upper
CHARACTER(LEN=default_string_length), &
DIMENSION(:), POINTER :: keyword
INTEGER :: handle, i, i_rep, iatom, ielem_found, &
iw, n_rep, natom
LOGICAL :: user_defined
REAL(KIND=dp), DIMENSION(:), POINTER :: mass
TYPE(atom_info_type), POINTER :: atom_info
TYPE(cp_logger_type), POINTER :: logger
TYPE(section_vals_type), POINTER :: kind_section
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)
atom_info => topology%atom_info
natom = topology%natoms
! Available external info
kind_section => section_vals_get_subs_vals(subsys_section, "KIND")
CALL section_vals_get(kind_section, n_repetition=n_rep)
ALLOCATE (keyword(n_rep))
ALLOCATE (mass(n_rep))
mass = HUGE(0.0_dp)
DO i_rep = 1, n_rep
CALL section_vals_val_get(kind_section, "_SECTION_PARAMETERS_", &
c_val=keyword(i_rep), i_rep_section=i_rep)
CALL uppercase(keyword(i_rep))
CALL section_vals_val_get(kind_section, i_rep_section=i_rep, &
keyword_name="MASS", n_rep_val=i)
IF (i > 0) CALL section_vals_val_get(kind_section, i_rep_section=i_rep, &
keyword_name="MASS", r_val=mass(i_rep))
END DO
!
DO iatom = 1, natom
!If we reach this point then we've definitely identified the element..
!Let's look if an external mass has been defined..
user_defined = .FALSE.
DO i = 1, SIZE(keyword)
atmname_upper = id2str(atom_info%id_atmname(iatom))
CALL uppercase(atmname_upper)
IF (TRIM(atmname_upper) == TRIM(keyword(i)) .AND. mass(i) /= HUGE(0.0_dp)) THEN
atom_info%atm_mass(iatom) = mass(i)
user_defined = .TRUE.
EXIT
END IF
END DO
! If name didn't match let's try with the element
IF (.NOT. user_defined) THEN
upper_sym_1 = id2str(atom_info%id_element(iatom))
CALL get_ptable_info(symbol=upper_sym_1, ielement=ielem_found, amass=atom_info%atm_mass(iatom))
END IF
IF (iw > 0) WRITE (iw, '(7X,A,A5,A,F12.5)') "In topology_set_atm_mass :: element = ", &
id2str(atom_info%id_element(iatom)), " a_mass ", atom_info%atm_mass(iatom)
END DO
DEALLOCATE (keyword)
DEALLOCATE (mass)
CALL timestop(handle)
CALL cp_print_key_finished_output(iw, logger, subsys_section, &
"PRINT%TOPOLOGY_INFO/UTIL_INFO")
END SUBROUTINE topology_set_atm_mass
! **************************************************************************************************
!> \brief Check and verify that all molecules of the same kind are bonded the same
!> \param topology ...
!> \param subsys_section ...
! **************************************************************************************************
SUBROUTINE topology_molecules_check(topology, subsys_section)
TYPE(topology_parameters_type), INTENT(INOUT) :: topology
TYPE(section_vals_type), POINTER :: subsys_section