ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
(Generate patch)

Comparing trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90 (file contents):
Revision 2279 by chrisfen, Tue Aug 30 18:23:50 2005 UTC vs.
Revision 2301 by gezelter, Thu Sep 15 22:05:21 2005 UTC

# Line 54 | Line 54 | module electrostatic_module
54  
55    PRIVATE
56  
57 + #define __FORTRAN90
58 + #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59 +
60    !! these prefactors convert the multipole interactions into kcal / mol
61    !! all were computed assuming distances are measured in angstroms
62    !! Charge-Charge, assuming charges are measured in electrons
# Line 68 | Line 71 | module electrostatic_module
71    !! This unit is also known affectionately as an esu centi-barn.
72    real(kind=dp), parameter :: pre14 = 69.13373_dp
73  
74 +  !! variables to handle different summation methods for long-range electrostatics:
75 +  integer, save :: summationMethod = NONE
76 +  real(kind=DP), save :: defaultCutoff = 0.0_DP
77 +  logical, save :: haveDefaultCutoff = .false.
78 +  real(kind=DP), save :: dampingAlpha = 0.0_DP
79 +  logical, save :: haveDampingAlpha = .false.
80 +  real(kind=DP), save :: dielectric = 0.0_DP
81 +  logical, save :: haveDielectric = .false.
82 +  real(kind=DP), save :: constERFC = 0.0_DP
83 +  real(kind=DP), save :: constEXP = 0.0_DP
84 +  logical, save :: haveDWAconstants = .false.
85 +
86 +
87 +  public :: setElectrostaticSummationMethod
88 +  public :: setElectrostaticCutoffRadius
89 +  public :: setDampedWolfAlpha
90 +  public :: setReactionFieldDielectric
91    public :: newElectrostaticType
92    public :: setCharge
93    public :: setDipoleMoment
# Line 95 | Line 115 | contains
115    type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
116  
117   contains
118 +
119 +  subroutine setElectrostaticSummationMethod(the_ESM)
120 +
121 +    integer, intent(in) :: the_ESM    
122 +
123 +    if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then
124 +       call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method")
125 +    endif
126 +
127 +  end subroutine setElectrostaticSummationMethod
128 +
129 +  subroutine setElectrostaticCutoffRadius(thisRcut)
130 +    real(kind=dp), intent(in) :: thisRcut
131 +    defaultCutoff = thisRcut
132 +    haveDefaultCutoff = .true.
133 +  end subroutine setElectrostaticCutoffRadius
134 +
135 +  subroutine setDampedWolfAlpha(thisAlpha)
136 +    real(kind=dp), intent(in) :: thisAlpha
137 +    dampingAlpha = thisAlpha
138 +    haveDampingAlpha = .true.
139 +  end subroutine setDampedWolfAlpha
140 +  
141 +  subroutine setReactionFieldDielectric(thisDielectric)
142 +    real(kind=dp), intent(in) :: thisDielectric
143 +    dielectric = thisDielectric
144 +    haveDielectric = .true.
145 +  end subroutine setReactionFieldDielectric
146  
147    subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
148         is_SplitDipole, is_Quadrupole, is_Tap, status)
# Line 306 | Line 354 | contains
354      dm = ElectrostaticMap(atid)%dipole_moment
355    end function getDipoleMoment
356  
357 +  subroutine checkSummationMethod()
358 +
359 +    if (summationMethod .eq. DAMPED_WOLF) then
360 +       if (.not.haveDWAconstants) then
361 +          
362 +          if (.not.haveDampingAlpha) then
363 +             call handleError("checkSummationMethod", "no Damping Alpha set!")
364 +          endif
365 +          
366 +          if (.not.have....)
367 +          constEXP =
368 +          constERFC =
369 +          
370 +          haveDWAconstants = .true.
371 +       endif
372 +    endif
373 +
374 +  end subroutine checkSummationMethod
375 +
376 +
377 +
378    subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, &
379 <       vpair, fpair, pot, eFrame, f, t, do_pot, corrMethod)
379 >       vpair, fpair, pot, eFrame, f, t, do_pot, corrMethod, rcuti)
380  
381      logical, intent(in) :: do_pot
382  
# Line 315 | Line 384 | contains
384      integer :: localError
385      integer, intent(in) :: corrMethod
386  
387 <    real(kind=dp), intent(in) :: rij, r2, sw
387 >    real(kind=dp), intent(in) :: rij, r2, sw, rcuti
388      real(kind=dp), intent(in), dimension(3) :: d
389      real(kind=dp), intent(inout) :: vpair
390      real(kind=dp), intent(inout), dimension(3) :: fpair
391  
392 <    real( kind = dp ) :: pot
392 >    real( kind = dp ) :: pot, swi
393      real( kind = dp ), dimension(9,nLocal) :: eFrame
394      real( kind = dp ), dimension(3,nLocal) :: f
395      real( kind = dp ), dimension(3,nLocal) :: t
# Line 333 | Line 402 | contains
402      logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
403      logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
404      logical :: i_is_Tap, j_is_Tap
336    logical :: use_damped_wolf, use_undamped_wolf
405      integer :: me1, me2, id1, id2
406      real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
407      real (kind=dp) :: qxx_i, qyy_i, qzz_i
# Line 343 | Line 411 | contains
411      real (kind=dp) :: cx2, cy2, cz2
412      real (kind=dp) :: ct_i, ct_j, ct_ij, a1
413      real (kind=dp) :: riji, ri, ri2, ri3, ri4
414 <    real (kind=dp) :: pref, vterm, epot, dudr    
414 >    real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2
415      real (kind=dp) :: xhat, yhat, zhat
416      real (kind=dp) :: dudx, dudy, dudz
417      real (kind=dp) :: scale, sc2, bigR, switcher, dswitcher
418 +    real (kind=dp) :: rcuti2, rcuti3, rcuti4
419  
351    use_damped_wolf = .false.
352    use_undamped_wolf = .false.
353    if (corrMethod .eq. 1) then
354       use_undamped_wolf = .true.
355    elseif (corrMethod .eq. 2) then
356       use_damped_wolf = .true.
357    endif
358
420      if (.not.allocated(ElectrostaticMap)) then
421         call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
422         return
423      end if
424 +
425 +    if (.not.summationMethodChecked) then
426 +       call checkSummationMethod()
427 +    endif
428  
429 +
430   #ifdef IS_MPI
431      me1 = atid_Row(atom1)
432      me2 = atid_Col(atom2)
# Line 377 | Line 443 | contains
443      yhat = d(2) * riji
444      zhat = d(3) * riji
445  
446 +    rcuti2 = rcuti*rcuti
447 +    rcuti3 = rcuti2*rcuti
448 +    rcuti4 = rcuti2*rcuti2
449 +
450 +    swi = 1.0d0 / sw
451 +
452      !! logicals
453      i_is_Charge = ElectrostaticMap(me1)%is_Charge
454      i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
# Line 520 | Line 592 | contains
592  
593         if (j_is_Charge) then
594  
595 <          vterm = pre11 * q_i * q_j * riji
596 <          vpair = vpair + vterm
525 <          epot = epot + sw*vterm
595 >          if (corrMethod .eq. 1) then
596 >             vterm = pre11 * q_i * q_j * (riji - rcuti)
597  
598 <          dudr  = - sw * vterm * riji
598 >             vpair = vpair + vterm
599 >             epot = epot + sw * vterm
600 >            
601 >             dudr  = - sw * pre11 * q_i * q_j * (riji*riji*riji - rcuti2*rcuti)
602 >            
603 >             dudx = dudx + dudr * d(1)
604 >             dudy = dudy + dudr * d(2)
605 >             dudz = dudz + dudr * d(3)
606  
607 <          dudx = dudx + dudr * xhat
608 <          dudy = dudy + dudr * yhat
531 <          dudz = dudz + dudr * zhat
607 >          else
608 >             vterm = pre11 * q_i * q_j * riji
609  
610 <       endif
610 >             vpair = vpair + vterm
611 >             epot = epot + sw * vterm
612 >            
613 >             dudr  = - sw * vterm * riji
614 >            
615 >             dudx = dudx + dudr * xhat
616 >             dudy = dudy + dudr * yhat
617 >             dudz = dudz + dudr * zhat
618  
535       if (j_is_Dipole) then
536
537          if (j_is_SplitDipole) then
538             BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
539             ri = 1.0_dp / BigR
540             scale = rij * ri
541          else
542             ri = riji
543             scale = 1.0_dp
619            endif
620  
621 <          ri2 = ri * ri
547 <          ri3 = ri2 * ri
548 <          sc2 = scale * scale
621 >       endif
622  
623 <          pref = pre12 * q_i * mu_j
551 <          vterm = - pref * ct_j * ri2 * scale
552 <          vpair = vpair + vterm
553 <          epot = epot + sw * vterm
623 >       if (j_is_Dipole) then
624  
625 <          !! this has a + sign in the () because the rij vector is
556 <          !! r_j - r_i and the charge-dipole potential takes the origin
557 <          !! as the point dipole, which is atom j in this case.
625 >          pref = sw * pre12 * q_i * mu_j
626  
627 <          dudx = dudx - pref * sw * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2)
628 <          dudy = dudy - pref * sw * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2)
629 <          dudz = dudz - pref * sw * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2)
627 >          if (corrMethod .eq. 1) then
628 >             ri2 = riji * riji
629 >             ri3 = ri2 * riji
630  
631 <          duduz_j(1) = duduz_j(1) - pref * sw * ri2 * xhat * scale
632 <          duduz_j(2) = duduz_j(2) - pref * sw * ri2 * yhat * scale
633 <          duduz_j(3) = duduz_j(3) - pref * sw * ri2 * zhat * scale
631 >             vterm = - pref * ct_j * (ri2 - rcuti2)
632 >             vpair = vpair + swi*vterm
633 >             epot = epot + vterm
634 >            
635 >             !! this has a + sign in the () because the rij vector is
636 >             !! r_j - r_i and the charge-dipole potential takes the origin
637 >             !! as the point dipole, which is atom j in this case.
638 >            
639 >             dudx = dudx - pref * ( ri3*( uz_j(1) - 3.0d0*ct_j*xhat) &
640 >                  - rcuti3*( uz_j(1) - 3.0d0*ct_j*d(1)*rcuti ) )
641 >             dudy = dudy - pref * ( ri3*( uz_j(2) - 3.0d0*ct_j*yhat) &
642 >                  - rcuti3*( uz_j(2) - 3.0d0*ct_j*d(2)*rcuti ) )
643 >             dudz = dudz - pref * ( ri3*( uz_j(3) - 3.0d0*ct_j*zhat) &
644 >                  - rcuti3*( uz_j(3) - 3.0d0*ct_j*d(3)*rcuti ) )
645 >            
646 >             duduz_j(1) = duduz_j(1) - pref*( ri2*xhat - d(1)*rcuti3 )
647 >             duduz_j(2) = duduz_j(2) - pref*( ri2*yhat - d(2)*rcuti3 )
648 >             duduz_j(3) = duduz_j(3) - pref*( ri2*zhat - d(3)*rcuti3 )
649  
650 +          else
651 +             if (j_is_SplitDipole) then
652 +                BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
653 +                ri = 1.0_dp / BigR
654 +                scale = rij * ri
655 +             else
656 +                ri = riji
657 +                scale = 1.0_dp
658 +             endif
659 +            
660 +             ri2 = ri * ri
661 +             ri3 = ri2 * ri
662 +             sc2 = scale * scale
663 +            
664 +             vterm = - pref * ct_j * ri2 * scale
665 +             vpair = vpair + swi * vterm
666 +             epot = epot + vterm
667 +            
668 +             !! this has a + sign in the () because the rij vector is
669 +             !! r_j - r_i and the charge-dipole potential takes the origin
670 +             !! as the point dipole, which is atom j in this case.
671 +            
672 +             dudx = dudx - pref * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2)
673 +             dudy = dudy - pref * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2)
674 +             dudz = dudz - pref * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2)
675 +            
676 +             duduz_j(1) = duduz_j(1) - pref * ri2 * xhat * scale
677 +             duduz_j(2) = duduz_j(2) - pref * ri2 * yhat * scale
678 +             duduz_j(3) = duduz_j(3) - pref * ri2 * zhat * scale
679 +
680 +          endif
681         endif
682  
683         if (j_is_Quadrupole) then
# Line 575 | Line 689 | contains
689            cz2 = cz_j * cz_j
690  
691  
692 <          pref =  pre14 * q_i / 3.0_dp
579 <          vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
580 <               qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
581 <               qzz_j * (3.0_dp*cz2 - 1.0_dp))
582 <          vpair = vpair + vterm
583 <          epot = epot + sw * vterm
692 >          pref =  sw * pre14 * q_i / 3.0_dp
693  
694 <          dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
695 <               qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
696 <               qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
697 <               qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
698 <          dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
699 <               qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
700 <               qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
701 <               qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
702 <          dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
703 <               qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
704 <               qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
705 <               qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
706 <
707 <          dudux_j(1) = dudux_j(1) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*xhat)
708 <          dudux_j(2) = dudux_j(2) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*yhat)
709 <          dudux_j(3) = dudux_j(3) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*zhat)
710 <
711 <          duduy_j(1) = duduy_j(1) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*xhat)
712 <          duduy_j(2) = duduy_j(2) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*yhat)
713 <          duduy_j(3) = duduy_j(3) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*zhat)
714 <
715 <          duduz_j(1) = duduz_j(1) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*xhat)
716 <          duduz_j(2) = duduz_j(2) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*yhat)
717 <          duduz_j(3) = duduz_j(3) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*zhat)
694 >          if (corrMethod .eq. 1) then
695 >             vterm1 = pref * ri3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
696 >                  qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
697 >                  qzz_j * (3.0_dp*cz2 - 1.0_dp) )
698 >             vterm2 = pref * rcuti3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
699 >                  qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
700 >                  qzz_j * (3.0_dp*cz2 - 1.0_dp) )
701 >             vpair = vpair + swi*( vterm1 - vterm2 )
702 >             epot = epot + ( vterm1 - vterm2 )
703 >            
704 >             dudx = dudx - (5.0_dp * &
705 >                  (vterm1*riji*xhat - vterm2*rcuti2*d(1))) + pref * ( &
706 >                  (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(1)) - &
707 >                  qxx_j*2.0_dp*(xhat - rcuti*d(1))) + &
708 >                  (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(1)) - &
709 >                  qyy_j*2.0_dp*(xhat - rcuti*d(1))) + &
710 >                  (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(1)) - &
711 >                  qzz_j*2.0_dp*(xhat - rcuti*d(1))) )
712 >             dudy = dudy - (5.0_dp * &
713 >                  (vterm1*riji*yhat - vterm2*rcuti2*d(2))) + pref * ( &
714 >                  (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(2)) - &
715 >                  qxx_j*2.0_dp*(yhat - rcuti*d(2))) + &
716 >                  (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(2)) - &
717 >                  qyy_j*2.0_dp*(yhat - rcuti*d(2))) + &
718 >                  (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(2)) - &
719 >                  qzz_j*2.0_dp*(yhat - rcuti*d(2))) )
720 >             dudz = dudz - (5.0_dp * &
721 >                  (vterm1*riji*zhat - vterm2*rcuti2*d(3))) + pref * ( &
722 >                  (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(3)) - &
723 >                  qxx_j*2.0_dp*(zhat - rcuti*d(3))) + &
724 >                  (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(3)) - &
725 >                  qyy_j*2.0_dp*(zhat - rcuti*d(3))) + &
726 >                  (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(3)) - &
727 >                  qzz_j*2.0_dp*(zhat - rcuti*d(3))) )
728 >            
729 >             dudux_j(1) = dudux_j(1) + pref * (ri3*(qxx_j*6.0_dp*cx_j*xhat) - &
730 >                  rcuti4*(qxx_j*6.0_dp*cx_j*d(1)))
731 >             dudux_j(2) = dudux_j(2) + pref * (ri3*(qxx_j*6.0_dp*cx_j*yhat) - &
732 >                  rcuti4*(qxx_j*6.0_dp*cx_j*d(2)))
733 >             dudux_j(3) = dudux_j(3) + pref * (ri3*(qxx_j*6.0_dp*cx_j*zhat) - &
734 >                  rcuti4*(qxx_j*6.0_dp*cx_j*d(3)))
735 >            
736 >             duduy_j(1) = duduy_j(1) + pref * (ri3*(qyy_j*6.0_dp*cy_j*xhat) - &
737 >                  rcuti4*(qyy_j*6.0_dp*cx_j*d(1)))
738 >             duduy_j(2) = duduy_j(2) + pref * (ri3*(qyy_j*6.0_dp*cy_j*yhat) - &
739 >                  rcuti4*(qyy_j*6.0_dp*cx_j*d(2)))
740 >             duduy_j(3) = duduy_j(3) + pref * (ri3*(qyy_j*6.0_dp*cy_j*zhat) - &
741 >                  rcuti4*(qyy_j*6.0_dp*cx_j*d(3)))
742 >            
743 >             duduz_j(1) = duduz_j(1) + pref * (ri3*(qzz_j*6.0_dp*cz_j*xhat) - &
744 >                  rcuti4*(qzz_j*6.0_dp*cx_j*d(1)))
745 >             duduz_j(2) = duduz_j(2) + pref * (ri3*(qzz_j*6.0_dp*cz_j*yhat) - &
746 >                  rcuti4*(qzz_j*6.0_dp*cx_j*d(2)))
747 >             duduz_j(3) = duduz_j(3) + pref * (ri3*(qzz_j*6.0_dp*cz_j*zhat) - &
748 >                  rcuti4*(qzz_j*6.0_dp*cx_j*d(3)))
749 >        
750 >          else
751 >             vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
752 >                  qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
753 >                  qzz_j * (3.0_dp*cz2 - 1.0_dp))
754 >             vpair = vpair + swi * vterm
755 >             epot = epot + vterm
756 >            
757 >             dudx = dudx - 5.0_dp*vterm*riji*xhat + pref * ri4 * ( &
758 >                  qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
759 >                  qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
760 >                  qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
761 >             dudy = dudy - 5.0_dp*vterm*riji*yhat + pref * ri4 * ( &
762 >                  qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
763 >                  qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
764 >                  qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
765 >             dudz = dudz - 5.0_dp*vterm*riji*zhat + pref * ri4 * ( &
766 >                  qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
767 >                  qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
768 >                  qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
769 >            
770 >             dudux_j(1) = dudux_j(1) + pref * ri3*(qxx_j*6.0_dp*cx_j*xhat)
771 >             dudux_j(2) = dudux_j(2) + pref * ri3*(qxx_j*6.0_dp*cx_j*yhat)
772 >             dudux_j(3) = dudux_j(3) + pref * ri3*(qxx_j*6.0_dp*cx_j*zhat)
773 >            
774 >             duduy_j(1) = duduy_j(1) + pref * ri3*(qyy_j*6.0_dp*cy_j*xhat)
775 >             duduy_j(2) = duduy_j(2) + pref * ri3*(qyy_j*6.0_dp*cy_j*yhat)
776 >             duduy_j(3) = duduy_j(3) + pref * ri3*(qyy_j*6.0_dp*cy_j*zhat)
777 >            
778 >             duduz_j(1) = duduz_j(1) + pref * ri3*(qzz_j*6.0_dp*cz_j*xhat)
779 >             duduz_j(2) = duduz_j(2) + pref * ri3*(qzz_j*6.0_dp*cz_j*yhat)
780 >             duduz_j(3) = duduz_j(3) + pref * ri3*(qzz_j*6.0_dp*cz_j*zhat)
781 >          
782 >          endif
783         endif
610
784      endif
785  
786      if (i_is_Dipole) then
787  
788         if (j_is_Charge) then
789  
790 <          if (i_is_SplitDipole) then
618 <             BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
619 <             ri = 1.0_dp / BigR
620 <             scale = rij * ri
621 <          else
622 <             ri = riji
623 <             scale = 1.0_dp
624 <          endif
790 >          pref = sw * pre12 * q_j * mu_i
791  
792 <          ri2 = ri * ri
793 <          ri3 = ri2 * ri
794 <          sc2 = scale * scale
792 >          if (corrMethod .eq. 1) then
793 >             ri2 = riji * riji
794 >             ri3 = ri2 * riji
795  
796 <          pref = pre12 * q_j * mu_i
797 <          vterm = pref * ct_i * ri2 * scale
798 <          vpair = vpair + vterm
799 <          epot = epot + sw * vterm
796 >             vterm = pref * ct_i * (ri2 - rcuti2)
797 >             vpair = vpair + swi * vterm
798 >             epot = epot + vterm
799 >            
800 >             !! this has a + sign in the () because the rij vector is
801 >             !! r_j - r_i and the charge-dipole potential takes the origin
802 >             !! as the point dipole, which is atom j in this case.
803 >            
804 >             dudx = dudx + pref * ( ri3*( uz_i(1) - 3.0d0*ct_i*xhat) &
805 >                  - rcuti3*( uz_i(1) - 3.0d0*ct_i*d(1)*rcuti ) )
806 >             dudy = dudy + pref * ( ri3*( uz_i(2) - 3.0d0*ct_i*yhat) &
807 >                  - rcuti3*( uz_i(2) - 3.0d0*ct_i*d(2)*rcuti ) )
808 >             dudz = dudz + pref * ( ri3*( uz_i(3) - 3.0d0*ct_i*zhat) &
809 >                  - rcuti3*( uz_i(3) - 3.0d0*ct_i*d(3)*rcuti ) )
810 >            
811 >             duduz_i(1) = duduz_i(1) - pref*( ri2*xhat - d(1)*rcuti3 )
812 >             duduz_i(2) = duduz_i(2) - pref*( ri2*yhat - d(2)*rcuti3 )
813 >             duduz_i(3) = duduz_i(3) - pref*( ri2*zhat - d(3)*rcuti3 )
814  
635          dudx = dudx + pref * sw * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
636          dudy = dudy + pref * sw * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
637          dudz = dudz + pref * sw * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
638
639          duduz_i(1) = duduz_i(1) + pref * sw * ri2 * xhat * scale
640          duduz_i(2) = duduz_i(2) + pref * sw * ri2 * yhat * scale
641          duduz_i(3) = duduz_i(3) + pref * sw * ri2 * zhat * scale
642       endif
643
644       if (j_is_Dipole) then
645
646          if (i_is_SplitDipole) then
647             if (j_is_SplitDipole) then
648                BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
649             else
650                BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
651             endif
652             ri = 1.0_dp / BigR
653             scale = rij * ri                
815            else
816 <             if (j_is_SplitDipole) then
817 <                BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
816 >             if (i_is_SplitDipole) then
817 >                BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
818                  ri = 1.0_dp / BigR
819 <                scale = rij * ri                            
820 <             else                
819 >                scale = rij * ri
820 >             else
821                  ri = riji
822                  scale = 1.0_dp
823               endif
824 +            
825 +             ri2 = ri * ri
826 +             ri3 = ri2 * ri
827 +             sc2 = scale * scale
828 +            
829 +             vterm = pref * ct_i * ri2 * scale
830 +             vpair = vpair + swi * vterm
831 +             epot = epot + vterm
832 +            
833 +             dudx = dudx + pref * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
834 +             dudy = dudy + pref * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
835 +             dudz = dudz + pref * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
836 +            
837 +             duduz_i(1) = duduz_i(1) + pref * ri2 * xhat * scale
838 +             duduz_i(2) = duduz_i(2) + pref * ri2 * yhat * scale
839 +             duduz_i(3) = duduz_i(3) + pref * ri2 * zhat * scale
840            endif
841 +       endif
842  
843 <          ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
843 >       if (j_is_Dipole) then
844  
845 <          ri2 = ri * ri
668 <          ri3 = ri2 * ri
669 <          ri4 = ri2 * ri2
670 <          sc2 = scale * scale
845 >          pref = sw * pre22 * mu_i * mu_j
846  
847 <          pref = pre22 * mu_i * mu_j
848 <          vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
849 <          vpair = vpair + vterm
850 <          epot = epot + sw * vterm
847 >          if (corrMethod .eq. 1) then
848 >             ri2 = riji * riji
849 >             ri3 = ri2 * riji
850 >             ri4 = ri2 * ri2
851  
852 <          a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
853 <
854 <          dudx = dudx + pref*sw*3.0d0*ri4*scale &
855 <                         *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
856 <          dudy = dudy + pref*sw*3.0d0*ri4*scale &
857 <                         *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
858 <          dudz = dudz + pref*sw*3.0d0*ri4*scale &
859 <                         *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
860 <
861 <          duduz_i(1) = duduz_i(1) + pref*sw*ri3 &
862 <                                     *(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
863 <          duduz_i(2) = duduz_i(2) + pref*sw*ri3 &
864 <                                     *(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
865 <          duduz_i(3) = duduz_i(3) + pref*sw*ri3 &
866 <                                     *(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
867 <
868 <          duduz_j(1) = duduz_j(1) + pref*sw*ri3 &
869 <                                     *(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
870 <          duduz_j(2) = duduz_j(2) + pref*sw*ri3 &
871 <                                     *(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
872 <          duduz_j(3) = duduz_j(3) + pref*sw*ri3 &
873 <                                     *(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
852 >             vterm = pref * (ri3 - rcuti3) * (ct_ij - 3.0d0 * ct_i * ct_j)
853 >             vpair = vpair + swi * vterm
854 >             epot = epot + vterm
855 >            
856 >             a1 = 5.0d0 * ct_i * ct_j - ct_ij
857 >            
858 >             dudx = dudx + pref*3.0d0*ri4 &
859 >                  *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) - &
860 >                  pref*3.0d0*rcuti4*(a1*rcuti*d(1)-ct_i*uz_j(1)-ct_j*uz_i(1))
861 >             dudy = dudy + pref*3.0d0*ri4 &
862 >                  *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) - &
863 >                  pref*3.0d0*rcuti4*(a1*rcuti*d(2)-ct_i*uz_j(2)-ct_j*uz_i(2))
864 >             dudz = dudz + pref*3.0d0*ri4 &
865 >                  *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) - &
866 >                  pref*3.0d0*rcuti4*(a1*rcuti*d(3)-ct_i*uz_j(3)-ct_j*uz_i(3))
867 >            
868 >             duduz_i(1) = duduz_i(1) + pref*(ri3*(uz_j(1) - 3.0d0*ct_j*xhat) &
869 >                  - rcuti3*(uz_j(1) - 3.0d0*ct_j*d(1)*rcuti))
870 >             duduz_i(2) = duduz_i(2) + pref*(ri3*(uz_j(2) - 3.0d0*ct_j*yhat) &
871 >                  - rcuti3*(uz_j(2) - 3.0d0*ct_j*d(2)*rcuti))
872 >             duduz_i(3) = duduz_i(3) + pref*(ri3*(uz_j(3) - 3.0d0*ct_j*zhat) &
873 >                  - rcuti3*(uz_j(3) - 3.0d0*ct_j*d(3)*rcuti))
874 >             duduz_j(1) = duduz_j(1) + pref*(ri3*(uz_i(1) - 3.0d0*ct_i*xhat) &
875 >                  - rcuti3*(uz_i(1) - 3.0d0*ct_i*d(1)*rcuti))
876 >             duduz_j(2) = duduz_j(2) + pref*(ri3*(uz_i(2) - 3.0d0*ct_i*yhat) &
877 >                  - rcuti3*(uz_i(2) - 3.0d0*ct_i*d(2)*rcuti))
878 >             duduz_j(3) = duduz_j(3) + pref*(ri3*(uz_i(3) - 3.0d0*ct_i*zhat) &
879 >                  - rcuti3*(uz_i(3) - 3.0d0*ct_i*d(3)*rcuti))
880 >          else
881 >            
882 >             if (i_is_SplitDipole) then
883 >                if (j_is_SplitDipole) then
884 >                   BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
885 >                else
886 >                   BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
887 >                endif
888 >                ri = 1.0_dp / BigR
889 >                scale = rij * ri                
890 >             else
891 >                if (j_is_SplitDipole) then
892 >                   BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
893 >                   ri = 1.0_dp / BigR
894 >                   scale = rij * ri                            
895 >                else                
896 >                   ri = riji
897 >                   scale = 1.0_dp
898 >                endif
899 >             endif
900 >            
901 >             ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
902 >            
903 >             ri2 = ri * ri
904 >             ri3 = ri2 * ri
905 >             ri4 = ri2 * ri2
906 >             sc2 = scale * scale
907 >            
908 >             vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
909 >             vpair = vpair + swi * vterm
910 >             epot = epot + vterm
911 >            
912 >             a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
913 >            
914 >             dudx = dudx + pref*3.0d0*ri4*scale &
915 >                  *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
916 >             dudy = dudy + pref*3.0d0*ri4*scale &
917 >                  *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
918 >             dudz = dudz + pref*3.0d0*ri4*scale &
919 >                  *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
920 >            
921 >             duduz_i(1) = duduz_i(1) + pref*ri3 &
922 >                  *(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
923 >             duduz_i(2) = duduz_i(2) + pref*ri3 &
924 >                  *(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
925 >             duduz_i(3) = duduz_i(3) + pref*ri3 &
926 >                  *(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
927 >            
928 >             duduz_j(1) = duduz_j(1) + pref*ri3 &
929 >                  *(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
930 >             duduz_j(2) = duduz_j(2) + pref*ri3 &
931 >                  *(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
932 >             duduz_j(3) = duduz_j(3) + pref*ri3 &
933 >                  *(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
934 >          endif
935         endif
700
936      endif
937  
938      if (i_is_Quadrupole) then
# Line 710 | Line 945 | contains
945            cy2 = cy_i * cy_i
946            cz2 = cz_i * cz_i
947  
948 <          pref = pre14 * q_j / 3.0_dp
714 <          vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
715 <               qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
716 <               qzz_i * (3.0_dp*cz2 - 1.0_dp))
717 <          vpair = vpair + vterm
718 <          epot = epot + sw * vterm
948 >          pref = sw * pre14 * q_j / 3.0_dp
949  
950 <          dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
951 <               qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
952 <               qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
953 <               qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
954 <          dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
955 <               qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
956 <               qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
957 <               qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
958 <          dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
959 <               qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
960 <               qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
961 <               qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
950 >          if (corrMethod .eq. 1) then
951 >             vterm1 = pref * ri3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
952 >                  qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
953 >                  qzz_i * (3.0_dp*cz2 - 1.0_dp) )
954 >             vterm2 = pref * rcuti3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
955 >                  qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
956 >                  qzz_i * (3.0_dp*cz2 - 1.0_dp) )
957 >             vpair = vpair + swi * ( vterm1 - vterm2 )
958 >             epot = epot + ( vterm1 - vterm2 )
959 >            
960 >             dudx = dudx - (5.0_dp*(vterm1*riji*xhat - vterm2*rcuti2*d(1))) + &
961 >                  pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(1)) - &
962 >                  qxx_i*2.0_dp*(xhat - rcuti*d(1))) + &
963 >                  (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(1)) - &
964 >                  qyy_i*2.0_dp*(xhat - rcuti*d(1))) + &
965 >                  (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(1)) - &
966 >                  qzz_i*2.0_dp*(xhat - rcuti*d(1))) )
967 >             dudy = dudy - (5.0_dp*(vterm1*riji*yhat - vterm2*rcuti2*d(2))) + &
968 >                  pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(2)) - &
969 >                  qxx_i*2.0_dp*(yhat - rcuti*d(2))) + &
970 >                  (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(2)) - &
971 >                  qyy_i*2.0_dp*(yhat - rcuti*d(2))) + &
972 >                  (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(2)) - &
973 >                  qzz_i*2.0_dp*(yhat - rcuti*d(2))) )
974 >             dudz = dudz - (5.0_dp*(vterm1*riji*zhat - vterm2*rcuti2*d(3))) + &
975 >                  pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(3)) - &
976 >                  qxx_i*2.0_dp*(zhat - rcuti*d(3))) + &
977 >                  (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(3)) - &
978 >                  qyy_i*2.0_dp*(zhat - rcuti*d(3))) + &
979 >                  (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(3)) - &
980 >                  qzz_i*2.0_dp*(zhat - rcuti*d(3))) )
981 >            
982 >             dudux_i(1) = dudux_i(1) + pref * (ri3*(qxx_i*6.0_dp*cx_i*xhat) - &
983 >                  rcuti4*(qxx_i*6.0_dp*cx_i*d(1)))
984 >             dudux_i(2) = dudux_i(2) + pref * (ri3*(qxx_i*6.0_dp*cx_i*yhat) - &
985 >                  rcuti4*(qxx_i*6.0_dp*cx_i*d(2)))
986 >             dudux_i(3) = dudux_i(3) + pref * (ri3*(qxx_i*6.0_dp*cx_i*zhat) - &
987 >                  rcuti4*(qxx_i*6.0_dp*cx_i*d(3)))
988 >            
989 >             duduy_i(1) = duduy_i(1) + pref * (ri3*(qyy_i*6.0_dp*cy_i*xhat) - &
990 >                  rcuti4*(qyy_i*6.0_dp*cx_i*d(1)))
991 >             duduy_i(2) = duduy_i(2) + pref * (ri3*(qyy_i*6.0_dp*cy_i*yhat) - &
992 >                  rcuti4*(qyy_i*6.0_dp*cx_i*d(2)))
993 >             duduy_i(3) = duduy_i(3) + pref * (ri3*(qyy_i*6.0_dp*cy_i*zhat) - &
994 >                  rcuti4*(qyy_i*6.0_dp*cx_i*d(3)))
995 >            
996 >             duduz_i(1) = duduz_i(1) + pref * (ri3*(qzz_i*6.0_dp*cz_i*xhat) - &
997 >                  rcuti4*(qzz_i*6.0_dp*cx_i*d(1)))
998 >             duduz_i(2) = duduz_i(2) + pref * (ri3*(qzz_i*6.0_dp*cz_i*yhat) - &
999 >                  rcuti4*(qzz_i*6.0_dp*cx_i*d(2)))
1000 >             duduz_i(3) = duduz_i(3) + pref * (ri3*(qzz_i*6.0_dp*cz_i*zhat) - &
1001 >                  rcuti4*(qzz_i*6.0_dp*cx_i*d(3)))
1002  
1003 <          dudux_i(1) = dudux_i(1) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*xhat)
1004 <          dudux_i(2) = dudux_i(2) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*yhat)
1005 <          dudux_i(3) = dudux_i(3) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*zhat)
1006 <
1007 <          duduy_i(1) = duduy_i(1) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*xhat)
1008 <          duduy_i(2) = duduy_i(2) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*yhat)
1009 <          duduy_i(3) = duduy_i(3) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*zhat)
1010 <
1011 <          duduz_i(1) = duduz_i(1) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*xhat)
1012 <          duduz_i(2) = duduz_i(2) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*yhat)
1013 <          duduz_i(3) = duduz_i(3) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*zhat)
1003 >          else
1004 >             vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
1005 >                  qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
1006 >                  qzz_i * (3.0_dp*cz2 - 1.0_dp))
1007 >             vpair = vpair + swi * vterm
1008 >             epot = epot + vterm
1009 >            
1010 >             dudx = dudx - 5.0_dp*vterm*riji*xhat + pref * ri4 * ( &
1011 >                  qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
1012 >                  qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
1013 >                  qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
1014 >             dudy = dudy - 5.0_dp*vterm*riji*yhat + pref * ri4 * ( &
1015 >                  qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
1016 >                  qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
1017 >                  qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
1018 >             dudz = dudz - 5.0_dp*vterm*riji*zhat + pref * ri4 * ( &
1019 >                  qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
1020 >                  qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
1021 >                  qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
1022 >            
1023 >             dudux_i(1) = dudux_i(1) + pref * ri3*(qxx_i*6.0_dp*cx_i*xhat)
1024 >             dudux_i(2) = dudux_i(2) + pref * ri3*(qxx_i*6.0_dp*cx_i*yhat)
1025 >             dudux_i(3) = dudux_i(3) + pref * ri3*(qxx_i*6.0_dp*cx_i*zhat)
1026 >            
1027 >             duduy_i(1) = duduy_i(1) + pref * ri3*(qyy_i*6.0_dp*cy_i*xhat)
1028 >             duduy_i(2) = duduy_i(2) + pref * ri3*(qyy_i*6.0_dp*cy_i*yhat)
1029 >             duduy_i(3) = duduy_i(3) + pref * ri3*(qyy_i*6.0_dp*cy_i*zhat)
1030 >            
1031 >             duduz_i(1) = duduz_i(1) + pref * ri3*(qzz_i*6.0_dp*cz_i*xhat)
1032 >             duduz_i(2) = duduz_i(2) + pref * ri3*(qzz_i*6.0_dp*cz_i*yhat)
1033 >             duduz_i(3) = duduz_i(3) + pref * ri3*(qzz_i*6.0_dp*cz_i*zhat)
1034 >          endif
1035         endif
1036      endif
1037  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines