-
Notifications
You must be signed in to change notification settings - Fork 4
/
Copy pathcalc_3c_colc_setup.f90
114 lines (99 loc) · 4.31 KB
/
calc_3c_colc_setup.f90
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
!
! ParaGauss, a program package for high-performance computations of
! molecular systems
!
! Copyright (C) 2014 T. Belling, T. Grauschopf, S. Krüger,
! F. Nörtemann, M. Staufer, M. Mayer, V. A. Nasluzov, U. Birkenheuer,
! A. Hu, A. V. Matveev, A. V. Shor, M. S. K. Fuchs-Rohr, K. M. Neyman,
! D. I. Ganyushin, T. Kerdcharoen, A. Woiterski, A. B. Gordienko,
! S. Majumder, M. H. i Rotllant, R. Ramakrishnan, G. Dixit,
! A. Nikodem, T. Soini, M. Roderus, N. Rösch
!
! This program is free software; you can redistribute it and/or modify
! it under the terms of the GNU General Public License version 2 as
! published by the Free Software Foundation [1].
!
! This program is distributed in the hope that it will be useful, but
! WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
! General Public License for more details.
!
! [1] http://www.gnu.org/licenses/gpl-2.0.html
!
! Please see the accompanying LICENSE file for further information.
!
subroutine calc_3c_colc_setup(la,lb)
#include "def.h"
use calc_3center_module
use cpksdervs_matrices, only: cpks_coul_grad_size
USE_MEMLOG
implicit none
integer(kind=i4_kind),intent(in):: la,lb
! *** end of interface ***
! local variables
integer(kind=i4_kind):: k_gr,k2dr,i_grad
allocate(coul_int(lc)%m(num, ncexps, n_independent_fcts, 2*lb+1, 2*la+1), &
stat=alloc_stat(14)) ! 1
MEMLOG(+size(coul_int(lc)%m))
ASSERT(alloc_stat(14).eq.0)
alloc_stat(14)=1
pointer_coul => coul_int(lc)%m
pointer_coul=0.0_r8_kind
if(integralpar_3cob_grad) then
do k_gr=1,6
allocate(coul_int_grad(k_gr)%l(lc)%m &
(num,ncexps,n_independent_fcts,2*lb+1,2*la+1),stat=alloc_stat(56))
ASSERT(alloc_stat(56).eq.0)
alloc_stat(56)=1 !2
#if 1
cpks_coul_grad_size=cpks_coul_grad_size+size(coul_int_grad(k_gr)%l(lc)%m)
#endif
MEMLOG(+size(coul_int_grad(k_gr)%l(lc)%m))
coul_int_grad(k_gr)%l(lc)%m=0.0_r8_kind
if(integralpar_dervs) then
do k2dr=1,6
allocate(coul_int_dervs(k_gr,k2dr)%l(lc)%m &
(num,ncexps,n_independent_fcts,2*lb+1,2*la+1),stat=alloc_stat(68))
! print*,'coul_int_dervs alloc',k_gr,k2dr,lc
MEMLOG(+size(coul_int_dervs(k_gr,k2dr)%l(lc)%m))
ASSERT(alloc_stat(68).eq.0)
alloc_stat(68)=1 !2
coul_int_dervs(k_gr,k2dr)%l(lc)%m=0.0_r8_kind
enddo
endif
enddo
lc_gr_tots: do i_grad=1,grad_dim ! only if moving_c
allocate(coul_int_grad_totsym(i_grad)%l(lc)%m(num,&
ncexps,n_independent_fcts,2*lb+1,2*la+1),&
stat=alloc_stat(60))
ASSERT(alloc_stat(60).eq.0)
alloc_stat(60)=1
MEMLOG(+size(coul_int_grad_totsym(i_grad)%l(lc)%m))
#if 1
cpks_coul_grad_size=cpks_coul_grad_size+&
size(coul_int_grad_totsym(i_grad)%l(lc)%m)
#endif
coul_int_grad_totsym(i_grad)%l(lc)%m=0.0_r8_kind
if(integralpar_dervs) then
do k2dr=1,grad_dim
allocate(coul_int_dervs_totsym(i_grad,k2dr)%l(lc)%m(num&
,ncexps,n_independent_fcts,2*lb+1,2*la+1),&
stat=alloc_stat(71))
MEMLOG(+size(coul_int_dervs_totsym(i_grad,k2dr)%l(lc)%m))
ASSERT(alloc_stat(71).eq.0)
alloc_stat(71)=1
coul_int_dervs_totsym(i_grad,k2dr)%l(lc)%m=0.0_r8_kind
enddo
do k2dr=1,6
allocate(coul_int_ca_dervs(i_grad,k2dr)%l(lc)%m(num,&
ncexps,n_independent_fcts,2*lb+1,2*la+1),&
stat=alloc_stat(75))
MEMLOG(+size(coul_int_ca_dervs(i_grad,k2dr)%l(lc)%m))
ASSERT(alloc_stat(75).eq.0)
alloc_stat(75)=1
coul_int_ca_dervs(i_grad,k2dr)%l(lc)%m=0.0_r8_kind
enddo
endif
enddo lc_gr_tots
endif
end subroutine calc_3c_colc_setup