kim-api  2.1.0+v2.1.0.GNU
An Application Programming Interface (API) for the Knowledgebase of Interatomic Models (KIM).
utility_forces_numer_deriv.f90
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--2019, Regents of the University of Minnesota.
23 ! All rights reserved.
24 !
25 ! Contributors:
26 ! Ellad B. Tadmor
27 ! Ryan S. Elliott
28 ! Stephen M. Whalen
29 !
30 
31 
32 module error
33  use, intrinsic :: iso_c_binding
34  implicit none
35  public
36 
37 contains
38  recursive subroutine my_error(message)
39  implicit none
40  character(len=*, kind=c_char), intent(in) :: message
41 
42  print *,"* Error : ", trim(message)
43  stop
44  end subroutine my_error
45 
46  recursive subroutine my_warning(message)
47  implicit none
48  character(len=*, kind=c_char), intent(in) :: message
49 
50  print *,"* Warning : ", trim(message)
51  end subroutine my_warning
52 end module error
53 
54 
55 !-------------------------------------------------------------------------------
56 !
57 ! module mod_neighborlist :
58 !
59 ! Module contains type and routines related to neighbor list calculation
60 !
61 !-------------------------------------------------------------------------------
62 
63 module mod_neighborlist
64 
65  use, intrinsic :: iso_c_binding
66 
67  public get_neigh
68 
69  type, bind(c) :: neighobject_type
70  real(c_double) :: cutoff
71  integer(c_int) :: number_of_particles
72  type(c_ptr) :: neighbor_list_pointer
73  end type neighobject_type
74 contains
75 
76 !-------------------------------------------------------------------------------
77 !
78 ! get_neigh neighbor list access function
79 !
80 ! This function implements Locator and Iterator mode
81 !
82 !-------------------------------------------------------------------------------
83 !-------------------------------------------------------------------------------
84 !
85 ! get_neigh neighbor list access function
86 !
87 ! This function implements Locator and Iterator mode
88 !
89 !-------------------------------------------------------------------------------
90 recursive subroutine get_neigh(data_object, number_of_neighbor_lists, cutoffs, &
91  neighbor_list_index, request, numnei, pnei1part, ierr) bind(c)
92  use error
93  implicit none
94 
95  !-- Transferred variables
96  type(c_ptr), value, intent(in) :: data_object
97  integer(c_int), value, intent(in) :: number_of_neighbor_lists
98  real(c_double), intent(in) :: cutoffs(*)
99  integer(c_int), value, intent(in) :: neighbor_list_index
100  integer(c_int), value, intent(in) :: request
101  integer(c_int), intent(out) :: numnei
102  type(c_ptr), intent(out) :: pnei1part
103  integer(c_int), intent(out) :: ierr
104 
105  !-- Local variables
106  integer(c_int) numberofparticles
107  type(neighobject_type), pointer :: neighobject
108  integer(c_int), pointer :: neighborlist(:,:)
109 
110  if (number_of_neighbor_lists > 1) then
111  call my_warning("Model requires too many neighbor lists")
112  ierr = 1
113  return
114  endif
115 
116  call c_f_pointer(data_object, neighobject)
117  numberofparticles = neighobject%number_of_particles
118  call c_f_pointer(neighobject%neighbor_list_pointer, neighborlist, &
119  [numberofparticles+1, numberofparticles])
120 
121  if (cutoffs(neighbor_list_index) > neighobject%cutoff) then
122  call my_warning("neighbor list cutoff too small for model cutoff")
123  ierr = 1
124  return
125  endif
126 
127  if ( (request.gt.numberofparticles) .or. (request.lt.1)) then
128  print *, request
129  call my_warning("Invalid part ID in get_neigh")
130  ierr = 1
131  return
132  endif
133 
134  ! set the returned number of neighbors for the returned part
135  numnei = neighborlist(1,request)
136 
137  ! set the location for the returned neighbor list
138  pnei1part = c_loc(neighborlist(2,request))
139 
140  ierr = 0
141  return
142 end subroutine get_neigh
143 
144 end module mod_neighborlist
145 
148  implicit none
149  public
150 
151 contains
152 
153 !-------------------------------------------------------------------------------
154 !
155 ! Check if we are compatible with the model
156 !
157 !-------------------------------------------------------------------------------
158 recursive subroutine check_model_compatibility(compute_arguments_handle, &
159  forces_optional, model_is_compatible, ierr)
160  use, intrinsic :: iso_c_binding
161  use error
162  implicit none
163 
164  !-- Transferred variables
165  type(kim_compute_arguments_handle_type), intent(in) :: compute_arguments_handle
166  logical, intent(out) :: forces_optional
167  logical, intent(out) :: model_is_compatible
168  integer(c_int), intent(out) :: ierr
169 
170  !-- Local variables
171  integer(c_int) i
172  integer(c_int) number_of_argument_names
173  integer(c_int) number_of_callback_names
174  type(kim_compute_argument_name_type) argument_name
175  type(kim_support_status_type) support_status
176  type(kim_compute_callback_name_type) callback_name
177 
178 
179  ! assume fail
180  model_is_compatible = .false.
181  ierr = 0
182 
183  ! check arguments
184  call kim_get_number_of_compute_argument_names(&
185  number_of_argument_names)
186  do i=1,number_of_argument_names
187  call kim_get_compute_argument_name(i, argument_name, &
188  ierr)
189  if (ierr /= 0) then
190  call my_warning("can't get argument name")
191  return
192  end if
193  call kim_get_argument_support_status( &
194  compute_arguments_handle, argument_name, support_status, ierr)
195  if (ierr /= 0) then
196  call my_warning("can't get argument support_status")
197  return
198  end if
199 
200  ! can only handle energy and forces as required args
201  if (support_status == kim_support_status_required) then
202  if (.not. ( &
203  (argument_name == kim_compute_argument_name_partial_energy) .or. &
204  (argument_name == kim_compute_argument_name_partial_forces))) then
205  call my_warning("unsupported required argument")
206  ierr = 0
207  return
208  end if
209  end if
210 
211  ! need both energy and forces not "notSupported"
212  if ((argument_name == kim_compute_argument_name_partial_energy) .and. &
213  (support_status == kim_support_status_not_supported)) then
214  call my_warning("model does not support energy")
215  ierr = 0
216  return
217  end if
218  if (argument_name == kim_compute_argument_name_partial_forces) then
219  if (support_status == kim_support_status_not_supported) then
220  call my_warning("model does not support forces")
221  ierr = 0
222  return
223  else if (support_status == kim_support_status_required) then
224  forces_optional = .false.
225  else if (support_status == kim_support_status_optional) then
226  forces_optional = .true.
227  else
228  call my_warning("unknown support_status for forces")
229  ierr = 0
230  return
231  end if
232  end if
233  end do
234 
235  ! check call backs
236  call kim_get_number_of_compute_callback_names( &
237  number_of_callback_names)
238  do i=1,number_of_callback_names
239  call kim_get_compute_callback_name(i, callback_name, &
240  ierr)
241  if (ierr /= 0) then
242  call my_warning("can't get call back name")
243  return
244  end if
245  call kim_get_callback_support_status( &
246  compute_arguments_handle, callback_name, support_status, ierr)
247  if (ierr /= 0) then
248  call my_warning("can't get call back support_status")
249  return
250  end if
251 
252  ! cannot handle any "required" call backs
253  if (support_status == kim_support_status_required) then
254  call my_warning("unsupported required call back")
255  ierr = 0
256  return
257  end if
258  end do
259 
260  ! got to here, then everything must be OK
261  model_is_compatible = .true.
262  ierr = 0
263  return
264 end subroutine check_model_compatibility
265 
266 !-------------------------------------------------------------------------------
267 !
268 ! Get number and identities of particle species supported by
269 ! KIM Model `modelname'
270 !
271 !-------------------------------------------------------------------------------
272 recursive subroutine get_model_supported_species(model_handle, max_species, &
273  model_species, num_species, ier)
274 use, intrinsic :: iso_c_binding
275 implicit none
276 
277 !-- Transferred variables
278 type(kim_model_handle_type), intent(in) :: model_handle
279 integer(c_int), intent(in) :: max_species
280 type(kim_species_name_type), intent(out) :: model_species(max_species)
281 integer(c_int), intent(out) :: num_species
282 integer(c_int), intent(out) :: ier
283 
284 !-- Local variables
285 integer(c_int) i
286 integer(c_int) total_num_species
287 type(kim_species_name_type) :: species_name
288 integer(c_int) species_is_supported
289 integer(c_int) code
290 
291 ! Initialize error flag
292 ier = 1
293 
294 call kim_get_number_of_species_names(total_num_species)
295 
296 if (total_num_species .gt. max_species) return
297 
298 num_species = 0
299 do i=1,total_num_species
300  call kim_get_species_name(i, species_name, ier)
301  call kim_get_species_support_and_code(model_handle, species_name, &
302  species_is_supported, code, ier)
303  if ((ier == 0) .and. (species_is_supported .ne. 0)) then
304  num_species = num_species + 1
305  model_species(num_species) = species_name
306  end if
307 end do
308 
309 ier = 0
310 return
311 
312 end subroutine get_model_supported_species
313 
314 recursive subroutine update_neighborlist(DIM,N,coords,cutoff,cutpad, &
315  do_update_list,coordsave, &
316  neighObject,ierr)
317 use, intrinsic :: iso_c_binding
319 implicit none
320 integer(c_int), parameter :: cd = c_double ! used for literal constants
321 
322 !-- Transferred variables
323 integer(c_int), intent(in) :: DIM
324 integer(c_int), intent(in) :: N
325 real(c_double), intent(in) :: coords(dim,n)
326 real(c_double), intent(in) :: cutoff
327 real(c_double), intent(in) :: cutpad
328 logical, intent(inout) :: do_update_list
329 real(c_double), intent(inout) :: coordsave(dim,n)
330 type(neighobject_type), intent(inout) :: neighObject
331 integer(c_int), intent(out) :: ierr
332 
333 !-- Local variables
334 real(c_double) disp, disp1, disp2, cutrange, dispvec(dim)
335 integer(c_int) i
336 
337 ! Initialize error code
338 !
339 ierr = 0
340 
341 ! Update neighbor lists if necessary
342 !
343 if (.not.do_update_list) then ! if update not requested
344 
345  ! check whether a neighbor list update is necessary even if it hasn't been
346  ! requested using the "two max sum" criterion
347  disp1 = 0.0_cd
348  disp2 = 0.0_cd
349  do i=1,n
350  dispvec(1:dim) = coords(1:dim,i) - coordsave(1:dim,i)
351  disp = sqrt( dot_product(dispvec,dispvec) )
352  if (disp >= disp1) then ! 1st position taken
353  disp2 = disp1 ! push current 1st into 2nd place
354  disp1 = disp ! and put this one into current 1st
355  else if (disp >= disp2) then ! 2nd position taken
356  disp2 = disp
357  endif
358  enddo
359  do_update_list = ( disp1 + disp2 > cutpad )
360 
361 endif
362 
363 if (do_update_list) then
364 
365  ! save current coordinates
366  coordsave(1:dim,1:n) = coords(1:dim,1:n)
367 
368  ! compute neighbor lists
369  cutrange = cutoff + cutpad
370  call neigh_pure_cluster_neighborlist(.false., n, coords, cutrange, &
371  neighobject)
372 
373  ! neighbor list uptodate, no need to compute again for now
374  do_update_list = .false.
375 endif
376 
377 return
378 
379 end subroutine update_neighborlist
380 
381 
382 
383 !-------------------------------------------------------------------------------
384 !
385 ! NEIGH_PURE_cluster_neighborlist
386 !
387 !-------------------------------------------------------------------------------
388 recursive subroutine neigh_pure_cluster_neighborlist(half, numberOfParticles, &
389  coords, cutoff, neighObject)
390  use, intrinsic :: iso_c_binding
391  use mod_neighborlist
392  implicit none
393 
394  !-- Transferred variables
395  logical, intent(in) :: half
396  integer(c_int), intent(in) :: numberOfParticles
397  real(c_double), dimension(3,numberOfParticles), &
398  intent(in) :: coords
399  real(c_double), intent(in) :: cutoff
400  type(neighobject_type), intent(inout) :: neighObject
401 
402  !-- Local variables
403  integer(c_int), pointer :: neighborList(:,:)
404  integer(c_int) i, j, a
405  real(c_double) dx(3)
406  real(c_double) r2
407  real(c_double) cutoff2
408 
409  call c_f_pointer(neighobject%neighbor_list_pointer, neighborlist, &
410  [numberofparticles+1, numberofparticles])
411 
412  neighobject%cutoff = cutoff
413 
414  cutoff2 = cutoff**2
415 
416  do i=1,numberofparticles
417  a = 1
418  do j=1,numberofparticles
419  dx(:) = coords(:, j) - coords(:, i)
420  r2 = dot_product(dx, dx)
421  if (r2.le.cutoff2) then
422  ! part j is a neighbor of part i
423  if ( (j .gt. i) .OR. ((.not. half) .AND. (i.ne.j)) ) then
424  a = a+1
425  neighborlist(a,i) = j
426  endif
427  endif
428  enddo
429  ! part i has a-1 neighbors
430  neighborlist(1,i) = a-1
431  enddo
432 
433  return
434 
435 end subroutine neigh_pure_cluster_neighborlist
436 
437 !-------------------------------------------------------------------------------
438 !
439 ! create_FCC_configuration subroutine
440 !
441 ! creates a cubic configuration of FCC particles with lattice spacing
442 ! `FCCspacing' and `nCellsPerSide' cells along each direction.
443 !
444 ! With periodic==.true. this will result in a total number of particles equal
445 ! to 4*(nCellsPerSide)**3 + 6*(nCellsPerSide)**2 + 3*(nCellsPerSide) + 1
446 !
447 ! With periodic==.false. this will result in a total number of particles equal
448 ! to 4*(nCellsPerSide)**3
449 !
450 ! Returns the Id of the particle situated in the middle of the configuration
451 ! (this particle will have the most neighbors.)
452 !
453 !-------------------------------------------------------------------------------
454 recursive subroutine create_fcc_configuration(FCCspacing, nCellsPerSide, &
455  periodic, coords, MiddlePartId)
456  use, intrinsic :: iso_c_binding
457  implicit none
458  integer(c_int), parameter :: cd = c_double ! used for literal constants
459 
460  !-- Transferred variables
461  real(c_double), intent(in) :: FCCspacing
462  integer(c_int), intent(in) :: nCellsPerSide
463  logical, intent(in) :: periodic
464  real(c_double), intent(out) :: coords(3,*)
465  integer(c_int), intent(out) :: MiddlePartId
466  !
467  ! cluster setup variables
468  !
469  real(c_double) FCCshifts(3,4)
470  real(c_double) latVec(3)
471  integer(c_int) a, i, j, k, m
472 
473  ! Create a cubic FCC cluster
474  !
475  fccshifts(1,1) = 0.0_cd
476  fccshifts(2,1) = 0.0_cd
477  fccshifts(3,1) = 0.0_cd
478  fccshifts(1,2) = 0.5_cd*fccspacing
479  fccshifts(2,2) = 0.5_cd*fccspacing
480  fccshifts(3,2) = 0.0_cd
481  fccshifts(1,3) = 0.5_cd*fccspacing
482  fccshifts(2,3) = 0.0_cd
483  fccshifts(3,3) = 0.5_cd*fccspacing
484  fccshifts(1,4) = 0.0_cd
485  fccshifts(2,4) = 0.5_cd*fccspacing
486  fccshifts(3,4) = 0.5_cd*fccspacing
487 
488  middlepartid = 1 ! Always put middle particle as #1
489  a = 1 ! leave space for middle particle as particle #1
490  do i=1,ncellsperside
491  latvec(1) = (i-1)*fccspacing
492  do j=1,ncellsperside
493  latvec(2) = (j-1)*fccspacing
494  do k=1,ncellsperside
495  latvec(3) = (k-1)*fccspacing
496  do m=1,4
497  a = a+1
498  coords(:,a) = latvec + fccshifts(:,m)
499  if ((i.eq.ncellsperside/2+1).and.(j.eq.ncellsperside/2+1) .and. &
500  (k.eq.ncellsperside/2+1) .and. (m.eq.1)) then
501  coords(:,1) = latvec + fccshifts(:,m) ! put middle particle as #1
502  a = a - 1
503  endif
504  enddo
505  enddo
506  if (.not. periodic) then
507  ! Add in the remaining three faces
508  ! pos-x face
509  latvec(1) = ncellsperside*fccspacing
510  latvec(2) = (i-1)*fccspacing
511  latvec(3) = (j-1)*fccspacing
512  a = a+1; coords(:,a) = latvec
513  a = a+1; coords(:,a) = latvec + fccshifts(:,4)
514  ! pos-y face
515  latvec(1) = (i-1)*fccspacing
516  latvec(2) = ncellsperside*fccspacing
517  latvec(3) = (j-1)*fccspacing
518  a = a+1; coords(:,a) = latvec
519  a = a+1; coords(:,a) = latvec + fccshifts(:,3)
520  ! pos-z face
521  latvec(1) = (i-1)*fccspacing
522  latvec(2) = (j-1)*fccspacing
523  latvec(3) = ncellsperside*fccspacing
524  a = a+1; coords(:,a) = latvec
525  a = a+1; coords(:,a) = latvec + fccshifts(:,2)
526  endif
527  enddo
528  if (.not. periodic) then
529  ! Add in the remaining three edges
530  latvec(1) = (i-1)*fccspacing
531  latvec(2) = ncellsperside*fccspacing
532  latvec(3) = ncellsperside*fccspacing
533  a = a+1; coords(:,a) = latvec
534  latvec(1) = ncellsperside*fccspacing
535  latvec(2) = (i-1)*fccspacing
536  latvec(3) = ncellsperside*fccspacing
537  a = a+1; coords(:,a) = latvec
538  latvec(1) = ncellsperside*fccspacing
539  latvec(2) = ncellsperside*fccspacing
540  latvec(3) = (i-1)*fccspacing
541  a = a+1; coords(:,a) = latvec
542  endif
543  enddo
544  if (.not. periodic) then
545  ! Add in the remaining corner
546  a = a+1; coords(:,a) = ncellsperside*fccspacing
547  endif
548 
549  return
550 
551 end subroutine create_fcc_configuration
552 
553 recursive subroutine compute_numer_deriv(partnum,dir,model_handle, &
554  compute_arguments_handle,DIM,N,coords,cutoff,cutpad,energy,do_update_list, &
555  coordsave,neighObject,deriv,deriv_err,ierr)
556 use, intrinsic :: iso_c_binding
557 use error
559 implicit none
560 integer(c_int), parameter :: cd = c_double ! used for literal constants
561 
562 !--Transferred variables
563 integer(c_int), intent(in) :: partnum
564 integer(c_int), intent(in) :: dir
565 type(kim_model_handle_type), intent(in) :: model_handle
566 type(kim_compute_arguments_handle_type), intent(in) :: compute_arguments_handle
567 integer(c_int), intent(in) :: DIM
568 integer(c_int), intent(in) :: N
569 real(c_double), intent(inout) :: coords(dim,n)
570 real(c_double), intent(in) :: cutoff
571 real(c_double), intent(in) :: cutpad
572 real(c_double), intent(inout) :: energy
573 logical, intent(inout) :: do_update_list
574 real(c_double), intent(inout) :: coordsave(dim,n)
575 type(neighobject_type), intent(inout) :: neighObject
576 real(c_double), intent(out) :: deriv
577 real(c_double), intent(out) :: deriv_err
578 integer(c_int), intent(out) :: ierr
579 
580 !-- Local variables
581 real(c_double), parameter :: eps_init = 1.e-6_cd
582 integer(c_int), parameter :: number_eps_levels = 15
583 real(c_double) eps, deriv_last, deriv_err_last
584 integer(c_int) i
585 
586 ! Initialize error flag
587 ierr = 0
588 
589 deriv_last = 0.0_cd ! initialize
590 
591 ! Outer loop of Ridders' method for computing numerical derivative
592 !
593 eps = eps_init
594 deriv_err_last = huge(1.0_cd)
595 do i=1,number_eps_levels
596  deriv = dfridr(eps,deriv_err)
597  if (ierr /= 0) then
598  call my_error("compute_numer_deriv")
599  endif
600  if (deriv_err>deriv_err_last) then
601  deriv = deriv_last
602  deriv_err = deriv_err_last
603  exit
604  endif
605  eps = eps*10.0_cd
606  deriv_last = deriv
607  deriv_err_last = deriv_err
608 enddo
609 
610 return
611 
612 contains
613 
614  !----------------------------------------------------------------------------
615  !
616  ! Compute numerical derivative using Ridders' method
617  !
618  ! Based on code from Numerical Recipes, Press et al., Second Ed., Cambridge,
619  ! 1992
620  !
621  ! Ref: Ridders, C. J. F., "Two algorithms for the calculation of F'(x)=D",
622  ! Advances in Engineering Software, Vol. 4, no. 2, pp. 75-76, 1982.
623  !
624  !
625  ! Returns the gradient grad() of a KIM-compliant interatomic model at the
626  ! current configuration by Ridders' method of polynomial extrapolation.
627  ! An estimate for the error in each component of the gradient is returned in
628  ! grad_err().
629  !
630  !----------------------------------------------------------------------------
631  real(c_double) recursive function dfridr(h,err)
632  implicit none
633 
634  !-- Transferred variables
635  real(c_double), intent(inout) :: h
636  real(c_double), intent(out) :: err
637 
638  !-- Local variables
639  integer(c_int), parameter :: NTAB=10 ! Maximum size of tableau
640  real(c_double), parameter :: CON=1.4_cd ! Stepsize incr. by CON at each iter
641  real(c_double), parameter :: CON2=con*con
642  real(c_double), parameter :: BIG=huge(1.0_cd)
643  real(c_double), parameter :: SAFE=2.0_cd ! Returns when error is SAFE worse
644  ! than the best so far
645  integer(c_int) i,j
646  real(c_double) errt,fac,hh,a(ntab,ntab),fp,fm,coordorig
647 
648  dfridr = 0.0_cd ! initialize
649 
650  if (abs(h).le.tiny(0.0_cd)) then ! avoid division by zero
651  ierr = 1
652  return
653  endif
654 
655  hh = h
656  coordorig = coords(dir,partnum)
657  coords(dir,partnum) = coordorig + hh
658  call update_neighborlist(dim,n,coords,cutoff,cutpad, &
659  do_update_list,coordsave, &
660  neighobject,ierr)
661  call kim_compute(model_handle, compute_arguments_handle, ierr)
662  if (ierr /= 0) then
663  call my_error("kim_api_model_compute")
664  endif
665  fp = energy
666  coords(dir,partnum) = coordorig - hh
667  call update_neighborlist(dim,n,coords,cutoff,cutpad,&
668  do_update_list,coordsave, &
669  neighobject,ierr)
670  call kim_compute(model_handle, compute_arguments_handle, ierr)
671  if (ierr /= 0) then
672  call my_error("kim_api_model_compute")
673  endif
674  fm = energy
675  coords(dir,partnum) = coordorig
676  call update_neighborlist(dim,n,coords,cutoff,cutpad,&
677  do_update_list,coordsave, &
678  neighobject,ierr)
679  a(1,1)=(fp-fm)/(2.0_cd*hh)
680  err=big
681  ! successive columns in the Neville tableau will go to smaller step sizes
682  ! and higher orders of extrapolation
683  do i=2,ntab
684  ! try new, smaller step size
685  hh=hh/con
686  coords(dir,partnum) = coordorig + hh
687  call update_neighborlist(dim,n,coords,cutoff,cutpad, &
688  do_update_list,coordsave, &
689  neighobject,ierr)
690  call kim_compute(model_handle, compute_arguments_handle, ierr)
691  if (ierr /= 0) then
692  call my_error("kim_api_model_compute")
693  endif
694  fp = energy
695  coords(dir,partnum) = coordorig - hh
696  call update_neighborlist(dim,n,coords,cutoff,cutpad, &
697  do_update_list,coordsave, &
698  neighobject,ierr)
699  call kim_compute(model_handle, compute_arguments_handle, ierr)
700  if (ierr /= 0) then
701  call my_error("kim_api_model_compute")
702  endif
703  fm = energy
704  coords(dir,partnum) = coordorig
705  call update_neighborlist(dim,n,coords,cutoff,cutpad, &
706  do_update_list,coordsave, &
707  neighobject,ierr)
708  a(1,i)=(fp-fm)/(2.0_cd*hh)
709  fac=con2
710  ! compute extrapolations of various orders, requiring no new function
711  ! evaluations
712  do j=2,i
713  a(j,i)=(a(j-1,i)*fac-a(j-1,i-1))/(fac-1.0_cd)
714  fac=con2*fac
715  ! The error strategy is to compute each new extrapolation to one order
716  ! lower, both at the present step size and the previous one.
717  errt = max(abs(a(j,i)-a(j-1,i)),abs(a(j,i)-a(j-1,i-1)))
718  if (errt.le.err) then ! if error is decreased, save the improved answer
719  err = errt
720  dfridr=a(j,i)
721  endif
722  enddo
723  if (abs(a(i,i)-a(i-1,i-1)).ge.safe*err) return ! if higher order is worse
724  ! by significant factor
725  ! `SAFE', then quit early.
726  enddo
727  return
728  end function dfridr
729 
730  end subroutine compute_numer_deriv
731 
732 end module mod_utilities
733 
734 !*******************************************************************************
735 !**
736 !** PROGRAM vc_forces_numer_deriv
737 !**
738 !** KIM compliant program to perform numerical derivative check on a model
739 !**
740 !*******************************************************************************
741 
742 !-------------------------------------------------------------------------------
743 !
744 ! Main program
745 !
746 !-------------------------------------------------------------------------------
748  use, intrinsic :: iso_c_binding
749  use error
750  use mod_neighborlist
751  use mod_utilities
752  implicit none
753  integer(c_int), parameter :: cd = c_double ! used for literal constants
754 
755  integer(c_int), parameter :: nCellsPerSide = 2
756  integer(c_int), parameter :: DIM = 3
757  real(c_double), parameter :: cutpad = 0.75_cd
758  integer(c_int), parameter :: max_species = 200 ! most species supported
759  real(c_double), parameter :: eps_prec = epsilon(1.0_cd)
760  real(c_double) FCCspacing
761 
762  integer(c_int), parameter :: &
763  N = 4*(ncellsperside)**3 + 6*(ncellsperside)**2 + 3*(ncellsperside) + 1
764  real(c_double), allocatable :: forces_num(:,:)
765  real(c_double), allocatable :: forces_num_err(:,:)
766  type(kim_species_name_type) :: model_species(max_species)
767  integer(c_int), target :: num_species
768  character(len=5, kind=c_char) :: passfail
769  real(c_double) :: forcediff
770  real(c_double) :: forcediff_sumsq
771  real(c_double) :: weight
772  real(c_double) :: weight_sum
773  real(c_double) :: alpha
774  real(c_double) :: term
775  real(c_double) :: term_max
776  real(c_double), allocatable :: cluster_coords(:,:)
777  real(c_double), allocatable :: cluster_disps(:,:)
778  type(kim_species_name_type), allocatable :: cluster_species(:)
779  integer(c_int) I,J,Imax,Jmax,species_code
780  integer(c_int) seed_size
781  integer(c_int), allocatable :: seed(:)
782 
783  !
784  ! neighbor list
785  !
786  type(neighobject_type), target :: neighObject
787  integer(c_int), allocatable, target :: neighborList(:,:)
788  real(c_double), allocatable :: coordsave(:,:)
789  logical do_update_list
790 
791  !
792  ! KIM variables
793  !
794  character(len=256, kind=c_char) :: testname = "vc_forces_numer_deriv"
795  character(len=256, kind=c_char) :: modelname
796 
797  type(kim_model_handle_type) :: model_handle
798  type(kim_compute_arguments_handle_type) :: compute_arguments_handle
799  integer(c_int) ierr, ierr2
800  integer(c_int) species_is_supported
801  integer(c_int), target :: numberOfParticles
802  integer(c_int), target :: particleSpeciesCodes(n)
803  integer(c_int), target :: particleContributing(n)
804  real(c_double) :: influence_distance
805  integer(c_int) :: number_of_neighbor_lists
806  real(c_double), allocatable :: cutoffs(:)
807  integer(c_int), allocatable :: &
808  model_will_not_request_neighbors_of_noncontributing_particles(:)
809  real(c_double) :: cutoff
810  real(c_double), target :: energy
811  real(c_double), target :: coords(3,n)
812  real(c_double), target :: forces_kim(3,n)
813  real(c_double) :: forces(3,n)
814  integer(c_int) middleDum
815  logical forces_optional
816  logical model_is_compatible
817  integer(c_int) number_of_model_routine_names
818  type(kim_model_routine_name_type) model_routine_name
819  integer(c_int) present
820  integer(c_int) required
821  integer(c_int) requested_units_accepted
822  real(c_double) rnd, deriv, deriv_err
823  real(c_double), pointer :: null_pointer
824 
825  nullify(null_pointer)
826 
827  numberofparticles = n
828 
829  term_max = 0.0_cd ! initialize
830 
831  ! Initialize error flag
832  ierr = 0
833 
834  ! Initialize seed for random number generator
835  !
836  ! NOTE: Here, we set the seed to a fixed value for reproducibility
837  call random_seed(size=seed_size) ! Get seed size for this CPU
838  allocate(seed(seed_size))
839  seed(:) = 13
840  call random_seed(put=seed)
841  deallocate(seed)
842 
843  ! Get KIM Model name to use
844  print '("Please enter a valid KIM model name: ")'
845  read(*,*) modelname
846 
847  ! Print output header
848  !
849  print *
850  print *,'VERIFICATION CHECK: NUMERICAL DERIVATIVE VERIFICATION OF FORCES'
851  print *
852  print '(120(''-''))'
853  print '("This is Test : ",A)', trim(testname)
854  print '("Results for KIM Model : ",A)', trim(modelname)
855 
856  ! Create empty KIM object
857  !
858  call kim_model_create(kim_numbering_one_based, &
859  kim_length_unit_a, &
860  kim_energy_unit_ev, &
861  kim_charge_unit_e, &
862  kim_temperature_unit_k, &
863  kim_time_unit_ps, &
864  trim(modelname), &
865  requested_units_accepted, &
866  model_handle, ierr)
867  if (ierr /= 0) then
868  call my_error("kim_api_create")
869  endif
870 
871  ! check that we are compatible
872  if (requested_units_accepted == 0) then
873  call my_error("Must adapt to model units")
874  end if
875 
876  ! check that we know about all required routines
877  call kim_get_number_of_model_routine_names(number_of_model_routine_names)
878  do i=1,number_of_model_routine_names
879  call kim_get_model_routine_name(i, model_routine_name, ierr)
880  if (ierr /= 0) then
881  call my_error("kim_get_model_routine_name")
882  endif
883  call kim_is_routine_present(model_handle, model_routine_name, present, &
884  required, ierr)
885  if (ierr /= 0) then
886  call my_error("kim_is_routine_present")
887  endif
888 
889  if ((present == 1) .and. (required == 1)) then
890  if (.not. ((model_routine_name == kim_model_routine_name_create) &
891  .or. (model_routine_name == &
892  kim_model_routine_name_compute_arguments_create) &
893  .or. (model_routine_name == kim_model_routine_name_compute) &
894  .or. (model_routine_name == kim_model_routine_name_refresh) &
895  .or. (model_routine_name == &
896  kim_model_routine_name_compute_arguments_destroy) &
897  .or. (model_routine_name == kim_model_routine_name_destroy))) then
898 
899  call my_error("Unknown required ModelRoutineName found.")
900  endif
901  endif
902  enddo
903 
904  ! create compute_arguments object
905  call kim_compute_arguments_create(model_handle, &
906  compute_arguments_handle, ierr)
907  if (ierr /= 0) then
908  call my_error("kim_model_compute_arguments_create")
909  endif
910 
911  call check_model_compatibility(compute_arguments_handle, forces_optional, &
912  model_is_compatible, ierr)
913  if (ierr /= 0) then
914  call my_error("error checking compatibility")
915  end if
916  if (.not. model_is_compatible) then
917  call my_error("incompatibility reported by check_model_compatibility")
918  end if
919 
920  ! Get list of particle species supported by the model
921  !
922  call get_model_supported_species(model_handle, max_species, model_species, &
923  num_species, ierr)
924  if (ierr /= 0) then
925  call my_error("Get_Model_Supported_Species")
926  endif
927  ! Setup random cluster
928  !
929  allocate(cluster_coords(3,n),cluster_disps(3,n),cluster_species(n))
930  do i=1,n
931  call random_number(rnd) ! return random number between 0 and 1
932  species_code = 1 + int(rnd*num_species)
933  cluster_species(i) = model_species(species_code)
934  enddo
935  fccspacing = 1.0_cd ! initially generate an FCC cluster with lattice
936  ! spacing equal to one. This is scaled below based
937  ! on the cutoff radius of the model.
938  call create_fcc_configuration(fccspacing, ncellsperside, .false., &
939  cluster_coords, middledum)
940  ! Generate random displacements for all particles
941  !
942  do i=1,n
943  do j=1,dim
944  call random_number(rnd) ! return random number between 0 and 1
945  cluster_disps(j,i) = 0.1_cd*(rnd-0.5_cd)
946  enddo
947  enddo
948 
949  ! register memory with the KIM system
950  ierr = 0
951  call kim_set_argument_pointer(compute_arguments_handle, &
952  kim_compute_argument_name_number_of_particles, numberofparticles, ierr2)
953  ierr = ierr + ierr2
954  call kim_set_argument_pointer(compute_arguments_handle, &
955  kim_compute_argument_name_particle_species_codes, particlespeciescodes, &
956  ierr2)
957  ierr = ierr + ierr2
958  call kim_set_argument_pointer(compute_arguments_handle, &
959  kim_compute_argument_name_particle_contributing, particlecontributing, &
960  ierr2)
961  ierr = ierr + ierr2
962  call kim_set_argument_pointer(compute_arguments_handle, &
963  kim_compute_argument_name_coordinates, coords, ierr2)
964  ierr = ierr + ierr2
965  call kim_set_argument_pointer(compute_arguments_handle, &
966  kim_compute_argument_name_partial_energy, energy, ierr2)
967  ierr = ierr + ierr2
968  call kim_set_argument_pointer(compute_arguments_handle, &
969  kim_compute_argument_name_partial_forces, forces_kim, ierr2)
970  ierr = ierr + ierr2
971  if (ierr /= 0) then
972  call my_error("set_argument_pointer")
973  endif
974 
975  ! Allocate storage for neighbor lists and
976  ! store pointers to neighbor list object and access function
977  !
978  allocate(neighborlist(n+1,n))
979  neighobject%neighbor_list_pointer = c_loc(neighborlist)
980  neighobject%number_of_particles = n
981 
982  ! Set pointer in KIM object to neighbor list routine and object
983  !
984  call kim_set_callback_pointer(compute_arguments_handle, &
985  kim_compute_callback_name_get_neighbor_list, kim_language_name_fortran, &
986  c_funloc(get_neigh), c_loc(neighobject), ierr)
987  if (ierr /= 0) then
988  call my_error("set_callback_pointer")
989  end if
990 
991  call kim_get_influence_distance(model_handle, influence_distance)
992  call kim_get_number_of_neighbor_lists(model_handle, &
993  number_of_neighbor_lists)
994  allocate(cutoffs(number_of_neighbor_lists), &
995  model_will_not_request_neighbors_of_noncontributing_particles( &
996  number_of_neighbor_lists))
997  call kim_get_neighbor_list_values(model_handle, cutoffs, &
998  model_will_not_request_neighbors_of_noncontributing_particles, ierr)
999  if (ierr /= 0) then
1000  call my_error("get_neighbor_list_values")
1001  end if
1002  cutoff = maxval(cutoffs)
1003 
1004  ! Scale reference FCC configuration based on cutoff radius.
1005  fccspacing = 0.75_cd*cutoff ! set the FCC spacing to a fraction
1006  ! of the cutoff radius
1007  do i=1,n
1008  cluster_coords(:,i) = fccspacing*cluster_coords(:,i)
1009  enddo
1010  print '("Using FCC lattice parameter: ",f12.5)', fccspacing
1011 
1012  do i=1,n
1013  call kim_get_species_support_and_code(model_handle, &
1014  cluster_species(i), species_is_supported, particlespeciescodes(i), ierr)
1015  enddo
1016  if (ierr /= 0) then
1017  call my_error("kim_api_get_species_code")
1018  endif
1019  do i=1,n
1020  particlecontributing(i) = 1 ! every particle contributes
1021  enddo
1022  do i=1,n
1023  coords(:,i) = cluster_coords(:,i) + cluster_disps(:,i)
1024  enddo
1025 
1026  ! Compute neighbor lists
1027  !
1028  do_update_list = .true.
1029  allocate(coordsave(dim,n))
1030  call update_neighborlist(dim,n,coords,cutoff,cutpad, &
1031  do_update_list,coordsave, &
1032  neighobject,ierr)
1033  if (ierr /= 0) then
1034  call my_error("update_neighborlist")
1035  endif
1036 
1037  ! Call model compute to get forces (gradient)
1038  !
1039  call kim_compute(model_handle, compute_arguments_handle, ierr)
1040  if (ierr /= 0) then
1041  call my_error("kim_api_model_compute")
1042  endif
1043 
1044  ! Copy forces in case model will overwrite forces_kim below
1045  !
1046  forces = forces_kim
1047 
1048  ! Print results to screen
1049  !
1050  print '(41(''=''))'
1051  print '("Energy = ",ES25.15)', energy
1052  print '(41(''=''))'
1053  print *
1054 
1055  ! Turn off force computation, if possible
1056  !
1057  if (forces_optional) then
1058  call kim_set_argument_pointer(compute_arguments_handle, &
1059  kim_compute_argument_name_partial_forces, null_pointer, ierr)
1060  if (ierr /= 0) then
1061  call my_error("set_argument_pointer")
1062  endif
1063  endif
1064 
1065  ! Compute gradient using numerical differentiation
1066  !
1067  allocate(forces_num(dim,n),forces_num_err(dim,n))
1068  do i=1,n
1069  do j=1,dim
1070  call compute_numer_deriv(i,j,model_handle,compute_arguments_handle,&
1071  dim,n,coords,cutoff, &
1072  cutpad, energy, do_update_list, &
1073  coordsave,neighobject,deriv,deriv_err,ierr)
1074  if (ierr /= 0) then
1075  call my_error("compute_numer_deriv")
1076  endif
1077  forces_num(j,i) = -deriv
1078  forces_num_err(j,i) = deriv_err
1079  enddo
1080  enddo
1081 
1082  ! Continue printing results to screen
1083  !
1084  print '(A6,2X,A4,2X,A3,2X,2A25,3A15,2X,A4)',"Part","Spec","Dir", &
1085  "Force_model", "Force_numer", "Force diff", "pred error", "weight", &
1086  "stat"
1087  forcediff_sumsq = 0.0_cd
1088  weight_sum = 0.0_cd
1089  do i=1,n
1090  do j=1,dim
1091  forcediff = abs(forces(j,i)-forces_num(j,i))
1092  if (forcediff<forces_num_err(j,i)) then
1093  passfail = "ideal"
1094  else
1095  passfail = " "
1096  endif
1097  weight = max(abs(forces_num(j,i)),eps_prec)/ &
1098  max(abs(forces_num_err(j,i)),eps_prec)
1099  term = weight*forcediff**2
1100  if (term.gt.term_max) then
1101  term_max = term
1102  imax = i
1103  jmax = j
1104  endif
1105  forcediff_sumsq = forcediff_sumsq + term
1106  weight_sum = weight_sum + weight
1107  if (j.eq.1) then
1108  print '(I6,2X,I4,2X,I3,2X,2ES25.15,3ES15.5,2X,A5)', &
1109  i,particlespeciescodes(i),j,forces(j,i),forces_num(j,i), &
1110  forcediff,forces_num_err(j,i),weight,passfail
1111  else
1112  print '(14X,I3,2X,2ES25.15,3ES15.5,2X,A5)', &
1113  j,forces(j,i),forces_num(j,i), &
1114  forcediff,forces_num_err(j,i),weight,passfail
1115  endif
1116  enddo
1117  print *
1118  enddo
1119  alpha = sqrt(forcediff_sumsq/weight_sum)/dble(dim*n)
1120  print *
1121  print '("alpha = |Force_model - Force_numer|_w/(DIM*N) = ",ES15.5," (units of force)")', &
1122  alpha
1123  print *
1124  print '(''Maximum term obtained for Part = '',I6,'', Dir = '',I1,' // &
1125  ''', forcediff = '',ES15.5, '', forcediff/force_model = '',ES15.5)', &
1126  imax,jmax,abs(forces(jmax,imax)-forces_num(jmax,imax)), &
1127  abs(forces(jmax,imax)-forces_num(jmax,imax))/abs(forces(jmax,imax))
1128 
1129  ! Free temporary storage
1130  !
1131  deallocate(forces_num)
1132  deallocate(forces_num_err)
1133  deallocate(neighborlist)
1134  deallocate(coordsave)
1135  deallocate(cutoffs, model_will_not_request_neighbors_of_noncontributing_particles)
1136 
1137  call kim_compute_arguments_destroy(model_handle, &
1138  compute_arguments_handle, ierr)
1139  if (ierr /= 0) then
1140  call my_error("kim_model_compute_arguments_destroy")
1141  endif
1142  call kim_model_destroy(model_handle)
1143 
1144  ! Print output footer
1145  !
1146  print *
1147  print '(120(''-''))'
1148 
1149  ! Free cluster storage
1150  !
1151  deallocate(cluster_coords,cluster_disps,cluster_species)
1152 
1153  stop
1154 
1155 end program vc_forces_numer_deriv
recursive subroutine check_model_compatibility(compute_arguments_handle, forces_optional, model_is_compatible, ierr)
recursive subroutine, public get_neigh(data_object, number_of_neighbor_lists, cutoffs, neighbor_list_index, request, numnei, pnei1part, ierr)
real(c_double) recursive function dfridr(h, err)
recursive subroutine my_warning(message)
recursive subroutine update_neighborlist(DIM, N, coords, cutoff, cutpad, do_update_list, coordsave, neighObject, ierr)
recursive subroutine get_model_supported_species(model_handle, max_species, model_species, num_species, ier)
recursive subroutine my_error(message)
program vc_forces_numer_deriv
recursive subroutine compute_numer_deriv(partnum, dir, model_handle, compute_arguments_handle, DIM, N, coords, cutoff, cutpad, energy, do_update_list, coordsave, neighObject, deriv, deriv_err, ierr)
recursive subroutine neigh_pure_cluster_neighborlist(half, numberOfParticles, coords, cutoff, neighObject)
recursive subroutine create_fcc_configuration(FCCspacing, nCellsPerSide, periodic, coords, MiddlePartId)