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