Map fragment gradient to system gradient coordinates with hydrogen cap redistribution
This function rebuilds the fragment to get local→global mappings and cap information, then redistributes gradients including hydrogen caps to their original atoms. Bond connectivity is accessed via sys_geom%bonds.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=dp), | intent(in) | :: | frag_grad(:,:) |
(3, natoms_frag) |
||
| integer, | intent(in) | :: | monomers(:) |
Monomer indices in fragment |
||
| type(system_geometry_t), | intent(in) | :: | sys_geom | |||
| real(kind=dp), | intent(inout) | :: | sys_grad(:,:) |
(3, total_atoms) |
||
| type(comm_t), | intent(in), | optional | :: | world_comm |
MPI communicator for abort |
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| integer, | private | :: | current_log_level | ||||
| type(error_t), | private | :: | error | ||||
| integer, | private | :: | frag_atom_idx | ||||
| type(physical_fragment_t), | private | :: | fragment | ||||
| integer, | private | :: | i_atom | ||||
| integer, | private | :: | i_mon | ||||
| integer, | private | :: | sys_atom_idx |
subroutine map_fragment_to_system_gradient(frag_grad, monomers, sys_geom, sys_grad, world_comm) !! Map fragment gradient to system gradient coordinates with hydrogen cap redistribution !! !! This function rebuilds the fragment to get local→global mappings and cap information, !! then redistributes gradients including hydrogen caps to their original atoms. !! Bond connectivity is accessed via sys_geom%bonds. use mqc_physical_fragment, only: build_fragment_from_indices, redistribute_cap_gradients use mqc_error, only: error_t use pic_logger, only: verbose_level real(dp), intent(in) :: frag_grad(:, :) !! (3, natoms_frag) integer, intent(in) :: monomers(:) !! Monomer indices in fragment type(system_geometry_t), intent(in) :: sys_geom real(dp), intent(inout) :: sys_grad(:, :) !! (3, total_atoms) type(comm_t), intent(in), optional :: world_comm !! MPI communicator for abort type(physical_fragment_t) :: fragment type(error_t) :: error integer :: i_mon, i_atom, frag_atom_idx, sys_atom_idx integer :: current_log_level ! Explicitly zero out the entire sys_grad array sys_grad = 0.0_dp ! Debug output call logger%configuration(level=current_log_level) if (current_log_level >= debug_level) then block character(len=256) :: debug_msg write (debug_msg, '(a,i0,a,*(i0,1x))') " Mapping fragment with ", size(monomers), " monomers: ", monomers call logger%debug(trim(debug_msg)) write (debug_msg, '(a,i0,a)') " Fragment has ", size(frag_grad, 2), " atoms" call logger%debug(trim(debug_msg)) end block end if 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) if (error%has_error()) then call logger%error(error%get_full_trace()) if (present(world_comm)) then call abort_comm(world_comm, 1) else error stop "Failed to build fragment in gradient mapping" end if end if ! Use new gradient redistribution with cap handling call redistribute_cap_gradients(fragment, frag_grad, sys_grad) ! Clean up call fragment%destroy() else ! Old code path for fragments without hydrogen caps ! Map fragment gradient to system positions (fixed-size monomers only) frag_atom_idx = 0 do i_mon = 1, size(monomers) do i_atom = 1, sys_geom%atoms_per_monomer frag_atom_idx = frag_atom_idx + 1 sys_atom_idx = (monomers(i_mon) - 1)*sys_geom%atoms_per_monomer + i_atom if (current_log_level >= debug_level .and. i_atom == 1) then block character(len=256) :: debug_msg write (debug_msg, '(a,i0,a,i0,a,i0)') & " Monomer ", monomers(i_mon), ": frag atoms ", & frag_atom_idx, " -> sys atom ", sys_atom_idx call logger%debug(trim(debug_msg)) end block end if sys_grad(:, sys_atom_idx) = frag_grad(:, frag_atom_idx) end do end do end if end subroutine map_fragment_to_system_gradient