diff --git a/src/deepmd_wrapper.F b/src/deepmd_wrapper.F index 25cd8d4c8e..a1a89b267f 100644 --- a/src/deepmd_wrapper.F +++ b/src/deepmd_wrapper.F @@ -15,80 +15,62 @@ ! ************************************************************************************************** MODULE deepmd_wrapper - USE ISO_C_BINDING, ONLY: C_CHAR,& C_DOUBLE,& C_INT,& - C_LOC,& C_NULL_CHAR,& C_NULL_PTR,& C_PTR + USE kinds, ONLY: dp #include "./base/base_uses.f90" IMPLICIT NONE PRIVATE - PUBLIC :: dp_deep_pot, dp_deep_pot_compute, deepmd_model_release - -#if(__DEEPMD) - INTERFACE - FUNCTION dp_deep_pot_c(c_model) BIND(C, name="DP_NewDeepPot") - USE ISO_C_BINDING, ONLY: C_PTR, & - C_CHAR, & - C_DOUBLE, & - C_INT, & - C_NULL_CHAR, & - C_LOC - CHARACTER(LEN=1, KIND=C_CHAR) :: c_model(*) - TYPE(C_PTR) :: dp_deep_pot_c - - END FUNCTION - - SUBROUTINE dp_deep_pot_compute_c(dp, natom, & - coord, atype, cell, & - energy, force, virial, & - atomic_energy, atomic_virial) BIND(C, name="DP_DeepPotCompute") - USE ISO_C_BINDING, ONLY: C_PTR, & - C_INT, & - C_CHAR, & - C_DOUBLE - TYPE(C_PTR), VALUE :: dp - INTEGER(C_INT), VALUE :: natom - REAL(C_DOUBLE), DIMENSION(natom, 3), INTENT(IN) :: coord - INTEGER(C_INT), DIMENSION(natom), INTENT(IN) :: atype - REAL(C_DOUBLE), DIMENSION(9), INTENT(IN) :: cell - REAL(C_DOUBLE), INTENT(OUT) :: energy - REAL(C_DOUBLE), DIMENSION(natom, 3), INTENT(OUT) :: force - REAL(C_DOUBLE), DIMENSION(9), INTENT(OUT) :: virial - REAL(C_DOUBLE), DIMENSION(natom), INTENT(OUT) :: atomic_energy - REAL(C_DOUBLE), DIMENSION(natom, 9), INTENT(OUT) :: atomic_virial - - END SUBROUTINE - END INTERFACE -#endif + PUBLIC :: deepmd_model_type, deepmd_model_load, deepmd_model_compute, deepmd_model_release + + TYPE deepmd_model_type + PRIVATE + TYPE(C_PTR) :: c_ptr = C_NULL_PTR + END TYPE deepmd_model_type CONTAINS ! ************************************************************************************************** !> \brief Load DP from a model file. -!> \param model Path to the model file. +!> \param filename Path to the model file. !> \return Pointer to the DP model. ! ************************************************************************************************** + FUNCTION deepmd_model_load(filename) RESULT(model) + CHARACTER(len=*), INTENT(INOUT) :: filename + TYPE(deepmd_model_type) :: model - FUNCTION dp_deep_pot(model) RESULT(pot) - CHARACTER(len=*), INTENT(INOUT) :: model - TYPE(C_PTR) :: pot + CHARACTER(LEN=*), PARAMETER :: routineN = 'deepmd_model_load' -#if(__DEEPMD) - pot = dp_deep_pot_c(c_model=TRIM(model)//C_NULL_CHAR) + INTEGER :: handle + INTERFACE + FUNCTION NewDeepPot(filename) BIND(C, name="DP_NewDeepPot") + IMPORT :: C_PTR, C_CHAR + CHARACTER(kind=C_CHAR), DIMENSION(*) :: filename + TYPE(C_PTR) :: NewDeepPot + END FUNCTION + END INTERFACE + + CALL timeset(routineN, handle) + +#if defined(__DEEPMD) + model%c_ptr = NewDeepPot(filename=TRIM(filename)//C_NULL_CHAR) #else + CPABORT("CP2K was compiled without libdeepmd_c library.") + MARK_USED(filename) MARK_USED(model) - pot = C_NULL_PTR #endif - END FUNCTION dp_deep_pot + + CALL timestop(handle) + END FUNCTION deepmd_model_load ! ************************************************************************************************** !> \brief Compute energy, force and virial from DP. -!> \param dp Pointer to the DP model. +!> \param model Pointer to the DP model. !> \param natom Number of atoms. !> \param coord Coordinates of the atoms. !> \param atype Atom types. @@ -99,25 +81,55 @@ END FUNCTION dp_deep_pot !> \param atomic_energy Atomic energies. !> \param atomic_virial Atomic virial tensors. ! ************************************************************************************************** - - SUBROUTINE dp_deep_pot_compute(dp, natom, coord, atype, cell, energy, force, virial, atomic_energy, atomic_virial) - USE, INTRINSIC :: ISO_C_BINDING - TYPE(C_PTR), VALUE :: dp - INTEGER(C_INT), VALUE :: natom - REAL(C_DOUBLE), DIMENSION(natom, 3), INTENT(IN) :: coord - INTEGER(C_INT), DIMENSION(natom), INTENT(IN) :: atype - REAL(C_DOUBLE), DIMENSION(9), INTENT(IN), OPTIONAL :: cell - REAL(C_DOUBLE), INTENT(OUT), OPTIONAL :: energy - REAL(C_DOUBLE), DIMENSION(natom, 3), INTENT(OUT), OPTIONAL :: force - REAL(C_DOUBLE), DIMENSION(9), INTENT(OUT), OPTIONAL :: virial - REAL(C_DOUBLE), DIMENSION(natom), INTENT(OUT), OPTIONAL :: atomic_energy - REAL(C_DOUBLE), DIMENSION(natom, 9), INTENT(OUT), OPTIONAL :: atomic_virial - -#if(__DEEPMD) - CALL dp_deep_pot_compute_c(dp, natom, coord, atype, cell, energy, force, virial, atomic_energy, atomic_virial) + SUBROUTINE deepmd_model_compute(model, natom, coord, atype, cell, energy, force, virial, & + atomic_energy, atomic_virial) + TYPE(deepmd_model_type) :: model + INTEGER :: natom + REAL(kind=dp), DIMENSION(natom, 3), INTENT(IN) :: coord + INTEGER, DIMENSION(natom), INTENT(IN) :: atype + REAL(kind=dp), DIMENSION(9), INTENT(IN) :: cell + REAL(kind=dp), INTENT(OUT) :: energy + REAL(kind=dp), DIMENSION(natom, 3), INTENT(OUT) :: force + REAL(kind=dp), DIMENSION(9), INTENT(OUT) :: virial + REAL(kind=dp), DIMENSION(natom), INTENT(OUT) :: atomic_energy + REAL(kind=dp), DIMENSION(natom, 9), INTENT(OUT) :: atomic_virial + + CHARACTER(LEN=*), PARAMETER :: routineN = 'deepmd_model_compute' + + INTEGER :: handle + INTERFACE + SUBROUTINE DeepPotCompute(model, natom, coord, atype, cell, energy, force, virial, & + atomic_energy, atomic_virial) BIND(C, name="DP_DeepPotCompute") + IMPORT :: C_PTR, C_INT, C_DOUBLE + TYPE(C_PTR), VALUE :: model + INTEGER(C_INT), VALUE :: natom + REAL(C_DOUBLE), DIMENSION(natom, 3) :: coord + INTEGER(C_INT), DIMENSION(natom) :: atype + REAL(C_DOUBLE), DIMENSION(9) :: cell + REAL(C_DOUBLE) :: energy + REAL(C_DOUBLE), DIMENSION(natom, 3) :: force + REAL(C_DOUBLE), DIMENSION(9) :: virial + REAL(C_DOUBLE), DIMENSION(natom) :: atomic_energy + REAL(C_DOUBLE), DIMENSION(natom, 9) :: atomic_virial + END SUBROUTINE + END INTERFACE + + CALL timeset(routineN, handle) + +#if defined(__DEEPMD) + CALL DeepPotCompute(model=model%c_ptr, & + natom=natom, & + coord=coord, & + atype=atype, & + cell=cell, & + energy=energy, & + force=force, & + virial=virial, & + atomic_energy=atomic_energy, & + atomic_virial=atomic_virial) #else CPABORT("CP2K was compiled without libdeepmd_c library.") - MARK_USED(dp) + MARK_USED(model) MARK_USED(natom) MARK_USED(coord) MARK_USED(atype) @@ -128,6 +140,8 @@ SUBROUTINE dp_deep_pot_compute(dp, natom, coord, atype, cell, energy, force, vir MARK_USED(atomic_energy) MARK_USED(atomic_virial) #endif + + CALL timestop(handle) END SUBROUTINE ! ************************************************************************************************** @@ -135,11 +149,9 @@ SUBROUTINE dp_deep_pot_compute(dp, natom, coord, atype, cell, energy, force, vir !> \param model Pointer to the DP model. ! ************************************************************************************************** SUBROUTINE deepmd_model_release(model) - USE ISO_C_BINDING, ONLY: C_PTR, & - C_NULL_PTR - TYPE(C_PTR), INTENT(INOUT) :: model + TYPE(deepmd_model_type) :: model - model = C_NULL_PTR + model%c_ptr = C_NULL_PTR END SUBROUTINE deepmd_model_release END MODULE deepmd_wrapper diff --git a/src/fist_nonbond_env_types.F b/src/fist_nonbond_env_types.F index a187bec62e..aa0b39156c 100644 --- a/src/fist_nonbond_env_types.F +++ b/src/fist_nonbond_env_types.F @@ -11,11 +11,11 @@ !> \author HAF ! ************************************************************************************************** MODULE fist_nonbond_env_types - USE ISO_C_BINDING, ONLY: C_PTR USE atomic_kind_types, ONLY: atomic_kind_type USE cell_types, ONLY: cell_release,& cell_type - USE deepmd_wrapper, ONLY: deepmd_model_release + USE deepmd_wrapper, ONLY: deepmd_model_release,& + deepmd_model_type USE fist_neighbor_list_types, ONLY: fist_neighbor_deallocate,& fist_neighbor_type USE kinds, ONLY: default_string_length,& @@ -74,7 +74,7 @@ MODULE fist_nonbond_env_types INTEGER, POINTER :: use_indices(:) => NULL() REAL(KIND=dp), POINTER :: force(:, :) => NULL() REAL(KIND=dp) :: virial(3, 3) = 0.0_dp - TYPE(C_PTR) :: model + TYPE(deepmd_model_type) :: model END TYPE ! ************************************************************************************************** diff --git a/src/manybody_deepmd.F b/src/manybody_deepmd.F index 650dd8c558..fae16bf890 100644 --- a/src/manybody_deepmd.F +++ b/src/manybody_deepmd.F @@ -7,31 +7,26 @@ MODULE manybody_deepmd - USE cp_log_handling, ONLY: cp_logger_get_default_io_unit - USE atomic_kind_types, ONLY: atomic_kind_type - USE bibliography, ONLY: Wang2018, & - Zeng2023, & - cite_reference - USE cell_types, ONLY: cell_type - USE fist_nonbond_env_types, ONLY: fist_nonbond_env_get, & - fist_nonbond_env_set, & - fist_nonbond_env_type, & - deepmd_data_type - USE kinds, ONLY: dp - USE message_passing, ONLY: mp_para_env_type - USE pair_potential_types, ONLY: pair_potential_pp_type, & - pair_potential_single_type, & - deepmd_pot_type, & - deepmd_type - USE particle_types, ONLY: particle_type - USE physcon, ONLY: angstrom, & - evolt -#ifdef __DEEPMD - USE ISO_C_BINDING, ONLY: C_PTR, C_LOC - USE deepmd_wrapper, ONLY: dp_deep_pot, & - dp_deep_pot_compute -#endif - + USE atomic_kind_types, ONLY: atomic_kind_type + USE bibliography, ONLY: Wang2018,& + Zeng2023,& + cite_reference + USE cell_types, ONLY: cell_type + USE cp_log_handling, ONLY: cp_logger_get_default_io_unit + USE deepmd_wrapper, ONLY: deepmd_model_compute,& + deepmd_model_load + USE fist_nonbond_env_types, ONLY: deepmd_data_type,& + fist_nonbond_env_get,& + fist_nonbond_env_set,& + fist_nonbond_env_type + USE kinds, ONLY: dp + USE message_passing, ONLY: mp_para_env_type + USE pair_potential_types, ONLY: deepmd_type,& + pair_potential_pp_type,& + pair_potential_single_type + USE particle_types, ONLY: particle_type + USE physcon, ONLY: angstrom,& + evolt #include "./base/base_uses.f90" IMPLICIT NONE @@ -65,44 +60,26 @@ SUBROUTINE deepmd_energy_store_force_virial(particle_set, cell, atomic_kind_set, CHARACTER(LEN=*), PARAMETER :: routineN = 'deepmd_energy_store_force_virial' -#ifdef __DEEPMD - INTEGER :: i, iat, iat_use, ikind, & - jkind, n_atoms, n_atoms_use, & - output_unit - LOGICAL, ALLOCATABLE :: use_atom(:) - REAL(kind=dp) :: lattice(3, 3) - TYPE(pair_potential_single_type), & - POINTER :: pot - TYPE(deepmd_data_type), POINTER :: deepmd_data - TYPE(deepmd_pot_type), POINTER :: deepmd - REAL(kind=dp), ALLOCATABLE :: dpmd_force(:), dpmd_virial(:), & - dpmd_atomic_energy(:), dpmd_atomic_virial(:), & - dpmd_cell(:), dpmd_coord(:) - INTEGER, ALLOCATABLE :: dpmd_atype(:), use_atom_type(:) - INTEGER :: dpmd_natoms - -#endif - INTEGER :: handle + INTEGER :: dpmd_natoms, handle, i, iat, iat_use, & + ikind, jkind, n_atoms, n_atoms_use, & + output_unit + INTEGER, ALLOCATABLE :: dpmd_atype(:), use_atom_type(:) + LOGICAL, ALLOCATABLE :: use_atom(:) + REAL(kind=dp) :: lattice(3, 3) + REAL(kind=dp), ALLOCATABLE :: dpmd_atomic_energy(:), & + dpmd_atomic_virial(:), dpmd_cell(:), & + dpmd_coord(:), dpmd_force(:), & + dpmd_virial(:) + TYPE(deepmd_data_type), POINTER :: deepmd_data + TYPE(pair_potential_single_type), POINTER :: pot CALL timeset(routineN, handle) -#ifndef __DEEPMD - CPABORT("In order to use DeePMD-kit you need to download and install libdeepmd_c") - MARK_USED(particle_set) - MARK_USED(cell) - MARK_USED(atomic_kind_set) - MARK_USED(potparm) - MARK_USED(fist_nonbond_env) - MARK_USED(pot_deepmd) - MARK_USED(para_env) -#else n_atoms = SIZE(particle_set) ALLOCATE (use_atom(n_atoms)) ALLOCATE (use_atom_type(n_atoms)) use_atom = .FALSE. use_atom_type = 0 - NULLIFY (deepmd) - DO ikind = 1, SIZE(atomic_kind_set) DO jkind = 1, SIZE(atomic_kind_set) pot => potparm%pot(ikind, jkind)%pot @@ -152,7 +129,7 @@ SUBROUTINE deepmd_energy_store_force_virial(particle_set, cell, atomic_kind_set, ALLOCATE (deepmd_data) CALL fist_nonbond_env_set(fist_nonbond_env, deepmd_data=deepmd_data) NULLIFY (deepmd_data%use_indices, deepmd_data%force) - deepmd_data%model = dp_deep_pot(pot%set(1)%deepmd%deepmd_file_name) + deepmd_data%model = deepmd_model_load(filename=pot%set(1)%deepmd%deepmd_file_name) END IF IF (ASSOCIATED(deepmd_data%force)) THEN IF (SIZE(deepmd_data%force, 2) /= n_atoms_use) THEN @@ -164,16 +141,16 @@ SUBROUTINE deepmd_energy_store_force_virial(particle_set, cell, atomic_kind_set, ALLOCATE (deepmd_data%use_indices(n_atoms_use)) END IF - CALL dp_deep_pot_compute(dp=deepmd_data%model, & - natom=dpmd_natoms, & - coord=dpmd_coord, & - atype=dpmd_atype, & - cell=dpmd_cell, & - energy=pot_deepmd, & - force=dpmd_force, & - virial=dpmd_virial, & - atomic_energy=dpmd_atomic_energy, & - atomic_virial=dpmd_atomic_virial) + CALL deepmd_model_compute(model=deepmd_data%model, & + natom=dpmd_natoms, & + coord=dpmd_coord, & + atype=dpmd_atype, & + cell=dpmd_cell, & + energy=pot_deepmd, & + force=dpmd_force, & + virial=dpmd_virial, & + atomic_energy=dpmd_atomic_energy, & + atomic_virial=dpmd_atomic_virial) ! convert units pot_deepmd = pot_deepmd/evolt @@ -202,7 +179,7 @@ SUBROUTINE deepmd_energy_store_force_virial(particle_set, cell, atomic_kind_set, DEALLOCATE (use_atom, use_atom_type, dpmd_coord, dpmd_force, & dpmd_virial, dpmd_atomic_energy, dpmd_atomic_virial, & dpmd_cell, dpmd_atype) -#endif + CALL timestop(handle) END SUBROUTINE deepmd_energy_store_force_virial @@ -214,27 +191,17 @@ END SUBROUTINE deepmd_energy_store_force_virial !> \param use_virial ... ! ************************************************************************************************** SUBROUTINE deepmd_add_force_virial(fist_nonbond_env, force, pv_nonbond, use_virial) - TYPE(fist_nonbond_env_type), POINTER :: fist_nonbond_env - REAL(KIND=dp) :: force(:, :), pv_nonbond(3, 3) - LOGICAL, OPTIONAL :: use_virial + TYPE(fist_nonbond_env_type), POINTER :: fist_nonbond_env + REAL(KIND=dp) :: force(:, :), pv_nonbond(3, 3) + LOGICAL, OPTIONAL :: use_virial CHARACTER(LEN=*), PARAMETER :: routineN = 'deepmd_add_force_virial' - INTEGER :: handle -#ifdef __DEEPMD - INTEGER :: iat, iat_use - TYPE(deepmd_data_type), POINTER :: deepmd_data -#endif + INTEGER :: handle, iat, iat_use + TYPE(deepmd_data_type), POINTER :: deepmd_data CALL timeset(routineN, handle) -#ifndef __DEEPMD - CPABORT("In order to use DeePMD-kit you need to download and install libdeepmd_c") - MARK_USED(fist_nonbond_env) - MARK_USED(force) - MARK_USED(pv_nonbond) - MARK_USED(use_virial) -#else CALL fist_nonbond_env_get(fist_nonbond_env, deepmd_data=deepmd_data) IF (.NOT. ASSOCIATED(deepmd_data)) RETURN @@ -249,7 +216,6 @@ SUBROUTINE deepmd_add_force_virial(fist_nonbond_env, force, pv_nonbond, use_viri pv_nonbond = pv_nonbond + deepmd_data%virial END IF -#endif CALL timestop(handle) END SUBROUTINE deepmd_add_force_virial