forked from lgcrego/Dynemol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
film_STO.f
433 lines (326 loc) · 13.7 KB
/
film_STO.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
module Psi_squared_cube_format
use type_m
use constants_m
use parameters_m , only : electron_state
use Babel_m , only : System_Characteristics
use Semi_Empirical_Parms , only : atom
use Structure_Builder , only : Extended_Cell
use Slater_Type_Orbitals , only : s_orb , p_orb , d_x2y2 , d_z2 , d_xyz
public :: Gaussian_Cube_Format
private
real*8 , parameter :: fringe = 8.d0
interface Gaussian_Cube_Format
module procedure Gaussian_Cube_Format_Real
module procedure Gaussian_Cube_Format_Cmplx
end interface Gaussian_Cube_Format
contains
!
!
!
!===================================================================
subroutine gaussian_cube_format_Real( bra , ket , it , t , el_hl )
!===================================================================
implicit none
real*8 , intent(in) :: bra(:), ket(:)
real*8 , intent(in) :: t
integer , intent(in) :: it
character(*) , optional , intent(in) :: el_hl
! local variables ...
integer :: many_steps = 131
integer :: medium_steps = 89
integer :: few_steps = 59
integer :: n_xyz_steps(3)
complex*16 :: TotalPsiKet
real*8 , allocatable :: xyz(:,:) , Psi(:,:,:)
real*8 :: x , y , z , x0 , y0 , z0 , dx , dy , dz , a , b , c , r , SlaterOrbital, dumb
integer :: AtNo , i , j , ix , iy , iz , k , l
character(len=2) :: prefix
character(len=5) :: string
character(len=31) :: f_name
! bra is never used, so to avoid compiler warnings ...
dumb = bra(1)
allocate(xyz(extended_cell%atoms,3))
write(string,'(i5.5)') it
prefix = merge( "el" , el_hl , .NOT. present(el_hl) )
f_name = 'MO.trunk/'//prefix//'_MO_shot'//string//'.cube'
OPEN(unit=4,file=f_name,status='unknown')
! bounding box for isosurfaces ...
CALL BoundingBox( extended_cell )
! fix number of steps for each direction according to aspect ratio ...
do i = 1 , 3
! default case
n_xyz_steps(i) = medium_steps
If( extended_cell%BoxSides(i) >= 2.5d1 ) n_xyz_steps(i) = many_steps
If( extended_cell%BoxSides(i) <= 1.2d1 ) n_xyz_steps(i) = few_steps
end do
! voxel dimensions
dx = extended_cell%BoxSides(1) / n_xyz_steps(1)
dy = extended_cell%BoxSides(2) / n_xyz_steps(2)
dz = extended_cell%BoxSides(3) / n_xyz_steps(3)
! translation to the center of mass
forall(i=1:extended_cell%atoms,j=1:3) xyz(i,j) = extended_cell%coord(i,j) - extended_cell%Center_of_Mass(j)
! initial corner of the volume Box
a = minval(xyz(:,1)) - fringe / two
b = minval(xyz(:,2)) - fringe / two
c = minval(xyz(:,3)) - fringe / two
! start writing Gaussian cube files
write(4,*) System_Characteristics
write(4,*) 'electron_state = ',electron_state,' / time = ', t
write(4,111) extended_cell%atoms , a/a_Bohr , b/a_Bohr , c/a_Bohr
write(4,111) n_xyz_steps(1) + 1 , dx/a_Bohr , 0.d0 , 0.d0
write(4,111) n_xyz_steps(2) + 1 , 0.d0 , dy/a_Bohr , 0.d0
write(4,111) n_xyz_steps(3) + 1 , 0.d0 , 0.d0 , dz/a_Bohr
! coordinates in a.u., because zeta is in units of [a_0^{-1}] ...
xyz = xyz / a_Bohr
DO i = 1 , extended_cell%atoms
write(4,113) extended_cell%AtNo(i), 0.0 , xyz(i,1) , xyz(i,2) , xyz(i,3)
END DO
!---------------------------------------------------------
! drawing the wavefunction denssity in cube format
!========================================================
! the order orbitals are stored
!
! S --> 1 --> l = 0 , m = 0
! Py --> 2 --> l = 1 , m = -1
! Pz --> 3 --> l = 1 , m = 0
! Px --> 4 --> l = 1 , m = +1
! Dxy --> 5 --> l = 2 , m = -2
! Dyz --> 6 --> l = 2 , m = -1
! Dz2 --> 7 --> l = 2 , m = 0
! Dxz --> 8 --> l = 2 , m = +1
! Dx2y2 --> 9 --> l = 2 , m = +2
!========================================================
allocate( Psi(n_xyz_steps(1)+1,n_xyz_steps(2)+1,n_xyz_steps(3)+1) , source = D_zero )
!$OMP parallel private(ix,iy,iz,x,y,z,x0,y0,z0,SlaterOrbital,r,AtNo,i,j,k,l,TotalPsiKet)
!$OMP single
DO ix = 0 , n_xyz_steps(1)
x = (a + ix * dx) / a_Bohr
!$OMP task untied
DO iy = 0 , n_xyz_steps(2)
y = (b + iy * dy) / a_Bohr
DO iz = 0 , n_xyz_steps(3)
z = (c + iz * dz) / a_Bohr
i = 0
TotalPsiKet = C_zero
DO k = 1 , extended_cell%atoms
if( Extended_Cell% QMMM(k) /= "QM" ) cycle
AtNo = extended_cell%AtNo(k)
! distance from the center of the nuclei, in a.u.
r = dsqrt((x-xyz(k,1))*(x-xyz(k,1)) + (y-xyz(k,2))*(y-xyz(k,2)) + (z-xyz(k,3))*(z-xyz(k,3)))
! coordinates centered on the nuclei, in a.u.
x0 = x - xyz(k,1)
y0 = y - xyz(k,2)
z0 = z - xyz(k,3)
do j = 1 , atom(AtNo)%DOS
i = i + 1
l = extended_cell%BasisPointer(k) + j
select case (j)
case( 1 )
SlaterOrbital = s_orb(r,AtNo,l)
case( 2 )
SlaterOrbital = p_orb(r,y0,AtNo,l)
case( 3 )
SlaterOrbital = p_orb(r,z0,AtNo,l)
case( 4 )
SlaterOrbital = p_orb(r,x0,AtNo,l)
case( 5 )
SlaterOrbital = d_xyz(r,x0,y0,AtNo,l)
case( 6 )
SlaterOrbital = d_xyz(r,y0,z0,AtNo,l)
case( 7 )
SlaterOrbital = d_z2(r,x0,y0,z0,AtNo,l)
case( 8 )
SlaterOrbital = d_xyz(r,x0,z0,AtNo,l)
case( 9 )
SlaterOrbital = d_x2y2(r,x0,y0,AtNo,l)
end select
TotalPsiKet = TotalPsiKet + ket(i) * SlaterOrbital
end do ! <== DOS
end do ! <== atoms
Psi( ix+1 , iy+1 , iz+1 ) = TotalPsiKet
END DO ! <== Z coord
END DO ! <== Y coord
!$OMP end task
END DO ! <== X coord
!$OMP end single
!$OMP end parallel
!::::::::::::::::::::::::::::::::::::::::::::::::::::::::
DO ix = 0 , n_xyz_steps(1)
DO iy = 0 , n_xyz_steps(2)
DO iz = 0 , n_xyz_steps(3)
write(4,112,advance='no') Psi( ix+1 , iy+1 , iz+1 )
If( (mod(iz,6) == 5) ) write(4,'(a)',advance='yes')
END DO
write(4,'(a)',advance='yes')
END DO
END DO
!::::::::::::::::::::::::::::::::::::::::::::::::::::::::
close(4)
deallocate(xyz , Psi)
include 'formats.h'
end subroutine Gaussian_Cube_Format_Real
!
!
!
!====================================================================
subroutine gaussian_cube_format_Cmplx( bra , ket , it , t , el_hl )
!====================================================================
implicit none
complex*16 , intent(in) :: bra(:), ket(:)
real*8 , intent(in) :: t
integer , intent(in) :: it
character(*) , optional , intent(in) :: el_hl
! local variables ...
integer :: many_steps = 131
integer :: medium_steps = 89
integer :: few_steps = 59
integer :: n_xyz_steps(3)
complex*16 :: TotalPsiBra , TotalPsiKet
real*8 , allocatable :: xyz(:,:) , Psi_2(:,:,:)
real*8 :: x , y , z , x0 , y0 , z0 , dx , dy , dz , a , b , c , r , SlaterOrbital
integer :: AtNo , i , j , ix , iy , iz , k , l
character(len=2) :: prefix
character(len=5) :: string
character(len=22) :: f_name
allocate(xyz(extended_cell%atoms,3))
write(string,'(i5.5)') it
prefix = merge( "el" , el_hl , .NOT. present(el_hl) )
f_name = prefix//'_dens_shot'//string//'.cube'
OPEN(unit=4,file=f_name,status='unknown')
! bounding box for isosurfaces ...
CALL BoundingBox( extended_cell )
! fix number of steps for each direction according to aspect ratio ...
do i = 1 , 3
! default case
n_xyz_steps(i) = medium_steps
If( extended_cell%BoxSides(i) >= 2.5d1 ) n_xyz_steps(i) = many_steps
If( extended_cell%BoxSides(i) <= 1.2d1 ) n_xyz_steps(i) = few_steps
end do
! voxel dimensions
dx = extended_cell%BoxSides(1) / n_xyz_steps(1)
dy = extended_cell%BoxSides(2) / n_xyz_steps(2)
dz = extended_cell%BoxSides(3) / n_xyz_steps(3)
! translation to the center of mass
forall(i=1:extended_cell%atoms,j=1:3) xyz(i,j) = extended_cell%coord(i,j) - extended_cell%Center_of_Mass(j)
! initial corner of the volume Box
a = minval(xyz(:,1)) - fringe / two
b = minval(xyz(:,2)) - fringe / two
c = minval(xyz(:,3)) - fringe / two
! start writing Gaussian cube files
write(4,*) System_Characteristics
write(4,*) 'electron_state = ',electron_state,' / time = ', t
write(4,111) extended_cell%atoms , a/a_Bohr , b/a_Bohr , c/a_Bohr
write(4,111) n_xyz_steps(1) + 1 , dx/a_Bohr , 0.d0 , 0.d0
write(4,111) n_xyz_steps(2) + 1 , 0.d0 , dy/a_Bohr , 0.d0
write(4,111) n_xyz_steps(3) + 1 , 0.d0 , 0.d0 , dz/a_Bohr
! coordinates in a.u., because zeta is in units of [a_0^{-1}] ...
xyz = xyz / a_Bohr
DO i = 1 , extended_cell%atoms
write(4,113) extended_cell%AtNo(i), 0.0 , xyz(i,1) , xyz(i,2) , xyz(i,3)
END DO
!---------------------------------------------------------
! drawing the wavefunction denssity in cube format
!========================================================
! the order orbitals are stored
!
! S --> 1 --> l = 0 , m = 0
! Py --> 2 --> l = 1 , m = -1
! Pz --> 3 --> l = 1 , m = 0
! Px --> 4 --> l = 1 , m = +1
! Dxy --> 5 --> l = 2 , m = -2
! Dyz --> 6 --> l = 2 , m = -1
! Dz2 --> 7 --> l = 2 , m = 0
! Dxz --> 8 --> l = 2 , m = +1
! Dx2y2 --> 9 --> l = 2 , m = +2
!========================================================
allocate( Psi_2(n_xyz_steps(1)+1,n_xyz_steps(2)+1,n_xyz_steps(3)+1) , source = D_zero )
!$OMP parallel private(ix,iy,iz,x,y,z,x0,y0,z0,SlaterOrbital,r,AtNo,i,j,k,l,TotalPsiBra,TotalPsiKet)
!$OMP single
DO ix = 0 , n_xyz_steps(1)
x = (a + ix * dx) / a_Bohr
!$OMP task untied
DO iy = 0 , n_xyz_steps(2)
y = (b + iy * dy) / a_Bohr
DO iz = 0 , n_xyz_steps(3)
z = (c + iz * dz) / a_Bohr
i = 0
TotalPsiBra = C_zero
TotalPsiKet = C_zero
DO k = 1 , extended_cell%atoms
if( Extended_Cell% QMMM(k) /= "QM" ) cycle
AtNo = extended_cell%AtNo(k)
! distance from the center of the nuclei
r = dsqrt((x-xyz(k,1))*(x-xyz(k,1)) + (y-xyz(k,2))*(y-xyz(k,2)) + (z-xyz(k,3))*(z-xyz(k,3)))
! coordinates centered on the nuclei
x0 = x - xyz(k,1)
y0 = y - xyz(k,2)
z0 = z - xyz(k,3)
do j = 1 , atom(AtNo)%DOS
i = i + 1
l = extended_cell%BasisPointer(k) + j
select case (j)
case( 1 )
SlaterOrbital = s_orb(r,AtNo,l)
case( 2 )
SlaterOrbital = p_orb(r,y0,AtNo,l)
case( 3 )
SlaterOrbital = p_orb(r,z0,AtNo,l)
case( 4 )
SlaterOrbital = p_orb(r,x0,AtNo,l)
case( 5 )
SlaterOrbital = d_xyz(r,x0,y0,AtNo,l)
case( 6 )
SlaterOrbital = d_xyz(r,y0,z0,AtNo,l)
case( 7 )
SlaterOrbital = d_z2(r,x0,y0,z0,AtNo,l)
case( 8 )
SlaterOrbital = d_xyz(r,x0,z0,AtNo,l)
case( 9 )
SlaterOrbital = d_x2y2(r,x0,y0,AtNo,l)
end select
TotalPsiBra = TotalPsiBra + bra(i) * SlaterOrbital
TotalPsiKet = TotalPsiKet + ket(i) * SlaterOrbital
end do ! <== DOS
end do ! <== atoms
Psi_2( ix+1 , iy+1 , iz+1 ) = cdabs( cdsqrt( TotalPsiBra * TotalPsiKet ) )
END DO ! <== Z coord
END DO ! <== Y coord
!$OMP end task
END DO ! <== X coord
!$OMP end single
!$OMP end parallel
!::::::::::::::::::::::::::::::::::::::::::::::::::::::::
DO ix = 0 , n_xyz_steps(1)
DO iy = 0 , n_xyz_steps(2)
DO iz = 0 , n_xyz_steps(3)
write(4,112,advance='no') Psi_2( ix+1 , iy+1 , iz+1 )
If( (mod(iz,6) == 5) ) write(4,'(a)',advance='yes')
END DO
write(4,'(a)',advance='yes')
END DO
END DO
!::::::::::::::::::::::::::::::::::::::::::::::::::::::::
close(4)
deallocate(xyz , Psi_2)
include 'formats.h'
end subroutine Gaussian_Cube_Format_Cmplx
!
!
!
!===========================
subroutine BoundingBox( a )
!===========================
implicit none
type(structure) :: a
! local variables ...
integer :: i
! size of the box
forall(i=1:3) &
a%BoxSides(i) = maxval(a%coord(:,i)) - minval(a%coord(:,i)) + fringe
! find the center of mass
forall(i=1:3) &
a%Center_of_Mass(i) = sum(a%coord(:,i)) / a%atoms
end subroutine BoundingBox
!
!
end module Psi_squared_cube_format