KIM API V2
ex_model_Ar_P_MLJ_F03.F03
Go to the documentation of this file.
1 !
2 ! CDDL HEADER START
3 !
4 ! The contents of this file are subject to the terms of the Common Development
5 ! and Distribution License Version 1.0 (the "License").
6 !
7 ! You can obtain a copy of the license at
8 ! http://www.opensource.org/licenses/CDDL-1.0. See the License for the
9 ! specific language governing permissions and limitations under the License.
10 !
11 ! When distributing Covered Code, include this CDDL HEADER in each file and
12 ! include the License file in a prominent location with the name LICENSE.CDDL.
13 ! If applicable, add the following below this CDDL HEADER, with the fields
14 ! enclosed by brackets "[]" replaced with your own identifying information:
15 !
16 ! Portions Copyright (c) [yyyy] [name of copyright owner]. All rights reserved.
17 !
18 ! CDDL HEADER END
19 !
20 
21 !
22 ! Copyright (c) 2013--2018, Regents of the University of Minnesota.
23 ! All rights reserved.
24 !
25 ! Contributors:
26 ! Ryan S. Elliott
27 ! Ellad B. Tadmor
28 ! Valeriu Smirichinski
29 ! Stephen M. Whalen
30 !
31 
32 !****************************************************************************
33 !**
34 !** MODULE ex_model_Ar_P_MLJ_F03
35 !**
36 !** Modified Lennard-Jones pair potential (with smooth cutoff) model for Ar
37 !**
38 !** Reference: Ashcroft and Mermin
39 !**
40 !** Language: Fortran 2003
41 !**
42 !****************************************************************************
43 
44 
45 #include "KIM_API_status.h"
46 #define THIS_FILE_NAME __FILE__
47 #define TRUEFALSE(TRUTH) merge(1,0,(TRUTH))
48 
50 
51 use, intrinsic :: iso_c_binding
52 use kim_api_f03
53 implicit none
54 
55 save
56 private
57 public compute_energy_forces, &
59 
60 ! Below are the definitions and values of all Model parameters
61 integer(c_int), parameter :: cd = c_double ! used for literal constants
62 integer(c_int), parameter :: DIM = 3 ! dimensionality of space
63 integer(c_int), parameter :: speccode = 1 ! internal species code
64 real(c_double), parameter :: model_cutoff = 8.15_cd ! cutoff radius
65  ! in angstroms
66 real(c_double), parameter :: model_cutsq = model_cutoff**2
67 
68 !-------------------------------------------------------------------------------
69 ! Below are the definitions and values of all additional model parameters
70 !
71 ! Recall that the Fortran 2003 format for declaring parameters is as follows:
72 !
73 ! integer(c_int), parameter :: parname = value ! This defines an integer
74 ! ! parameter called `parname'
75 ! ! with a value equal to
76 ! ! `value' (a number)
77 !
78 ! real(c_double), parameter :: parname = value ! This defines a real(c_double)
79 ! ! parameter called `parname'
80 ! ! with a value equal to
81 ! ! `value' (a number)
82 !-------------------------------------------------------------------------------
83 real(c_double), parameter :: lj_epsilon = 0.0104_cd
84 real(c_double), parameter :: lj_sigma = 3.40_cd
85 real(c_double), parameter :: lj_cutnorm = model_cutoff/lj_sigma
86 real(c_double), parameter :: lj_A = 12.0_cd*lj_epsilon*(-26.0_cd &
87  + 7.0_cd*lj_cutnorm**6)/(lj_cutnorm**14 &
88  *lj_sigma**2)
89 real(c_double), parameter :: lj_B = 96.0_cd*lj_epsilon*(7.0_cd &
90  - 2.0_cd*lj_cutnorm**6)/(lj_cutnorm**13*lj_sigma)
91 real(c_double), parameter :: lj_C = 28.0_cd*lj_epsilon*(-13.0_cd &
92  + 4.0_cd*lj_cutnorm**6)/(lj_cutnorm**12)
93 
94 contains
95 
96 !-------------------------------------------------------------------------------
97 !
98 ! Calculate pair potential phi(r)
99 !
100 !-------------------------------------------------------------------------------
101 subroutine calc_phi(r,phi)
102 implicit none
103 
104 !-- Transferred variables
105 real(c_double), intent(in) :: r
106 real(c_double), intent(out) :: phi
107 
108 !-- Local variables
109 real(c_double) rsq,sor,sor6,sor12
110 
111 rsq = r*r ! r^2
112 sor = lj_sigma/r ! (sig/r)
113 sor6 = sor*sor*sor !
114 sor6 = sor6*sor6 ! (sig/r)^6
115 sor12= sor6*sor6 ! (sig/r)^12
116 if (r .gt. model_cutoff) then
117  ! Argument exceeds cutoff radius
118  phi = 0.0_cd
119 else
120  phi = 4.0_cd*lj_epsilon*(sor12-sor6) + lj_a*rsq + lj_b*r + lj_c
121 endif
122 
123 end subroutine calc_phi
124 
125 !-------------------------------------------------------------------------------
126 !
127 ! Calculate pair potential phi(r) and its derivative dphi(r)
128 !
129 !-------------------------------------------------------------------------------
130 subroutine calc_phi_dphi(r,phi,dphi)
131 implicit none
132 
133 !-- Transferred variables
134 real(c_double), intent(in) :: r
135 real(c_double), intent(out) :: phi,dphi
136 
137 !-- Local variables
138 real(c_double) rsq,sor,sor6,sor12
139 
140 rsq = r*r ! r^2
141 sor = lj_sigma/r ! (sig/r)
142 sor6 = sor*sor*sor !
143 sor6 = sor6*sor6 ! (sig/r)^6
144 sor12= sor6*sor6 ! (sig/r)^12
145 if (r .gt. model_cutoff) then
146  ! Argument exceeds cutoff radius
147  phi = 0.0_cd
148  dphi = 0.0_cd
149 else
150  phi = 4.0_cd*lj_epsilon*(sor12-sor6) + lj_a*rsq + lj_b*r + lj_c
151  dphi = 24.0_cd*lj_epsilon*(-2.0_cd*sor12+sor6)/r + 2.0_cd*lj_a*r + lj_b
152 endif
153 
154 end subroutine calc_phi_dphi
155 
156 !-------------------------------------------------------------------------------
157 !
158 ! Compute energy and forces on particles from the positions.
159 !
160 !-------------------------------------------------------------------------------
161 integer(c_int) function compute_energy_forces(pkim) bind(c)
162 implicit none
163 
164 !-- Transferred variables
165 type(c_ptr), intent(in) :: pkim
166 
167 !-- Local variables
168 real(c_double) :: Rij(DIM)
169 real(c_double) :: r,Rsqij,phi,dphi,dEidr = 0.0_cd
170 integer(c_int) :: i,j,jj,numnei,part_ret,comp_force,comp_enepot,comp_virial, &
171  comp_energy
172 character (len=80) :: error_message
173 
174 !-- KIM variables
175 integer(c_int), pointer :: N; type(c_ptr) :: pN
176 real(c_double), pointer :: energy; type(c_ptr) :: penergy
177 real(c_double), pointer :: coor(:,:); type(c_ptr) :: pcoor
178 real(c_double), pointer :: force(:,:); type(c_ptr) :: pforce
179 real(c_double), pointer :: enepot(:); type(c_ptr) :: penepot
180 real(c_double), pointer :: Rij_list(:,:); type(c_ptr) :: pRij_list
181 integer(c_int), pointer :: nei1part(:); type(c_ptr) :: pnei1part
182 integer(c_int), pointer :: particleSpecies(:);type(c_ptr) :: pparticleSpecies
183 real(c_double), pointer :: virial(:); type(c_ptr) :: pvirial
184 integer(c_int) idum
185 
186 
187 ! Check to see if we have been asked to compute the forces, energyperpart,
188 ! energy and virial
189 !
190 call kim_api_getm_compute(pkim, compute_energy_forces, &
191  "energy", comp_energy, 1, &
192  "forces", comp_force, 1, &
193  "particleEnergy", comp_enepot, 1, &
194  "virial", comp_virial, 1)
195 if (compute_energy_forces.lt.kim_status_ok) then
196  idum = kim_api_report_error(__line__, this_file_name, &
197  "kim_api_getm_compute", compute_energy_forces)
198  return
199 endif
200 
201 ! Unpack data from KIM object
202 !
203 call kim_api_getm_data(pkim, compute_energy_forces, &
204  "numberOfParticles", pn, 1, &
205  "particleSpecies", pparticlespecies,1, &
206  "coordinates", pcoor, 1, &
207  "energy", penergy, truefalse(comp_energy.eq.1), &
208  "forces", pforce, truefalse(comp_force.eq.1), &
209  "particleEnergy", penepot, truefalse(comp_enepot.eq.1), &
210  "virial", pvirial, truefalse(comp_virial.eq.1))
211 if (compute_energy_forces.lt.kim_status_ok) then
212  idum = kim_api_report_error(__line__, this_file_name, &
213  "kim_api_getm_data_f", compute_energy_forces)
214  return
215 endif
216 
217 call c_f_pointer(pn, n)
218 call c_f_pointer(pparticlespecies, particlespecies, [n])
219 call c_f_pointer(pcoor, coor, [dim,n])
220 if (comp_energy.eq.1) call c_f_pointer(penergy, energy)
221 if (comp_force.eq.1) call c_f_pointer(pforce, force, [dim,n])
222 if (comp_enepot.eq.1) call c_f_pointer(penepot, enepot, [n])
223 if (comp_virial.eq.1) call c_f_pointer(pvirial, virial, [6])
224 
225 
226 ! Check to be sure that the species are correct
227 !
228 compute_energy_forces = kim_status_fail ! assume an error
229 do i = 1,n
230  if (particlespecies(i).ne.speccode) then
231  idum = kim_api_report_error(__line__, this_file_name, &
232  "Unexpected species detected", &
234  return
235  endif
236 enddo
237 compute_energy_forces = kim_status_ok ! everything is ok
238 
239 ! Initialize potential energies, forces, virial term
240 !
241 if (comp_enepot.eq.1) enepot = 0.0_cd
242 if (comp_energy.eq.1) energy = 0.0_cd
243 if (comp_force.eq.1) force = 0.0_cd
244 if (comp_virial.eq.1) virial = 0.0_cd
245 
246 
247 !
248 ! Compute energy and forces
249 !
250 
251 ! Loop over particles and compute energy and forces
252 !
253 do i=1,n
254  compute_energy_forces = kim_api_get_neigh(pkim,1,i,part_ret,numnei, &
255  pnei1part,prij_list)
256  if (compute_energy_forces.ne.kim_status_ok) then
257  ! some sort of problem, exit
258  idum = kim_api_report_error(__line__, this_file_name, &
259  "kim_api_get_neigh", &
261  compute_energy_forces = kim_status_fail
262  return
263  endif
264 
265  call c_f_pointer(pnei1part, nei1part, [numnei])
266 
267  ! Loop over the neighbors of particle i
268  !
269  do jj = 1, numnei
270 
271  j = nei1part(jj) ! get neighbor ID
272 
273  ! compute relative position vector
274  !
275  rij(:) = coor(:,j) - coor(:,i) ! distance vector between i j
276 
277  ! compute energy and forces
278  !
279  rsqij = dot_product(rij,rij) ! compute square distance
280  if ( rsqij .lt. model_cutsq ) then ! particles are interacting?
281 
282  r = sqrt(rsqij) ! compute distance
283  if (comp_force.eq.1.or.comp_virial.eq.1) then
284  call calc_phi_dphi(r,phi,dphi) ! compute pair potential
285  ! and it derivative
286  deidr = 0.5_cd*dphi
287  else
288  call calc_phi(r,phi) ! compute just pair potential
289  endif
290 
291  ! contribution to energy
292  !
293  if (comp_enepot.eq.1) then
294  enepot(i) = enepot(i) + 0.5_cd*phi ! accumulate energy
295  endif
296  if (comp_energy.eq.1) then
297  energy = energy + 0.5_cd*phi
298  endif
299 
300  ! contribution to virial tensor, virial(i,j)=r(i)*r(j)*(dV/dr)/r
301  !
302  if (comp_virial.eq.1) then
303  virial(1) = virial(1) + rij(1)*rij(1)*deidr/r
304  virial(2) = virial(2) + rij(2)*rij(2)*deidr/r
305  virial(3) = virial(3) + rij(3)*rij(3)*deidr/r
306  virial(4) = virial(4) + rij(2)*rij(3)*deidr/r
307  virial(5) = virial(5) + rij(1)*rij(3)*deidr/r
308  virial(6) = virial(6) + rij(1)*rij(2)*deidr/r
309  endif
310 
311  ! contribution to forces
312  !
313  if (comp_force.eq.1) then
314  force(:,i) = force(:,i) + deidr*rij/r ! accumulate force on i
315  force(:,j) = force(:,j) - deidr*rij/r ! accumulate force on j
316  endif
317 
318  endif
319 
320  enddo ! loop on jj
321 
322 enddo
323 
324 ! Everything is great
325 !
326 compute_energy_forces = kim_status_ok
327 return
328 
329 end function compute_energy_forces
330 
331 end module ex_model_ar_p_mlj_f03
332 
333 !-------------------------------------------------------------------------------
334 !
335 ! Model initialization routine (REQUIRED)
336 !
337 !-------------------------------------------------------------------------------
338 integer(c_int) function model_init(pkim) bind(c)
339 use, intrinsic :: iso_c_binding
341 use kim_api_f03
342 implicit none
343 
344 !-- Transferred variables
345 type(c_ptr), intent(in) :: pkim
346 
347 !-- Local variables
348 integer(c_int), parameter :: one=1
349 integer(c_int) ier, idum
350 
351 !-- KIM variables
352 real(c_double), pointer :: cutoff; type(c_ptr) :: pcutoff
353 
354 ! store pointer to compute function in KIM object
355 ier = kim_api_set_method(pkim,"compute",one,c_funloc(compute_energy_forces))
356 if (ier.lt.kim_status_ok) then
357  idum = kim_api_report_error(__line__, this_file_name, &
358  "kim_api_set_method", ier)
359  goto 42
360 endif
361 
362 ! store model cutoff in KIM object
363 pcutoff = kim_api_get_data(pkim,"cutoff",ier)
364 if (ier.lt.kim_status_ok) then
365  idum = kim_api_report_error(__line__, this_file_name, &
366  "kim_api_get_data", ier)
367  goto 42
368 endif
369 call c_f_pointer(pcutoff, cutoff)
370 cutoff = model_cutoff
371 
372 ier = kim_status_ok
373 42 continue
374 model_init = ier
375 return
376 
377 end function model_init
real(c_double), parameter, public model_cutoff
integer(c_int) function, public compute_energy_forces(pkim)