Map fragment Hessian to system Hessian coordinates with hydrogen cap redistribution Bond connectivity is accessed via sys_geom%bonds
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=dp), | intent(in) | :: | frag_hess(:,:) |
(3natoms_frag, 3natoms_frag) |
||
| integer, | intent(in) | :: | monomers(:) | |||
| type(system_geometry_t), | intent(in) | :: | sys_geom | |||
| real(kind=dp), | intent(inout) | :: | sys_hess(:,:) |
(3total_atoms, 3total_atoms) |
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| type(error_t), | private | :: | error | ||||
| type(physical_fragment_t), | private | :: | fragment |
subroutine map_fragment_to_system_hessian(frag_hess, monomers, sys_geom, sys_hess) !! Map fragment Hessian to system Hessian coordinates with hydrogen cap redistribution !! Bond connectivity is accessed via sys_geom%bonds use mqc_physical_fragment, only: build_fragment_from_indices, redistribute_cap_hessian use mqc_error, only: error_t real(dp), intent(in) :: frag_hess(:, :) !! (3*natoms_frag, 3*natoms_frag) integer, intent(in) :: monomers(:) type(system_geometry_t), intent(in) :: sys_geom real(dp), intent(inout) :: sys_hess(:, :) !! (3*total_atoms, 3*total_atoms) type(physical_fragment_t) :: fragment type(error_t) :: error ! Zero out sys_hess = 0.0_dp if (allocated(sys_geom%bonds)) then ! Rebuild fragment to get local→global mapping and cap information call build_fragment_from_indices(sys_geom, monomers, fragment, error, sys_geom%bonds) call redistribute_cap_hessian(fragment, frag_hess, sys_hess) call fragment%destroy() else ! Old code path for fragments without hydrogen caps ! Map fragment Hessian to system positions (fixed-size monomers only) block integer :: i_mon, j_mon, i_atom, j_atom integer :: frag_atom_i, frag_atom_j, sys_atom_i, sys_atom_j integer :: frag_row_start, frag_col_start, sys_row_start, sys_col_start integer :: n_monomers n_monomers = size(monomers) frag_atom_i = 0 ! Map each monomer's atoms do i_mon = 1, n_monomers do i_atom = 1, sys_geom%atoms_per_monomer frag_atom_i = frag_atom_i + 1 sys_atom_i = (monomers(i_mon) - 1)*sys_geom%atoms_per_monomer + i_atom frag_row_start = (frag_atom_i - 1)*3 + 1 sys_row_start = (sys_atom_i - 1)*3 + 1 ! Map this atom's Hessian blocks with all other atoms in fragment frag_atom_j = 0 do j_mon = 1, n_monomers do j_atom = 1, sys_geom%atoms_per_monomer frag_atom_j = frag_atom_j + 1 sys_atom_j = (monomers(j_mon) - 1)*sys_geom%atoms_per_monomer + j_atom frag_col_start = (frag_atom_j - 1)*3 + 1 sys_col_start = (sys_atom_j - 1)*3 + 1 ! Copy the 3×3 block for this atom pair sys_hess(sys_row_start:sys_row_start + 2, sys_col_start:sys_col_start + 2) = & frag_hess(frag_row_start:frag_row_start + 2, frag_col_start:frag_col_start + 2) end do end do end do end do end block end if end subroutine map_fragment_to_system_hessian