ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/nano_mpi/src/force_lj_module.F90
Revision: 4
Committed: Mon Jun 10 17:18:36 2002 UTC (22 years, 1 month ago) by chuckv
File size: 14089 byte(s)
Log Message:
Import Root

File Contents

# User Rev Content
1 chuckv 4 module lj_module
2     use definitions, ONLY : DP,ndim
3     use parameter
4     use simulation
5     use second_deriv
6     use status, ONLY: error,info,warning
7     use force_utilities
8     #ifdef MPI
9     use mpi_module
10     #endif
11    
12     integer, parameter :: n_ljatypes = 12
13     real( kind = DP ),allocatable, dimension(:) :: lj_eps
14     real( kind = DP ),allocatable, dimension(:) :: lj_sigma
15     integer, allocatable, dimension(:) :: ljatype ! to be fixed in a module
16    
17    
18    
19    
20     public :: lj_eps,lj_sigma
21     public :: calc_lj_forces,initialize_lj
22     private :: mass_weight
23    
24     contains
25    
26     subroutine allocate_lj_module(n_size_atype)
27     integer, intent(in) :: n_size_atype
28    
29     allocate(ljatype(n_size_atype))
30     allocate(lj_eps(n_size_atype))
31     allocate(lj_sigma(n_size_atype))
32     end subroutine allocate_lj_module
33    
34     subroutine deallocate_lj_module()
35     deallocate(lj_eps)
36     deallocate(lj_sigma)
37     deallocate(ljatype)
38     end subroutine deallocate_lj_module
39    
40     subroutine calc_lj_forces(update_nlist, nmflag,pe)
41    
42     ! include 'headers/sizes.h'
43     ! include 'headers/fileio.h'
44    
45    
46     #ifdef MPI
47     real( kind = DP ), dimension(3,ncol) :: efr
48     real( kind = DP ) :: pot_local
49     #else
50     real( kind = DP ), dimension(3,natoms) :: efr
51     #endif
52    
53     real( kind = DP ), intent(out), optional :: pe
54     logical, intent(in) :: nmflag
55     logical, intent(in) :: update_nlist
56     logical :: do_pot
57    
58     integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
59     integer :: nlist
60     integer :: j_start
61     integer :: tag_i,tag_j
62     real( kind = DP ) :: r, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
63     real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq
64    
65    
66    
67    
68     #ifndef MPI
69     integer :: nrow
70     integer :: ncol
71    
72     nrow = natoms - 1
73     ncol = natoms
74     #else
75     j_start = 1
76     #endif
77    
78    
79     do_pot = .false.
80     if (present(pe)) do_pot = .true.
81    
82     #ifndef MPI
83     if (do_pot) pot = 0.0E0_DP
84     f = 0.0E0_DP
85     e = 0.0E0_DP
86     #else
87     f_row = 0.0E0_DP
88     f_col = 0.0E0_DP
89    
90     pot_local = 0.0E0_DP
91    
92     e_row = 0.0E0_DP
93     e_col = 0.0E0_DP
94     e_tmp = 0.0E0_DP
95     #endif
96     efr = 0.0E0_DP
97    
98     ! communicate MPI positions
99     #ifdef MPI
100     call gather(q,q_row,plan_row3)
101     call gather(q,q_col,plan_col3)
102     #endif
103    
104     if (update_nlist) then
105    
106     ! save current configuration, contruct neighbor list,
107     ! and calculate forces
108     call save_nlist()
109    
110     nlist = 0
111    
112     do i = 1, nrow
113     point(i) = nlist + 1
114     #ifdef MPI
115     tag_i = tag_row(i)
116     rxi = q_row(1,i)
117     ryi = q_row(2,i)
118     rzi = q_row(3,i)
119     #else
120     j_start = i + 1
121     rxi = q(1,i)
122     ryi = q(2,i)
123     rzi = q(3,i)
124     #endif
125    
126     inner: do j = j_start, ncol
127     #ifdef MPI
128     tag_j = tag_col(j)
129     if (newtons_thrd) then
130     if (tag_i <= tag_j) then
131     if (mod(tag_i + tag_j,2) == 0) cycle inner
132     else
133     if (mod(tag_i + tag_j,2) == 1) cycle inner
134     endif
135     endif
136    
137     rxij = wrap(rxi - q_col(1,j), 1)
138     ryij = wrap(ryi - q_col(2,j), 2)
139     rzij = wrap(rzi - q_col(3,j), 3)
140     #else
141     rxij = wrap(rxi - q(1,j), 1)
142     ryij = wrap(ryi - q(2,j), 2)
143     rzij = wrap(rzi - q(3,j), 3)
144     #endif
145     rijsq = rxij*rxij + ryij*ryij + rzij*rzij
146    
147     #ifdef MPI
148     if (rijsq <= rlstsq .AND. &
149     tag_j /= tag_i) then
150     #else
151     if (rijsq < rlstsq) then
152     #endif
153    
154     nlist = nlist + 1
155     list(nlist) = j
156    
157    
158     if (rijsq < rcutsq) then
159    
160     r = dsqrt(rijsq)
161    
162     call LJ_mix(r,pot,dudr,d2,i,j)
163     #ifdef MPI
164     e_row(i) = e_row(i) + pot*0.5
165     e_col(i) = e_col(i) + pot*0.5
166     #else
167     if (do_pot) pe = pe + pot
168     #endif
169     efr(1,j) = -rxij
170     efr(2,j) = -ryij
171     efr(3,j) = -rzij
172    
173     do dim = 1, 3
174    
175     drdx1 = efr(j,dim) / r
176     ftmp = dudr * drdx1
177    
178     #ifdef MPI
179     f_col(dim,j) = f_col(dim,j) - ftmp
180     f_row(dim,i) = f_row(dim,i) + ftmp
181     #else
182     f(dim,j) = f(dim,j) - ftmp
183     f(dim,i) = f(dim,i) + ftmp
184     #endif
185     if (nmflag) then
186     idim = 3 * (i-1) + dim
187     jdim = 3 * (j-1) + dim
188    
189     do dim2 = 1, 3
190    
191     kt1 = d2 * efr(dim,j) * efr(dim2,j)/r/r
192     kt2 = - dudr * efr(dim,j) * efr(dim2,j)/r/r/r
193    
194     if (dim.eq.dim2) then
195     kt3 = dudr / r
196     else
197     kt3 = 0.0E0_DP
198     endif
199    
200     ! The factor of 2 below is to compensate for
201     ! overcounting.
202     ! Mass weighting is done separately...
203    
204     ktmp = (kt1+kt2+kt3)/2.0E0_DP
205     idim2 = 3 * (i-1) + dim2
206     jdim2 = 3 * (j-1) + dim2
207    
208     d(idim, idim2) = d(idim,idim2) + ktmp
209     d(idim2, idim) = d(idim2,idim) + ktmp
210    
211     d(idim, jdim2) = d(idim,jdim2) - ktmp
212     d(idim2, jdim) = d(idim2,jdim) - ktmp
213    
214     d(jdim, idim2) = d(jdim,idim2) - ktmp
215     d(jdim2, idim) = d(jdim2,idim) - ktmp
216    
217     d(jdim, jdim2) = d(jdim,jdim2) + ktmp
218     d(jdim2, jdim) = d(jdim2,jdim) + ktmp
219    
220     enddo
221     endif
222     enddo
223     endif
224     endif
225     enddo inner
226     enddo
227     #ifdef MPI
228     point(nrow + 1) = nlist + 1
229     #else
230     point(natoms) = nlist + 1
231     #endif
232    
233     else
234     ! use the list to find the neighbors
235     do i = 1, nrow
236     JBEG = POINT(i)
237     JEND = POINT(i+1) - 1
238     ! check thiat molecule i has neighbors
239     if (jbeg .le. jend) then
240     #ifdef MPI
241     rxi = q_row(1,i)
242     ryi = q_row(2,i)
243     rzi = q_row(3,i)
244     #else
245     rxi = q(1,i)
246     ryi = q(2,i)
247     rzi = q(3,i)
248     #endif
249     do jnab = jbeg, jend
250     j = list(jnab)
251     #ifdef MPI
252     rxij = wrap(rxi - q_col(1,j), 1)
253     ryij = wrap(ryi - q_col(2,j), 2)
254     rzij = wrap(rzi - q_col(3,j), 3)
255     #else
256     rxij = wrap(rxi - q(1,j), 1)
257     ryij = wrap(ryi - q(2,j), 2)
258     rzij = wrap(rzi - q(3,j), 3)
259     #endif
260     rijsq = rxij*rxij + ryij*ryij + rzij*rzij
261    
262     if (rijsq < rcutsq) then
263    
264     r = dsqrt(rijsq)
265    
266     call LJ_mix(r,pot,dudr,d2,i,j)
267     #ifdef MPI
268     e_row(i) = e_row(i) + pot*0.5
269     e_col(i) = e_col(i) + pot*0.5
270     #else
271     if (do_pot) pe = pe + pot
272     #endif
273    
274    
275     efr(1,j) = -rxij
276     efr(2,j) = -ryij
277     efr(3,j) = -rzij
278    
279     do dim = 1, 3
280    
281     drdx1 = efr(j,dim) / r
282     ftmp = dudr * drdx1
283     #ifdef MPI
284     f_col(dim,j) = f_col(dim,j) - ftmp
285     f_row(dim,i) = f_row(dim,i) + ftmp
286     #else
287     f(dim,j) = f(dim,j) - ftmp
288     f(dim,i) = f(dim,i) + ftmp
289     #endif
290     if (nmflag) then
291     idim = 3 * (i-1) + dim
292     jdim = 3 * (j-1) + dim
293    
294     do dim2 = 1, 3
295    
296     kt1 = d2 * efr(dim,j) * efr(dim2,j)/r/r
297     kt2 = - dudr * efr(dim,j) * efr(dim2,j)/r/r/r
298    
299     if (dim.eq.dim2) then
300     kt3 = dudr / r
301     else
302     kt3 = 0.0E0_DP
303     endif
304    
305     ! The factor of 2 below is to compensate for
306     ! overcounting.
307     ! Mass weighting is done separately...
308    
309     ktmp = (kt1+kt2+kt3)/2.0E0_DP
310     idim2 = 3 * (i-1) + dim2
311     jdim2 = 3 * (j-1) + dim2
312    
313     d(idim, idim2) = d(idim,idim2) + ktmp
314     d(idim2, idim) = d(idim2,idim) + ktmp
315    
316     d(idim, jdim2) = d(idim,jdim2) - ktmp
317     d(idim2, jdim) = d(idim2,jdim) - ktmp
318    
319     d(jdim, idim2) = d(jdim,idim2) - ktmp
320     d(jdim2, idim) = d(jdim2,idim) - ktmp
321    
322     d(jdim, jdim2) = d(jdim,jdim2) + ktmp
323     d(jdim2, jdim) = d(jdim2,jdim) + ktmp
324    
325     enddo
326     endif
327     enddo
328    
329     endif
330     enddo
331     endif
332     enddo
333     endif
334    
335     #ifdef MPI
336     !!distribute forces
337     call scatter(f_row,f,plan_row3)
338     if (newtons_thrd) then
339     call scatter(f_col,f_tmp,plan_col3)
340     do i = 1,nlocal
341     do dim = 1,3
342     f(dim,i) = f(dim,i) + f_tmp(dim,i)
343     end do
344     end do
345     endif
346    
347    
348     if (do_pot) then
349     ! scatter/gather pot_row into the members of my column
350     call scatter(e_row,e_tmp,plan_row)
351    
352     ! scatter/gather pot_local into all other procs
353     ! add resultant to get total pot
354     do i = 1, nlocal
355     pot_local = pot_local + e_tmp(i)
356     enddo
357     if (newtons_thrd) then
358     e_tmp = 0.0E0_DP
359     call scatter(e_col,e_tmp,plan_col)
360     do i = 1, nlocal
361     pot_local = pot_local + e_tmp(i)
362     enddo
363     endif
364    
365    
366     call mpi_reduce(pot_local,pe,1,mpi_double_precision, &
367     mpi_sum,0,mpi_comm_world,mpi_err)
368     endif
369     #endif
370    
371    
372    
373     if (nmflag) then
374     call mass_weight()
375     endif
376    
377     return
378     end subroutine calc_lj_forces
379    
380     subroutine LJ_mix(r,pot,dudr,d2,atom1,atom2)
381    
382     ! include 'headers/sizes.h'
383     ! include 'headers/atom.h'
384    
385    
386     integer :: atom1, atom2, id1, id2
387     real( kind = DP ) :: r, pot, dudr, d2
388     real( kind = DP ) :: u1, dudr1, d21
389     real( kind = DP ) :: this_sigma, this_eps
390    
391     #ifdef MPI
392     id1 = ident_row(atom1)
393     id2 = ident_col(atom2)
394     #else
395     id1 = ident(atom1)
396     id2 = ident(atom2)
397     #endif
398     if (id1.eq.id2) then
399     this_sigma = lj_sigma(id1)
400     this_eps = lj_eps(id1)
401     else
402     ! use Lorentz-Berthelot mixing rules:
403     this_sigma = 0.5E0_DP * (lj_sigma(id1) + lj_sigma(id2))
404     this_eps = dsqrt(lj_eps(id1)*lj_eps(id2))
405     endif
406    
407     call LJ_pot(r, this_sigma, this_eps, u1, dudr1, d21)
408    
409     pot = u1
410     dudr = dudr1
411     d2 = d21
412    
413     return
414     end subroutine LJ_mix
415    
416     subroutine LJ_pot(r, sigma, epsilon, u, dudr, d2)
417    
418    
419     real( kind = DP ) :: r,sigma, epsilon, u, dudr, d2
420     real( kind = DP ) :: sigma2, sigma6, r2, r6, rcut2, rcut6
421     real( kind = DP ) :: t6, t12, tp6, tp12, delta
422    
423     sigma2 = sigma*sigma
424     sigma6 = sigma2*sigma2*sigma2
425    
426     r2 = r*r
427     r6 = r2*r2*r2
428    
429    
430     rcut2 = rcut*rcut
431     rcut6 = rcut2*rcut2*rcut2
432    
433     t6 = sigma6 / r6
434     t12 = t6*t6
435    
436     tp6 = sigma6 / rcut6
437     tp12 = tp6*tp6
438    
439     delta = -4.0E0_DP*epsilon * (tp12 - tp6)
440    
441     if (r.le.rcut) then
442     u = 4.0E0_DP * epsilon * (t12 - t6) + delta
443     dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / r
444     d2 = 24.0E0_DP * epsilon * (26.0E0_DP*t12 - 7.0E0_DP*t6)/r/r
445     else
446     u = 0.0E0_DP
447     dudr = 0.0E0_DP
448     d2 = 0.0E0_DP
449     endif
450    
451     return
452     end subroutine LJ_pot
453    
454     subroutine initialize_lj()
455     use model_module
456     include 'headers/atom.h'
457    
458     integer :: n_atypes
459    
460     n_atypes = get_max_atype()
461     call allocate_lj_module(n_atypes)
462    
463    
464     ljatype(1) = H_atom
465     ljatype(2) = He_atom
466     ljatype(3) = C_atom
467     ljatype(4) = N_atom
468     ljatype(5) = O_atom
469     ljatype(6) = F_atom
470     ljatype(7) = Ne_atom
471     ljatype(8) = S_atom
472     ljatype(9) = Cl_atom
473     ljatype(10) = Ar_atom
474     ljatype(11) = Br_atom
475     ljatype(12) = Kr_atom
476    
477     lj_sigma(H_atom) = 2.81E0_DP
478     lj_sigma(He_atom) = 2.28E0_DP
479     lj_sigma(C_atom) = 3.35E0_DP
480     lj_sigma(N_atom) = 3.31E0_DP
481     lj_sigma(O_atom) = 2.95E0_DP
482     lj_sigma(F_atom) = 2.83E0_DP
483     lj_sigma(Ne_atom) = 2.72E0_DP
484     lj_sigma(S_atom) = 3.52E0_DP
485     lj_sigma(Cl_atom) = 3.35E0_DP
486     lj_sigma(Ar_atom) = 3.41E0_DP
487     lj_sigma(Br_atom) = 3.54E0_DP
488     lj_sigma(Kr_atom) = 3.83E0_DP
489    
490     lj_eps(H_atom) = 0.01708992E0_DP
491     lj_eps(He_atom) = 0.02026944E0_DP
492     lj_eps(C_atom) = 0.10174464E0_DP
493     lj_eps(N_atom) = 0.07412256E0_DP
494     lj_eps(O_atom) = 0.12241152E0_DP
495     lj_eps(F_atom) = 0.10492416E0_DP
496     lj_eps(Ne_atom) = 0.0933984E0_DP
497     lj_eps(S_atom) = 0.3636576E0_DP
498     lj_eps(Cl_atom) = 0.3447792E0_DP
499     lj_eps(Ar_atom) = 0.23806656E0_DP
500     lj_eps(Br_atom) = 0.51110784E0_DP
501     lj_eps(Kr_atom) = 0.3259008E0_DP
502    
503     end subroutine initialize_lj
504    
505     subroutine mass_weight()
506     integer ia, ja, dim, dim2, idim, idim2
507     real( kind = DP ) :: mt, m1, m2, wt
508    
509    
510     do ia = 1, natoms
511     m1 = mass(ia)
512     do ja = 1, natoms
513     m2 = mass(ja)
514     wt = 1.0E0_DP/dsqrt(m1*m2)
515     do dim = 1, 3
516     idim = 3 * (ia-1) + dim
517     do dim2 = 1, 3
518     idim2 = 3 * (ja-1) + dim2
519     d(idim,idim2) = d(idim,idim2)*wt
520     enddo
521     enddo
522     enddo
523     enddo
524    
525     end subroutine mass_weight
526    
527    
528    
529    
530     end module lj_module