ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 302
Committed: Mon Mar 10 14:53:36 2003 UTC (21 years, 5 months ago) by gezelter
File size: 16478 byte(s)
Log Message:
Added stuff for Dipole/Dipole and ReactionField

File Contents

# User Rev Content
1 chuckv 292 !! Calculates Long Range forces.
2     !! @author Charles F. Vardeman II
3     !! @author Matthew Meineke
4 gezelter 302 !! @version $Id: do_Forces.F90,v 1.6 2003-03-10 14:53:36 gezelter Exp $, $Date: 2003-03-10 14:53:36 $, $Name: not supported by cvs2svn $, $Revision: 1.6 $
5 chuckv 292
6    
7    
8     module do_Forces
9     use simulation
10     use definitions
11 chuckv 298 use forceGlobals
12     use atype_typedefs
13 chuckv 292 use neighborLists
14 chuckv 298
15 chuckv 292
16     use lj_FF
17     use sticky_FF
18     use dp_FF
19     use gb_FF
20    
21     #ifdef IS_MPI
22     use mpiSimulation
23     #endif
24     implicit none
25     PRIVATE
26    
27 chuckv 298 public :: do_force_loop
28 chuckv 292
29     contains
30    
31     !! FORCE routine Calculates Lennard Jones forces.
32     !------------------------------------------------------------->
33 gezelter 297 subroutine do_force_loop(q,A,mu,u_l,f,t,tau,potE,do_pot,FFerror)
34 chuckv 292 !! Position array provided by C, dimensioned by getNlocal
35     real ( kind = dp ), dimension(3,getNlocal()) :: q
36     !! Rotation Matrix for each long range particle in simulation.
37     real( kind = dp), dimension(9,getNlocal()) :: A
38    
39     !! Magnitude dipole moment
40     real( kind = dp ), dimension(3,getNlocal()) :: mu
41     !! Unit vectors for dipoles (lab frame)
42     real( kind = dp ), dimension(3,getNlocal()) :: u_l
43     !! Force array provided by C, dimensioned by getNlocal
44     real ( kind = dp ), dimension(3,getNlocal()) :: f
45     !! Torsion array provided by C, dimensioned by getNlocal
46     real( kind = dp ), dimension(3,getNlocal()) :: t
47    
48     !! Stress Tensor
49     real( kind = dp), dimension(9) :: tau
50 chuckv 295
51 chuckv 292 real ( kind = dp ) :: potE
52     logical ( kind = 2) :: do_pot
53     integer :: FFerror
54    
55    
56     type(atype), pointer :: Atype_i
57     type(atype), pointer :: Atype_j
58    
59    
60    
61    
62    
63    
64     #ifdef IS_MPI
65     real( kind = DP ) :: pot_local
66    
67     !! Local arrays needed for MPI
68    
69     #endif
70    
71    
72    
73     real( kind = DP ) :: pe
74     logical :: update_nlist
75    
76    
77     integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
78     integer :: nlist
79     integer :: j_start
80 chuckv 295
81     real( kind = DP ) :: r_ij, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
82    
83     real( kind = DP ) :: rx_ij, ry_ij, rz_ij, rijsq
84 chuckv 292 real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut
85    
86     real( kind = DP ) :: dielectric = 0.0_dp
87    
88     ! a rig that need to be fixed.
89     #ifdef IS_MPI
90     real( kind = dp ) :: pe_local
91     integer :: nlocal
92     #endif
93     integer :: nrow
94     integer :: ncol
95     integer :: natoms
96     integer :: neighborListSize
97     integer :: listerror
98     !! should we calculate the stress tensor
99     logical :: do_stress = .false.
100    
101    
102     FFerror = 0
103    
104     ! Make sure we are properly initialized.
105     if (.not. isFFInit) then
106     write(default_error,*) "ERROR: lj_FF has not been properly initialized"
107     FFerror = -1
108     return
109     endif
110     #ifdef IS_MPI
111     if (.not. isMPISimSet()) then
112     write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
113     FFerror = -1
114     return
115     endif
116     #endif
117    
118     !! initialize local variables
119     natoms = getNlocal()
120     call getRcut(rcut,rcut2=rcutsq)
121     call getRlist(rlist,rlistsq)
122 chuckv 295
123 chuckv 292 !! Find ensemble
124     if (isEnsemble("NPT")) do_stress = .true.
125 chuckv 295 !! set to wrap
126     if (isPBC()) wrap = .true.
127 chuckv 292
128 chuckv 295
129 chuckv 292
130    
131     !! See if we need to update neighbor lists
132     call check(q,update_nlist)
133    
134     !--------------WARNING...........................
135     ! Zero variables, NOTE:::: Forces are zeroed in C
136     ! Zeroing them here could delete previously computed
137     ! Forces.
138     !------------------------------------------------
139 chuckv 295 call zero_module_variables()
140 chuckv 292
141    
142     ! communicate MPI positions
143     #ifdef IS_MPI
144     call gather(q,qRow,plan_row3d)
145     call gather(q,qCol,plan_col3d)
146    
147     call gather(mu,muRow,plan_row3d)
148     call gather(mu,muCol,plan_col3d)
149    
150     call gather(u_l,u_lRow,plan_row3d)
151     call gather(u_l,u_lCol,plan_col3d)
152    
153     call gather(A,ARow,plan_row_rotation)
154     call gather(A,ACol,plan_col_rotation)
155     #endif
156    
157    
158 gezelter 297 #ifdef IS_MPI
159 chuckv 292
160 gezelter 297 if (update_nlist) then
161    
162     ! save current configuration, contruct neighbor list,
163     ! and calculate forces
164     call save_neighborList(q)
165    
166     neighborListSize = getNeighborListSize()
167     nlist = 0
168    
169     nrow = getNrow(plan_row)
170     ncol = getNcol(plan_col)
171     nlocal = getNlocal()
172    
173     do i = 1, nrow
174     point(i) = nlist + 1
175     Atype_i => identPtrListRow(i)%this
176    
177     inner: do j = 1, ncol
178     Atype_j => identPtrListColumn(j)%this
179    
180     call get_interatomic_vector(i,j,qRow(:,i),qCol(:,j),&
181     rxij,ryij,rzij,rijsq,r)
182    
183     ! skip the loop if the atoms are identical
184     if (mpi_cycle_jLoop(i,j)) cycle inner:
185    
186     if (rijsq < rlistsq) then
187    
188     nlist = nlist + 1
189    
190     if (nlist > neighborListSize) then
191     call expandList(listerror)
192     if (listerror /= 0) then
193     FFerror = -1
194     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
195     return
196     end if
197     endif
198    
199     list(nlist) = j
200    
201    
202     if (rijsq < rcutsq) then
203     call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
204     endif
205     endif
206     enddo inner
207     enddo
208 chuckv 292
209 gezelter 297 point(nrow + 1) = nlist + 1
210    
211     else !! (update)
212 chuckv 292
213 gezelter 297 ! use the list to find the neighbors
214     do i = 1, nrow
215     JBEG = POINT(i)
216     JEND = POINT(i+1) - 1
217     ! check thiat molecule i has neighbors
218     if (jbeg .le. jend) then
219 chuckv 292
220 gezelter 297 Atype_i => identPtrListRow(i)%this
221     do jnab = jbeg, jend
222     j = list(jnab)
223     Atype_j = identPtrListColumn(j)%this
224     call get_interatomic_vector(i,j,qRow(:,i),qCol(:,j),&
225     rxij,ryij,rzij,rijsq,r)
226    
227     call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
228     enddo
229     endif
230     enddo
231     endif
232    
233     #else
234    
235     if (update_nlist) then
236 chuckv 292
237 gezelter 297 ! save current configuration, contruct neighbor list,
238     ! and calculate forces
239     call save_neighborList(q)
240    
241     neighborListSize = getNeighborListSize()
242     nlist = 0
243    
244    
245     do i = 1, natoms-1
246     point(i) = nlist + 1
247     Atype_i => identPtrList(i)%this
248    
249     inner: do j = i+1, natoms
250     Atype_j => identPtrList(j)%this
251     call get_interatomic_vector(i,j,q(:,i),q(:,j),&
252     rxij,ryij,rzij,rijsq,r)
253 chuckv 295
254 gezelter 297 if (rijsq < rlistsq) then
255 chuckv 292
256 gezelter 297 nlist = nlist + 1
257    
258     if (nlist > neighborListSize) then
259     call expandList(listerror)
260     if (listerror /= 0) then
261     FFerror = -1
262     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
263     return
264     end if
265     endif
266    
267     list(nlist) = j
268 chuckv 295
269 chuckv 292
270 gezelter 297 if (rijsq < rcutsq) then
271     call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
272     endif
273     endif
274 chuckv 292 enddo inner
275 gezelter 297 enddo
276    
277     point(natoms) = nlist + 1
278    
279     else !! (update)
280    
281     ! use the list to find the neighbors
282     do i = 1, nrow
283     JBEG = POINT(i)
284     JEND = POINT(i+1) - 1
285     ! check thiat molecule i has neighbors
286     if (jbeg .le. jend) then
287    
288     Atype_i => identPtrList(i)%this
289     do jnab = jbeg, jend
290     j = list(jnab)
291     Atype_j = identPtrList(j)%this
292     call get_interatomic_vector(i,j,q(:,i),q(:,j),&
293     rxij,ryij,rzij,rijsq,r)
294     call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
295     enddo
296     endif
297     enddo
298     endif
299 chuckv 292
300 gezelter 297 #endif
301    
302    
303 chuckv 292 #ifdef IS_MPI
304 gezelter 297 !! distribute all reaction field stuff (or anything for post-pair):
305     call scatter(rflRow,rflTemp1,plan_row3d)
306     call scatter(rflCol,rflTemp2,plan_col3d)
307     do i = 1,nlocal
308     rflTemp(1:3,i) = rflTemp1(1:3,i) + rflTemp2(1:3,i)
309     end do
310 chuckv 292 #endif
311    
312 gezelter 297 ! This is the post-pair loop:
313     #ifdef IS_MPI
314 chuckv 292
315 gezelter 297 if (system_has_postpair_atoms) then
316     do i = 1, nlocal
317     Atype_i => identPtrListRow(i)%this
318     call do_postpair(i, Atype_i)
319     enddo
320     endif
321 chuckv 295
322 chuckv 292 #else
323 gezelter 297
324     if (system_has_postpair_atoms) then
325     do i = 1, natoms
326     Atype_i => identPtr(i)%this
327     call do_postpair(i, Atype_i)
328     enddo
329     endif
330    
331 chuckv 292 #endif
332    
333 chuckv 295
334 gezelter 297
335    
336 chuckv 292 #ifdef IS_MPI
337 chuckv 295 !!distribute forces
338    
339 gezelter 297 call scatter(fRow,fTemp1,plan_row3d)
340     call scatter(fCol,fTemp2,plan_col3d)
341 chuckv 295
342 gezelter 297
343 chuckv 295 do i = 1,nlocal
344 gezelter 297 fTemp(1:3,i) = fTemp1(1:3,i) + fTemp2(1:3,i)
345 chuckv 295 end do
346    
347     if (do_torque) then
348 gezelter 297 call scatter(tRow,tTemp1,plan_row3d)
349     call scatter(tCol,tTemp2,plan_col3d)
350 chuckv 295
351     do i = 1,nlocal
352 gezelter 297 tTemp(1:3,i) = tTemp1(1:3,i) + tTemp2(1:3,i)
353 chuckv 295 end do
354     endif
355    
356     if (do_pot) then
357     ! scatter/gather pot_row into the members of my column
358     call scatter(eRow,eTemp,plan_row)
359    
360     ! scatter/gather pot_local into all other procs
361     ! add resultant to get total pot
362     do i = 1, nlocal
363     pe_local = pe_local + eTemp(i)
364     enddo
365    
366     eTemp = 0.0E0_DP
367     call scatter(eCol,eTemp,plan_col)
368     do i = 1, nlocal
369     pe_local = pe_local + eTemp(i)
370     enddo
371    
372     pe = pe_local
373     endif
374     #else
375     ! Copy local array into return array for c
376 gezelter 297 f = f+fTemp
377     t = t+tTemp
378 chuckv 295 #endif
379    
380     potE = pe
381    
382    
383     if (do_stress) then
384     #ifdef IS_MPI
385     mpi_allreduce = (tau,tauTemp,9,mpi_double_precision,mpi_sum, &
386     mpi_comm_world,mpi_err)
387     #else
388     tau = tauTemp
389     #endif
390     endif
391    
392     end subroutine do_force_loop
393    
394    
395    
396    
397    
398    
399    
400    
401    
402    
403     !! Calculate any pre-force loop components and update nlist if necessary.
404     subroutine do_preForce(updateNlist)
405     logical, intent(inout) :: updateNlist
406    
407    
408    
409     end subroutine do_preForce
410    
411    
412    
413    
414    
415    
416    
417    
418    
419    
420    
421    
422    
423     !! Calculate any post force loop components, i.e. reaction field, etc.
424     subroutine do_postForce()
425    
426    
427    
428     end subroutine do_postForce
429    
430    
431    
432    
433    
434    
435    
436    
437    
438    
439    
440    
441    
442    
443    
444    
445    
446     subroutine do_pair(atype_i,atype_j,i,j,r_ij,rx_ij,ry_ij,rz_ij)
447     type (atype ), pointer, intent(inout) :: atype_i
448     type (atype ), pointer, intent(inout) :: atype_j
449     integer :: i
450     integer :: j
451     real ( kind = dp ), intent(inout) :: rx_ij
452     real ( kind = dp ), intent(inout) :: ry_ij
453     real ( kind = dp ), intent(inout) :: rz_ij
454    
455    
456     real( kind = dp ) :: fx = 0.0_dp
457     real( kind = dp ) :: fy = 0.0_dp
458     real( kind = dp ) :: fz = 0.0_dp
459 gezelter 297
460 chuckv 295 real( kind = dp ) :: drdx = 0.0_dp
461     real( kind = dp ) :: drdy = 0.0_dp
462     real( kind = dp ) :: drdz = 0.0_dp
463    
464    
465 gezelter 302 #ifdef IS_MPI
466    
467 gezelter 297 if (Atype_i%is_LJ .and. Atype_j%is_LJ) then
468     call getLJForce(r,pot,dudr,ljAtype_i,ljAtype_j,fx,fy,fz)
469     endif
470 chuckv 295
471 gezelter 297 if (Atype_i%is_dp .and. Atype_j%is_dp) then
472 chuckv 295
473 gezelter 297 call dipole_dipole(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
474     ulRow(:,i), ulCol(:,j), rt, rrf, pot)
475 gezelter 302
476     if (do_reaction_field) then
477     call accumulate_rf(i, j, r_ij, rflRow(:,i), rflCol(:j), &
478     ulRow(:i), ulCol(:,j), rt, rrf)
479     endif
480    
481     endif
482    
483     if (Atype_i%is_sticky .and. Atype_j%is_sticky) then
484     call getstickyforce(r,pot,dudr,ljAtype_i,ljAtype_j)
485     endif
486    
487 gezelter 297 #else
488 gezelter 302
489     if (Atype_i%is_LJ .and. Atype_j%is_LJ) then
490     call getLJForce(r,pot,dudr,ljAtype_i,ljAtype_j,fx,fy,fz)
491     endif
492    
493     if (Atype_i%is_dp .and. Atype_j%is_dp) then
494 gezelter 297 call dipole_dipole(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
495     ul(:,i), ul(:,j), rt, rrf, pot)
496 chuckv 295
497 gezelter 297 if (do_reaction_field) then
498     call accumulate_rf(i, j, r_ij, rfl(:,i), rfl(:j), &
499     ul(:,i), ul(:,j), rt, rrf)
500     endif
501 chuckv 295
502 gezelter 297 endif
503    
504     if (Atype_i%is_sticky .and. Atype_j%is_sticky) then
505     call getstickyforce(r,pot,dudr,ljAtype_i,ljAtype_j)
506     endif
507    
508 gezelter 302 #endif
509    
510 chuckv 295
511     #ifdef IS_MPI
512 chuckv 292 eRow(i) = eRow(i) + pot*0.5
513     eCol(i) = eCol(i) + pot*0.5
514     #else
515 chuckv 295 pe = pe + pot
516 chuckv 292 #endif
517 chuckv 295
518 chuckv 292 drdx = -rxij / r
519     drdy = -ryij / r
520     drdz = -rzij / r
521    
522     fx = dudr * drdx
523     fy = dudr * drdy
524     fz = dudr * drdz
525    
526     #ifdef IS_MPI
527     fCol(1,j) = fCol(1,j) - fx
528     fCol(2,j) = fCol(2,j) - fy
529     fCol(3,j) = fCol(3,j) - fz
530    
531     fRow(1,j) = fRow(1,j) + fx
532     fRow(2,j) = fRow(2,j) + fy
533     fRow(3,j) = fRow(3,j) + fz
534     #else
535 chuckv 295 fTemp(1,j) = fTemp(1,j) - fx
536     fTemp(2,j) = fTemp(2,j) - fy
537     fTemp(3,j) = fTemp(3,j) - fz
538     fTemp(1,i) = fTemp(1,i) + fx
539     fTemp(2,i) = fTemp(2,i) + fy
540     fTemp(3,i) = fTemp(3,i) + fz
541 chuckv 292 #endif
542    
543     if (do_stress) then
544     tauTemp(1) = tauTemp(1) + fx * rxij
545     tauTemp(2) = tauTemp(2) + fx * ryij
546     tauTemp(3) = tauTemp(3) + fx * rzij
547     tauTemp(4) = tauTemp(4) + fy * rxij
548     tauTemp(5) = tauTemp(5) + fy * ryij
549     tauTemp(6) = tauTemp(6) + fy * rzij
550     tauTemp(7) = tauTemp(7) + fz * rxij
551     tauTemp(8) = tauTemp(8) + fz * ryij
552     tauTemp(9) = tauTemp(9) + fz * rzij
553     endif
554    
555    
556    
557 chuckv 295 end subroutine do_pair
558 chuckv 292
559    
560    
561 chuckv 295
562    
563    
564    
565    
566    
567    
568    
569    
570    
571    
572    
573    
574     subroutine get_interatomic_vector(q_i,q_j,rx_ij,ry_ij,rz_ij,r_sq,r_ij)
575     !---------------- Arguments-------------------------------
576     !! index i
577    
578     !! Position array
579     real (kind = dp), dimension(3) :: q_i
580     real (kind = dp), dimension(3) :: q_j
581     !! x component of vector between i and j
582     real ( kind = dp ), intent(out) :: rx_ij
583     !! y component of vector between i and j
584     real ( kind = dp ), intent(out) :: ry_ij
585     !! z component of vector between i and j
586     real ( kind = dp ), intent(out) :: rz_ij
587     !! magnitude of r squared
588     real ( kind = dp ), intent(out) :: r_sq
589     !! magnitude of vector r between atoms i and j.
590     real ( kind = dp ), intent(out) :: r_ij
591     !! wrap into periodic box.
592     logical, intent(in) :: wrap
593    
594     !--------------- Local Variables---------------------------
595     !! Distance between i and j
596     real( kind = dp ) :: d(3)
597     !---------------- END DECLARATIONS-------------------------
598    
599    
600     ! Find distance between i and j
601     d(1:3) = q_i(1:3) - q_j(1:3)
602    
603     ! Wrap back into periodic box if necessary
604     if ( wrap ) then
605     d(1:3) = d(1:3) - thisSim%box(1:3) * sign(1.0_dp,thisSim%box(1:3)) * &
606     int(abs(d(1:3)/thisSim%box(1:3) + 0.5_dp)
607     end if
608 chuckv 292
609 chuckv 295 ! Find Magnitude of the vector
610     r_sq = dot_product(d,d)
611     r_ij = sqrt(r_sq)
612 chuckv 292
613 chuckv 295 ! Set each component for force calculation
614     rx_ij = d(1)
615     ry_ij = d(2)
616     rz_ij = d(3)
617 chuckv 292
618    
619 chuckv 295 end subroutine get_interatomic_vector
620 chuckv 292
621 chuckv 295 subroutine zero_module_variables()
622 chuckv 292
623 chuckv 295 #ifndef IS_MPI
624 chuckv 292
625 chuckv 295 pe = 0.0E0_DP
626     tauTemp = 0.0_dp
627     fTemp = 0.0_dp
628     tTemp = 0.0_dp
629 chuckv 292 #else
630 chuckv 295 qRow = 0.0_dp
631     qCol = 0.0_dp
632    
633     muRow = 0.0_dp
634     muCol = 0.0_dp
635    
636     u_lRow = 0.0_dp
637     u_lCol = 0.0_dp
638    
639     ARow = 0.0_dp
640     ACol = 0.0_dp
641    
642     fRow = 0.0_dp
643     fCol = 0.0_dp
644    
645    
646     tRow = 0.0_dp
647     tCol = 0.0_dp
648    
649    
650 chuckv 292
651 chuckv 295 eRow = 0.0_dp
652     eCol = 0.0_dp
653     eTemp = 0.0_dp
654     #endif
655 chuckv 292
656 chuckv 295 end subroutine zero_module_variables
657 chuckv 292
658 chuckv 295 #ifdef IS_MPI
659     !! Function to properly build neighbor lists in MPI using newtons 3rd law.
660     !! We don't want 2 processors doing the same i j pair twice.
661     !! Also checks to see if i and j are the same particle.
662     function mpi_cycle_jLoop(i,j) result(do_cycle)
663     !--------------- Arguments--------------------------
664     ! Index i
665     integer,intent(in) :: i
666     ! Index j
667     integer,intent(in) :: j
668     ! Result do_cycle
669     logical :: do_cycle
670     !--------------- Local variables--------------------
671     integer :: tag_i
672     integer :: tag_j
673     !--------------- END DECLARATIONS------------------
674     tag_i = tagRow(i)
675     tag_j = tagColumn(j)
676 chuckv 292
677 chuckv 295 do_cycle = .false.
678 chuckv 292
679 chuckv 295 if (tag_i == tag_j) then
680     do_cycle = .true.
681     return
682     end if
683    
684     if (tag_i < tag_j) then
685     if (mod(tag_i + tag_j,2) == 0) do_cycle = .true.
686     return
687     else
688     if (mod(tag_i + tag_j,2) == 1) do_cycle = .true.
689     endif
690     end function mpi_cycle_jLoop
691     #endif
692    
693 chuckv 292 end module do_Forces