!**************************************************************************** ! ! CDDL HEADER START ! ! The contents of this file are subject to the terms of the Common Development ! and Distribution License Version 1.0 (the "License"). ! ! You can obtain a copy of the license at ! http://www.opensource.org/licenses/CDDL-1.0. See the License for the ! specific language governing permissions and limitations under the License. ! ! When distributing Covered Code, include this CDDL HEADER in each file and ! include the License file in a prominent location with the name LICENSE.CDDL. ! If applicable, add the following below this CDDL HEADER, with the fields ! enclosed by brackets "[]" replaced with your own identifying information: ! ! Portions Copyright (c) [yyyy] [name of copyright owner]. All rights reserved. ! ! CDDL HEADER END ! ! ! Copyright (c) 2020, Regents of the University of Minnesota. ! All rights reserved. ! ! Contributors: ! Ellad B. Tadmor ! !**************************************************************************** !** !** MODULE EAM_Mendelev_2019_CuZr !** !** EAM potential for Cu-Zr by Mendelev (2019) !** !** Language: Fortran 2003 !** !** !**************************************************************************** module EAM_Mendelev_2019_CuZr use, intrinsic :: iso_c_binding use kim_model_headers_module implicit none save private public Compute_Energy_Forces, & destroy, & compute_arguments_create, & compute_arguments_destroy, & model_cutoff, & speccodeCu, & speccodeZr, & buffer_type ! Below are the definitions and values of all Model parameters integer(c_int), parameter :: cd = c_double ! used for literal constants integer(c_int), parameter :: DIM = 3 ! dimensionality of space integer(c_int), parameter :: speccodeCu = 1 ! internal species code for Cu integer(c_int), parameter :: speccodeZr = 2 ! internal species code for Zr real(c_double), parameter :: model_cutoff = 7.60_cd ! cutoff radius in angstroms real(c_double), parameter :: pair_cutoff(3) = (/6.00_cd, 7.60_cd, 7.60_cd/) ! pair cutoff for 1-1, 1-2, 2-2 real(c_double), parameter :: dens_cutoff(3) = (/6.00_cd, 6.00_cd, 5.60_cd/) ! density func cutoff for 1-1, 1-2, 2-2 real(c_double), parameter :: model_cutsq = model_cutoff**2 type, bind(c) :: buffer_type real(c_double) :: influence_distance real(c_double) :: cutoff(1) integer(c_int) :: & model_will_not_request_neighbors_of_noncontributing_particles(1) end type buffer_type contains !------------------------------------------------------------------------------- ! ! Calculate pair potential phi(r) and its derivative dphi(r) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine calc_phi_dphi(r,phi,dphi,spec1,spec2,yes_phi,yes_dphi) implicit none !-- Transferred variables real(c_double), intent(in) :: r real(c_double), intent(inout) :: phi,dphi integer(c_int), intent(in) :: spec1,spec2 logical, intent(in) :: yes_phi, yes_dphi !-- Local variables integer(c_int) typ ! Identify interaction type ! typ=1 for 1-1; typ=2 for 1-2 or 2-1; typ=3 for 2-2 typ = spec1+spec2-1 if (r .gt. pair_cutoff(typ)) then ! Argument exceeds cutoff radius if (yes_phi) phi = 0.0_cd if (yes_dphi) dphi = 0.0_cd return endif select case (typ) ! 1-1 interactions case(1) if (yes_phi) & phi = (1.2110196703541E+04/r)*(0.1818*exp(-2.9677556220903E+01*r)+0.5099* & exp(-8.7391128834241E+00*r)+0.2802*exp(-3.7365898129381E+00*r)+0.02817* & exp(-1.8696860419169E+00*r))*H0(1.0000-r)+exp(1.1026565103477E+01- & 1.0167211017722E+01*r+6.0017702915006E+00*r**2-1.9598299733506E+00*r**3)* & H(r-1.0000)*H0(1.8000-r)+(3.3519281301971E+00*(H(2.8-r)*(2.8-r)**4)- & 4.7447602323833E+01*(H(2.8-r)*(2.8-r)**5)+1.1106454537813E+02*(H(2.8-r)* & (2.8-r)**6)-1.2256379390195E+02*(H(2.8-r)*(2.8-r)**7)+4.9145722026502E+01*& (H(2.8-r)*(2.8-r)**8)+4.0605833179061E+00*(H(4.8-r)*(4.8-r)**4)+ & 2.5958091214976E+00*(H(4.8-r)*(4.8-r)**5)+5.5656640545299E+00*(H(4.8-r)* & (4.8-r)**6)+1.5184323060743E+00*(H(4.8-r)*(4.8-r)**7)+3.9696001635415E-01*& (H(4.8-r)*(4.8-r)**8)-2.1402913758299E-01*(H(6.0-r)*(6.0-r)**4)+ & 1.1714811538458E+00*(H(6.0-r)*(6.0-r)**5)-1.9913969426765E+00*(H(6.0-r)* & (6.0-r)**6)+1.3862043035438E+00*(H(6.0-r)*(6.0-r)**7)-3.4520315264743E-01*& (H(6.0-r)*(6.0-r)**8))*H(r-1.8000) if (yes_dphi) & dphi = -( & (1.2110196703541E+04/r)*((1/r+2.9677556220903E+01)*0.1818* & exp(-2.9677556220903E+01*r)+(1/r+8.7391128834241E+00)*0.5099* & exp(-8.7391128834241E+00*r)+(1/r+3.7365898129381E+00)*0.2802* & exp(-3.7365898129381E+00*r)+(1/r+1.8696860419169E+00)*0.02817* & exp(-1.8696860419169E+00*r))*H0(1.0000-r)-(-1.0167211017722E+01+ & 1.2003540583001E+01*r-5.8794899200518E+00*r**2)* & exp(1.1026565103477E+01-1.0167211017722E+01*r+6.0017702915006E+00* & r**2-1.9598299733506E+00*r**3)*H(r-1.0000)*H0(1.8000-r)+ & (3.3519281301971E+00*(H(2.8-r)*4*(2.8-r)**3)-4.7447602323833E+01* & (H(2.8-r)*5*(2.8-r)**4)+1.1106454537813E+02*(H(2.8-r)*6*(2.8-r)**5)- & 1.2256379390195E+02*(H(2.8-r)*7*(2.8-r)**6)+4.9145722026502E+01* & (H(2.8-r)*8*(2.8-r)**7)+4.0605833179061E+00*(H(4.8-r)*4*(4.8-r)**3)+ & 2.5958091214976E+00*(H(4.8-r)*5*(4.8-r)**4)+5.5656640545299E+00* & (H(4.8-r)*6*(4.8-r)**5)+1.5184323060743E+00*(H(4.8-r)*7*(4.8-r)**6)+ & 3.9696001635415E-01*(H(4.8-r)*8*(4.8-r)**7)-2.1402913758299E-01* & (H(6.0-r)*4*(6.0-r)**3)+1.1714811538458E+00*(H(6.0-r)*5*(6.0-r)**4)- & 1.9913969426765E+00*(H(6.0-r)*6*(6.0-r)**5)+1.3862043035438E+00* & (H(6.0-r)*7*(6.0-r)**6)-3.4520315264743E-01*(H(6.0-r)*8*(6.0-r)**7))* & H(r-1.8000) & ) ! 1-2 interactions case(2) if (yes_phi) & phi = (1.6703719591091E+04/r)*(0.1818*exp(-3.1401495864361E+01*r)+0.5099* & exp(-9.2467592353084E+00*r)+0.2802*exp(-3.9536445886722E+00*r)+0.02817* & exp(-1.9782942394547E+00*r))*H0(1.0000-r)+exp(9.1019930660972E+00- & 3.7511598234520E+00*r-5.0415926201351E-01*r**2+2.0120128716553E-01*r**3)* & H(r-1.0000)*H0(2.2000-r)+(6.0060951714678E+00*(H(4.0-r)*(4.0-r)**4)+ & 5.1271330965401E+00*(H(4.0-r)*(4.0-r)**5)+4.1248957113155E+00*(H(4.0-r)* & (4.0-r)**6)+1.8597579947697E+00*(H(4.0-r)*(4.0-r)**7)+1.3200992217759E-01*& (H(4.0-r)*(4.0-r)**8)-6.9673315120402E+00*(H(5.8-r)*(5.8-r)**4)+ & 1.4395210270079E+00*(H(5.8-r)*(5.8-r)**5)-6.5366946764648E+00*(H(5.8-r)* & (5.8-r)**6)+7.2008171787267E-01*(H(5.8-r)*(5.8-r)**7)-2.9778417035783E-01*& (H(5.8-r)*(5.8-r)**8)+3.9042809819224E-01*(H(7.6-r)*(7.6-r)**4)- & 1.1693902326951E+00*(H(7.6-r)*(7.6-r)**5)+1.2316955134740E+00*(H(7.6-r)* & (7.6-r)**6)-5.5333525392045E-01*(H(7.6-r)*(7.6-r)**7)+9.0516165873998E-02*& (H(7.6-r)*(7.6-r)**8))*H(r-2.2000) if (yes_dphi) & dphi = -( & (1.6703719591091E+04/r)*((1/r+3.1401495864361E+01)*0.1818* & exp(-3.1401495864361E+01*r)+(1/r+9.2467592353084E+00)*0.5099* & exp(-9.2467592353084E+00*r)+(1/r+3.9536445886722E+00)*0.2802* & exp(-3.9536445886722E+00*r)+(1/r+1.9782942394547E+00)*0.02817* & exp(-1.9782942394547E+00*r))*H0(1.0000-r)-(-3.7511598234520E+00- & 1.0083185240270E+00*r+6.0360386149660E-01*r**2)* & exp(9.1019930660972E+00-3.7511598234520E+00*r-5.0415926201351E-01* & r**2+2.0120128716553E-01*r**3)*H(r-1.0000)*H0(2.2000-r)+ & (6.0060951714678E+00*(H(4.0-r)*4*(4.0-r)**3)+5.1271330965401E+00* & (H(4.0-r)*5*(4.0-r)**4)+4.1248957113155E+00*(H(4.0-r)*6*(4.0-r)**5)+ & 1.8597579947697E+00*(H(4.0-r)*7*(4.0-r)**6)+1.3200992217759E-01* & (H(4.0-r)*8*(4.0-r)**7)-6.9673315120402E+00*(H(5.8-r)*4*(5.8-r)**3)+ & 1.4395210270079E+00*(H(5.8-r)*5*(5.8-r)**4)-6.5366946764648E+00* & (H(5.8-r)*6*(5.8-r)**5)+7.2008171787267E-01*(H(5.8-r)*7*(5.8-r)**6)- & 2.9778417035783E-01*(H(5.8-r)*8*(5.8-r)**7)+3.9042809819224E-01* & (H(7.6-r)*4*(7.6-r)**3)-1.1693902326951E+00*(H(7.6-r)*5*(7.6-r)**4)+ & 1.2316955134740E+00*(H(7.6-r)*6*(7.6-r)**5)-5.5333525392045E-01* & (H(7.6-r)*7*(7.6-r)**6)+9.0516165873998E-02*(H(7.6-r)*8*(7.6-r)**7))* & H(r-2.2000) & ) ! 2-2 interactions case(3) if (yes_phi) & phi = (2.3039613229091E+04/r)*(0.1818*exp(-3.3035595072498E+01*r)+0.5099* & exp(-9.7279503865046E+00*r)+0.2802*exp(-4.1593878920967E+00*r)+0.02817* & exp(-2.0812424895674E+00*r))*H0(1.0000-r)+exp(1.2333392311315E+01- & 1.0847321975247E+01*r+4.5733524457570E+00*r**2-8.5266291503854E-01*r**3)* & H(r-1.0000)*H0(2.3000-r)+(-1.4261501929757E+01*(H(3.5-r)*(3.5-r)**4)+ & 1.5850036758176E+01*(H(3.5-r)*(3.5-r)**5)-1.1325102264291E+01*(H(3.5-r)* & (3.5-r)**6)-4.0971114831366E+00*(H(3.5-r)*(3.5-r)**7)+3.6739378016909E+00* & (H(3.5-r)*(3.5-r)**8)+1.3066813393823E+00*(H(6.0-r)*(6.0-r)**4)- & 6.0542710718094E-01*(H(6.0-r)*(6.0-r)**5)+1.0055527194350E+00*(H(6.0-r)* & (6.0-r)**6)-1.4918186777562E-01*(H(6.0-r)*(6.0-r)**7)+3.2773112059590E-02* & (H(6.0-r)*(6.0-r)**8)+1.1433120304691E-02*(H(7.6-r)*(7.6-r)**4)- & 2.1982172508973E-02*(H(7.6-r)*(7.6-r)**5)-1.2542439692607E-02*(H(7.6-r)* & (7.6-r)**6)+2.5062673874258E-02*(H(7.6-r)*(7.6-r)**7)-7.5442887837418E-03* & (H(7.6-r)*(7.6-r)**8))*H(r-2.3000) if (yes_dphi) & dphi = -( & (2.3039613229091E+04/r)*((1/r+3.3035595072498E+01)*0.1818* & exp(-3.3035595072498E+01*r)+(1/r+9.7279503865046E+00)*0.5099* & exp(-9.7279503865046E+00*r)+(1/r+4.1593878920967E+00)*0.2802* & exp(-4.1593878920967E+00*r)+(1/r+2.0812424895674E+00)*0.02817* & exp(-2.0812424895674E+00*r))*H0(1.0000-r)-(-1.0847321975247E+01+ & 9.1467048915139E+00*r-2.5579887451156E+00*r**2)* & exp(1.2333392311315E+01-1.0847321975247E+01*r+4.5733524457570E+00* & r**2-8.5266291503854E-01*r**3)*H(r-1.0000)*H0(2.3000-r)+ & (-1.4261501929757E+01*(H(3.5-r)*4*(3.5-r)**3)+1.5850036758176E+01* & (H(3.5-r)*5*(3.5-r)**4)-1.1325102264291E+01*(H(3.5-r)*6*(3.5-r)**5)- & 4.0971114831366E+00*(H(3.5-r)*7*(3.5-r)**6)+3.6739378016909E+00* & (H(3.5-r)*8*(3.5-r)**7)+1.3066813393823E+00*(H(6.0-r)*4*(6.0-r)**3)- & 6.0542710718094E-01*(H(6.0-r)*5*(6.0-r)**4)+1.0055527194350E+00* & (H(6.0-r)*6*(6.0-r)**5)-1.4918186777562E-01*(H(6.0-r)*7*(6.0-r)**6)+ & 3.2773112059590E-02*(H(6.0-r)*8*(6.0-r)**7)+1.1433120304691E-02* & (H(7.6-r)*4*(7.6-r)**3)-2.1982172508973E-02*(H(7.6-r)*5*(7.6-r)**4)- & 1.2542439692607E-02*(H(7.6-r)*6*(7.6-r)**5)+2.5062673874258E-02* & (H(7.6-r)*7*(7.6-r)**6)-7.5442887837418E-03*(H(7.6-r)*8*(7.6-r)**7))* & H(r-2.3000) & ) end select return end subroutine calc_phi_dphi !------------------------------------------------------------------------------- ! ! Calculate electron density g(r) and its derivative dg(r) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine calc_g_dg(r,g,dg,spec1,spec2,yes_g,yes_dg) implicit none !-- Transferred variables real(c_double), intent(in) :: r real(c_double), intent(inout) :: g,dg integer(c_int), intent(in) :: spec1,spec2 logical, intent(in) :: yes_g,yes_dg !-- Local variables integer(c_int) typ ! Identify interaction type ! typ=1 for 1-1; typ=2 for 1-2 or 2-1; typ=3 for 2-2 typ = spec1+spec2-1 if (r .gt. dens_cutoff(typ)) then ! Argument exceeds cutoff radius if (yes_g) g = 0.0_cd if (yes_dg) dg = 0.0_cd return endif select case (typ) ! 1-1 interactions case(1) if (yes_g) & g = 1.539907037884*(1.9999999875362E-02*(H(2.4-r)*(2.4-r)**4)+ & 1.9987533420669E-02*(H(3.2-r)*(3.2-r)**4)+1.8861676713565E-02* & (H(4.5-r)*(4.5-r)**4)+6.6082982694659E-03*(H(6.0-r)*(6.0-r)**4)) if (yes_dg) & dg = -( & 1.539907037884*(1.9999999875362E-02*(H(2.4-r)*4*(2.4-r)**3)+ & 1.9987533420669E-02*(H(3.2-r)*4*(3.2-r)**3)+1.8861676713565E-02* & (H(4.5-r)*4*(4.5-r)**3)+6.6082982694659E-03*(H(6.0-r)*4*(6.0-r)**3)) & ) ! 1-2 interactions case(2) if (yes_g) & g = 1.0779076111400E-01*(H(2.8-r)*(2.8-r)**4)+1.0882724199700E-01* & (H(3.2-r)*(3.2-r)**4)+1.6470644009500E-01*(H(4.4-r)*(4.4-r)**4)+ & 2.3681839912700E-02*(H(6.0-r)*(6.0-r)**4) if (yes_dg) & dg = -( & 1.0779076111400E-01*(H(2.8-r)*4*(2.8-r)**3)+1.0882724199700E-01* & (H(3.2-r)*4*(3.2-r)**3)+1.6470644009500E-01*(H(4.4-r)*4*(4.4-r)**3)+ & 2.3681839912700E-02*(H(6.0-r)*4*(6.0-r)**3) & ) ! 2-2 interactions case(3) if (yes_g) & g = 0.901364818051*(7.7718711248373E-01*(H(5.6-r)*(5.6-r)**4)- & 4.8102928454986E-01*(H(5.6-r)*(5.6-r)**5)+1.4501312593993E-01*(H(5.6-r)* & (5.6-r)**6)-2.1292226813959E-02*(H(5.6-r)*(5.6-r)**7)+1.2209217625670E-03*& (H(5.6-r)*(5.6-r)**8)) if (yes_dg) & dg = -( & 0.901364818051*(7.7718711248373E-01*(H(5.6-r)*4*(5.6-r)**3)- & 4.8102928454986E-01*(H(5.6-r)*5*(5.6-r)**4)+1.4501312593993E-01* & (H(5.6-r)*6*(5.6-r)**5)-2.1292226813959E-02*(H(5.6-r)*7*(5.6-r)**6)+ & 1.2209217625670E-03*(H(5.6-r)*8*(5.6-r)**7)) & ) end select return end subroutine calc_g_dg !------------------------------------------------------------------------------- ! ! Calculate embedding function U(rho) and first derivative dU(rho) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine calc_U_dU(ro,U,dU,spec,yes_U,yes_dU) implicit none !-- Transferred variables real(c_double), intent(in) :: ro real(c_double), intent(inout) :: U,dU integer(c_int), intent(in) :: spec logical, intent(in) :: yes_U,yes_dU select case (spec) ! species 1 case(1) if (yes_U) & U = 1.0000000000000E+00*(-(ro/1.539907037884)**0.5)-5.7112865649408E-05* & (H((ro/1.539907037884)-9)*((ro/1.539907037884)-9)**4)+ & 3.0303487333648E-04*(H((ro/1.539907037884)-11)*((ro/1.539907037884)-11)**4)- & 5.4720795296134E-04*(H((ro/1.539907037884)-13)*((ro/1.539907037884)-13)**4)+ & 4.6278681464721E-04*(H((ro/1.539907037884)-15)*((ro/1.539907037884)-15)**4)- & 1.0310712451906E-03*(H((ro/1.539907037884)-16)*((ro/1.539907037884)-16)**4)+ & 3.0634000239833E-03*(H((ro/1.539907037884)-16.5)* & ((ro/1.539907037884)-16.5)**4)-2.8308102136994E-03*(H((ro/1.539907037884)-17)* & ((ro/1.539907037884)-17)**4)+6.4044567482688E-04*(H((ro/1.539907037884)-18)* & ((ro/1.539907037884)-18)**4)+1.9890362320891E-06*(H((ro/1.539907037884)-20)* & ((ro/1.539907037884)-20)**4) if (yes_dU) & dU = (1.0000000000000E+00*(-0.5/(ro/1.539907037884)**0.5)-5.7112865649408E-05* & (H((ro/1.539907037884)-9)*4*((ro/1.539907037884)-9)**3)+3.0303487333648E-04* & (H((ro/1.539907037884)-11)*4*((ro/1.539907037884)-11)**3)-5.4720795296134E-04* & (H((ro/1.539907037884)-13)*4*((ro/1.539907037884)-13)**3)+4.6278681464721E-04* & (H((ro/1.539907037884)-15)*4*((ro/1.539907037884)-15)**3)-1.0310712451906E-03* & (H((ro/1.539907037884)-16)*4*((ro/1.539907037884)-16)**3)+3.0634000239833E-03* & (H((ro/1.539907037884)-16.5)*4*((ro/1.539907037884)-16.5)**3)-2.8308102136994E-03* & (H((ro/1.539907037884)-17)*4*((ro/1.539907037884)-17)**3)+6.4044567482688E-04* & (H((ro/1.539907037884)-18)*4*((ro/1.539907037884)-18)**3))/1.539907037884+ & 1.9890362320891E-06*((H((ro/1.539907037884)-20)*4*((ro/1.539907037884)-20)**3)/ & 1.539907037884) ! species 2 case(2) if (yes_U) & U = 1.0000000000000E+00*(-(ro/0.901364818051)**0.5)-1.9162462126235E-07* & (H((ro/0.901364818051)-60)*((ro/0.901364818051)-60)**4)+4.6418727035037E-07* & (H((ro/0.901364818051)-70)*((ro/0.901364818051)-70)**4)+6.6448294272955E-07* & (H((ro/0.901364818051)-80)*((ro/0.901364818051)-80)**4)-2.0680252960229E-06* & (H((ro/0.901364818051)-85)*((ro/0.901364818051)-85)**4)+1.1387131464983E-06* & (H((ro/0.901364818051)-90)*((ro/0.901364818051)-90)**4)-1.8281354613667E-09* & (H((ro/0.901364818051)-100)*((ro/0.901364818051)-100)**4) if (yes_dU) & dU = (1.0000000000000E+00*(-0.5/(ro/0.901364818051)**0.5)-1.9162462126235E-07* & (H((ro/0.901364818051)-60)*4*((ro/0.901364818051)-60)**3)+4.6418727035037E-07* & (H((ro/0.901364818051)-70)*4*((ro/0.901364818051)-70)**3)+6.6448294272955E-07* & (H((ro/0.901364818051)-80)*4*((ro/0.901364818051)-80)**3)-2.0680252960229E-06* & (H((ro/0.901364818051)-85)*4*((ro/0.901364818051)-85)**3)+1.1387131464983E-06* & (H((ro/0.901364818051)-90)*4*((ro/0.901364818051)-90)**3))/0.901364818051- & 1.8281354613667E-09*((H((ro/0.901364818051)-100)*4*((ro/0.901364818051)-100)**3)/ & 0.901364818051) end select return end subroutine calc_U_dU !------------------------------------------------------------------------------- ! ! Heaviside function H (1 at x=0) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive real(c_double) function H(x) implicit none !-- Transferred variables real(c_double), intent(in) :: x if (x<0.0_cd) then H = 0.0_cd else H = 1.0_cd endif end function H !------------------------------------------------------------------------------- ! ! Heaviside function H0 (0 at x=0) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive real(c_double) function H0(x) implicit none !-- Transferred variables real(c_double), intent(in) :: x if (x<=0.0_cd) then H0 = 0.0_cd else H0 = 1.0_cd endif end function H0 !------------------------------------------------------------------------------- ! ! Compute energy and forces on atoms from the positions. ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine Compute_Energy_Forces(model_compute_handle, & model_compute_arguments_handle, ierr) bind(c) implicit none !-- Transferred variables type(kim_model_compute_handle_type), intent(in) :: model_compute_handle type(kim_model_compute_arguments_handle_type), intent(in) :: & model_compute_arguments_handle integer(c_int), intent(out) :: ierr !-- Local variables real(c_double) :: Rij(DIM) real(c_double) :: r,Rsqij,phi,dphi,g,dg,dU,U,dphieff real(c_double) :: dphii,dUi,dphij,dUj integer(c_int) :: i,j,jj,numnei,comp_force,comp_particleEnergy,comp_virial,comp_energy integer(c_int) :: ierr2 real(c_double), allocatable :: rho(:), derU(:) integer(c_int) :: speci,specj real(c_double) :: dum !-- KIM variables integer(c_int), pointer :: N real(c_double), pointer :: energy real(c_double), pointer :: coor(:,:) real(c_double), pointer :: force(:,:) real(c_double), pointer :: particleEnergy(:) integer(c_int), pointer :: nei1part(:) integer(c_int), pointer :: particleSpeciesCodes(:) integer(c_int), pointer :: particleContributing(:) real(c_double), pointer :: virial(:) ! Avoid unused argument warning if (model_compute_handle .eq. KIM_MODEL_COMPUTE_NULL_HANDLE) continue ! Unpack data from KIM object ! ierr = 0 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_NUMBER_OF_PARTICLES, N, ierr2) ierr = ierr + ierr2 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTICLE_SPECIES_CODES, N, particleSpeciesCodes, & ierr2) ierr = ierr + ierr2 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTICLE_CONTRIBUTING, N, particleContributing, & ierr2) ierr = ierr + ierr2 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_COORDINATES, dim, N, coor, ierr2) ierr = ierr + ierr2 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_ENERGY, energy, ierr2) ierr = ierr + ierr2 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_FORCES, dim, N, force, ierr2) ierr = ierr + ierr2 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_PARTICLE_ENERGY, n, particleEnergy, ierr2) ierr = ierr + ierr2 call kim_get_argument_pointer( & model_compute_arguments_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_VIRIAL, 6, virial, ierr2) ierr = ierr + ierr2 if (ierr /= 0) then call kim_log_entry(model_compute_arguments_handle, & KIM_LOG_VERBOSITY_ERROR, & "Failed to retrieve data from KIM API compute-arguments object") return endif ! Check to see if we have been asked to compute the forces, energyperpart, ! energy and virial ! if (associated(energy)) then comp_energy = 1 else comp_energy = 0 end if if (associated(force)) then comp_force = 1 else comp_force = 0 end if if (associated(particleEnergy)) then comp_particleEnergy = 1 else comp_particleEnergy = 0 end if if (associated(virial)) then comp_virial = 1 else comp_virial = 0 end if ! Check to be sure that the species are correct ! ierr = 1 ! assume an error do i = 1,N if (particleSpeciesCodes(i).ne.speccodeCu .and. particleSpeciesCodes(i).ne.speccodeZr) then call kim_log_entry(model_compute_arguments_handle, & KIM_LOG_VERBOSITY_ERROR, "Unexpected species code detected") return endif enddo ierr = 0 ! everything is ok ! Initialize potential energies, forces, virial term, and electron density ! ! Note: that the variable `particleEnergy' does not need to be initialized ! because it's initial value is set during the embedding energy ! calculation. if (comp_energy.eq.1) energy = 0.0_cd if (comp_force.eq.1) force = 0.0_cd if (comp_virial.eq.1) virial = 0.0_cd dphieff = 0.0_cd allocate( rho(N) ) ! pair functional electron density rho = 0.0_cd ! EAM embedded energy deriv if (comp_force.eq.1.or.comp_virial.eq.1) allocate( derU(N) ) ! Loop over particles in the neighbor list a first time, ! to compute electron density (=coordination) ! do i = 1,N if(particleContributing(i).eq.1) then ! Get particle species speci = particleSpeciesCodes(i) ! Get neighbor list of current atom call kim_get_neighbor_list( & model_compute_arguments_handle, 1, i, numnei, nei1part, ierr) if (ierr /= 0) then ! some sort of problem, exit call kim_log_entry(model_compute_arguments_handle, & KIM_LOG_VERBOSITY_ERROR, "GetNeighborList failed") ierr = 1 return endif ! Loop over the neighbors of atom i and compute charge density ! do jj = 1, numnei j = nei1part(jj) ! get neighbor ID if ( .not.( (particleContributing(j).eq.1) .and. & j.lt.i) ) then ! effective half list ! compute relative position vector ! Rij(:) = coor(:,j) - coor(:,i) ! distance vector between i j ! compute contribution to electron density ! Rsqij = dot_product(Rij,Rij) ! compute square distance if ( Rsqij .lt. model_cutsq ) then ! particles are interacting? ! Get particle species specj = particleSpeciesCodes(j) r = sqrt(Rsqij) ! compute distance call calc_g_dg(r,g,dum,speci,specj,.true.,.false.) ! compute electron density rho(i) = rho(i) + g ! accumulate electron density if (particleContributing(j).eq.1) then rho(j) = rho(j) + g ! WARNING: This model assumes that the charge density ! is symmetric between species, i.e. ! g_{\alpha\beta}(r) = g_{\beta\alpha}(r) ! where \alpha and \beta are species. ! If this is not the case, calc_g_dg would need ! to be called again with the species order ! reversed to get the contribution for atom j. endif endif endif ! if ( i.lt.j ) enddo ! loop on jj endif ! Check on whether particle is contributing enddo ! infinite do loop (terminated by exit statements above) ! Now that we know the electron densities, calculate embedding part of energy ! U and its derivative U' (derU) ! do i = 1,N if(particleContributing(i).eq.1) then ! Get particle species speci = particleSpeciesCodes(i) if (comp_force.eq.1.or.comp_virial.eq.1) then call calc_U_dU(rho(i),U,dU,speci,.true.,.true.) ! compute embedding energy ! and its derivative derU(i) = dU ! store dU for later use else call calc_U_dU(rho(i),U,dum,speci,.true.,.false.) ! compute just embedding energy endif ! accumulate the embedding energy contribution ! ! Assuming U(rho=0) = 0.0_cd ! if (comp_particleEnergy.eq.1) then ! accumulate embedding energy contribution particleEnergy(i) = U endif if (comp_energy.eq.1) then energy = energy + U endif endif ! Check on whether particle is contributing enddo ! Loop over particles in the neighbor list a second time, to compute ! the forces and complete energy calculation ! do i = 1,N if(particleContributing(i).eq.1) then ! Get particle species speci = particleSpeciesCodes(i) ! Get neighbor list of current atom call kim_get_neighbor_list( & model_compute_arguments_handle, 1, i, numnei, nei1part, ierr) if (ierr /= 0) then ! some sort of problem, exit call kim_log_entry(model_compute_arguments_handle, & KIM_LOG_VERBOSITY_ERROR, "GetNeighborList failed") ierr = 1 return endif ! Loop over the neighbors of atom i ! do jj = 1, numnei j = nei1part(jj) ! get neighbor ID if ( .not.( (particleContributing(j).eq.1) .and. & j.lt.i) ) then ! effective half list ! compute relative position vector ! Rij(:) = coor(:,j) - coor(:,i) ! distance vector between i j ! compute energy and forces ! Rsqij = dot_product(Rij,Rij) ! compute square distance if ( Rsqij .lt. model_cutsq ) then ! particles are interacting? ! Get particle species specj = particleSpeciesCodes(j) r = sqrt(Rsqij) ! compute distance if (comp_force.eq.1.or.comp_virial.eq.1) then call calc_phi_dphi(r,phi,dphi,speci,specj,.true.,.true.) ! compute pair potential ! and it derivative call calc_g_dg(r,dum,dg,speci,specj,.false.,.true.) ! compute elect dens first deriv if (particleContributing(j).eq.1) then dphii = 0.5_cd*dphi dphij = 0.5_cd*dphi dUi = derU(i)*dg dUj = derU(j)*dg ! WARNING: This model assumes that the else ! charge density is symmetric dphii = 0.5_cd*dphi ! between species. See comment dphij = 0.0_cd ! above. dUi = derU(i)*dg dUj = 0.0_cd endif dphieff = dphii + dphij + dUi + dUj else call calc_phi_dphi(r,phi,dum,speci,specj,.true.,.false.) ! compute just pair potential endif ! contribution to energy ! if (comp_particleEnergy.eq.1) then particleEnergy(i) = particleEnergy(i) + 0.5_cd*phi ! accumulate energy if (particleContributing(j).eq.1) then particleEnergy(j) = particleEnergy(j) + 0.5_cd*phi ! accumulate energy endif endif if (comp_energy.eq.1) then if (particleContributing(j).eq.1) then energy = energy + phi ! accumulate energy else energy = energy + 0.5_cd*phi endif endif ! contribution to virial tensor ! if (comp_virial.eq.1) then virial(1) = virial(1) + Rij(1)*Rij(1)*dphieff/r virial(2) = virial(2) + Rij(2)*Rij(2)*dphieff/r virial(3) = virial(3) + Rij(3)*Rij(3)*dphieff/r virial(4) = virial(4) + Rij(2)*Rij(3)*dphieff/r virial(5) = virial(5) + Rij(1)*Rij(3)*dphieff/r virial(6) = virial(6) + Rij(1)*Rij(2)*dphieff/r endif ! contribution to forces ! if (comp_force.eq.1) then force(:,i) = force(:,i) + dphieff*Rij/r ! accumulate force on atom i force(:,j) = force(:,j) - dphieff*Rij/r ! accumulate force on atom j endif endif endif ! if ( i.lt.j ) enddo ! loop on jj endif ! Check on whether particle is contributing enddo ! infinite do loop (terminated by exit statements above) ! Free temporary storage ! deallocate( rho ) if (comp_force.eq.1.or.comp_virial.eq.1) deallocate( derU ) ! Everything is great ierr = 0 return end subroutine Compute_Energy_Forces !------------------------------------------------------------------------------- ! ! Model destroy routine ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine destroy(model_destroy_handle, ierr) bind(c) use, intrinsic :: iso_c_binding implicit none !-- Transferred variables type(kim_model_destroy_handle_type), intent(inout) :: model_destroy_handle integer(c_int), intent(out) :: ierr type(buffer_type), pointer :: buf; type(c_ptr) :: pbuf call kim_get_model_buffer_pointer(model_destroy_handle, pbuf) call c_f_pointer(pbuf, buf) call kim_log_entry(model_destroy_handle, KIM_LOG_VERBOSITY_INFORMATION, & "deallocating model buffer") deallocate(buf) ierr = 0 ! everything is good return end subroutine destroy !------------------------------------------------------------------------------- ! ! Model compute arguments create routine (REQUIRED) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine compute_arguments_create(model_compute_handle, & model_compute_arguments_create_handle, ierr) bind(c) use, intrinsic :: iso_c_binding use kim_model_compute_arguments_create_module implicit none !-- Transferred variables type(kim_model_compute_handle_type), intent(in) :: model_compute_handle type(kim_model_compute_arguments_create_handle_type), intent(inout) :: & model_compute_arguments_create_handle integer(c_int), intent(out) :: ierr integer(c_int) :: ierr2 ! Avoid unused argument warning if (model_compute_handle .eq. KIM_MODEL_COMPUTE_NULL_HANDLE) continue ierr = 0 ierr2 = 0 ! register arguments call kim_set_argument_support_status( & model_compute_arguments_create_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_ENERGY, & KIM_SUPPORT_STATUS_OPTIONAL, ierr2) ierr = ierr + ierr2 call kim_set_argument_support_status( & model_compute_arguments_create_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_FORCES, & KIM_SUPPORT_STATUS_OPTIONAL, ierr2) ierr = ierr + ierr2 call kim_set_argument_support_status( & model_compute_arguments_create_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_PARTICLE_ENERGY, & KIM_SUPPORT_STATUS_OPTIONAL, ierr2) ierr = ierr + ierr2 call kim_set_argument_support_status( & model_compute_arguments_create_handle, & KIM_COMPUTE_ARGUMENT_NAME_PARTIAL_VIRIAL, & KIM_SUPPORT_STATUS_OPTIONAL, ierr2) ierr = ierr + ierr2 ! register call backs (ProcessDEDrTerm, ProcessD2EDr2Term ! NONE if (ierr /= 0) then ierr = 1 call kim_log_entry(model_compute_arguments_create_handle, & KIM_LOG_VERBOSITY_ERROR, & "Unable to successfully create compute-arguments object") endif return end subroutine compute_arguments_create !------------------------------------------------------------------------------- ! ! Model compute arguments destroy routine (REQUIRED) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine compute_arguments_destroy(model_compute_handle, & model_compute_arguments_destroy_handle, ierr) bind(c) use, intrinsic :: iso_c_binding implicit none !-- Transferred variables type(kim_model_compute_handle_type), intent(in) :: model_compute_handle type(kim_model_compute_arguments_destroy_handle_type), intent(inout) :: & model_compute_arguments_destroy_handle integer(c_int), intent(out) :: ierr ! avoid unused dummy argument warnings if (model_compute_handle .eq. KIM_MODEL_COMPUTE_NULL_HANDLE) continue if (model_compute_arguments_destroy_handle .eq. & KIM_MODEL_COMPUTE_ARGUMENTS_DESTROY_NULL_HANDLE) continue ierr = 0 return end subroutine compute_arguments_destroy end module EAM_Mendelev_2019_CuZr !------------------------------------------------------------------------------- ! ! Model initialization routine (REQUIRED) ! !------------------------------------------------------------------------------- ! The "recursive" keyword is added below make this routine thread safe recursive subroutine create(model_create_handle, requested_length_unit, & requested_energy_unit, requested_charge_unit, requested_temperature_unit, & requested_time_unit, ierr) bind(c) use, intrinsic :: iso_c_binding use EAM_Mendelev_2019_CuZr use kim_model_headers_module implicit none !-- Transferred variables type(kim_model_create_handle_type), intent(inout) :: model_create_handle type(kim_length_unit_type), intent(in) :: requested_length_unit type(kim_energy_unit_type), intent(in) :: requested_energy_unit type(kim_charge_unit_type), intent(in) :: requested_charge_unit type(kim_temperature_unit_type), intent(in) :: requested_temperature_unit type(kim_time_unit_type), intent(in) :: requested_time_unit integer(c_int), intent(out) :: ierr !-- KIM variables integer(c_int) :: ierr2 type(buffer_type), pointer :: buf ierr = 0 ierr2 = 0 call kim_set_units(model_create_handle, & KIM_LENGTH_UNIT_A, & KIM_ENERGY_UNIT_EV, & KIM_CHARGE_UNIT_UNUSED, & KIM_TEMPERATURE_UNIT_UNUSED, & KIM_TIME_UNIT_UNUSED, & ierr2) ierr = ierr + ierr2 ! register species call kim_set_species_code(model_create_handle, & KIM_SPECIES_NAME_CU, speccodeCu, ierr2) ierr = ierr + ierr2 call kim_set_species_code(model_create_handle, & KIM_SPECIES_NAME_ZR, speccodeZr, ierr2) ierr = ierr + ierr2 ! register numbering call kim_set_model_numbering(model_create_handle, & KIM_NUMBERING_ONE_BASED, ierr2); ierr = ierr + ierr2 ! register function pointers call kim_set_routine_pointer(model_create_handle, & KIM_MODEL_ROUTINE_NAME_COMPUTE, KIM_LANGUAGE_NAME_FORTRAN, & 1, c_funloc(Compute_Energy_Forces), ierr2) ierr = ierr + ierr2 call kim_set_routine_pointer( & model_create_handle, KIM_MODEL_ROUTINE_NAME_COMPUTE_ARGUMENTS_CREATE, & KIM_LANGUAGE_NAME_FORTRAN, 1, c_funloc(compute_arguments_create), ierr2) ierr = ierr + ierr2 call kim_set_routine_pointer( & model_create_handle, KIM_MODEL_ROUTINE_NAME_COMPUTE_ARGUMENTS_DESTROY, & KIM_LANGUAGE_NAME_FORTRAN, 1, c_funloc(compute_arguments_destroy), ierr2) ierr = ierr + ierr2 call kim_set_routine_pointer(model_create_handle, & KIM_MODEL_ROUTINE_NAME_DESTROY, KIM_LANGUAGE_NAME_FORTRAN, & 1, c_funloc(destroy), ierr2) ierr = ierr + ierr2 ! allocate buffer allocate(buf) ! store model buffer in KIM object call kim_set_model_buffer_pointer(model_create_handle, & c_loc(buf)) ! set buffer values buf%influence_distance = model_cutoff buf%cutoff(1) = model_cutoff buf%model_will_not_request_neighbors_of_noncontributing_particles(1) = 1 ! register influence distance call kim_set_influence_distance_pointer( & model_create_handle, buf%influence_distance) ! register cutoff call kim_set_neighbor_list_pointers(model_create_handle, & 1, buf%cutoff, & buf%model_will_not_request_neighbors_of_noncontributing_particles) if (ierr /= 0) then ierr = 1 call kim_log_entry(model_create_handle, KIM_LOG_VERBOSITY_ERROR, & "Unable to successfully initialize model") endif return end subroutine create