Build a fragment from explicit atom list (for GMBE intersection fragments)
Similar to build_fragment_from_indices but takes atom indices directly instead of monomer indices. Used for building intersection fragments in GMBE calculations. Intersection fragments are ALWAYS NEUTRAL (charge=0, multiplicity=1).
Example: atom_indices = [3, 4, 5] builds fragment from atoms 3, 4, 5 of the system
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(system_geometry_t), | intent(in) | :: | sys_geom | |||
| integer, | intent(in) | :: | atom_indices(:) |
0-indexed atom indices |
||
| integer, | intent(in) | :: | n_atoms |
Number of atoms in list |
||
| type(physical_fragment_t), | intent(out) | :: | fragment | |||
| type(error_t), | intent(out) | :: | error | |||
| type(bond_t), | intent(in), | optional | :: | bonds(:) |
Connectivity for capping |
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| integer, | private | :: | atom_global_idx | ||||
| integer, | private | :: | frag_atom_idx | ||||
| integer, | private | :: | i | ||||
| integer, | private | :: | n_caps |
subroutine build_fragment_from_atom_list(sys_geom, atom_indices, n_atoms, fragment, error, bonds) !! Build a fragment from explicit atom list (for GMBE intersection fragments) !! !! Similar to build_fragment_from_indices but takes atom indices directly instead of !! monomer indices. Used for building intersection fragments in GMBE calculations. !! Intersection fragments are ALWAYS NEUTRAL (charge=0, multiplicity=1). !! !! Example: atom_indices = [3, 4, 5] builds fragment from atoms 3, 4, 5 of the system type(system_geometry_t), intent(in) :: sys_geom integer, intent(in) :: atom_indices(:) !! 0-indexed atom indices integer, intent(in) :: n_atoms !! Number of atoms in list type(physical_fragment_t), intent(out) :: fragment type(error_t), intent(out) :: error type(bond_t), intent(in), optional :: bonds(:) !! Connectivity for capping integer :: i, frag_atom_idx, atom_global_idx integer :: n_caps ! Count how many caps we need call count_hydrogen_caps(atom_indices(1:n_atoms), bonds, n_caps) ! Allocate arrays with space for original atoms + caps fragment%n_atoms = n_atoms + n_caps fragment%n_caps = n_caps allocate (fragment%element_numbers(fragment%n_atoms)) allocate (fragment%coordinates(3, fragment%n_atoms)) if (n_caps > 0) allocate (fragment%cap_replaces_atom(n_caps)) allocate (fragment%local_to_global(n_atoms)) ! Only non-cap atoms ! Copy original atoms and build local→global mapping (atom_indices are 0-indexed, add 1 for Fortran arrays) do i = 1, n_atoms atom_global_idx = atom_indices(i) + 1 ! Convert to 1-indexed fragment%element_numbers(i) = sys_geom%element_numbers(atom_global_idx) fragment%coordinates(:, i) = sys_geom%coordinates(:, atom_global_idx) fragment%local_to_global(i) = atom_global_idx ! Store 1-indexed global position end do ! Add hydrogen caps at end (if any) if (present(bonds) .and. n_caps > 0) then call add_hydrogen_caps(atom_indices(1:n_atoms), bonds, sys_geom, fragment, n_atoms) end if ! Intersection fragments are ALWAYS NEUTRAL ! Rationale: For polypeptides, intersections are backbone atoms; ! charged side chains are in non-overlapping regions fragment%charge = 0 fragment%multiplicity = 1 call fragment%compute_nelec() ! Validate: check for spatially overlapping atoms call check_duplicate_atoms(fragment, error) if (error%has_error()) then call error%add_context("mqc_physical_fragment:build_fragment_from_atom_list") return end if end subroutine build_fragment_from_atom_list