LCOV - code coverage report
Current view: top level - src - topology_constraint_util.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:06f838d) Lines: 87.3 % 726 634
Test Date: 2026-06-05 07:04:50 Functions: 100.0 % 10 10

            Line data    Source code
       1              : !--------------------------------------------------------------------------------------------------!
       2              : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3              : !   Copyright 2000-2026 CP2K developers group <https://cp2k.org>                                   !
       4              : !                                                                                                  !
       5              : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6              : !--------------------------------------------------------------------------------------------------!
       7              : 
       8              : ! **************************************************************************************************
       9              : !> \brief Collection of subroutine needed for topology related things
      10              : !> \par History
      11              : !>     jgh (23-05-2004) Last atom of molecule information added
      12              : ! **************************************************************************************************
      13              : MODULE topology_constraint_util
      14              :    USE atomic_kind_types,               ONLY: atomic_kind_type,&
      15              :                                               get_atomic_kind,&
      16              :                                               is_hydrogen
      17              :    USE cell_types,                      ONLY: cell_transform_input_cartesian,&
      18              :                                               use_perd_x,&
      19              :                                               use_perd_xy,&
      20              :                                               use_perd_xyz,&
      21              :                                               use_perd_xz,&
      22              :                                               use_perd_y,&
      23              :                                               use_perd_yz,&
      24              :                                               use_perd_z
      25              :    USE colvar_methods,                  ONLY: colvar_eval_mol_f
      26              :    USE colvar_types,                    ONLY: &
      27              :         colvar_clone, colvar_counters, colvar_create, colvar_p_reallocate, colvar_release, &
      28              :         colvar_setup, colvar_type, dist_colvar_id, torsion_colvar_id, xyz_diag_colvar_id, &
      29              :         xyz_outerdiag_colvar_id
      30              :    USE colvar_utils,                    ONLY: post_process_colvar
      31              :    USE cp_log_handling,                 ONLY: cp_get_default_logger,&
      32              :                                               cp_logger_type,&
      33              :                                               cp_to_string
      34              :    USE cp_output_handling,              ONLY: cp_print_key_finished_output,&
      35              :                                               cp_print_key_unit_nr
      36              :    USE input_constants,                 ONLY: do_constr_atomic,&
      37              :                                               do_constr_molec
      38              :    USE input_section_types,             ONLY: section_vals_get,&
      39              :                                               section_vals_get_subs_vals,&
      40              :                                               section_vals_type,&
      41              :                                               section_vals_val_get,&
      42              :                                               section_vals_val_set
      43              :    USE kinds,                           ONLY: default_string_length,&
      44              :                                               dp
      45              :    USE memory_utilities,                ONLY: reallocate
      46              :    USE molecule_kind_types,             ONLY: &
      47              :         atom_type, bond_type, colvar_constraint_type, fixd_constraint_type, g3x3_constraint_type, &
      48              :         g4x6_constraint_type, get_molecule_kind, molecule_kind_type, set_molecule_kind, &
      49              :         setup_colvar_counters, vsite_constraint_type
      50              :    USE molecule_types,                  ONLY: get_molecule,&
      51              :                                               global_constraint_type,&
      52              :                                               local_colvar_constraint_type,&
      53              :                                               local_constraint_type,&
      54              :                                               local_g3x3_constraint_type,&
      55              :                                               local_g4x6_constraint_type,&
      56              :                                               molecule_type,&
      57              :                                               set_molecule
      58              :    USE particle_types,                  ONLY: particle_type
      59              :    USE qmmm_ff_fist,                    ONLY: qmmm_ff_precond_only_qm
      60              :    USE qmmm_types_low,                  ONLY: qmmm_env_mm_type
      61              :    USE topology_types,                  ONLY: constr_list_type,&
      62              :                                               constraint_info_type,&
      63              :                                               topology_parameters_type
      64              : #include "./base/base_uses.f90"
      65              : 
      66              :    IMPLICIT NONE
      67              : 
      68              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'topology_constraint_util'
      69              : 
      70              :    PRIVATE
      71              :    PUBLIC :: topology_constraint_pack
      72              : 
      73              : CONTAINS
      74              : 
      75              : ! **************************************************************************************************
      76              : !> \brief Pack in all the information needed for the constraints
      77              : !> \param molecule_kind_set ...
      78              : !> \param molecule_set ...
      79              : !> \param topology ...
      80              : !> \param qmmm_env ...
      81              : !> \param particle_set ...
      82              : !> \param input_file ...
      83              : !> \param subsys_section ...
      84              : !> \param gci ...
      85              : ! **************************************************************************************************
      86        10726 :    SUBROUTINE topology_constraint_pack(molecule_kind_set, molecule_set, &
      87              :                                        topology, qmmm_env, particle_set, input_file, subsys_section, gci)
      88              :       TYPE(molecule_kind_type), DIMENSION(:), POINTER    :: molecule_kind_set
      89              :       TYPE(molecule_type), DIMENSION(:), POINTER         :: molecule_set
      90              :       TYPE(topology_parameters_type), INTENT(INOUT)      :: topology
      91              :       TYPE(qmmm_env_mm_type), OPTIONAL, POINTER          :: qmmm_env
      92              :       TYPE(particle_type), DIMENSION(:), POINTER         :: particle_set
      93              :       TYPE(section_vals_type), POINTER                   :: input_file, subsys_section
      94              :       TYPE(global_constraint_type), POINTER              :: gci
      95              : 
      96              :       CHARACTER(len=*), PARAMETER :: routineN = 'topology_constraint_pack'
      97              : 
      98              :       CHARACTER(LEN=2)                                   :: element_symbol
      99              :       CHARACTER(LEN=default_string_length)               :: molname, name
     100              :       CHARACTER(LEN=default_string_length), &
     101        10726 :          DIMENSION(:), POINTER                           :: atom_typeh, cnds
     102              :       INTEGER :: cind, first, first_atom, gind, handle, handle2, i, ii, itype, iw, j, k, k1loc, &
     103              :          k2loc, kk, last, last_atom, m, n_start_colv, natom, nbond, ncolv_glob, ncolv_mol, &
     104              :          nfixd_list_gci, nfixd_restart, nfixd_restraint, nfixed_atoms, ng3x3, ng3x3_restraint, &
     105              :          ng4x6, ng4x6_restraint, nhdist, nmolecule, nrep, nvsite, nvsite_restraint, offset
     106        10726 :       INTEGER, DIMENSION(:), POINTER                     :: constr_x_glob, inds, molecule_list
     107              :       LOGICAL :: exclude_mm, exclude_qm, fix_atom_mm, fix_atom_molname, fix_atom_qm, &
     108              :          fix_atom_qmmm, fix_fixed_atom, found_molname, is_qm, ishbond, ldummy, &
     109              :          restart_restraint_clv, restart_restraint_pos, use_clv_info
     110        10726 :       LOGICAL, ALLOCATABLE, DIMENSION(:)                 :: missed_molname
     111              :       REAL(KIND=dp)                                      :: rmod, rvec(3)
     112        10726 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: hdist, r
     113        10726 :       TYPE(atom_type), DIMENSION(:), POINTER             :: atom_list
     114              :       TYPE(atomic_kind_type), POINTER                    :: atomic_kind
     115        10726 :       TYPE(bond_type), DIMENSION(:), POINTER             :: bond_list
     116              :       TYPE(colvar_constraint_type), DIMENSION(:), &
     117        10726 :          POINTER                                         :: colv_list
     118              :       TYPE(colvar_counters)                              :: ncolv
     119        10726 :       TYPE(constr_list_type), DIMENSION(:), POINTER      :: constr_x_mol
     120              :       TYPE(constraint_info_type), POINTER                :: cons_info
     121              :       TYPE(cp_logger_type), POINTER                      :: logger
     122        10726 :       TYPE(fixd_constraint_type), DIMENSION(:), POINTER  :: fixd_list, fixd_list_gci
     123        10726 :       TYPE(g3x3_constraint_type), DIMENSION(:), POINTER  :: g3x3_list
     124        10726 :       TYPE(g4x6_constraint_type), DIMENSION(:), POINTER  :: g4x6_list
     125              :       TYPE(local_colvar_constraint_type), DIMENSION(:), &
     126        10726 :          POINTER                                         :: lcolv
     127              :       TYPE(local_constraint_type), POINTER               :: lci
     128              :       TYPE(local_g3x3_constraint_type), DIMENSION(:), &
     129        10726 :          POINTER                                         :: lg3x3
     130              :       TYPE(local_g4x6_constraint_type), DIMENSION(:), &
     131        10726 :          POINTER                                         :: lg4x6
     132              :       TYPE(molecule_kind_type), POINTER                  :: molecule_kind
     133              :       TYPE(molecule_type), POINTER                       :: molecule
     134              :       TYPE(section_vals_type), POINTER                   :: colvar_func_info, colvar_rest, &
     135              :                                                             fixd_restr_rest, hbonds_section
     136        10726 :       TYPE(vsite_constraint_type), DIMENSION(:), POINTER :: vsite_list
     137              : 
     138        10726 :       NULLIFY (logger, constr_x_mol, constr_x_glob)
     139        21452 :       logger => cp_get_default_logger()
     140              :       iw = cp_print_key_unit_nr(logger, subsys_section, "PRINT%TOPOLOGY_INFO/UTIL_INFO", &
     141        10726 :                                 extension=".subsysLog")
     142        10726 :       CALL timeset(routineN, handle)
     143        10726 :       CALL timeset(routineN//"_1", handle2)
     144              : 
     145        10726 :       cons_info => topology%cons_info
     146              :       hbonds_section => section_vals_get_subs_vals(input_file, &
     147        10726 :                                                    "MOTION%CONSTRAINT%HBONDS")
     148              :       fixd_restr_rest => section_vals_get_subs_vals(input_file, &
     149        10726 :                                                     "MOTION%CONSTRAINT%FIX_ATOM_RESTART")
     150        10726 :       CALL section_vals_get(fixd_restr_rest, explicit=restart_restraint_pos)
     151              :       colvar_rest => section_vals_get_subs_vals(input_file, &
     152        10726 :                                                 "MOTION%CONSTRAINT%COLVAR_RESTART")
     153        10726 :       CALL section_vals_get(colvar_rest, explicit=restart_restraint_clv)
     154              :       colvar_func_info => section_vals_get_subs_vals(subsys_section, &
     155        10726 :                                                      "COLVAR%COLVAR_FUNC_INFO")
     156        10726 :       CALL section_vals_get(colvar_func_info, explicit=use_clv_info)
     157              :       !-----------------------------------------------------------------------------
     158              :       !-----------------------------------------------------------------------------
     159              :       ! 1. NULLIFY the molecule_set(imol)%lci via set_molecule_set
     160              :       !-----------------------------------------------------------------------------
     161       321606 :       DO i = 1, topology%nmol
     162       310880 :          molecule => molecule_set(i)
     163       310880 :          NULLIFY (lci)
     164              :          ! only allocate the lci if constraints are active. Can this stuff be distributed ?
     165              :          IF (topology%const_atom .OR. topology%const_hydr .OR. &
     166              :              topology%const_33 .OR. topology%const_46 .OR. &
     167       310880 :              topology%const_colv .OR. topology%const_vsite) THEN
     168        43692 :             ALLOCATE (lci)
     169        43692 :             NULLIFY (lci%lcolv)
     170        43692 :             NULLIFY (lci%lg3x3)
     171        43692 :             NULLIFY (lci%lg4x6)
     172              :          END IF
     173       321606 :          CALL set_molecule(molecule, lci=lci)
     174              :       END DO
     175        10726 :       ALLOCATE (gci)
     176              :       NULLIFY (gci%lcolv, &
     177        10726 :                gci%lg3x3, &
     178        10726 :                gci%lg4x6, &
     179        10726 :                gci%fixd_list, &
     180        10726 :                gci%colv_list, &
     181        10726 :                gci%g3x3_list, &
     182        10726 :                gci%g4x6_list, &
     183        10726 :                gci%vsite_list)
     184        10726 :       gci%ntot = 0
     185        10726 :       gci%ng3x3 = 0
     186        10726 :       gci%ng4x6 = 0
     187        10726 :       gci%nvsite = 0
     188        10726 :       gci%ng3x3_restraint = 0
     189        10726 :       gci%ng4x6_restraint = 0
     190        10726 :       gci%nvsite_restraint = 0
     191        10726 :       CALL setup_colvar_counters(gci%colv_list, gci%ncolv)
     192              :       gci%nrestraint = gci%ng3x3_restraint + &
     193              :                        gci%ng4x6_restraint + &
     194              :                        gci%nvsite_restraint + &
     195        10726 :                        gci%ncolv%nrestraint
     196        10726 :       CALL timestop(handle2)
     197        10726 :       CALL timeset(routineN//"_2", handle2)
     198              :       !-----------------------------------------------------------------------------
     199              :       !-----------------------------------------------------------------------------
     200              :       ! 2. Add more stuff to COLVAR constraint if constraint hydrogen is on
     201              :       !-----------------------------------------------------------------------------
     202        10726 :       IF (topology%const_hydr) THEN
     203           16 :          topology%const_colv = .TRUE.
     204           16 :          NULLIFY (atom_typeh, hdist)
     205           98 :          ALLOCATE (constr_x_mol(SIZE(molecule_kind_set)))
     206           66 :          DO i = 1, SIZE(molecule_kind_set)
     207           50 :             ALLOCATE (constr_x_mol(i)%constr(1))
     208           66 :             constr_x_mol(i)%constr(1) = 1
     209              :          END DO
     210           16 :          CALL section_vals_val_get(hbonds_section, "MOLECULE", n_rep_val=nrep)
     211           16 :          IF (nrep /= 0) THEN
     212            4 :             NULLIFY (inds)
     213           36 :             DO i = 1, SIZE(molecule_kind_set)
     214           36 :                constr_x_mol(i)%constr(1) = 0
     215              :             END DO
     216            4 :             CALL section_vals_val_get(hbonds_section, "MOLECULE", i_vals=inds)
     217           32 :             DO i = 1, SIZE(inds)
     218           32 :                constr_x_mol(inds(i))%constr(1) = 1
     219              :             END DO
     220              :          ELSE
     221           12 :             CALL section_vals_val_get(hbonds_section, "MOLNAME", n_rep_val=nrep)
     222           12 :             IF (nrep /= 0) THEN
     223            2 :                NULLIFY (cnds)
     224           10 :                DO i = 1, SIZE(molecule_kind_set)
     225           10 :                   constr_x_mol(i)%constr(1) = 0
     226              :                END DO
     227            2 :                CALL section_vals_val_get(hbonds_section, "MOLNAME", c_vals=cnds)
     228            4 :                DO i = 1, SIZE(cnds)
     229            2 :                   found_molname = .FALSE.
     230           10 :                   DO k = 1, SIZE(molecule_kind_set)
     231            8 :                      molecule_kind => molecule_kind_set(k)
     232            8 :                      name = molecule_kind%name
     233            8 :                      ldummy = qmmm_ff_precond_only_qm(id1=name)
     234           10 :                      IF (cnds(i) == name) THEN
     235            4 :                         constr_x_mol(k)%constr(1) = 1
     236            4 :                         found_molname = .TRUE.
     237              :                      END IF
     238              :                   END DO
     239            4 :                   CALL print_warning_molname(found_molname, cnds(i))
     240              :                END DO
     241              :             END IF
     242              :          END IF
     243           16 :          CALL section_vals_val_get(hbonds_section, "ATOM_TYPE", n_rep_val=nrep)
     244           16 :          IF (nrep /= 0) &
     245            8 :             CALL section_vals_val_get(hbonds_section, "ATOM_TYPE", c_vals=atom_typeh)
     246           16 :          CALL section_vals_val_get(hbonds_section, "TARGETS", n_rep_val=nrep)
     247           16 :          IF (nrep /= 0) &
     248            8 :             CALL section_vals_val_get(hbonds_section, "TARGETS", r_vals=hdist)
     249           16 :          IF (ASSOCIATED(hdist)) THEN
     250            8 :             CPASSERT(SIZE(hdist) == SIZE(atom_typeh))
     251              :          END IF
     252           16 :          CALL section_vals_val_get(hbonds_section, "exclude_qm", l_val=exclude_qm)
     253           16 :          CALL section_vals_val_get(hbonds_section, "exclude_mm", l_val=exclude_mm)
     254           16 :          nhdist = 0
     255           66 :          DO i = 1, SIZE(molecule_kind_set)
     256           50 :             molecule_kind => molecule_kind_set(i)
     257           50 :             IF (constr_x_mol(i)%constr(1) == 0) CYCLE
     258              :             CALL get_molecule_kind(molecule_kind=molecule_kind, &
     259              :                                    bond_list=bond_list, nbond=nbond, atom_list=atom_list, &
     260           42 :                                    molecule_list=molecule_list)
     261              :             ! Let's tag all requested atoms involving Hydrogen
     262              :             ! on the first molecule of this kind
     263           42 :             molecule => molecule_set(molecule_list(1))
     264           42 :             CALL get_molecule(molecule, first_atom=first_atom, last_atom=last_atom)
     265           42 :             natom = last_atom - first_atom + 1
     266          464 :             DO k = 1, nbond
     267          364 :                ishbond = .FALSE.
     268          364 :                j = bond_list(k)%a
     269          364 :                IF (j < 1 .OR. j > natom) CYCLE
     270          364 :                atomic_kind => atom_list(j)%atomic_kind
     271          364 :                CALL get_atomic_kind(atomic_kind=atomic_kind, name=name)
     272          364 :                is_qm = qmmm_ff_precond_only_qm(id1=name)
     273          364 :                IF ((name(1:1) == "H") .OR. is_hydrogen(atomic_kind)) ishbond = .TRUE.
     274          364 :                IF (is_qm .AND. exclude_qm) ishbond = .FALSE.
     275          292 :                IF (.NOT. (is_qm) .AND. exclude_mm) ishbond = .FALSE.
     276          332 :                IF (.NOT. ishbond) THEN
     277          364 :                   j = bond_list(k)%b
     278          364 :                   IF (j < 1 .OR. j > natom) CYCLE
     279          344 :                   atomic_kind => atom_list(j)%atomic_kind
     280          344 :                   CALL get_atomic_kind(atomic_kind=atomic_kind, name=name)
     281          344 :                   is_qm = qmmm_ff_precond_only_qm(id1=name)
     282          344 :                   IF ((name(1:1) == "H") .OR. is_hydrogen(atomic_kind)) ishbond = .TRUE.
     283          344 :                   IF (is_qm .AND. exclude_qm) ishbond = .FALSE.
     284          288 :                   IF (.NOT. (is_qm) .AND. exclude_mm) ishbond = .FALSE.
     285              :                END IF
     286          354 :                IF (ishbond) THEN
     287          180 :                   nhdist = nhdist + 1
     288              :                END IF
     289              :             END DO
     290              :          END DO
     291           16 :          n_start_colv = cons_info%nconst_colv
     292           16 :          cons_info%nconst_colv = nhdist + n_start_colv
     293           16 :          CALL reallocate(cons_info%const_colv_mol, 1, cons_info%nconst_colv)
     294           16 :          CALL reallocate(cons_info%const_colv_molname, 1, cons_info%nconst_colv)
     295           16 :          CALL reallocate(cons_info%const_colv_target, 1, cons_info%nconst_colv)
     296           16 :          CALL reallocate(cons_info%const_colv_target_growth, 1, cons_info%nconst_colv)
     297           16 :          CALL colvar_p_reallocate(cons_info%colvar_set, 1, cons_info%nconst_colv)
     298              :          ! Fill in Restraints info
     299           16 :          CALL reallocate(cons_info%colv_intermolecular, 1, cons_info%nconst_colv)
     300           16 :          CALL reallocate(cons_info%colv_restraint, 1, cons_info%nconst_colv)
     301           16 :          CALL reallocate(cons_info%colv_k0, 1, cons_info%nconst_colv)
     302           16 :          CALL reallocate(cons_info%colv_exclude_qm, 1, cons_info%nconst_colv)
     303           16 :          CALL reallocate(cons_info%colv_exclude_mm, 1, cons_info%nconst_colv)
     304              :          ! Bonds involving hydrogens are by their nature only intramolecular
     305          196 :          cons_info%colv_intermolecular(n_start_colv + 1:cons_info%nconst_colv) = .FALSE.
     306          196 :          cons_info%colv_exclude_qm(n_start_colv + 1:cons_info%nconst_colv) = .FALSE.
     307          196 :          cons_info%colv_exclude_mm(n_start_colv + 1:cons_info%nconst_colv) = .FALSE.
     308          196 :          cons_info%colv_restraint(n_start_colv + 1:cons_info%nconst_colv) = cons_info%hbonds_restraint
     309          196 :          cons_info%colv_k0(n_start_colv + 1:cons_info%nconst_colv) = cons_info%hbonds_k0
     310              :          !
     311           16 :          nhdist = 0
     312           66 :          DO i = 1, SIZE(molecule_kind_set)
     313           50 :             IF (constr_x_mol(i)%constr(1) == 0) CYCLE
     314           42 :             molecule_kind => molecule_kind_set(i)
     315              :             CALL get_molecule_kind(molecule_kind=molecule_kind, &
     316              :                                    bond_list=bond_list, nbond=nbond, atom_list=atom_list, &
     317           42 :                                    molecule_list=molecule_list)
     318           42 :             molecule => molecule_set(molecule_list(1))
     319           42 :             CALL get_molecule(molecule, first_atom=first_atom, last_atom=last_atom)
     320           42 :             natom = last_atom - first_atom + 1
     321           42 :             offset = first_atom - 1
     322          464 :             DO k = 1, nbond
     323          364 :                ishbond = .FALSE.
     324          364 :                j = bond_list(k)%a
     325          364 :                IF (j < 1 .OR. j > natom) CYCLE
     326          364 :                atomic_kind => atom_list(j)%atomic_kind
     327          364 :                CALL get_atomic_kind(atomic_kind=atomic_kind, name=name)
     328          364 :                is_qm = qmmm_ff_precond_only_qm(id1=name)
     329          364 :                IF ((name(1:1) == "H") .OR. is_hydrogen(atomic_kind)) ishbond = .TRUE.
     330          364 :                IF (is_qm .AND. exclude_qm) ishbond = .FALSE.
     331          292 :                IF (.NOT. (is_qm) .AND. exclude_mm) ishbond = .FALSE.
     332          332 :                IF (.NOT. ishbond) THEN
     333          364 :                   j = bond_list(k)%b
     334          364 :                   IF (j < 1 .OR. j > natom) CYCLE
     335          344 :                   atomic_kind => atom_list(j)%atomic_kind
     336          344 :                   CALL get_atomic_kind(atomic_kind=atomic_kind, name=name)
     337          344 :                   is_qm = qmmm_ff_precond_only_qm(id1=name)
     338          344 :                   IF ((name(1:1) == "H") .OR. is_hydrogen(atomic_kind)) ishbond = .TRUE.
     339          344 :                   IF (is_qm .AND. exclude_qm) ishbond = .FALSE.
     340          288 :                   IF (.NOT. (is_qm) .AND. exclude_mm) ishbond = .FALSE.
     341              :                END IF
     342          354 :                IF (ishbond) THEN
     343          180 :                   nhdist = nhdist + 1
     344          720 :                   rvec = particle_set(offset + bond_list(k)%a)%r - particle_set(offset + bond_list(k)%b)%r
     345          720 :                   rmod = SQRT(DOT_PRODUCT(rvec, rvec))
     346          180 :                   IF (ASSOCIATED(hdist)) THEN
     347           32 :                      IF (SIZE(hdist) > 0) THEN
     348           32 :                         IF (bond_list(k)%a == j) atomic_kind => atom_list(bond_list(k)%b)%atomic_kind
     349           32 :                         IF (bond_list(k)%b == j) atomic_kind => atom_list(bond_list(k)%a)%atomic_kind
     350              :                         CALL get_atomic_kind(atomic_kind=atomic_kind, &
     351           32 :                                              name=name, element_symbol=element_symbol)
     352           32 :                         ldummy = qmmm_ff_precond_only_qm(id1=name)
     353           32 :                         DO m = 1, SIZE(hdist)
     354           32 :                            IF (TRIM(name) == TRIM(atom_typeh(m))) EXIT
     355           32 :                            IF (TRIM(element_symbol) == TRIM(atom_typeh(m))) EXIT
     356              :                         END DO
     357           32 :                         IF (m <= SIZE(hdist)) THEN
     358           32 :                            rmod = hdist(m)
     359              :                         END IF
     360              :                      END IF
     361              :                   END IF
     362          180 :                   cons_info%const_colv_mol(nhdist + n_start_colv) = i
     363          180 :                   cons_info%const_colv_molname(nhdist + n_start_colv) = "UNDEF"
     364          180 :                   cons_info%const_colv_target(nhdist + n_start_colv) = rmod
     365          180 :                   cons_info%const_colv_target_growth(nhdist + n_start_colv) = 0.0_dp
     366              :                   CALL colvar_create(cons_info%colvar_set(nhdist + n_start_colv)%colvar, &
     367          180 :                                      dist_colvar_id)
     368          180 :                   cons_info%colvar_set(nhdist + n_start_colv)%colvar%dist_param%i_at = bond_list(k)%a
     369          180 :                   cons_info%colvar_set(nhdist + n_start_colv)%colvar%dist_param%j_at = bond_list(k)%b
     370          180 :                   CALL colvar_setup(cons_info%colvar_set(nhdist + n_start_colv)%colvar)
     371              :                END IF
     372              :             END DO
     373              :          END DO
     374           66 :          DO j = 1, SIZE(constr_x_mol)
     375           66 :             DEALLOCATE (constr_x_mol(j)%constr)
     376              :          END DO
     377           80 :          DEALLOCATE (constr_x_mol)
     378              :       END IF
     379              : 
     380        10726 :       CALL timestop(handle2)
     381        10726 :       CALL timeset(routineN//"_3", handle2)
     382              :       !-----------------------------------------------------------------------------
     383              :       !-----------------------------------------------------------------------------
     384              :       ! 3. Set the COLVAR constraint molecule_kind_set(ikind)%colv_list
     385              :       !-----------------------------------------------------------------------------
     386        10726 :       IF (topology%const_colv) THEN
     387              :          ! Post Process of COLVARS..
     388          586 :          DO ii = 1, SIZE(cons_info%colvar_set)
     389          586 :             CALL post_process_colvar(cons_info%colvar_set(ii)%colvar, particle_set)
     390              :          END DO
     391              :          ! Real constraint/restraint part..
     392              :          CALL give_constraint_array(cons_info%const_colv_mol, &
     393              :                                     cons_info%const_colv_molname, &
     394              :                                     cons_info%colv_intermolecular, &
     395              :                                     constr_x_mol, &
     396              :                                     constr_x_glob, &
     397              :                                     molecule_kind_set, &
     398              :                                     cons_info%colv_exclude_qm, &
     399          136 :                                     cons_info%colv_exclude_mm)
     400              :          ! Intramolecular constraints
     401          136 :          gind = 0
     402          136 :          cind = 0
     403          714 :          DO ii = 1, SIZE(molecule_kind_set)
     404          578 :             molecule_kind => molecule_kind_set(ii)
     405              :             CALL get_molecule_kind(molecule_kind=molecule_kind, &
     406          578 :                                    nmolecule=nmolecule, molecule_list=molecule_list)
     407          578 :             ncolv_mol = SIZE(constr_x_mol(ii)%constr)
     408         1660 :             ALLOCATE (colv_list(ncolv_mol))
     409              :             ! Starting index of the first molecule of this kind.
     410              :             ! We need the index if no target is provided in the input file
     411              :             ! for the collective variable.. The target will be computed on the
     412              :             ! first molecule of the kind...
     413          578 :             molecule => molecule_set(molecule_list(1))
     414          578 :             CALL get_molecule(molecule, first_atom=first_atom)
     415              :             CALL setup_colv_list(colv_list, constr_x_mol(ii)%constr, gind, &
     416              :                                  cons_info, topology, particle_set, restart_restraint_clv, &
     417          578 :                                  colvar_rest, first_atom)
     418          578 :             CALL setup_colvar_counters(colv_list, ncolv)
     419          578 :             CALL set_molecule_kind(molecule_kind, colv_list=colv_list, ncolv=ncolv)
     420         3978 :             DO j = 1, nmolecule
     421         2108 :                molecule => molecule_set(molecule_list(j))
     422         2108 :                CALL get_molecule(molecule, first_atom=first_atom, last_atom=last_atom)
     423         7014 :                ALLOCATE (lcolv(ncolv_mol))
     424              :                CALL setup_lcolv(lcolv, constr_x_mol(ii)%constr, first_atom, last_atom, &
     425         2108 :                                 cons_info, particle_set, colvar_func_info, use_clv_info, cind)
     426         2686 :                CALL set_molecule(molecule=molecule, lcolv=lcolv)
     427              :             END DO
     428              :          END DO
     429          714 :          DO j = 1, SIZE(constr_x_mol)
     430          714 :             DEALLOCATE (constr_x_mol(j)%constr)
     431              :          END DO
     432          136 :          DEALLOCATE (constr_x_mol)
     433              :          ! Intermolecular constraints
     434          136 :          ncolv_glob = 0
     435          136 :          IF (ASSOCIATED(constr_x_glob)) THEN
     436           44 :             ncolv_glob = SIZE(constr_x_glob)
     437          198 :             ALLOCATE (colv_list(ncolv_glob))
     438              :             CALL setup_colv_list(colv_list, constr_x_glob, gind, cons_info, &
     439              :                                  topology, particle_set, restart_restraint_clv, colvar_rest, &
     440           44 :                                  first_atom=1)
     441           44 :             CALL setup_colvar_counters(colv_list, ncolv)
     442          198 :             ALLOCATE (lcolv(ncolv_glob))
     443              :             CALL setup_lcolv(lcolv, constr_x_glob, 1, SIZE(particle_set), cons_info, &
     444           44 :                              particle_set, colvar_func_info, use_clv_info, cind)
     445           44 :             gci%colv_list => colv_list
     446           44 :             gci%lcolv => lcolv
     447           44 :             gci%ncolv = ncolv
     448              :             ! Total number of Intermolecular constraints
     449           44 :             gci%ntot = gci%ncolv%ntot + gci%ntot
     450           88 :             DEALLOCATE (constr_x_glob)
     451              :          END IF
     452              :       END IF
     453              : 
     454        10726 :       CALL timestop(handle2)
     455        10726 :       CALL timeset(routineN//"_4", handle2)
     456              :       !-----------------------------------------------------------------------------
     457              :       !-----------------------------------------------------------------------------
     458              :       ! 4. Set the group 3x3 constraint g3x3_list
     459              :       !-----------------------------------------------------------------------------
     460        10726 :       IF (topology%const_33) THEN
     461              :          CALL give_constraint_array(cons_info%const_g33_mol, &
     462              :                                     cons_info%const_g33_molname, &
     463              :                                     cons_info%g33_intermolecular, &
     464              :                                     constr_x_mol, &
     465              :                                     constr_x_glob, &
     466              :                                     molecule_kind_set, &
     467              :                                     cons_info%g33_exclude_qm, &
     468          156 :                                     cons_info%g33_exclude_mm)
     469              :          ! Intramolecular constraints
     470          426 :          DO ii = 1, SIZE(molecule_kind_set)
     471          270 :             molecule_kind => molecule_kind_set(ii)
     472              :             CALL get_molecule_kind(molecule_kind=molecule_kind, &
     473              :                                    nmolecule=nmolecule, &
     474          270 :                                    molecule_list=molecule_list)
     475          270 :             ng3x3 = SIZE(constr_x_mol(ii)%constr)
     476          852 :             ALLOCATE (g3x3_list(ng3x3))
     477          270 :             CALL setup_g3x3_list(g3x3_list, constr_x_mol(ii)%constr, cons_info, ng3x3_restraint)
     478          270 :             CALL set_molecule_kind(molecule_kind, ng3x3=ng3x3, ng3x3_restraint=ng3x3_restraint, g3x3_list=g3x3_list)
     479        37320 :             DO j = 1, nmolecule
     480        36354 :                molecule => molecule_set(molecule_list(j))
     481        36354 :                CALL get_molecule(molecule, first_atom=first_atom, last_atom=last_atom)
     482      2524548 :                ALLOCATE (lg3x3(ng3x3))
     483        36354 :                CALL setup_lg3x3(lg3x3, g3x3_list, first_atom, last_atom)
     484        36624 :                CALL set_molecule(molecule=molecule, lg3x3=lg3x3)
     485              :             END DO
     486              :          END DO
     487          426 :          DO j = 1, SIZE(constr_x_mol)
     488          426 :             DEALLOCATE (constr_x_mol(j)%constr)
     489              :          END DO
     490          156 :          DEALLOCATE (constr_x_mol)
     491              :          ! Intermolecular constraints
     492          156 :          IF (ASSOCIATED(constr_x_glob)) THEN
     493            4 :             ng3x3 = SIZE(constr_x_glob)
     494           16 :             ALLOCATE (g3x3_list(ng3x3))
     495            4 :             CALL setup_g3x3_list(g3x3_list, constr_x_glob, cons_info, ng3x3_restraint)
     496          280 :             ALLOCATE (lg3x3(ng3x3))
     497            4 :             CALL setup_lg3x3(lg3x3, g3x3_list, first_atom, last_atom)
     498            4 :             gci%g3x3_list => g3x3_list
     499            4 :             gci%lg3x3 => lg3x3
     500            4 :             gci%ng3x3 = ng3x3
     501            4 :             gci%ng3x3_restraint = ng3x3_restraint
     502              :             ! Total number of Intermolecular constraints
     503            4 :             gci%ntot = 3*gci%ng3x3 + gci%ntot
     504            8 :             DEALLOCATE (constr_x_glob)
     505              :          END IF
     506              :       END IF
     507              : 
     508        10726 :       CALL timestop(handle2)
     509        10726 :       CALL timeset(routineN//"_5", handle2)
     510              :       !-----------------------------------------------------------------------------
     511              :       !-----------------------------------------------------------------------------
     512              :       ! 5. Set the group 4x6 constraint g4x6_list
     513              :       !-----------------------------------------------------------------------------
     514        10726 :       IF (topology%const_46) THEN
     515              :          CALL give_constraint_array(cons_info%const_g46_mol, &
     516              :                                     cons_info%const_g46_molname, &
     517              :                                     cons_info%g46_intermolecular, &
     518              :                                     constr_x_mol, &
     519              :                                     constr_x_glob, &
     520              :                                     molecule_kind_set, &
     521              :                                     cons_info%g46_exclude_qm, &
     522           16 :                                     cons_info%g46_exclude_mm)
     523              :          ! Intramolecular constraints
     524           36 :          DO ii = 1, SIZE(molecule_kind_set)
     525           20 :             molecule_kind => molecule_kind_set(ii)
     526              :             CALL get_molecule_kind(molecule_kind=molecule_kind, &
     527           20 :                                    nmolecule=nmolecule, molecule_list=molecule_list)
     528           20 :             ng4x6 = SIZE(constr_x_mol(ii)%constr)
     529           64 :             ALLOCATE (g4x6_list(ng4x6))
     530           20 :             CALL setup_g4x6_list(g4x6_list, constr_x_mol(ii)%constr, cons_info, ng4x6_restraint)
     531           20 :             CALL set_molecule_kind(molecule_kind, ng4x6=ng4x6, ng4x6_restraint=ng4x6_restraint, g4x6_list=g4x6_list)
     532          726 :             DO j = 1, nmolecule
     533          650 :                molecule => molecule_set(molecule_list(j))
     534          650 :                CALL get_molecule(molecule, first_atom=first_atom, last_atom=last_atom)
     535        99580 :                ALLOCATE (lg4x6(ng4x6))
     536          650 :                CALL setup_lg4x6(lg4x6, g4x6_list, first_atom, last_atom)
     537          670 :                CALL set_molecule(molecule=molecule, lg4x6=lg4x6)
     538              :             END DO
     539              :          END DO
     540           36 :          DO j = 1, SIZE(constr_x_mol)
     541           36 :             DEALLOCATE (constr_x_mol(j)%constr)
     542              :          END DO
     543           16 :          DEALLOCATE (constr_x_mol)
     544              :          ! Intermolecular constraints
     545           16 :          IF (ASSOCIATED(constr_x_glob)) THEN
     546            4 :             ng4x6 = SIZE(constr_x_glob)
     547           16 :             ALLOCATE (g4x6_list(ng4x6))
     548            4 :             CALL setup_g4x6_list(g4x6_list, constr_x_glob, cons_info, ng4x6_restraint)
     549          616 :             ALLOCATE (lg4x6(ng4x6))
     550            4 :             CALL setup_lg4x6(lg4x6, g4x6_list, first_atom, last_atom)
     551            4 :             gci%g4x6_list => g4x6_list
     552            4 :             gci%lg4x6 => lg4x6
     553            4 :             gci%ng4x6 = ng4x6
     554            4 :             gci%ng4x6_restraint = ng4x6_restraint
     555              :             ! Total number of Intermolecular constraints
     556            4 :             gci%ntot = 6*gci%ng4x6 + gci%ntot
     557            8 :             DEALLOCATE (constr_x_glob)
     558              :          END IF
     559              :       END IF
     560              : 
     561        10726 :       CALL timestop(handle2)
     562        10726 :       CALL timeset(routineN//"_6", handle2)
     563              :       !-----------------------------------------------------------------------------
     564              :       !-----------------------------------------------------------------------------
     565              :       ! 6. Set the group vsite constraint vsite_list
     566              :       !-----------------------------------------------------------------------------
     567        10726 :       IF (topology%const_vsite) THEN
     568              :          CALL give_constraint_array(cons_info%const_vsite_mol, &
     569              :                                     cons_info%const_vsite_molname, &
     570              :                                     cons_info%vsite_intermolecular, &
     571              :                                     constr_x_mol, &
     572              :                                     constr_x_glob, &
     573              :                                     molecule_kind_set, &
     574              :                                     cons_info%vsite_exclude_qm, &
     575            8 :                                     cons_info%vsite_exclude_mm)
     576              :          ! Intramolecular constraints
     577           18 :          DO ii = 1, SIZE(molecule_kind_set)
     578           10 :             molecule_kind => molecule_kind_set(ii)
     579              :             CALL get_molecule_kind(molecule_kind=molecule_kind, &
     580           10 :                                    nmolecule=nmolecule, molecule_list=molecule_list)
     581           10 :             nvsite = SIZE(constr_x_mol(ii)%constr)
     582           36 :             ALLOCATE (vsite_list(nvsite))
     583           10 :             CALL setup_vsite_list(vsite_list, constr_x_mol(ii)%constr, cons_info, nvsite_restraint)
     584              :             CALL set_molecule_kind(molecule_kind, nvsite=nvsite, nvsite_restraint=nvsite_restraint, &
     585           28 :                                    vsite_list=vsite_list)
     586              :          END DO
     587           18 :          DO j = 1, SIZE(constr_x_mol)
     588           18 :             DEALLOCATE (constr_x_mol(j)%constr)
     589              :          END DO
     590            8 :          DEALLOCATE (constr_x_mol)
     591              :          ! Intermolecular constraints
     592            8 :          IF (ASSOCIATED(constr_x_glob)) THEN
     593            0 :             nvsite = SIZE(constr_x_glob)
     594            0 :             ALLOCATE (vsite_list(nvsite))
     595            0 :             CALL setup_vsite_list(vsite_list, constr_x_glob, cons_info, nvsite_restraint)
     596            0 :             gci%vsite_list => vsite_list
     597            0 :             gci%nvsite = nvsite
     598            0 :             gci%nvsite_restraint = nvsite_restraint
     599              :             ! Total number of Intermolecular constraints
     600            0 :             gci%ntot = gci%nvsite + gci%ntot
     601            0 :             DEALLOCATE (constr_x_glob)
     602              :          END IF
     603              :       END IF
     604        10726 :       CALL timestop(handle2)
     605        10726 :       CALL timeset(routineN//"_7", handle2)
     606              :       !-----------------------------------------------------------------------------
     607              :       !-----------------------------------------------------------------------------
     608              :       ! 7. Set the group fixed_atom constraint fixd_list
     609              :       !-----------------------------------------------------------------------------
     610        10726 :       IF (topology%const_atom) THEN
     611        30574 :          ALLOCATE (fixd_list_gci(SIZE(particle_set)))
     612          110 :          nfixd_list_gci = 0
     613          226 :          ALLOCATE (missed_molname(SIZE(cons_info%fixed_molnames, 1)))
     614          116 :          missed_molname = .TRUE.
     615          110 :          nfixd_restart = 0
     616         5036 :          DO i = 1, SIZE(molecule_kind_set)
     617         4926 :             molecule_kind => molecule_kind_set(i)
     618              :             CALL get_molecule_kind(molecule_kind=molecule_kind, &
     619         4926 :                                    nmolecule=nmolecule, molecule_list=molecule_list, name=molname)
     620         4926 :             is_qm = qmmm_ff_precond_only_qm(id1=molname)
     621         4938 :             WHERE (molname == cons_info%fixed_molnames)
     622              :                missed_molname = .FALSE.
     623              :             END WHERE
     624              :             ! Try to figure out how many atoms of the list belong to this molecule_kind
     625         4926 :             nfixed_atoms = 0
     626        17634 :             DO j = 1, nmolecule
     627        12708 :                molecule => molecule_set(molecule_list(j))
     628        12708 :                CALL get_molecule(molecule, first_atom=first, last_atom=last)
     629        12708 :                fix_atom_molname = .FALSE.
     630        12708 :                IF (ASSOCIATED(cons_info%fixed_molnames)) THEN
     631        14274 :                   DO k = 1, SIZE(cons_info%fixed_molnames)
     632        14274 :                      IF (cons_info%fixed_molnames(k) == molname) THEN
     633           48 :                         fix_atom_molname = .TRUE.
     634           48 :                         IF (is_qm .AND. cons_info%fixed_exclude_qm(k)) fix_atom_molname = .FALSE.
     635           44 :                         IF ((.NOT. is_qm) .AND. cons_info%fixed_exclude_mm(k)) fix_atom_molname = .FALSE.
     636              :                      END IF
     637              :                   END DO
     638              :                END IF
     639        47548 :                DO k = first, last
     640        29914 :                   fix_atom_qmmm = .FALSE.
     641        29914 :                   IF (PRESENT(qmmm_env)) THEN
     642          324 :                      SELECT CASE (cons_info%freeze_qm)
     643              :                      CASE (do_constr_atomic)
     644            0 :                         IF (ANY(qmmm_env%qm_atom_index == k)) fix_atom_qmmm = .TRUE.
     645              :                      CASE (do_constr_molec)
     646          336 :                         IF (ANY(qmmm_env%qm_molecule_index == molecule_list(j))) fix_atom_qmmm = .TRUE.
     647              :                      END SELECT
     648          394 :                      SELECT CASE (cons_info%freeze_mm)
     649              :                      CASE (do_constr_atomic)
     650          840 :                         IF (ALL(qmmm_env%qm_atom_index /= k)) fix_atom_qmmm = .TRUE.
     651              :                      CASE (do_constr_molec)
     652          408 :                         IF (ALL(qmmm_env%qm_molecule_index /= molecule_list(j))) fix_atom_qmmm = .TRUE.
     653              :                      END SELECT
     654              :                   END IF
     655      3861838 :                   IF (ANY(cons_info%fixed_atoms == k) .OR. fix_atom_qmmm .OR. fix_atom_molname) THEN
     656        10196 :                      nfixed_atoms = nfixed_atoms + 1
     657              :                   END IF
     658              :                END DO
     659              :             END DO
     660        39006 :             ALLOCATE (fixd_list(nfixed_atoms))
     661         4926 :             kk = 0
     662         4926 :             nfixd_restraint = 0
     663         4926 :             IF (nfixed_atoms /= 0) THEN
     664        10122 :                DO j = 1, nmolecule
     665         5942 :                   molecule => molecule_set(molecule_list(j))
     666         5942 :                   CALL get_molecule(molecule, first_atom=first, last_atom=last)
     667         5942 :                   fix_atom_molname = .FALSE.
     668         5942 :                   IF (ASSOCIATED(cons_info%fixed_molnames)) THEN
     669         5942 :                      DO k1loc = 1, SIZE(cons_info%fixed_molnames)
     670         5942 :                         IF (cons_info%fixed_molnames(k1loc) == molname) THEN
     671           44 :                            fix_atom_molname = .TRUE.
     672           44 :                            itype = cons_info%fixed_mol_type(k1loc)
     673           44 :                            EXIT
     674              :                         END IF
     675              :                      END DO
     676              :                   END IF
     677        21162 :                   DO k = first, last
     678              :                      ! FIXED LIST ATOMS
     679        11040 :                      fix_fixed_atom = .FALSE.
     680      2891634 :                      DO k2loc = 1, SIZE(cons_info%fixed_atoms)
     681      2891634 :                         IF (cons_info%fixed_atoms(k2loc) == k) THEN
     682        10012 :                            fix_fixed_atom = .TRUE.
     683        10012 :                            itype = cons_info%fixed_type(k2loc)
     684        10012 :                            EXIT
     685              :                         END IF
     686              :                      END DO
     687              :                      ! QMMM FIXED ATOMS (QM OR MM)
     688        11040 :                      fix_atom_qmmm = .FALSE.
     689        11040 :                      fix_atom_mm = .FALSE.
     690        11040 :                      fix_atom_qm = .FALSE.
     691        11040 :                      IF (PRESENT(qmmm_env)) THEN
     692          224 :                         SELECT CASE (cons_info%freeze_qm)
     693              :                         CASE (do_constr_atomic)
     694            0 :                            IF (ANY(qmmm_env%qm_atom_index == k)) THEN
     695            0 :                               fix_atom_qmmm = .TRUE.
     696            0 :                               fix_atom_qm = .TRUE.
     697            0 :                               itype = cons_info%freeze_qm_type
     698              :                            END IF
     699              :                         CASE (do_constr_molec)
     700          224 :                            IF (ANY(qmmm_env%qm_molecule_index == molecule_list(j))) THEN
     701            6 :                               fix_atom_qmmm = .TRUE.
     702            6 :                               fix_atom_qm = .TRUE.
     703            6 :                               itype = cons_info%freeze_qm_type
     704              :                            END IF
     705              :                         END SELECT
     706          294 :                         SELECT CASE (cons_info%freeze_mm)
     707              :                         CASE (do_constr_atomic)
     708          840 :                            IF (ALL(qmmm_env%qm_atom_index /= k)) THEN
     709           42 :                               fix_atom_qmmm = .TRUE.
     710           42 :                               fix_atom_mm = .TRUE.
     711           42 :                               itype = cons_info%freeze_mm_type
     712              :                            END IF
     713              :                         CASE (do_constr_molec)
     714          308 :                            IF (ALL(qmmm_env%qm_molecule_index /= molecule_list(j))) THEN
     715           84 :                               fix_atom_qmmm = .TRUE.
     716           84 :                               fix_atom_mm = .TRUE.
     717           84 :                               itype = cons_info%freeze_mm_type
     718              :                            END IF
     719              :                         END SELECT
     720              :                         ! We should never reach this point but let's check it anyway
     721          126 :                         IF (fix_atom_qm .AND. fix_atom_mm) THEN
     722              :                            CALL cp_abort(__LOCATION__, &
     723              :                                          "Atom number: "//cp_to_string(k)// &
     724            0 :                                          " has been defined both QM and MM. General Error!")
     725              :                         END IF
     726              :                      END IF
     727              :                      ! Check that the fixed atom constraint/restraint is unique
     728              :                      IF ((fix_fixed_atom .AND. fix_atom_qmmm) .OR. (fix_fixed_atom .AND. fix_atom_molname) &
     729        11040 :                          .OR. (fix_atom_qmmm .AND. fix_atom_molname)) THEN
     730              :                         CALL cp_abort(__LOCATION__, &
     731              :                                       "Atom number: "//cp_to_string(k)// &
     732              :                                       " has been constrained/restrained to be fixed in more than one"// &
     733            0 :                                       " input section. Check and correct your input file!")
     734              :                      END IF
     735              :                      ! Let's store the atom index
     736        16982 :                      IF (fix_fixed_atom .OR. fix_atom_qmmm .OR. fix_atom_molname) THEN
     737        10196 :                         IF (ASSOCIATED(topology%cell_muc)) THEN
     738        10196 :                            IF (topology%cell_muc%input_cell_canonicalized .AND. itype /= use_perd_xyz) THEN
     739              :                               CALL cp_abort(__LOCATION__, &
     740              :                                             "Partial FIXED_ATOMS components cannot be transformed "// &
     741              :                                             "after CELL%CANONICALIZE. Use COMPONENTS_TO_FIX XYZ or "// &
     742            0 :                                             "disable CELL%CANONICALIZE for this input.")
     743              :                            END IF
     744              :                         END IF
     745        10196 :                         kk = kk + 1
     746        10196 :                         fixd_list(kk)%fixd = k
     747        71372 :                         fixd_list(kk)%coord = particle_set(k)%r
     748        10196 :                         fixd_list(kk)%itype = itype
     749              :                         ! Possibly Restraint
     750        10196 :                         IF (fix_fixed_atom) THEN
     751        10012 :                            fixd_list(kk)%restraint%active = cons_info%fixed_restraint(k2loc)
     752        10012 :                            fixd_list(kk)%restraint%k0 = cons_info%fixed_k0(k2loc)
     753          184 :                         ELSEIF (fix_atom_qm) THEN
     754            6 :                            fixd_list(kk)%restraint%active = cons_info%fixed_qm_restraint
     755            6 :                            fixd_list(kk)%restraint%k0 = cons_info%fixed_qm_k0
     756          178 :                         ELSEIF (fix_atom_mm) THEN
     757          126 :                            fixd_list(kk)%restraint%active = cons_info%fixed_mm_restraint
     758          126 :                            fixd_list(kk)%restraint%k0 = cons_info%fixed_mm_k0
     759           52 :                         ELSEIF (fix_atom_molname) THEN
     760           52 :                            fixd_list(kk)%restraint%active = cons_info%fixed_mol_restraint(k1loc)
     761           52 :                            fixd_list(kk)%restraint%k0 = cons_info%fixed_mol_k0(k1loc)
     762              :                         ELSE
     763              :                            ! Should never reach this point
     764            0 :                            CPABORT("")
     765              :                         END IF
     766        10196 :                         IF (fixd_list(kk)%restraint%active) THEN
     767           38 :                            nfixd_restraint = nfixd_restraint + 1
     768           38 :                            nfixd_restart = nfixd_restart + 1
     769              :                            ! Check that we use the components that we really want..
     770            0 :                            SELECT CASE (itype)
     771              :                            CASE (use_perd_x)
     772            0 :                               fixd_list(kk)%coord(2) = HUGE(0.0_dp)
     773            0 :                               fixd_list(kk)%coord(3) = HUGE(0.0_dp)
     774              :                            CASE (use_perd_y)
     775            0 :                               fixd_list(kk)%coord(1) = HUGE(0.0_dp)
     776            0 :                               fixd_list(kk)%coord(3) = HUGE(0.0_dp)
     777              :                            CASE (use_perd_z)
     778            0 :                               fixd_list(kk)%coord(1) = HUGE(0.0_dp)
     779            0 :                               fixd_list(kk)%coord(2) = HUGE(0.0_dp)
     780              :                            CASE (use_perd_xy)
     781            0 :                               fixd_list(kk)%coord(3) = HUGE(0.0_dp)
     782              :                            CASE (use_perd_xz)
     783            0 :                               fixd_list(kk)%coord(2) = HUGE(0.0_dp)
     784              :                            CASE (use_perd_yz)
     785           38 :                               fixd_list(kk)%coord(1) = HUGE(0.0_dp)
     786              :                            END SELECT
     787           38 :                            IF (restart_restraint_pos) THEN
     788              :                               ! Read  coord0 value for restraint
     789              :                               CALL section_vals_val_get(fixd_restr_rest, "_DEFAULT_KEYWORD_", &
     790           14 :                                                         i_rep_val=nfixd_restart, r_vals=r)
     791            0 :                               SELECT CASE (itype)
     792              :                               CASE (use_perd_x)
     793            0 :                                  CPASSERT(SIZE(r) == 1)
     794            0 :                                  fixd_list(kk)%coord(1) = r(1)
     795              :                               CASE (use_perd_y)
     796            0 :                                  CPASSERT(SIZE(r) == 1)
     797            0 :                                  fixd_list(kk)%coord(2) = r(1)
     798              :                               CASE (use_perd_z)
     799            0 :                                  CPASSERT(SIZE(r) == 1)
     800            0 :                                  fixd_list(kk)%coord(3) = r(1)
     801              :                               CASE (use_perd_xy)
     802            0 :                                  CPASSERT(SIZE(r) == 2)
     803            0 :                                  fixd_list(kk)%coord(1) = r(1)
     804            0 :                                  fixd_list(kk)%coord(2) = r(2)
     805              :                               CASE (use_perd_xz)
     806            0 :                                  CPASSERT(SIZE(r) == 2)
     807            0 :                                  fixd_list(kk)%coord(1) = r(1)
     808            0 :                                  fixd_list(kk)%coord(3) = r(2)
     809              :                               CASE (use_perd_yz)
     810            0 :                                  CPASSERT(SIZE(r) == 2)
     811            0 :                                  fixd_list(kk)%coord(2) = r(1)
     812            0 :                                  fixd_list(kk)%coord(3) = r(2)
     813              :                               CASE (use_perd_xyz)
     814           14 :                                  CPASSERT(SIZE(r) == 3)
     815           98 :                                  fixd_list(kk)%coord(1:3) = r(1:3)
     816           14 :                                  IF (ASSOCIATED(topology%cell_muc)) &
     817           28 :                                     CALL cell_transform_input_cartesian(topology%cell_muc, fixd_list(kk)%coord)
     818              :                               END SELECT
     819              :                            ELSE
     820              :                               ! Write coord0 value for restraint
     821            0 :                               SELECT CASE (itype)
     822              :                               CASE (use_perd_x)
     823            0 :                                  ALLOCATE (r(1))
     824            0 :                                  r(1) = fixd_list(kk)%coord(1)
     825              :                               CASE (use_perd_y)
     826            0 :                                  ALLOCATE (r(1))
     827            0 :                                  r(1) = fixd_list(kk)%coord(2)
     828              :                               CASE (use_perd_z)
     829            0 :                                  ALLOCATE (r(1))
     830            0 :                                  r(1) = fixd_list(kk)%coord(3)
     831              :                               CASE (use_perd_xy)
     832            0 :                                  ALLOCATE (r(2))
     833            0 :                                  r(1) = fixd_list(kk)%coord(1)
     834            0 :                                  r(2) = fixd_list(kk)%coord(2)
     835              :                               CASE (use_perd_xz)
     836            0 :                                  ALLOCATE (r(2))
     837            0 :                                  r(1) = fixd_list(kk)%coord(1)
     838            0 :                                  r(2) = fixd_list(kk)%coord(3)
     839              :                               CASE (use_perd_yz)
     840            0 :                                  ALLOCATE (r(2))
     841            0 :                                  r(1) = fixd_list(kk)%coord(1)
     842            0 :                                  r(2) = fixd_list(kk)%coord(3)
     843              :                               CASE (use_perd_xyz)
     844           24 :                                  ALLOCATE (r(3))
     845          120 :                                  r(1:3) = fixd_list(kk)%coord(1:3)
     846              :                               END SELECT
     847              :                               CALL section_vals_val_set(fixd_restr_rest, "_DEFAULT_KEYWORD_", &
     848           24 :                                                         i_rep_val=nfixd_restart, r_vals_ptr=r)
     849              :                            END IF
     850              :                         END IF
     851              :                      END IF
     852              :                   END DO
     853              :                END DO
     854              :             END IF
     855         4926 :             IF (iw > 0) THEN
     856            0 :                WRITE (iw, *) "MOLECULE KIND:", i, " NR. FIXED ATOMS:", SIZE(fixd_list(:)%fixd), " LIST::", fixd_list(:)%fixd
     857              :             END IF
     858              :             CALL set_molecule_kind(molecule_kind, nfixd=nfixed_atoms, nfixd_restraint=nfixd_restraint, &
     859         4926 :                                    fixd_list=fixd_list)
     860        25318 :             fixd_list_gci(nfixd_list_gci + 1:nfixd_list_gci + nfixed_atoms) = fixd_list
     861         9962 :             nfixd_list_gci = nfixd_list_gci + nfixed_atoms
     862              :          END DO
     863          110 :          IF (iw > 0) THEN
     864            0 :             WRITE (iw, *) "TOTAL NUMBER OF FIXED ATOMS:", nfixd_list_gci
     865              :          END IF
     866          116 :          CPASSERT(COUNT(missed_molname) == 0)
     867          110 :          DEALLOCATE (missed_molname)
     868              :          ! Intermolecular constraints
     869          110 :          IF (gci%ntot /= 0) THEN
     870           16 :             ALLOCATE (fixd_list(nfixd_list_gci))
     871           10 :             fixd_list(1:nfixd_list_gci) = fixd_list_gci(1:nfixd_list_gci)
     872            2 :             gci%fixd_list => fixd_list
     873              :          END IF
     874          110 :          DEALLOCATE (fixd_list_gci)
     875              :       END IF
     876              :       ! Final setup of the number of possible restraints
     877              :       gci%nrestraint = gci%ng3x3_restraint + &
     878              :                        gci%ng4x6_restraint + &
     879              :                        gci%nvsite_restraint + &
     880        10726 :                        gci%ncolv%nrestraint
     881              :       CALL cp_print_key_finished_output(iw, logger, subsys_section, &
     882        10726 :                                         "PRINT%TOPOLOGY_INFO/UTIL_INFO")
     883        10726 :       CALL timestop(handle2)
     884        10726 :       CALL timestop(handle)
     885        10726 :    END SUBROUTINE topology_constraint_pack
     886              : 
     887              : ! **************************************************************************************************
     888              : !> \brief Setup the colv_list for the packing of constraints
     889              : !> \param colv_list ...
     890              : !> \param ilist ...
     891              : !> \param gind ...
     892              : !> \param cons_info ...
     893              : !> \param topology ...
     894              : !> \param particle_set ...
     895              : !> \param restart_restraint_clv ...
     896              : !> \param colvar_rest ...
     897              : !> \param first_atom ...
     898              : !> \par History
     899              : !>      Updated 2007 for intermolecular constraints
     900              : !> \author Teodoro Laino [2007]
     901              : ! **************************************************************************************************
     902          622 :    SUBROUTINE setup_colv_list(colv_list, ilist, gind, cons_info, topology, &
     903              :                               particle_set, restart_restraint_clv, colvar_rest, first_atom)
     904              : 
     905              :       TYPE(colvar_constraint_type), DIMENSION(:), &
     906              :          POINTER                                         :: colv_list
     907              :       INTEGER, DIMENSION(:), POINTER                     :: ilist
     908              :       INTEGER, INTENT(INOUT)                             :: gind
     909              :       TYPE(constraint_info_type), POINTER                :: cons_info
     910              :       TYPE(topology_parameters_type), INTENT(INOUT)      :: topology
     911              :       TYPE(particle_type), DIMENSION(:), OPTIONAL, &
     912              :          POINTER                                         :: particle_set
     913              :       LOGICAL, INTENT(IN)                                :: restart_restraint_clv
     914              :       TYPE(section_vals_type), POINTER                   :: colvar_rest
     915              :       INTEGER, INTENT(IN)                                :: first_atom
     916              : 
     917              :       INTEGER                                            :: j, kdim, kk, ncolv_mol
     918              :       REAL(KIND=dp)                                      :: rmod
     919              :       TYPE(colvar_type), POINTER                         :: local_colvar
     920              : 
     921          622 :       ncolv_mol = 0
     922         1070 :       DO kk = 1, SIZE(ilist)
     923          448 :          j = ilist(kk)
     924          448 :          ncolv_mol = ncolv_mol + 1
     925          448 :          kdim = SIZE(cons_info%colvar_set(j)%colvar%i_atom)
     926         1344 :          ALLOCATE (colv_list(ncolv_mol)%i_atoms(kdim))
     927          448 :          colv_list(ncolv_mol)%inp_seq_num = j
     928          448 :          colv_list(ncolv_mol)%type_id = cons_info%colvar_set(j)%colvar%type_id
     929         3020 :          colv_list(ncolv_mol)%i_atoms = cons_info%colvar_set(j)%colvar%i_atom
     930          448 :          colv_list(ncolv_mol)%use_points = cons_info%colvar_set(j)%colvar%use_points
     931              :          ! Restraint
     932          448 :          colv_list(ncolv_mol)%restraint%active = cons_info%colv_restraint(j)
     933          448 :          colv_list(ncolv_mol)%restraint%k0 = cons_info%colv_k0(j)
     934          448 :          IF (cons_info%const_colv_target(j) == -HUGE(0.0_dp)) THEN
     935              :             ! Let's compute the value..
     936          100 :             NULLIFY (local_colvar)
     937              :             CALL colvar_clone(local_colvar, cons_info%colvar_set(j)%colvar, &
     938          100 :                               i_atom_offset=first_atom - 1)
     939          100 :             CALL colvar_eval_mol_f(local_colvar, topology%cell, particle_set)
     940          100 :             colv_list(ncolv_mol)%expected_value = local_colvar%ss
     941          100 :             CALL colvar_release(local_colvar)
     942              :          ELSE
     943          348 :             colv_list(ncolv_mol)%expected_value = cons_info%const_colv_target(j)
     944              :          END IF
     945          448 :          colv_list(ncolv_mol)%expected_value_growth_speed = cons_info%const_colv_target_growth(j)
     946              :          ! In case of Restraint let's check for possible restart values
     947          448 :          IF (colv_list(ncolv_mol)%restraint%active .AND. &
     948              :              (colv_list(ncolv_mol)%expected_value_growth_speed == 0.0_dp)) THEN
     949           96 :             gind = gind + 1
     950           96 :             IF (restart_restraint_clv) THEN
     951              :                CALL section_vals_val_get(colvar_rest, "_DEFAULT_KEYWORD_", &
     952           14 :                                          i_rep_val=gind, r_val=rmod)
     953           14 :                colv_list(ncolv_mol)%expected_value = rmod
     954              :             ELSE
     955           82 :                rmod = colv_list(ncolv_mol)%expected_value
     956              :                CALL section_vals_val_set(colvar_rest, "_DEFAULT_KEYWORD_", &
     957           82 :                                          i_rep_val=gind, r_val=rmod)
     958              :             END IF
     959              :          END IF
     960              :          ! Only if torsion let's take into account the singularity in the definition
     961              :          ! of the dihedral
     962         1070 :          IF (cons_info%colvar_set(j)%colvar%type_id == torsion_colvar_id) THEN
     963           38 :             cons_info%colvar_set(j)%colvar%torsion_param%o0 = colv_list(ncolv_mol)%expected_value
     964              :          END IF
     965              :       END DO
     966          622 :    END SUBROUTINE setup_colv_list
     967              : 
     968              : ! **************************************************************************************************
     969              : !> \brief Setup the g3x3_list for the packing of constraints
     970              : !> \param g3x3_list ...
     971              : !> \param ilist ...
     972              : !> \param cons_info ...
     973              : !> \param ng3x3_restraint ...
     974              : !> \par History
     975              : !>      Updated 2007 for intermolecular constraints
     976              : !> \author Teodoro Laino [2007]
     977              : ! **************************************************************************************************
     978          274 :    SUBROUTINE setup_g3x3_list(g3x3_list, ilist, cons_info, ng3x3_restraint)
     979              :       TYPE(g3x3_constraint_type), DIMENSION(:), POINTER  :: g3x3_list
     980              :       INTEGER, DIMENSION(:), POINTER                     :: ilist
     981              :       TYPE(constraint_info_type), POINTER                :: cons_info
     982              :       INTEGER, INTENT(OUT)                               :: ng3x3_restraint
     983              : 
     984              :       INTEGER                                            :: j, ng3x3
     985              : 
     986          274 :       ng3x3_restraint = 0
     987          434 :       DO ng3x3 = 1, SIZE(ilist)
     988          160 :          j = ilist(ng3x3)
     989          160 :          g3x3_list(ng3x3)%a = cons_info%const_g33_a(j)
     990          160 :          g3x3_list(ng3x3)%b = cons_info%const_g33_b(j)
     991          160 :          g3x3_list(ng3x3)%c = cons_info%const_g33_c(j)
     992          160 :          g3x3_list(ng3x3)%dab = cons_info%const_g33_dab(j)
     993          160 :          g3x3_list(ng3x3)%dac = cons_info%const_g33_dac(j)
     994          160 :          g3x3_list(ng3x3)%dbc = cons_info%const_g33_dbc(j)
     995              :          ! Restraint
     996          160 :          g3x3_list(ng3x3)%restraint%active = cons_info%g33_restraint(j)
     997          160 :          g3x3_list(ng3x3)%restraint%k0 = cons_info%g33_k0(j)
     998          434 :          IF (g3x3_list(ng3x3)%restraint%active) ng3x3_restraint = ng3x3_restraint + 1
     999              :       END DO
    1000              : 
    1001          274 :    END SUBROUTINE setup_g3x3_list
    1002              : 
    1003              : ! **************************************************************************************************
    1004              : !> \brief Setup the g4x6_list for the packing of constraints
    1005              : !> \param g4x6_list ...
    1006              : !> \param ilist ...
    1007              : !> \param cons_info ...
    1008              : !> \param ng4x6_restraint ...
    1009              : !> \par History
    1010              : !>      Updated 2007 for intermolecular constraints
    1011              : !> \author Teodoro Laino [2007]
    1012              : ! **************************************************************************************************
    1013           24 :    SUBROUTINE setup_g4x6_list(g4x6_list, ilist, cons_info, ng4x6_restraint)
    1014              :       TYPE(g4x6_constraint_type), DIMENSION(:), POINTER  :: g4x6_list
    1015              :       INTEGER, DIMENSION(:), POINTER                     :: ilist
    1016              :       TYPE(constraint_info_type), POINTER                :: cons_info
    1017              :       INTEGER, INTENT(OUT)                               :: ng4x6_restraint
    1018              : 
    1019              :       INTEGER                                            :: j, ng4x6
    1020              : 
    1021           24 :       ng4x6 = 0
    1022           24 :       ng4x6_restraint = 0
    1023           40 :       DO ng4x6 = 1, SIZE(ilist)
    1024           16 :          j = ilist(ng4x6)
    1025           16 :          g4x6_list(ng4x6)%a = cons_info%const_g46_a(j)
    1026           16 :          g4x6_list(ng4x6)%b = cons_info%const_g46_b(j)
    1027           16 :          g4x6_list(ng4x6)%c = cons_info%const_g46_c(j)
    1028           16 :          g4x6_list(ng4x6)%d = cons_info%const_g46_d(j)
    1029           16 :          g4x6_list(ng4x6)%dab = cons_info%const_g46_dab(j)
    1030           16 :          g4x6_list(ng4x6)%dac = cons_info%const_g46_dac(j)
    1031           16 :          g4x6_list(ng4x6)%dbc = cons_info%const_g46_dbc(j)
    1032           16 :          g4x6_list(ng4x6)%dad = cons_info%const_g46_dad(j)
    1033           16 :          g4x6_list(ng4x6)%dbd = cons_info%const_g46_dbd(j)
    1034           16 :          g4x6_list(ng4x6)%dcd = cons_info%const_g46_dcd(j)
    1035              :          ! Restraint
    1036           16 :          g4x6_list(ng4x6)%restraint%active = cons_info%g46_restraint(j)
    1037           16 :          g4x6_list(ng4x6)%restraint%k0 = cons_info%g46_k0(j)
    1038           40 :          IF (g4x6_list(ng4x6)%restraint%active) ng4x6_restraint = ng4x6_restraint + 1
    1039              :       END DO
    1040              : 
    1041           24 :    END SUBROUTINE setup_g4x6_list
    1042              : 
    1043              : ! **************************************************************************************************
    1044              : !> \brief Setup the vsite_list for the packing of constraints
    1045              : !> \param vsite_list ...
    1046              : !> \param ilist ...
    1047              : !> \param cons_info ...
    1048              : !> \param nvsite_restraint ...
    1049              : !> \par History
    1050              : !> \author Marcel Baer [2008]
    1051              : ! **************************************************************************************************
    1052           10 :    SUBROUTINE setup_vsite_list(vsite_list, ilist, cons_info, nvsite_restraint)
    1053              :       TYPE(vsite_constraint_type), DIMENSION(:), POINTER :: vsite_list
    1054              :       INTEGER, DIMENSION(:), POINTER                     :: ilist
    1055              :       TYPE(constraint_info_type), POINTER                :: cons_info
    1056              :       INTEGER, INTENT(OUT)                               :: nvsite_restraint
    1057              : 
    1058              :       INTEGER                                            :: j, nvsite
    1059              : 
    1060           10 :       nvsite = 0
    1061           10 :       nvsite_restraint = 0
    1062           18 :       DO nvsite = 1, SIZE(ilist)
    1063            8 :          j = ilist(nvsite)
    1064            8 :          vsite_list(nvsite)%a = cons_info%const_vsite_a(j)
    1065            8 :          vsite_list(nvsite)%b = cons_info%const_vsite_b(j)
    1066            8 :          vsite_list(nvsite)%c = cons_info%const_vsite_c(j)
    1067            8 :          vsite_list(nvsite)%d = cons_info%const_vsite_d(j)
    1068            8 :          vsite_list(nvsite)%wbc = cons_info%const_vsite_wbc(j)
    1069            8 :          vsite_list(nvsite)%wdc = cons_info%const_vsite_wdc(j)
    1070              :          ! Restraint
    1071            8 :          vsite_list(nvsite)%restraint%active = cons_info%vsite_restraint(j)
    1072            8 :          vsite_list(nvsite)%restraint%k0 = cons_info%vsite_k0(j)
    1073           18 :          IF (vsite_list(nvsite)%restraint%active) nvsite_restraint = nvsite_restraint + 1
    1074              :       END DO
    1075              : 
    1076           10 :    END SUBROUTINE setup_vsite_list
    1077              : ! **************************************************************************************************
    1078              : !> \brief Setup the lcolv for the packing of constraints
    1079              : !> \param lcolv ...
    1080              : !> \param ilist ...
    1081              : !> \param first_atom ...
    1082              : !> \param last_atom ...
    1083              : !> \param cons_info ...
    1084              : !> \param particle_set ...
    1085              : !> \param colvar_func_info ...
    1086              : !> \param use_clv_info ...
    1087              : !> \param cind ...
    1088              : !> \par History
    1089              : !>      Updated 2007 for intermolecular constraints
    1090              : !> \author Teodoro Laino [2007]
    1091              : ! **************************************************************************************************
    1092         2152 :    SUBROUTINE setup_lcolv(lcolv, ilist, first_atom, last_atom, cons_info, &
    1093              :                           particle_set, colvar_func_info, use_clv_info, &
    1094              :                           cind)
    1095              :       TYPE(local_colvar_constraint_type), DIMENSION(:), &
    1096              :          POINTER                                         :: lcolv
    1097              :       INTEGER, DIMENSION(:), POINTER                     :: ilist
    1098              :       INTEGER, INTENT(IN)                                :: first_atom, last_atom
    1099              :       TYPE(constraint_info_type), POINTER                :: cons_info
    1100              :       TYPE(particle_type), DIMENSION(:), OPTIONAL, &
    1101              :          POINTER                                         :: particle_set
    1102              :       TYPE(section_vals_type), POINTER                   :: colvar_func_info
    1103              :       LOGICAL, INTENT(IN)                                :: use_clv_info
    1104              :       INTEGER, INTENT(INOUT)                             :: cind
    1105              : 
    1106              :       INTEGER                                            :: ind, k, kk
    1107         2152 :       REAL(KIND=dp), DIMENSION(:), POINTER               :: r_vals
    1108              : 
    1109         4446 :       DO kk = 1, SIZE(ilist)
    1110         2294 :          k = ilist(kk)
    1111         2294 :          lcolv(kk)%init = .FALSE.
    1112         2294 :          lcolv(kk)%lambda = 0.0_dp
    1113         2294 :          lcolv(kk)%sigma = 0.0_dp
    1114              : 
    1115              :          ! Set Up colvar variable
    1116         2294 :          NULLIFY (lcolv(kk)%colvar, lcolv(kk)%colvar_old)
    1117              :          ! Colvar
    1118              :          CALL colvar_clone(lcolv(kk)%colvar, cons_info%colvar_set(k)%colvar, &
    1119         2294 :                            i_atom_offset=first_atom - 1)
    1120              : 
    1121              :          ! Some COLVARS may need additional information for evaluating the
    1122              :          ! functional form: this is the case for COLVARS which depend on the
    1123              :          ! initial position of the atoms: This information is stored in a proper
    1124              :          ! container in the COLVAR_RESTART section..
    1125         2294 :          IF ((lcolv(kk)%colvar%type_id == xyz_diag_colvar_id) .OR. &
    1126              :              (lcolv(kk)%colvar%type_id == xyz_outerdiag_colvar_id)) THEN
    1127           12 :             cind = cind + 1
    1128           12 :             IF (use_clv_info) THEN
    1129              :                CALL section_vals_val_get(colvar_func_info, "_DEFAULT_KEYWORD_", &
    1130            0 :                                          i_rep_val=cind, r_vals=r_vals)
    1131            0 :                SELECT CASE (lcolv(kk)%colvar%type_id)
    1132              :                CASE (xyz_diag_colvar_id)
    1133            0 :                   CPASSERT(SIZE(r_vals) == 3)
    1134            0 :                   lcolv(kk)%colvar%xyz_diag_param%r0 = r_vals
    1135              :                CASE (xyz_outerdiag_colvar_id)
    1136            0 :                   CPASSERT(SIZE(r_vals) == 6)
    1137            0 :                   lcolv(kk)%colvar%xyz_outerdiag_param%r0(:, 1) = r_vals(1:3)
    1138            0 :                   lcolv(kk)%colvar%xyz_outerdiag_param%r0(:, 2) = r_vals(4:6)
    1139              :                END SELECT
    1140              :             ELSE
    1141            6 :                SELECT CASE (lcolv(kk)%colvar%type_id)
    1142              :                CASE (xyz_diag_colvar_id)
    1143            6 :                   ALLOCATE (r_vals(3))
    1144            6 :                   ind = first_atom - 1 + lcolv(kk)%colvar%xyz_diag_param%i_atom
    1145           24 :                   r_vals = particle_set(ind)%r
    1146           48 :                   lcolv(kk)%colvar%xyz_diag_param%r0 = r_vals
    1147              :                CASE (xyz_outerdiag_colvar_id)
    1148            6 :                   ALLOCATE (r_vals(6))
    1149            6 :                   ind = first_atom - 1 + lcolv(kk)%colvar%xyz_outerdiag_param%i_atoms(1)
    1150           24 :                   r_vals(1:3) = particle_set(ind)%r
    1151            6 :                   ind = first_atom - 1 + lcolv(kk)%colvar%xyz_outerdiag_param%i_atoms(2)
    1152           24 :                   r_vals(4:6) = particle_set(ind)%r
    1153           48 :                   lcolv(kk)%colvar%xyz_outerdiag_param%r0(:, 1) = r_vals(1:3)
    1154           60 :                   lcolv(kk)%colvar%xyz_outerdiag_param%r0(:, 2) = r_vals(4:6)
    1155              :                END SELECT
    1156              :                CALL section_vals_val_set(colvar_func_info, "_DEFAULT_KEYWORD_", &
    1157           12 :                                          i_rep_val=cind, r_vals_ptr=r_vals)
    1158              :             END IF
    1159              :          END IF
    1160              : 
    1161              :          ! Setup Colvar_old
    1162         2294 :          CALL colvar_clone(lcolv(kk)%colvar_old, lcolv(kk)%colvar)
    1163              : 
    1164              :          ! Check for consistency in the constraint definition
    1165        14100 :          IF (ANY(lcolv(kk)%colvar%i_atom > last_atom) .OR. &
    1166         2152 :              ANY(lcolv(kk)%colvar%i_atom < first_atom)) THEN
    1167            0 :             WRITE (*, '(T5,"|",T8,A)') "Error in constraints setup!"
    1168            0 :             WRITE (*, '(T5,"|",T8,A)') "A constraint has been defined for a molecule type", &
    1169            0 :                " but the atoms specified in the constraint and the atoms defined for", &
    1170            0 :                " the molecule DO NOT match!", &
    1171            0 :                "This could be very probable due to a wrong connectivity, or an error", &
    1172            0 :                " in the constraint specification in the input file.", &
    1173            0 :                " Please check it carefully!"
    1174            0 :             CPABORT("")
    1175              :          END IF
    1176              :       END DO
    1177         2152 :    END SUBROUTINE setup_lcolv
    1178              : 
    1179              : ! **************************************************************************************************
    1180              : !> \brief Setup the lg3x3 for the packing of constraints
    1181              : !> \param lg3x3 ...
    1182              : !> \param g3x3_list ...
    1183              : !> \param first_atom ...
    1184              : !> \param last_atom ...
    1185              : !> \par History
    1186              : !>      Updated 2007 for intermolecular constraints
    1187              : !> \author Teodoro Laino [2007]
    1188              : ! **************************************************************************************************
    1189        36358 :    SUBROUTINE setup_lg3x3(lg3x3, g3x3_list, first_atom, last_atom)
    1190              :       TYPE(local_g3x3_constraint_type), DIMENSION(:), &
    1191              :          POINTER                                         :: lg3x3
    1192              :       TYPE(g3x3_constraint_type), DIMENSION(:), POINTER  :: g3x3_list
    1193              :       INTEGER, INTENT(IN)                                :: first_atom, last_atom
    1194              : 
    1195              :       INTEGER                                            :: kk
    1196              : 
    1197        62600 :       DO kk = 1, SIZE(lg3x3)
    1198        26242 :          lg3x3(kk)%init = .FALSE.
    1199        26242 :          lg3x3(kk)%scale = 0.0_dp
    1200        26242 :          lg3x3(kk)%scale_old = 0.0_dp
    1201       104968 :          lg3x3(kk)%fa = 0.0_dp
    1202       104968 :          lg3x3(kk)%fb = 0.0_dp
    1203       104968 :          lg3x3(kk)%fc = 0.0_dp
    1204       104968 :          lg3x3(kk)%ra_old = 0.0_dp
    1205       104968 :          lg3x3(kk)%rb_old = 0.0_dp
    1206       104968 :          lg3x3(kk)%rc_old = 0.0_dp
    1207       104968 :          lg3x3(kk)%va = 0.0_dp
    1208       104968 :          lg3x3(kk)%vb = 0.0_dp
    1209       104968 :          lg3x3(kk)%vc = 0.0_dp
    1210       104968 :          lg3x3(kk)%lambda = 0.0_dp
    1211              :          IF ((g3x3_list(kk)%a + first_atom - 1 < first_atom) .OR. &
    1212              :              (g3x3_list(kk)%b + first_atom - 1 < first_atom) .OR. &
    1213              :              (g3x3_list(kk)%c + first_atom - 1 < first_atom) .OR. &
    1214              :              (g3x3_list(kk)%a + first_atom - 1 > last_atom) .OR. &
    1215        26242 :              (g3x3_list(kk)%b + first_atom - 1 > last_atom) .OR. &
    1216        36358 :              (g3x3_list(kk)%c + first_atom - 1 > last_atom)) THEN
    1217            0 :             WRITE (*, '(T5,"|",T8,A)') "Error in constraints setup!"
    1218            0 :             WRITE (*, '(T5,"|",T8,A)') "A constraint has been defined for a molecule type", &
    1219            0 :                " but the atoms specified in the constraint and the atoms defined for", &
    1220            0 :                " the molecule DO NOT match!", &
    1221            0 :                "This could be very probable due to a wrong connectivity, or an error", &
    1222            0 :                " in the constraint specification in the input file.", &
    1223            0 :                " Please check it carefully!"
    1224            0 :             CPABORT("")
    1225              :          END IF
    1226              :       END DO
    1227              : 
    1228        36358 :    END SUBROUTINE setup_lg3x3
    1229              : 
    1230              : ! **************************************************************************************************
    1231              : !> \brief Setup the lg4x6 for the packing of constraints
    1232              : !> \param lg4x6 ...
    1233              : !> \param g4x6_list ...
    1234              : !> \param first_atom ...
    1235              : !> \param last_atom ...
    1236              : !> \par History
    1237              : !>      Updated 2007 for intermolecular constraints
    1238              : !> \author Teodoro Laino [2007]
    1239              : ! **************************************************************************************************
    1240          654 :    SUBROUTINE setup_lg4x6(lg4x6, g4x6_list, first_atom, last_atom)
    1241              :       TYPE(local_g4x6_constraint_type), DIMENSION(:), &
    1242              :          POINTER                                         :: lg4x6
    1243              :       TYPE(g4x6_constraint_type), DIMENSION(:), POINTER  :: g4x6_list
    1244              :       INTEGER, INTENT(IN)                                :: first_atom, last_atom
    1245              : 
    1246              :       INTEGER                                            :: kk
    1247              : 
    1248         1048 :       DO kk = 1, SIZE(lg4x6)
    1249          394 :          lg4x6(kk)%init = .FALSE.
    1250          394 :          lg4x6(kk)%scale = 0.0_dp
    1251          394 :          lg4x6(kk)%scale_old = 0.0_dp
    1252         1576 :          lg4x6(kk)%fa = 0.0_dp
    1253         1576 :          lg4x6(kk)%fb = 0.0_dp
    1254         1576 :          lg4x6(kk)%fc = 0.0_dp
    1255         1576 :          lg4x6(kk)%fd = 0.0_dp
    1256         1576 :          lg4x6(kk)%fe = 0.0_dp
    1257         1576 :          lg4x6(kk)%ff = 0.0_dp
    1258         1576 :          lg4x6(kk)%ra_old = 0.0_dp
    1259         1576 :          lg4x6(kk)%rb_old = 0.0_dp
    1260         1576 :          lg4x6(kk)%rc_old = 0.0_dp
    1261         1576 :          lg4x6(kk)%rd_old = 0.0_dp
    1262         1576 :          lg4x6(kk)%re_old = 0.0_dp
    1263         1576 :          lg4x6(kk)%rf_old = 0.0_dp
    1264         1576 :          lg4x6(kk)%va = 0.0_dp
    1265         1576 :          lg4x6(kk)%vb = 0.0_dp
    1266         1576 :          lg4x6(kk)%vc = 0.0_dp
    1267         1576 :          lg4x6(kk)%vd = 0.0_dp
    1268         1576 :          lg4x6(kk)%ve = 0.0_dp
    1269         1576 :          lg4x6(kk)%vf = 0.0_dp
    1270         2758 :          lg4x6(kk)%lambda = 0.0_dp
    1271              :          IF ((g4x6_list(kk)%a + first_atom - 1 < first_atom) .OR. &
    1272              :              (g4x6_list(kk)%b + first_atom - 1 < first_atom) .OR. &
    1273              :              (g4x6_list(kk)%c + first_atom - 1 < first_atom) .OR. &
    1274              :              (g4x6_list(kk)%d + first_atom - 1 < first_atom) .OR. &
    1275              :              (g4x6_list(kk)%a + first_atom - 1 > last_atom) .OR. &
    1276              :              (g4x6_list(kk)%b + first_atom - 1 > last_atom) .OR. &
    1277          394 :              (g4x6_list(kk)%c + first_atom - 1 > last_atom) .OR. &
    1278          654 :              (g4x6_list(kk)%d + first_atom - 1 > last_atom)) THEN
    1279            0 :             WRITE (*, '(T5,"|",T8,A)') "Error in constraints setup!"
    1280            0 :             WRITE (*, '(T5,"|",T8,A)') "A constrained has been defined for a molecule type", &
    1281            0 :                " but the atoms specified in the constraint and the atoms defined for", &
    1282            0 :                " the molecule DO NOT match!", &
    1283            0 :                "This could be very probable due to a wrong connectivity, or an error", &
    1284            0 :                " in the constraint specification in the input file.", &
    1285            0 :                " Please check it carefully!"
    1286            0 :             CPABORT("")
    1287              :          END IF
    1288              :       END DO
    1289              : 
    1290          654 :    END SUBROUTINE setup_lg4x6
    1291              : 
    1292              : ! **************************************************************************************************
    1293              : !> \brief Gives back a list of molecule to which apply the constraint
    1294              : !> \param const_mol ...
    1295              : !> \param const_molname ...
    1296              : !> \param const_intermolecular ...
    1297              : !> \param constr_x_mol ...
    1298              : !> \param constr_x_glob ...
    1299              : !> \param molecule_kind_set ...
    1300              : !> \param exclude_qm ...
    1301              : !> \param exclude_mm ...
    1302              : !> \par History
    1303              : !>      Updated 2007 for intermolecular constraints
    1304              : !> \author Teodoro Laino [2006]
    1305              : ! **************************************************************************************************
    1306          316 :    SUBROUTINE give_constraint_array(const_mol, const_molname, const_intermolecular, &
    1307              :                                     constr_x_mol, constr_x_glob, molecule_kind_set, exclude_qm, exclude_mm)
    1308              : 
    1309              :       INTEGER, DIMENSION(:), POINTER                     :: const_mol
    1310              :       CHARACTER(LEN=default_string_length), &
    1311              :          DIMENSION(:), POINTER                           :: const_molname
    1312              :       LOGICAL, DIMENSION(:), POINTER                     :: const_intermolecular
    1313              :       TYPE(constr_list_type), DIMENSION(:), POINTER      :: constr_x_mol
    1314              :       INTEGER, DIMENSION(:), POINTER                     :: constr_x_glob
    1315              :       TYPE(molecule_kind_type), DIMENSION(:), POINTER    :: molecule_kind_set
    1316              :       LOGICAL, DIMENSION(:), POINTER                     :: exclude_qm, exclude_mm
    1317              : 
    1318              :       CHARACTER(len=*), PARAMETER :: routineN = 'give_constraint_array'
    1319              : 
    1320              :       CHARACTER(LEN=default_string_length)               :: myname, name
    1321              :       INTEGER                                            :: handle, i, iglob, isize, k
    1322              :       LOGICAL                                            :: found_molname, is_qm
    1323              :       TYPE(molecule_kind_type), POINTER                  :: molecule_kind
    1324              : 
    1325          316 :       CALL timeset(routineN, handle)
    1326          316 :       NULLIFY (molecule_kind)
    1327         1826 :       ALLOCATE (constr_x_mol(SIZE(molecule_kind_set)))
    1328         1194 :       DO i = 1, SIZE(constr_x_mol)
    1329          878 :          NULLIFY (constr_x_mol(i)%constr)
    1330         1194 :          ALLOCATE (constr_x_mol(i)%constr(0))
    1331              :       END DO
    1332          316 :       CPASSERT(SIZE(const_mol) == SIZE(const_molname))
    1333          316 :       iglob = 0
    1334          950 :       DO i = 1, SIZE(const_mol)
    1335          950 :          IF (const_intermolecular(i)) THEN
    1336              :             ! Intermolecular constraint
    1337           74 :             iglob = iglob + 1
    1338           74 :             CALL reallocate(constr_x_glob, 1, iglob)
    1339           74 :             constr_x_glob(iglob) = i
    1340              :          ELSE
    1341              :             ! Intramolecular constraint
    1342          560 :             IF (const_mol(i) /= 0) THEN
    1343          476 :                k = const_mol(i)
    1344          476 :                IF (k > SIZE(molecule_kind_set)) &
    1345              :                   CALL cp_abort(__LOCATION__, &
    1346              :                                 "A constraint has been specified providing the molecule index. But the"// &
    1347              :                                 " molecule index ("//cp_to_string(k)//") is out of range of the possible"// &
    1348            0 :                                 " molecule kinds ("//cp_to_string(SIZE(molecule_kind_set))//").")
    1349          476 :                isize = SIZE(constr_x_mol(k)%constr)
    1350          476 :                CALL reallocate(constr_x_mol(k)%constr, 1, isize + 1)
    1351          476 :                constr_x_mol(k)%constr(isize + 1) = i
    1352              :             ELSE
    1353           84 :                myname = const_molname(i)
    1354           84 :                found_molname = .FALSE.
    1355          304 :                DO k = 1, SIZE(molecule_kind_set)
    1356          220 :                   molecule_kind => molecule_kind_set(k)
    1357          220 :                   name = molecule_kind%name
    1358          220 :                   is_qm = qmmm_ff_precond_only_qm(id1=name)
    1359          220 :                   IF (is_qm .AND. exclude_qm(i)) CYCLE
    1360          152 :                   IF (.NOT. is_qm .AND. exclude_mm(i)) CYCLE
    1361          292 :                   IF (name == myname) THEN
    1362           82 :                      isize = SIZE(constr_x_mol(k)%constr)
    1363           82 :                      CALL reallocate(constr_x_mol(k)%constr, 1, isize + 1)
    1364           82 :                      constr_x_mol(k)%constr(isize + 1) = i
    1365           82 :                      found_molname = .TRUE.
    1366              :                   END IF
    1367              :                END DO
    1368           84 :                CALL print_warning_molname(found_molname, myname)
    1369              :             END IF
    1370              :          END IF
    1371              :       END DO
    1372          316 :       CALL timestop(handle)
    1373          316 :    END SUBROUTINE give_constraint_array
    1374              : 
    1375              : ! **************************************************************************************************
    1376              : !> \brief Prints a warning message if undefined molnames are used to define constraints
    1377              : !> \param found ...
    1378              : !> \param name ...
    1379              : !> \author Teodoro Laino [2007] - Zurich University
    1380              : ! **************************************************************************************************
    1381           86 :    SUBROUTINE print_warning_molname(found, name)
    1382              :       LOGICAL, INTENT(IN)                                :: found
    1383              :       CHARACTER(LEN=*), INTENT(IN)                       :: name
    1384              : 
    1385           86 :       IF (.NOT. found) &
    1386              :          CALL cp_warn(__LOCATION__, &
    1387              :                       " MOLNAME ("//TRIM(name)//") was defined for constraints, but this molecule name "// &
    1388              :                       "is not defined. Please check carefully your PDB, PSF (has priority over PDB) or "// &
    1389              :                       "input driven CP2K coordinates. In case you may not find the reason for this warning "// &
    1390              :                       "it may be a good idea to print all molecule information (including kind name) activating "// &
    1391            6 :                       "the print_key MOLECULES specific of the SUBSYS%PRINT section. ")
    1392              : 
    1393           86 :    END SUBROUTINE print_warning_molname
    1394              : 
    1395              : END MODULE topology_constraint_util
        

Generated by: LCOV version 2.0-1