~ubuntu-branches/ubuntu/saucy/nwchem/saucy

« back to all changes in this revision

Viewing changes to src/nwdft/zora/gridNuclearPotential.F

  • Committer: Package Import Robot
  • Author(s): Michael Banck, Michael Banck, Daniel Leidert
  • Date: 2012-02-09 20:02:41 UTC
  • mfrom: (1.1.1)
  • Revision ID: package-import@ubuntu.com-20120209200241-jgk03qfsphal4ug2
Tags: 6.1-1
* New upstream release.

[ Michael Banck ]
* debian/patches/02_makefile_flags.patch: Updated.
* debian/patches/02_makefile_flags.patch: Use internal blas and lapack code.
* debian/patches/02_makefile_flags.patch: Define GCC4 for LINUX and LINUX64
  (Closes: #632611 and LP: #791308).
* debian/control (Build-Depends): Added openssh-client.
* debian/rules (USE_SCALAPACK, SCALAPACK): Removed variables (Closes:
  #654658).
* debian/rules (LIBDIR, USE_MPIF4, ARMCI_NETWORK): New variables.
* debian/TODO: New file.
* debian/control (Build-Depends): Removed libblas-dev, liblapack-dev and
  libscalapack-mpi-dev.
* debian/patches/04_show_testsuite_diff_output.patch: New patch, shows the
  diff output for failed tests.
* debian/patches/series: Adjusted.
* debian/testsuite: Optionally run all tests if "all" is passed as option.
* debian/rules: Run debian/testsuite with "all" if DEB_BUILD_OPTIONS
  contains "checkall".

[ Daniel Leidert ]
* debian/control: Used wrap-and-sort. Added Vcs-Svn and Vcs-Browser fields.
  (Priority): Moved to extra according to policy section 2.5.
  (Standards-Version): Bumped to 3.9.2.
  (Description): Fixed a typo.
* debian/watch: Added.
* debian/patches/03_hurd-i386_define_path_max.patch: Added.
  - Define MAX_PATH if not defines to fix FTBFS on hurd.
* debian/patches/series: Adjusted.

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
c
 
2
c     == Calculate the point nuclear potential on the grid ==
 
3
      subroutine gridNuclearPotentialPoint(geom,
 
4
     &                                natoms,      ! IN : number of atoms
 
5
     &                                nqpts,       ! IN : number of grid points
 
6
     &                                qxyz,        ! IN : grid points
 
7
     &                                qwght,       ! IN : quadrature weightings
 
8
     &                                closegridpts,
 
9
     &                                amat_nucl)   ! OUT : nuclear potential
 
10
      implicit none
 
11
#include "geom.fh"
 
12
#include "stdio.fh"
 
13
#include "zora.fh"
 
14
#include "global.fh" ! Added for using ga_nodeid()
 
15
 
 
16
      integer iatom,natoms,geom,igrid,nqpts
 
17
      double precision nucCharge, nucCoords(3)
 
18
      double precision qxyz(3,nqpts),qwght(nqpts)
 
19
      character*16 tags(natoms)
 
20
      logical lSuccess
 
21
      double precision rx,ry,rz,dist
 
22
      double precision amat_nucl(nqpts)
 
23
      integer closegridpts(*)
 
24
      logical ofinite
 
25
c     == get the total nuclear potential on the grid ==
 
26
      do igrid = 1,nqpts   ! == loop over the grid points ==
 
27
       amat_nucl(igrid) = 0.d0  
 
28
       do iatom = 1,natoms ! == loop over the atoms ==
 
29
c
 
30
c     == get an atom (needs error handling) ==
 
31
        lSuccess = geom_cent_get(geom, iatom, tags, nucCoords, 
 
32
     &             nucCharge)
 
33
c
 
34
c     == distance from the grid points to the atom centers == 
 
35
        rx = nucCoords(1) - qxyz(1,igrid)
 
36
        ry = nucCoords(2) - qxyz(2,igrid)
 
37
        rz = nucCoords(3) - qxyz(3,igrid)
 
38
        dist = dsqrt(rx*rx + ry*ry + rz*rz)
 
39
        if (dist.gt.zoracutoff) then ! == check cutoff ==
 
40
           amat_nucl(igrid) = amat_nucl(igrid) - 
 
41
     &                        nucCharge/dist
 
42
        else
 
43
          closegridpts(igrid) = igrid
 
44
        end if
 
45
       end do
 
46
      end do
 
47
      return
 
48
      end 
 
49
c
 
50
c     == Finite nucleus 1 ==
 
51
      subroutine gridNuclearPotentialFinite(
 
52
     &               geom,
 
53
     &               natoms,      ! IN : number of atoms
 
54
     &               nqpts,       ! IN : number of grid points
 
55
     &               qxyz,        ! IN : grid points
 
56
     &               qwght,       ! IN : quadrature weightings
 
57
     &               zetanuc_arr, ! IN : sqrt(zetanuc) for Gaussian Nuclear Model
 
58
     &               closegridpts,
 
59
     &               amat_nucl)   ! OUT : nuclear potential
 
60
c
 
61
c Evaluates: V_L^{Gauss}(r) = \sum_L -Z_L/r_L erf(\tilde{r}_L)
 
62
c            \tilde{r}_L = srqt(zetanuc) r_L
 
63
c            r_L = |r_grid-r_N|
 
64
c            r_N   , nuclear position
 
65
c            r_grid, grid point
 
66
c            according to: Autschbach, J. ChemPhysChem 2009,V10,P 22774-2283
 
67
c            page 2276
 
68
c Author: Fredy W. Aquino 
 
69
c Date  : 04-20-11
 
70
c
 
71
      implicit none
 
72
#include "geom.fh"
 
73
#include "stdio.fh"
 
74
#include "zora.fh"
 
75
#include "global.fh" ! Added for using ga_nodeid()
 
76
 
 
77
      integer iatom,natoms,geom,igrid,nqpts
 
78
      double precision nucCharge, nucCoords(3)
 
79
      double precision qxyz(3,nqpts),qwght(nqpts),
 
80
     &                 rtemp
 
81
      character*16 tags(natoms)
 
82
      logical lSuccess
 
83
      double precision rx,ry,rz,dist
 
84
      double precision amat_nucl(nqpts),
 
85
     &                zetanuc_arr(natoms),erf
 
86
      integer closegridpts(*)
 
87
      logical ofinite
 
88
      external merf,get_zetanuc_arr
 
89
      do igrid = 1,nqpts   ! == loop over the grid points ==
 
90
       amat_nucl(igrid) = 0.d0  
 
91
       do iatom = 1,natoms ! == loop over the atoms ==
 
92
c     == get an atom (needs error handling) ==
 
93
        lSuccess = geom_cent_get(geom, iatom, tags, nucCoords, 
 
94
     &             nucCharge)
 
95
c     == distance from the grid points to the atom centers == 
 
96
        rx = nucCoords(1) - qxyz(1,igrid)
 
97
        ry = nucCoords(2) - qxyz(2,igrid)
 
98
        rz = nucCoords(3) - qxyz(3,igrid)
 
99
        dist = dsqrt(rx*rx + ry*ry + rz*rz)
 
100
        rtemp = zetanuc_arr(iatom)*dist
 
101
        call merf(rtemp,erf) ! Evaluate erf(sqrt(zetanuc)r_L)
 
102
        if (dist.gt.zoracutoff) then ! == check cutoff ==          
 
103
           amat_nucl(igrid) = amat_nucl(igrid) - 
 
104
     &                        nucCharge/dist*erf
 
105
        else
 
106
          closegridpts(igrid) = igrid
 
107
        end if
 
108
       end do ! end-loop-iatom
 
109
      end do ! end-loop-igrid
 
110
      return
 
111
      end
 
112
c
 
113
c     == Finite nucleus 2 ==
 
114
      subroutine gridNuclearPotentialFinite2(
 
115
     &               geom,
 
116
     &               natoms,      ! IN : number of atoms
 
117
     &               nqpts,       ! IN : number of grid points
 
118
     &               qxyz,        ! IN : grid points
 
119
     &               qwght,       ! IN : quadrature weightings
 
120
     &               closegridpts,
 
121
     &               amat_nucl)   ! OUT : nuclear potential
 
122
c Evaluates: V_L^{Gauss}(r) = - Z_L/r_L P(1/2,\tilde{r}_L^2)
 
123
c            \tilde{r}_L = srqt(zetanuc) r_L
 
124
c            r_L = |r_grid-r_N|
 
125
c            r_N   , nuclear position
 
126
c            r_grid, grid point
 
127
c            according to: Autschbach, J. ChemPhysChem 2009,V10,P 22774-2283
 
128
c            Eq. (5). 
 
129
c Author: Fredy W. Aquino 
 
130
c Date  : 04-20-11
 
131
      implicit none
 
132
#include "geom.fh"
 
133
#include "stdio.fh"
 
134
#include "zora.fh"
 
135
      integer iatom,natoms,geom,igrid,nqpts
 
136
      double precision nucCharge, nucCoords(3)
 
137
      double precision qxyz(3,nqpts),qwght(nqpts)
 
138
      character*16 tags(natoms)
 
139
      logical lSuccess
 
140
      double precision rx,ry,rz,dist
 
141
      double precision amat_nucl(nqpts),
 
142
     &                zetanuc_arr(natoms)
 
143
      double precision a_coeff,rtemp,gammap,gammaq
 
144
      integer closegridpts(*)
 
145
c
 
146
      external gratio,  ! Evaluates Incomplete Gamma Function P(a,x)
 
147
     &         get_zetanuc_arr
 
148
c
 
149
      call get_zetanuc_arr(geom,natoms,zetanuc_arr) ! zetanuc_arr(i) i=1,natoms
 
150
c
 
151
      a_coeff=0.5
 
152
      do igrid = 1,nqpts   ! == loop over the grid points ==
 
153
       amat_nucl(igrid) = 0.d0  
 
154
       do iatom = 1,natoms ! == loop over the atoms ==
 
155
c     == get an atom (needs error handling) ==
 
156
        lSuccess = geom_cent_get(geom, iatom, tags, nucCoords, 
 
157
     &             nucCharge)
 
158
c     == distance from the grid points to the atom centers == 
 
159
        rx = nucCoords(1) - qxyz(1,igrid)
 
160
        ry = nucCoords(2) - qxyz(2,igrid)
 
161
        rz = nucCoords(3) - qxyz(3,igrid)
 
162
        dist = dsqrt(rx*rx + ry*ry + rz*rz)
 
163
        rtemp = real(zetanuc_arr(iatom)*dist*dist)
 
164
        call gratio(a_coeff,rtemp,gammap,gammaq,0)! Evaluate P(1/2,zetanuc r_L^2)
 
165
        if (dist.gt.zoracutoff) then ! == check cutoff ==          
 
166
           amat_nucl(igrid) = amat_nucl(igrid) - 
 
167
     &                        nucCharge/dist*gammap
 
168
        else
 
169
          closegridpts(igrid) = igrid
 
170
        end if
 
171
       end do ! end-loop-iatom
 
172
      end do ! end-loop-igrid
 
173
      return
 
174
      end
 
175
c $Id: gridNuclearPotential.F 21488 2011-11-09 22:17:38Z niri $