ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/lj_FF.F90
Revision: 253
Committed: Thu Jan 30 15:20:21 2003 UTC (21 years, 5 months ago) by chuckv
File size: 20938 byte(s)
Log Message:
Added a generic util code directory and moved Linux_ifc_machdep to it.
MPI changes to compile MPI modules.

File Contents

# Content
1 !! Calculates Long Range forces Lennard-Jones interactions.
2 !! Corresponds to the force field defined in lj_FF.cpp
3 !! @author Charles F. Vardeman II
4 !! @author Matthew Meineke
5 !! @version $Id: lj_FF.F90,v 1.13 2003-01-30 15:20:21 chuckv Exp $, $Date: 2003-01-30 15:20:21 $, $Name: not supported by cvs2svn $, $Revision: 1.13 $
6
7
8
9 module lj_ff
10 use simulation
11 use definitions
12 use generic_atypes
13 #ifdef IS_MPI
14 use mpiSimulation
15 #endif
16 implicit none
17 PRIVATE
18
19 !! Number of lj_atypes in lj_atype_list
20 integer, save :: n_lj_atypes = 0
21
22 !! Global list of lj atypes in simulation
23 type (lj_atype), pointer :: ljListHead => null()
24 type (lj_atype), pointer :: ljListTail => null()
25
26
27 !! LJ mixing array
28 type (lj_atype), dimension(:,:), pointer :: ljMixed => null()
29
30
31 !! Neighbor list and commom arrays
32 !! This should eventally get moved to a neighbor list type
33 integer, allocatable, dimension(:) :: point
34 integer, allocatable, dimension(:) :: list
35 integer, parameter :: listMultiplier = 80
36 integer :: nListAllocs = 0
37 integer, parameter :: maxListAllocs = 5
38
39 logical, save :: firstTime = .True.
40
41 !! Atype identity pointer lists
42 #ifdef IS_MPI
43 !! Row lj_atype pointer list
44 type (lj_identPtrList), dimension(:), pointer :: identPtrListRow => null()
45 !! Column lj_atype pointer list
46 type (lj_identPtrList), dimension(:), pointer :: identPtrListColumn => null()
47 #else
48 type( lj_identPtrList ), dimension(:), pointer :: identPtrList => null()
49 #endif
50
51
52 !! Logical has lj force field module been initialized?
53 logical, save :: isljFFinit = .false.
54
55
56 !! Public methods and data
57 public :: new_lj_atype
58 public :: do_lj_ff
59 public :: getLjPot
60 public :: init_ljFF
61
62
63
64
65 contains
66
67 !! Adds a new lj_atype to the list.
68 subroutine new_lj_atype(ident,mass,epsilon,sigma,status)
69 real( kind = dp ), intent(in) :: mass
70 real( kind = dp ), intent(in) :: epsilon
71 real( kind = dp ), intent(in) :: sigma
72 integer, intent(in) :: ident
73 integer, intent(out) :: status
74
75 type (lj_atype), pointer :: newLJ_atype
76 integer :: alloc_error
77 integer :: atype_counter = 0
78 integer :: alloc_size
79 integer :: err_stat
80 status = 0
81
82
83
84 ! allocate a new atype
85 allocate(newLJ_atype,stat=alloc_error)
86 if (alloc_error /= 0 ) then
87 status = -1
88 return
89 end if
90
91 ! assign our new lj_atype information
92 newLJ_atype%mass = mass
93 newLJ_atype%epsilon = epsilon
94 newLJ_atype%sigma = sigma
95 newLJ_atype%sigma2 = sigma * sigma
96 newLJ_atype%sigma6 = newLJ_atype%sigma2 * newLJ_atype%sigma2 &
97 * newLJ_atype%sigma2
98 ! assume that this atype will be successfully added
99 newLJ_atype%atype_ident = ident
100 newLJ_atype%atype_number = n_lj_atypes + 1
101
102 call add_atype(newLJ_atype,ljListHead,ljListTail,err_stat)
103 if (err_stat /= 0 ) then
104 status = -1
105 return
106 endif
107
108 n_lj_atypes = n_lj_atypes + 1
109
110
111 end subroutine new_lj_atype
112
113
114 subroutine init_ljFF(nComponents,ident, status)
115 !! Number of components in ident array
116 integer, intent(inout) :: nComponents
117 !! Array of identities nComponents long corresponding to
118 !! ljatype ident.
119 integer, dimension(nComponents),intent(inout) :: ident
120 !! Result status, success = 0, error = -1
121 integer, intent(out) :: Status
122
123 integer :: alloc_stat
124
125 integer :: thisStat
126 integer :: i
127 #ifdef IS_MPI
128 integer, allocatable, dimension(:) :: identRow
129 integer, allocatable, dimension(:) :: identCol
130 integer :: nrow
131 integer :: ncol
132 #endif
133 status = 0
134 !! if were're not in MPI, we just update ljatypePtrList
135 #ifndef IS_MPI
136 call create_IdentPtrlst(ident,ljListHead,identPtrList,thisStat)
137 if ( thisStat /= 0 ) then
138 status = -1
139 return
140 endif
141
142 !! Allocate pointer lists
143 allocate(point(nComponents),stat=alloc_stat)
144 if (alloc_stat /=0) then
145 status = -1
146 return
147 endif
148
149 allocate(list(nComponents*listMultiplier),stat=alloc_stat)
150 if (alloc_stat /=0) then
151 status = -1
152 return
153 endif
154
155 ! if were're in MPI, we also have to worry about row and col lists
156 #else
157 ! We can only set up forces if mpiSimulation has been setup.
158 if (.not. isMPISimSet()) then
159 status = -1
160 return
161 endif
162 nrow = getNrow()
163 ncol = getNcol()
164 !! Allocate temperary arrays to hold gather information
165 allocate(identRow(nrow),stat=alloc_stat)
166 if (alloc_stat /= 0 ) then
167 status = -1
168 return
169 endif
170
171 allocate(identCol(ncol),stat=alloc_stat)
172 if (alloc_stat /= 0 ) then
173 status = -1
174 return
175 endif
176 !! Gather idents into row and column idents
177 call gather(ident,identRow,plan_row)
178 call gather(ident,identCol,plan_col)
179
180 !! Create row and col pointer lists
181 call create_IdentPtrlst(identRow,ljListTail,identPtrListRow,thisStat)
182 if (thisStat /= 0 ) then
183 status = -1
184 return
185 endif
186
187 call create_IdentPtrlst(identCol,ljListTail,identPtrListColumn,thisStat)
188 if (thisStat /= 0 ) then
189 status = -1
190 return
191 endif
192
193 !! free temporary ident arrays
194 deallocate(identCol)
195 deallocate(identRow)
196
197
198 !! Allocate neighbor lists for mpi simulations.
199 if (.not. allocated(point)) then
200 allocate(point(nrow),stat=alloc_stat)
201 if (alloc_stat /=0) then
202 status = -1
203 return
204 endif
205
206 allocate(list(ncol*listMultiplier),stat=alloc_stat)
207 if (alloc_stat /=0) then
208 status = -1
209 return
210 endif
211 else
212 deallocate(list)
213 deallocate(point)
214 allocate(point(nrow),stat=alloc_stat)
215 if (alloc_stat /=0) then
216 status = -1
217 return
218 endif
219
220 allocate(list(ncol*listMultiplier),stat=alloc_stat)
221 if (alloc_stat /=0) then
222 status = -1
223 return
224 endif
225 endif
226
227 #endif
228
229 call createMixingList(thisStat)
230 if (thisStat /= 0) then
231 status = -1
232 return
233 endif
234 isljFFinit = .true.
235
236 end subroutine init_ljFF
237
238
239
240
241
242
243 subroutine createMixingList(status)
244 integer :: listSize
245 integer :: status
246 integer :: i
247 integer :: j
248
249 integer :: outerCounter = 0
250 integer :: innerCounter = 0
251 type (lj_atype), pointer :: tmpPtrOuter => null()
252 type (lj_atype), pointer :: tmpPtrInner => null()
253 status = 0
254
255 listSize = getListLen(ljListHead)
256 if (listSize == 0) then
257 status = -1
258 return
259 end if
260
261
262 if (.not. associated(ljMixed)) then
263 allocate(ljMixed(listSize,listSize))
264 else
265 status = -1
266 return
267 end if
268
269
270
271 tmpPtrOuter => ljListHead
272 tmpPtrInner => tmpPtrOuter%next
273 do while (associated(tmpPtrOuter))
274 outerCounter = outerCounter + 1
275 ! do self mixing rule
276 ljMixed(outerCounter,outerCounter)%sigma = tmpPtrOuter%sigma
277
278 ljMixed(outerCounter,outerCounter)%sigma2 = ljMixed(outerCounter,outerCounter)%sigma &
279 * ljMixed(outerCounter,outerCounter)%sigma
280
281 ljMixed(outerCounter,outerCounter)%sigma6 = ljMixed(outerCounter,outerCounter)%sigma2 &
282 * ljMixed(outerCounter,outerCounter)%sigma2 &
283 * ljMixed(outerCounter,outerCounter)%sigma2
284
285 ljMixed(outerCounter,outerCounter)%epsilon = tmpPtrOuter%epsilon
286
287 innerCounter = outerCounter + 1
288 do while (associated(tmpPtrInner))
289
290 ljMixed(outerCounter,innerCounter)%sigma = &
291 calcLJMix("sigma",tmpPtrOuter%sigma, &
292 tmpPtrInner%sigma)
293
294 ljMixed(outerCounter,innerCounter)%sigma2 = &
295 ljMixed(outerCounter,innerCounter)%sigma &
296 * ljMixed(outerCounter,innerCounter)%sigma
297
298 ljMixed(outerCounter,innerCounter)%sigma6 = &
299 ljMixed(outerCounter,innerCounter)%sigma2 &
300 * ljMixed(outerCounter,innerCounter)%sigma2 &
301 * ljMixed(outerCounter,innerCounter)%sigma2
302
303 ljMixed(outerCounter,innerCounter)%epsilon = &
304 calcLJMix("epsilon",tmpPtrOuter%epsilon, &
305 tmpPtrInner%epsilon)
306 ljMixed(innerCounter,outerCounter)%sigma = ljMixed(outerCounter,innerCounter)%sigma
307 ljMixed(innerCounter,outerCounter)%sigma2 = ljMixed(outerCounter,innerCounter)%sigma2
308 ljMixed(innerCounter,outerCounter)%sigma6 = ljMixed(outerCounter,innerCounter)%sigma6
309 ljMixed(innerCounter,outerCounter)%epsilon = ljMixed(outerCounter,innerCounter)%epsilon
310
311
312 tmpPtrInner => tmpPtrInner%next
313 innerCounter = innerCounter + 1
314 end do
315 ! advance pointers
316 tmpPtrOuter => tmpPtrOuter%next
317 if (associated(tmpPtrOuter)) then
318 tmpPtrInner => tmpPtrOuter%next
319 endif
320
321 end do
322
323 end subroutine createMixingList
324
325
326
327
328
329
330
331
332 !! FORCE routine Calculates Lennard Jones forces.
333 !------------------------------------------------------------->
334 subroutine do_lj_ff(q,f,potE,do_pot)
335 !! Position array provided by C, dimensioned by getNlocal
336 real ( kind = dp ), dimension(3,getNlocal()) :: q
337 !! Force array provided by C, dimensioned by getNlocal
338 real ( kind = dp ), dimension(3,getNlocal()) :: f
339 real ( kind = dp ) :: potE
340 logical ( kind = 2) :: do_pot
341
342 type(lj_atype), pointer :: ljAtype_i
343 type(lj_atype), pointer :: ljAtype_j
344
345 #ifdef IS_MPI
346 real( kind = DP ), dimension(3,getNcol()) :: efr
347 real( kind = DP ) :: pot_local
348 #else
349 real( kind = DP ), dimension(3,getNlocal()) :: efr
350 #endif
351
352 !! Local arrays needed for MPI
353 #ifdef IS_MPI
354 real(kind = dp), dimension(3:getNrow()) :: qRow
355 real(kind = dp), dimension(3:getNcol()) :: qCol
356
357 real(kind = dp), dimension(3:getNrow()) :: fRow
358 real(kind = dp), dimension(3:getNcol()) :: fCol
359 real(kind = dp), dimension(3:getNlocal()) :: fMPITemp
360
361 real(kind = dp), dimension(getNrow()) :: eRow
362 real(kind = dp), dimension(getNcol()) :: eCol
363
364 real(kind = dp), dimension(getNlocal()) :: eTemp
365 #endif
366
367
368
369 real( kind = DP ) :: pe
370 logical :: update_nlist
371
372
373 integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
374 integer :: nlist
375 integer :: j_start
376 integer :: tag_i,tag_j
377 real( kind = DP ) :: r, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
378 real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq
379 real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut
380
381 ! a rig that need to be fixed.
382 #ifdef IS_MPI
383 logical :: newtons_thrd = .true.
384 real( kind = dp ) :: pe_local
385 integer :: nlocal
386 #endif
387 integer :: nrow
388 integer :: ncol
389 integer :: natoms
390
391
392
393
394 ! Make sure we are properly initialized.
395 if (.not. isljFFInit) then
396 write(default_error,*) "ERROR: lj_FF has not been properly initialized"
397 return
398 endif
399 #ifdef IS_MPI
400 if (.not. isMPISimSet()) then
401 write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
402 return
403 endif
404 #endif
405
406 !! initialize local variables
407 natoms = getNlocal()
408 call getRcut(rcut,rcut2=rcutsq)
409 call getRlist(rlist,rlistsq)
410
411 #ifndef IS_MPI
412 nrow = natoms - 1
413 ncol = natoms
414 #else
415 nrow = getNrow(plan_row)
416 ncol = getNcol(plan_col)
417 nlocal = natoms
418 j_start = 1
419 #endif
420
421
422 !! See if we need to update neighbor lists
423 call check(q,update_nlist)
424 if (firstTime) then
425 update_nlist = .true.
426 firstTime = .false.
427 endif
428
429 !--------------WARNING...........................
430 ! Zero variables, NOTE:::: Forces are zeroed in C
431 ! Zeroing them here could delete previously computed
432 ! Forces.
433 !------------------------------------------------
434 #ifndef IS_MPI
435 ! nloops = nloops + 1
436 pe = 0.0E0_DP
437
438 #else
439 fRow = 0.0E0_DP
440 fCol = 0.0E0_DP
441
442 pe_local = 0.0E0_DP
443
444 eRow = 0.0E0_DP
445 eCol = 0.0E0_DP
446 eTemp = 0.0E0_DP
447 #endif
448 efr = 0.0E0_DP
449
450 ! communicate MPI positions
451 #ifdef IS_MPI
452 call gather(q,qRow,plan_row3d)
453 call gather(q,qCol,plan_col3d)
454 #endif
455
456
457 if (update_nlist) then
458
459 ! save current configuration, contruct neighbor list,
460 ! and calculate forces
461 call save_nlist(q)
462
463 nlist = 0
464
465
466
467 do i = 1, nrow
468 point(i) = nlist + 1
469 #ifdef IS_MPI
470 ljAtype_i => identPtrListRow(i)%this
471 tag_i = tagRow(i)
472 rxi = qRow(1,i)
473 ryi = qRow(2,i)
474 rzi = qRow(3,i)
475 #else
476 ljAtype_i => identPtrList(i)%this
477 j_start = i + 1
478 rxi = q(1,i)
479 ryi = q(2,i)
480 rzi = q(3,i)
481 #endif
482
483 inner: do j = j_start, ncol
484 #ifdef IS_MPI
485 ! Assign identity pointers and tags
486 ljAtype_j => identPtrListColumn(j)%this
487 tag_j = tagCol(j)
488 if (newtons_thrd) then
489 if (tag_i <= tag_j) then
490 if (mod(tag_i + tag_j,2) == 0) cycle inner
491 else
492 if (mod(tag_i + tag_j,2) == 1) cycle inner
493 endif
494 endif
495
496 rxij = wrap(rxi - qCol(1,j), 1)
497 ryij = wrap(ryi - qCol(2,j), 2)
498 rzij = wrap(rzi - qCol(3,j), 3)
499 #else
500 ljAtype_j => identPtrList(j)%this
501 rxij = wrap(rxi - q(1,j), 1)
502 ryij = wrap(ryi - q(2,j), 2)
503 rzij = wrap(rzi - q(3,j), 3)
504
505 #endif
506 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
507
508 #ifdef IS_MPI
509 if (rijsq <= rlistsq .AND. &
510 tag_j /= tag_i) then
511 #else
512
513 if (rijsq < rlistsq) then
514 #endif
515
516 nlist = nlist + 1
517 if (nlist > size(list)) then
518 !! "Change how nlist size is done"
519 write(DEFAULT_ERROR,*) "ERROR: nlist > list size"
520 endif
521 list(nlist) = j
522
523
524 if (rijsq < rcutsq) then
525
526 r = dsqrt(rijsq)
527
528 call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j)
529
530 #ifdef IS_MPI
531 eRow(i) = eRow(i) + pot*0.5
532 eCol(i) = eCol(i) + pot*0.5
533 #else
534 pe = pe + pot
535 #endif
536
537 efr(1,j) = -rxij
538 efr(2,j) = -ryij
539 efr(3,j) = -rzij
540
541 do dim = 1, 3
542
543
544 drdx1 = efr(dim,j) / r
545 ftmp = dudr * drdx1
546
547
548 #ifdef IS_MPI
549 fCol(dim,j) = fCol(dim,j) - ftmp
550 fRow(dim,i) = fRow(dim,i) + ftmp
551 #else
552
553 f(dim,j) = f(dim,j) - ftmp
554 f(dim,i) = f(dim,i) + ftmp
555
556 #endif
557 enddo
558 endif
559 endif
560 enddo inner
561 enddo
562
563 #ifdef IS_MPI
564 point(nrow + 1) = nlist + 1
565 #else
566 point(natoms) = nlist + 1
567 #endif
568
569 else
570
571 ! use the list to find the neighbors
572 do i = 1, nrow
573 JBEG = POINT(i)
574 JEND = POINT(i+1) - 1
575 ! check thiat molecule i has neighbors
576 if (jbeg .le. jend) then
577 #ifdef IS_MPI
578 ljAtype_i => identPtrListRow(i)%this
579 rxi = qRow(1,i)
580 ryi = qRow(2,i)
581 rzi = qRow(3,i)
582 #else
583 ljAtype_i => identPtrList(i)%this
584 rxi = q(1,i)
585 ryi = q(2,i)
586 rzi = q(3,i)
587 #endif
588 do jnab = jbeg, jend
589 j = list(jnab)
590 #ifdef IS_MPI
591 ljAtype_j = identPtrListColumn(j)%this
592 rxij = wrap(rxi - q_col(1,j), 1)
593 ryij = wrap(ryi - q_col(2,j), 2)
594 rzij = wrap(rzi - q_col(3,j), 3)
595 #else
596 ljAtype_j = identPtrList(j)%this
597 rxij = wrap(rxi - q(1,j), 1)
598 ryij = wrap(ryi - q(2,j), 2)
599 rzij = wrap(rzi - q(3,j), 3)
600 #endif
601 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
602
603 if (rijsq < rcutsq) then
604
605 r = dsqrt(rijsq)
606
607 call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j)
608 #ifdef IS_MPI
609 eRow(i) = eRow(i) + pot*0.5
610 eCol(i) = eCol(i) + pot*0.5
611 #else
612 pe = pe + pot
613 #endif
614
615
616 efr(1,j) = -rxij
617 efr(2,j) = -ryij
618 efr(3,j) = -rzij
619
620 do dim = 1, 3
621
622 drdx1 = efr(dim,j) / r
623 ftmp = dudr * drdx1
624 #ifdef IS_MPI
625 fCol(dim,j) = fCol(dim,j) - ftmp
626 fRow(dim,i) = fRow(dim,i) + ftmp
627 #else
628 f(dim,j) = f(dim,j) - ftmp
629 f(dim,i) = f(dim,i) + ftmp
630 #endif
631 enddo
632 endif
633 enddo
634 endif
635 enddo
636 endif
637
638
639
640 #ifdef IS_MPI
641 !!distribute forces
642 call scatter(fRow,f,plan_row3d)
643
644 call scatter(fCol,fMPITemp,plan_col3d)
645
646 do i = 1,nlocal
647 f(1:3,i) = f(1:3,i) + fMPITemp(1:3,i)
648 end do
649
650
651
652 if (do_pot) then
653 ! scatter/gather pot_row into the members of my column
654 call scatter(eRow,eTemp,plan_row)
655
656 ! scatter/gather pot_local into all other procs
657 ! add resultant to get total pot
658 do i = 1, nlocal
659 pe_local = pe_local + eTemp(i)
660 enddo
661 if (newtons_thrd) then
662 eTemp = 0.0E0_DP
663 call scatter(eCol,eTemp,plan_col)
664 do i = 1, nlocal
665 pe_local = pe_local + eTemp(i)
666 enddo
667 endif
668 endif
669 potE = pe_local
670 #endif
671
672 potE = pe
673
674
675
676
677 end subroutine do_lj_ff
678
679 !! Calculates the potential between two lj particles based on two lj_atype pointers, optionally returns second
680 !! derivatives.
681 subroutine getLjPot(r,pot,dudr,atype1,atype2,d2,status)
682 ! arguments
683 !! Length of vector between particles
684 real( kind = dp ), intent(in) :: r
685 !! Potential Energy
686 real( kind = dp ), intent(out) :: pot
687 !! Derivatve wrt postion
688 real( kind = dp ), intent(out) :: dudr
689 !! Second Derivative, optional, used mainly for normal mode calculations.
690 real( kind = dp ), intent(out), optional :: d2
691
692 type (lj_atype), pointer :: atype1
693 type (lj_atype), pointer :: atype2
694
695 integer, intent(out), optional :: status
696
697 ! local Variables
698 real( kind = dp ) :: sigma
699 real( kind = dp ) :: sigma2
700 real( kind = dp ) :: sigma6
701 real( kind = dp ) :: epsilon
702
703 real( kind = dp ) :: rcut
704 real( kind = dp ) :: rcut2
705 real( kind = dp ) :: rcut6
706 real( kind = dp ) :: r2
707 real( kind = dp ) :: r6
708
709 real( kind = dp ) :: t6
710 real( kind = dp ) :: t12
711 real( kind = dp ) :: tp6
712 real( kind = dp ) :: tp12
713 real( kind = dp ) :: delta
714
715 logical :: doSec = .false.
716
717 integer :: errorStat
718
719 !! Optional Argument Checking
720 ! Check to see if we need to do second derivatives
721
722 if (present(d2)) doSec = .true.
723 if (present(status)) status = 0
724
725 ! Look up the correct parameters in the mixing matrix
726 sigma = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma
727 sigma2 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma2
728 sigma6 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma6
729 epsilon = ljMixed(atype1%atype_ident,atype2%atype_ident)%epsilon
730
731
732
733
734 call getRcut(rcut,rcut2=rcut2,rcut6=rcut6,status=errorStat)
735
736 r2 = r * r
737 r6 = r2 * r2 * r2
738
739 t6 = sigma6/ r6
740 t12 = t6 * t6
741
742
743
744 tp6 = sigma6 / rcut6
745 tp12 = tp6*tp6
746
747 delta = -4.0E0_DP*epsilon * (tp12 - tp6)
748
749 if (r.le.rcut) then
750 pot = 4.0E0_DP * epsilon * (t12 - t6) + delta
751 dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / r
752 if(doSec) d2 = 24.0E0_DP * epsilon * (26.0E0_DP*t12 - 7.0E0_DP*t6)/r/r
753 else
754 pot = 0.0E0_DP
755 dudr = 0.0E0_DP
756 if(doSec) d2 = 0.0E0_DP
757 endif
758
759 return
760
761
762
763 end subroutine getLjPot
764
765
766 !! Calculates the mixing for sigma or epslon based on character string for initialzition of mixing array.
767 function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
768 character(len=*) :: thisParam
769 real(kind = dp) :: param1
770 real(kind = dp) :: param2
771 real(kind = dp ) :: myMixParam
772 integer, optional :: status
773
774
775 myMixParam = 0.0_dp
776
777 if (present(status)) status = 0
778
779 select case (thisParam)
780
781 case ("sigma")
782 myMixParam = 0.5_dp * (param1 + param2)
783 case ("epsilon")
784 myMixParam = sqrt(param1 * param2)
785 case default
786 status = -1
787 end select
788
789 end function calcLJMix
790
791
792
793 end module lj_ff