mimic_request_system_data Subroutine

public subroutine mimic_request_system_data(sizes, system_data)

return data about atomic species, masses and multipoles in the MM part

Arguments

TypeIntentOptionalAttributesName
type(sizes_type), intent(inout) :: sizes

sizes data of the MM part

type(system_type), intent(inout) :: system_data

system data


Contents


Source Code

subroutine mimic_request_system_data(sizes, system_data)

    !> sizes data of the MM part
    type(sizes_type), intent(inout) :: sizes
    !> system data
    !> @todo intent(out)?
    type(system_type), intent(inout) :: system_data

    integer :: max_dim
    integer :: alloc_stat

    call timer_start("mimic_request_system_data")
    max_dim = maxval(sizes%atoms_pcode)

    allocate(system_data%multipole_values(maxval(sizes%multipoles_patom), &
                                          max_dim, communicator%num_clients), stat=alloc_stat)
    if (alloc_stat /= 0) then
        call handle_error(SEVERITY_FATAL, &
                          TYPE_MEM, &
                          "system_data%multipole_values", &
                          __FILE__, __LINE__)
    endif
    ! allocate(system_data%species(max_dim, num_clients))
    allocate(system_data%masses(sizes%num_species), stat=alloc_stat)
    if (alloc_stat /= 0) then
        call handle_error(SEVERITY_FATAL, &
                          TYPE_MEM, &
                          "system_data%masses", &
                          __FILE__, __LINE__)
    endif
    allocate(system_data%elements(sizes%num_species), stat=alloc_stat)
    if (alloc_stat /= 0) then
        call handle_error(SEVERITY_FATAL, &
                          TYPE_MEM, &
                          "system_data%elements", &
                          __FILE__, __LINE__)
    endif

    call communicator%gather_multipoles(system_data%species, &
                                        system_data%multipole_values, &
                                        sizes%atoms_pcode, sizes%multipoles_patom)

    call communicator%gather_masses(system_data%masses, &
                                    sizes%types_pcode, &
                                    system_data%species_map)

    allocate(system_data%atom_fragment_ids(maxval(sizes%atoms_pfragment),&
             maxval(sizes%frag_num),&
             communicator%num_clients), stat=alloc_stat)
    if (alloc_stat /= 0) then
        call handle_error(SEVERITY_FATAL, &
                          TYPE_MEM, &
                          "system_data%atom_fragment_ids", &
                          __FILE__, __LINE__)
    endif

    call communicator%gather_atom_fragment_ids(sizes%atoms_pcode,&
                            sizes%atoms_pfragment,&
                            system_data%atom_fragment_ids)

    call communicator%gather_elements(sizes%types_pcode, &
                                      system_data%elements, &
                                      system_data%species_map)
    call timer_stop

end subroutine mimic_request_system_data