ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 292
Committed: Thu Mar 6 14:52:44 2003 UTC (21 years, 4 months ago) by chuckv
File size: 18596 byte(s)
Log Message:
Changed code to use one generic force loop. Only one super atype that
has all information is now used.

File Contents

# Content
1 !! Calculates Long Range forces.
2 !! @author Charles F. Vardeman II
3 !! @author Matthew Meineke
4 !! @version $Id: do_Forces.F90,v 1.1 2003-03-06 14:52:44 chuckv Exp $, $Date: 2003-03-06 14:52:44 $, $Name: not supported by cvs2svn $, $Revision: 1.1 $
5
6
7
8 module do_Forces
9 use simulation
10 use definitions
11 use generic_atypes
12 use neighborLists
13
14 use lj_FF
15 use sticky_FF
16 use dp_FF
17 use gb_FF
18
19 #ifdef IS_MPI
20 use mpiSimulation
21 #endif
22 implicit none
23 PRIVATE
24
25 !! Number of lj_atypes in lj_atype_list
26 integer, save :: n_atypes = 0
27
28 !! Global list of lj atypes in simulation
29 type (atype), pointer :: ListHead => null()
30 type (atype), pointer :: ListTail => null()
31
32
33 !! LJ mixing array
34 ! type (lj_atype), dimension(:,:), pointer :: ljMixed => null()
35
36
37 logical, save :: firstTime = .True.
38
39 !! Atype identity pointer lists
40 #ifdef IS_MPI
41 !! Row lj_atype pointer list
42 type (identPtrList), dimension(:), pointer :: identPtrListRow => null()
43 !! Column lj_atype pointer list
44 type (identPtrList), dimension(:), pointer :: identPtrListColumn => null()
45 #else
46 type(identPtrList ), dimension(:), pointer :: identPtrList => null()
47 #endif
48
49
50 !! Logical has lj force field module been initialized?
51 logical, save :: isFFinit = .false.
52
53
54 !! Public methods and data
55 public :: new_atype
56 public :: do_lj_ff
57 public :: getLjPot
58 public :: init_FF
59
60
61
62
63 contains
64
65 !! Adds a new lj_atype to the list.
66 subroutine new_atype(ident,mass,epsilon,sigma, &
67 is_LJ,is_Sticky,is_DP,is_GB,w0,v0,dipoleMoment,status)
68 real( kind = dp ), intent(in) :: mass
69 real( kind = dp ), intent(in) :: epsilon
70 real( kind = dp ), intent(in) :: sigma
71 real( kind = dp ), intent(in) :: w0
72 real( kind = dp ), intent(in) :: v0
73 real( kind = dp ), intent(in) :: dipoleMoment
74
75 integer, intent(in) :: ident
76 integer, intent(out) :: status
77 integer, intent(in) :: is_Sticky
78 integer, intent(in) :: is_DP
79 integer, intent(in) :: is_GB
80 integer, intent(in) :: is_LJ
81
82
83 type (atype), pointer :: the_new_atype
84 integer :: alloc_error
85 integer :: atype_counter = 0
86 integer :: alloc_size
87 integer :: err_stat
88 status = 0
89
90
91
92 ! allocate a new atype
93 allocate(the_new_atype,stat=alloc_error)
94 if (alloc_error /= 0 ) then
95 status = -1
96 return
97 end if
98
99 ! assign our new atype information
100 the_new_atype%mass = mass
101 the_new_atype%epsilon = epsilon
102 the_new_atype%sigma = sigma
103 the_new_atype%sigma2 = sigma * sigma
104 the_new_atype%sigma6 = the_new_atype%sigma2 * the_new_atype%sigma2 &
105 * the_new_atype%sigma2
106 the_new_atype%w0 = w0
107 the_new_atype%v0 = v0
108 the_new_atype%dipoleMoment = dipoleMoment
109
110
111 ! assume that this atype will be successfully added
112 the_new_atype%atype_ident = ident
113 the_new_atype%atype_number = n_lj_atypes + 1
114
115 if ( is_Sticky /= 0 ) the_new_atype%is_Sticky = .true.
116 if ( is_GB /= 0 ) the_new_atype%is_GB = .true.
117 if ( is_LJ /= 0 ) the_new_atype%is_LJ = .true.
118 if ( is_DP /= 0 ) the_new_atype%is_DP = .true.
119
120 call add_atype(the_new_atype,ListHead,ListTail,err_stat)
121 if (err_stat /= 0 ) then
122 status = -1
123 return
124 endif
125
126 n_atypes = n_atypes + 1
127
128
129 end subroutine new_atype
130
131
132 subroutine init_FF(nComponents,ident, status)
133 !! Number of components in ident array
134 integer, intent(inout) :: nComponents
135 !! Array of identities nComponents long corresponding to
136 !! ljatype ident.
137 integer, dimension(nComponents),intent(inout) :: ident
138 !! Result status, success = 0, error = -1
139 integer, intent(out) :: Status
140
141 integer :: alloc_stat
142
143 integer :: thisStat
144 integer :: i
145
146 integer :: myNode
147 #ifdef IS_MPI
148 integer, allocatable, dimension(:) :: identRow
149 integer, allocatable, dimension(:) :: identCol
150 integer :: nrow
151 integer :: ncol
152 #endif
153 status = 0
154
155
156
157
158 !! if were're not in MPI, we just update ljatypePtrList
159 #ifndef IS_MPI
160 call create_IdentPtrlst(ident,ListHead,identPtrList,thisStat)
161 if ( thisStat /= 0 ) then
162 status = -1
163 return
164 endif
165
166
167 ! if were're in MPI, we also have to worry about row and col lists
168 #else
169
170 ! We can only set up forces if mpiSimulation has been setup.
171 if (.not. isMPISimSet()) then
172 write(default_error,*) "MPI is not set"
173 status = -1
174 return
175 endif
176 nrow = getNrow(plan_row)
177 ncol = getNcol(plan_col)
178 mynode = getMyNode()
179 !! Allocate temperary arrays to hold gather information
180 allocate(identRow(nrow),stat=alloc_stat)
181 if (alloc_stat /= 0 ) then
182 status = -1
183 return
184 endif
185
186 allocate(identCol(ncol),stat=alloc_stat)
187 if (alloc_stat /= 0 ) then
188 status = -1
189 return
190 endif
191
192 !! Gather idents into row and column idents
193
194 call gather(ident,identRow,plan_row)
195 call gather(ident,identCol,plan_col)
196
197
198 !! Create row and col pointer lists
199
200 call create_IdentPtrlst(identRow,ListHead,identPtrListRow,thisStat)
201 if (thisStat /= 0 ) then
202 status = -1
203 return
204 endif
205
206 call create_IdentPtrlst(identCol,ListHead,identPtrListColumn,thisStat)
207 if (thisStat /= 0 ) then
208 status = -1
209 return
210 endif
211
212 !! free temporary ident arrays
213 if (allocated(identCol)) then
214 deallocate(identCol)
215 end if
216 if (allocated(identCol)) then
217 deallocate(identRow)
218 endif
219
220 #endif
221
222 call initForce_Modules(thisStat)
223 if (thisStat /= 0) then
224 status = -1
225 return
226 endif
227
228 !! Create neighbor lists
229 call expandList(thisStat)
230 if (thisStat /= 0) then
231 status = -1
232 return
233 endif
234
235 isFFinit = .true.
236
237
238 end subroutine init_FF
239
240
241
242
243 subroutine initForce_Modules(thisStat)
244 integer, intent(out) :: thisStat
245 integer :: my_status
246
247 thisStat = 0
248 call init_lj_FF(ListHead,my_status)
249 if (my_status /= 0) then
250 thisStat = -1
251 return
252 end if
253
254 end subroutine initForce_Modules
255
256
257
258
259 !! FORCE routine Calculates Lennard Jones forces.
260 !------------------------------------------------------------->
261 subroutine do__force_loop(q,A,mu,u_l,f,t,tau,potE,do_pot,FFerror)
262 !! Position array provided by C, dimensioned by getNlocal
263 real ( kind = dp ), dimension(3,getNlocal()) :: q
264 !! Rotation Matrix for each long range particle in simulation.
265 real( kind = dp), dimension(9,getNlocal()) :: A
266
267 !! Magnitude dipole moment
268 real( kind = dp ), dimension(3,getNlocal()) :: mu
269 !! Unit vectors for dipoles (lab frame)
270 real( kind = dp ), dimension(3,getNlocal()) :: u_l
271 !! Force array provided by C, dimensioned by getNlocal
272 real ( kind = dp ), dimension(3,getNlocal()) :: f
273 !! Torsion array provided by C, dimensioned by getNlocal
274 real( kind = dp ), dimension(3,getNlocal()) :: t
275
276 !! Stress Tensor
277 real( kind = dp), dimension(9) :: tau
278 real( kind = dp), dimension(9) :: tauTemp
279 real ( kind = dp ) :: potE
280 logical ( kind = 2) :: do_pot
281 integer :: FFerror
282
283
284 type(atype), pointer :: Atype_i
285 type(atype), pointer :: Atype_j
286
287
288
289
290
291
292 #ifdef IS_MPI
293 real( kind = DP ) :: pot_local
294
295 !! Local arrays needed for MPI
296 real(kind = dp), dimension(3,getNrow(plan_row)) :: qRow = 0.0_dp
297 real(kind = dp), dimension(3,getNcol(plan_col)) :: qCol = 0.0_dp
298
299 real(kind = dp), dimension(3,getNrow(plan_row)) :: muRow = 0.0_dp
300 real(kind = dp), dimension(3,getNcol(plan_col)) :: muCol = 0.0_dp
301
302 real(kind = dp), dimension(3,getNrow(plan_row)) :: u_lRow = 0.0_dp
303 real(kind = dp), dimension(3,getNcol(plan_col)) :: u_lCol = 0.0_dp
304
305 real(kind = dp), dimension(3,getNrow(plan_row)) :: ARow = 0.0_dp
306 real(kind = dp), dimension(3,getNcol(plan_col)) :: ACol = 0.0_dp
307
308
309
310 real(kind = dp), dimension(3,getNrow(plan_row)) :: fRow = 0.0_dp
311 real(kind = dp), dimension(3,getNcol(plan_col)) :: fCol = 0.0_dp
312 real(kind = dp), dimension(3,getNlocal()) :: fMPITemp = 0.0_dp
313
314 real(kind = dp), dimension(3,getNrow(plan_row)) :: tRow = 0.0_dp
315 real(kind = dp), dimension(3,getNcol(plan_col)) :: tCol = 0.0_dp
316 real(kind = dp), dimension(3,getNlocal()) :: tMPITemp = 0.0_dp
317
318
319 real(kind = dp), dimension(getNrow(plan_row)) :: eRow = 0.0_dp
320 real(kind = dp), dimension(getNcol(plan_col)) :: eCol = 0.0_dp
321
322 real(kind = dp), dimension(getNlocal()) :: eTemp = 0.0_dp
323
324 #endif
325
326
327
328 real( kind = DP ) :: pe
329 logical :: update_nlist
330
331
332 integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
333 integer :: nlist
334 integer :: j_start
335 integer :: tag_i,tag_j
336 real( kind = DP ) :: r, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
337 real( kind = dp ) :: fx,fy,fz
338 real( kind = DP ) :: drdx, drdy, drdz
339 real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq
340 real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut
341
342 real( kind = DP ) :: dielectric = 0.0_dp
343
344 ! a rig that need to be fixed.
345 #ifdef IS_MPI
346 logical :: newtons_thrd = .true.
347 real( kind = dp ) :: pe_local
348 integer :: nlocal
349 #endif
350 integer :: nrow
351 integer :: ncol
352 integer :: natoms
353 integer :: neighborListSize
354 integer :: listerror
355 !! should we calculate the stress tensor
356 logical :: do_stress = .false.
357
358
359 FFerror = 0
360
361 ! Make sure we are properly initialized.
362 if (.not. isFFInit) then
363 write(default_error,*) "ERROR: lj_FF has not been properly initialized"
364 FFerror = -1
365 return
366 endif
367 #ifdef IS_MPI
368 if (.not. isMPISimSet()) then
369 write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
370 FFerror = -1
371 return
372 endif
373 #endif
374
375 !! initialize local variables
376 natoms = getNlocal()
377 call getRcut(rcut,rcut2=rcutsq)
378 call getRlist(rlist,rlistsq)
379 !! Find ensemble
380 if (isEnsemble("NPT")) do_stress = .true.
381
382 #ifndef IS_MPI
383 nrow = natoms - 1
384 ncol = natoms
385 #else
386 nrow = getNrow(plan_row)
387 ncol = getNcol(plan_col)
388 nlocal = natoms
389 j_start = 1
390 #endif
391
392
393 !! See if we need to update neighbor lists
394 call check(q,update_nlist)
395 ! if (firstTime) then
396 ! update_nlist = .true.
397 ! firstTime = .false.
398 ! endif
399
400 !--------------WARNING...........................
401 ! Zero variables, NOTE:::: Forces are zeroed in C
402 ! Zeroing them here could delete previously computed
403 ! Forces.
404 !------------------------------------------------
405 #ifndef IS_MPI
406 ! nloops = nloops + 1
407 pe = 0.0E0_DP
408
409 #else
410 fRow = 0.0E0_DP
411 fCol = 0.0E0_DP
412
413 pe_local = 0.0E0_DP
414
415 eRow = 0.0E0_DP
416 eCol = 0.0E0_DP
417 eTemp = 0.0E0_DP
418 #endif
419
420 ! communicate MPI positions
421 #ifdef IS_MPI
422 call gather(q,qRow,plan_row3d)
423 call gather(q,qCol,plan_col3d)
424
425 call gather(mu,muRow,plan_row3d)
426 call gather(mu,muCol,plan_col3d)
427
428 call gather(u_l,u_lRow,plan_row3d)
429 call gather(u_l,u_lCol,plan_col3d)
430
431 call gather(A,ARow,plan_row_rotation)
432 call gather(A,ACol,plan_col_rotation)
433
434 #endif
435
436
437 if (update_nlist) then
438
439 ! save current configuration, contruct neighbor list,
440 ! and calculate forces
441 call save_neighborList(q)
442
443 neighborListSize = getNeighborListSize()
444 nlist = 0
445
446
447
448 do i = 1, nrow
449 point(i) = nlist + 1
450 #ifdef IS_MPI
451 ljAtype_i => identPtrListRow(i)%this
452 tag_i = tagRow(i)
453 rxi = qRow(1,i)
454 ryi = qRow(2,i)
455 rzi = qRow(3,i)
456 #else
457 ljAtype_i => identPtrList(i)%this
458 j_start = i + 1
459 rxi = q(1,i)
460 ryi = q(2,i)
461 rzi = q(3,i)
462 #endif
463
464 inner: do j = j_start, ncol
465 #ifdef IS_MPI
466 ! Assign identity pointers and tags
467 ljAtype_j => identPtrListColumn(j)%this
468 tag_j = tagColumn(j)
469 if (newtons_thrd) then
470 if (tag_i <= tag_j) then
471 if (mod(tag_i + tag_j,2) == 0) cycle inner
472 else
473 if (mod(tag_i + tag_j,2) == 1) cycle inner
474 endif
475 endif
476
477 rxij = wrap(rxi - qCol(1,j), 1)
478 ryij = wrap(ryi - qCol(2,j), 2)
479 rzij = wrap(rzi - qCol(3,j), 3)
480 #else
481 ljAtype_j => identPtrList(j)%this
482 rxij = wrap(rxi - q(1,j), 1)
483 ryij = wrap(ryi - q(2,j), 2)
484 rzij = wrap(rzi - q(3,j), 3)
485
486 #endif
487 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
488
489 #ifdef IS_MPI
490 if (rijsq <= rlistsq .AND. &
491 tag_j /= tag_i) then
492 #else
493
494 if (rijsq < rlistsq) then
495 #endif
496
497 nlist = nlist + 1
498 if (nlist > neighborListSize) then
499 call expandList(listerror)
500 if (listerror /= 0) then
501 FFerror = -1
502 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
503 return
504 end if
505 endif
506 list(nlist) = j
507
508
509 if (rijsq < rcutsq) then
510
511 r = dsqrt(rijsq)
512
513 call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j)
514
515 #ifdef IS_MPI
516 eRow(i) = eRow(i) + pot*0.5
517 eCol(i) = eCol(i) + pot*0.5
518 #else
519 pe = pe + pot
520 #endif
521
522 drdx = -rxij / r
523 drdy = -ryij / r
524 drdz = -rzij / r
525
526 fx = dudr * drdx
527 fy = dudr * drdy
528 fz = dudr * drdz
529
530 #ifdef IS_MPI
531 fCol(1,j) = fCol(1,j) - fx
532 fCol(2,j) = fCol(2,j) - fy
533 fCol(3,j) = fCol(3,j) - fz
534
535 fRow(1,j) = fRow(1,j) + fx
536 fRow(2,j) = fRow(2,j) + fy
537 fRow(3,j) = fRow(3,j) + fz
538 #else
539 f(1,j) = f(1,j) - fx
540 f(2,j) = f(2,j) - fy
541 f(3,j) = f(3,j) - fz
542 f(1,i) = f(1,i) + fx
543 f(2,i) = f(2,i) + fy
544 f(3,i) = f(3,i) + fz
545 #endif
546
547 if (do_stress) then
548 tauTemp(1) = tauTemp(1) + fx * rxij
549 tauTemp(2) = tauTemp(2) + fx * ryij
550 tauTemp(3) = tauTemp(3) + fx * rzij
551 tauTemp(4) = tauTemp(4) + fy * rxij
552 tauTemp(5) = tauTemp(5) + fy * ryij
553 tauTemp(6) = tauTemp(6) + fy * rzij
554 tauTemp(7) = tauTemp(7) + fz * rxij
555 tauTemp(8) = tauTemp(8) + fz * ryij
556 tauTemp(9) = tauTemp(9) + fz * rzij
557 endif
558 endif
559 enddo inner
560 enddo
561
562 #ifdef IS_MPI
563 point(nrow + 1) = nlist + 1
564 #else
565 point(natoms) = nlist + 1
566 #endif
567
568 else
569
570 ! use the list to find the neighbors
571 do i = 1, nrow
572 JBEG = POINT(i)
573 JEND = POINT(i+1) - 1
574 ! check thiat molecule i has neighbors
575 if (jbeg .le. jend) then
576 #ifdef IS_MPI
577 ljAtype_i => identPtrListRow(i)%this
578 rxi = qRow(1,i)
579 ryi = qRow(2,i)
580 rzi = qRow(3,i)
581 #else
582 ljAtype_i => identPtrList(i)%this
583 rxi = q(1,i)
584 ryi = q(2,i)
585 rzi = q(3,i)
586 #endif
587 do jnab = jbeg, jend
588 j = list(jnab)
589 #ifdef IS_MPI
590 ljAtype_j = identPtrListColumn(j)%this
591 rxij = wrap(rxi - qCol(1,j), 1)
592 ryij = wrap(ryi - qCol(2,j), 2)
593 rzij = wrap(rzi - qCol(3,j), 3)
594 #else
595 ljAtype_j = identPtrList(j)%this
596 rxij = wrap(rxi - q(1,j), 1)
597 ryij = wrap(ryi - q(2,j), 2)
598 rzij = wrap(rzi - q(3,j), 3)
599 #endif
600 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
601
602 if (rijsq < rcutsq) then
603
604 r = dsqrt(rijsq)
605
606 call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j)
607 #ifdef IS_MPI
608 eRow(i) = eRow(i) + pot*0.5
609 eCol(i) = eCol(i) + pot*0.5
610 #else
611 pe = pe + pot
612 #endif
613
614 drdx = -rxij / r
615 drdy = -ryij / r
616 drdz = -rzij / r
617
618 fx = dudr * drdx
619 fy = dudr * drdy
620 fz = dudr * drdz
621
622 #ifdef IS_MPI
623 fCol(1,j) = fCol(1,j) - fx
624 fCol(2,j) = fCol(2,j) - fy
625 fCol(3,j) = fCol(3,j) - fz
626
627 fRow(1,j) = fRow(1,j) + fx
628 fRow(2,j) = fRow(2,j) + fy
629 fRow(3,j) = fRow(3,j) + fz
630 #else
631 f(1,j) = f(1,j) - fx
632 f(2,j) = f(2,j) - fy
633 f(3,j) = f(3,j) - fz
634 f(1,i) = f(1,i) + fx
635 f(2,i) = f(2,i) + fy
636 f(3,i) = f(3,i) + fz
637 #endif
638
639 if (do_stress) then
640 tauTemp(1) = tauTemp(1) + fx * rxij
641 tauTemp(2) = tauTemp(2) + fx * ryij
642 tauTemp(3) = tauTemp(3) + fx * rzij
643 tauTemp(4) = tauTemp(4) + fy * rxij
644 tauTemp(5) = tauTemp(5) + fy * ryij
645 tauTemp(6) = tauTemp(6) + fy * rzij
646 tauTemp(7) = tauTemp(7) + fz * rxij
647 tauTemp(8) = tauTemp(8) + fz * ryij
648 tauTemp(9) = tauTemp(9) + fz * rzij
649 endif
650
651
652 endif
653 enddo
654 endif
655 enddo
656 endif
657
658
659
660 #ifdef IS_MPI
661 !!distribute forces
662
663 call scatter(fRow,f,plan_row3d)
664 call scatter(fCol,fMPITemp,plan_col3d)
665
666 do i = 1,nlocal
667 f(1:3,i) = f(1:3,i) + fMPITemp(1:3,i)
668 end do
669
670
671 call scatter(tRow,t,plan_row3d)
672 call scatter(tCol,tMPITemp,plan_col3d)
673
674 do i = 1,nlocal
675 t(1:3,i) = t(1:3,i) + tMPITemp(1:3,i)
676 end do
677
678
679 if (do_pot) then
680 ! scatter/gather pot_row into the members of my column
681 call scatter(eRow,eTemp,plan_row)
682
683 ! scatter/gather pot_local into all other procs
684 ! add resultant to get total pot
685 do i = 1, nlocal
686 pe_local = pe_local + eTemp(i)
687 enddo
688 if (newtons_thrd) then
689 eTemp = 0.0E0_DP
690 call scatter(eCol,eTemp,plan_col)
691 do i = 1, nlocal
692 pe_local = pe_local + eTemp(i)
693 enddo
694 endif
695 pe = pe_local
696 endif
697
698 #endif
699
700 potE = pe
701
702
703 if (do_stress) then
704 #ifdef IS_MPI
705 mpi_allreduce = (tau,tauTemp,9,mpi_double_precision,mpi_sum, &
706 mpi_comm_world,mpi_err)
707 #else
708 tau = tauTemp
709 #endif
710 endif
711
712
713 end subroutine do_force_loop
714
715
716
717 end module do_Forces