1      SUBROUTINE new4ind_N5_osu(rtdb,d_v2,
2     1                                kax_v2_alpha_offset,
3     1                                size_2e)
4C     $Id: new_ga4ind_N5.F 25748 2014-06-08 07:53:05Z d3y133 $
5C     This is a Fortran77 program generated by Tensor Contraction Engine v.1.0
6C     Copyright (c) Battelle & Pacific Northwest National Laboratory (2002)
7C     t ( p1 p2 h3 h4 )_t
8      IMPLICIT NONE
9#include "rtdb.fh"
10#include "global.fh"
11#include "mafdecls.fh"
12#include "stdio.fh"
13#include "util.fh"
14#include "bas.fh"
15#include "schwarz.fh"
16#include "sym.fh"
17#include "sf.fh"
18#include "errquit.fh"
19#include "tce.fh"
20#include "tce_main.fh"
21c
22c written by K. Kowalski
23c
24c
25c     max. number of p2 groups =200
26c
27c
28      integer rtdb                 ! Run-time database
29      integer d_v2                 ! MO integrals
30      integer kax_v2_alpha_offset  ! MO integrals offset
31      integer size_2e              ! 2e file size
32      integer iter
33      INTEGER azone1,azone2,azone3,azone4
34      INTEGER g1b,g2b,g3b,g4b
35      INTEGER igi1,igi2,igi3,igi4
36      INTEGER ii,i,j,k,l,N,ipos1,ipos2,ipos3,ipos4
37      INTEGER size_4a,l_4a,k_4a
38      INTEGER size_aaaa, size_full_length
39      INTEGER d_agaa,size_agaa,l_offset_agaa,k_offset_agaa
40      INTEGER d_ggaa,size_ggaa,l_offset_ggaa,k_offset_ggaa
41      INTEGER d_ggga,size_ggga,l_offset_ggga,k_offset_ggga
42      INTEGER l_loc_aaaa,k_loc_aaaa,size_loc_aaaa
43      INTEGER l_loc_agaa,k_loc_agaa,size_loc_agaa
44      INTEGER d_agaa_f,size_agaa_f,l_offset_agaa_f,k_offset_agaa_f
45      INTEGER l_loc_agaa_f,k_loc_agaa_f,size_loc_agaa_f
46      INTEGER l_loc_ggaa,k_loc_ggaa,size_loc_ggaa
47      INTEGER l_loc_ggga,k_loc_ggga,size_loc_ggga
48      INTEGER l_loc_gggg,k_loc_gggg,size_loc_gggg
49c
50      integer key_aaaa,key_agaa,key_ggaa,key_ggga
51      integer offset_aaaa,offset_agaa,offset_ggaa,offset_ggga
52      integer key_gggg,offset_gggg
53c
54      integer max_size_temp,sumx
55c
56      integer tot_azone1_sh,tot_azone2_sh
57      integer tot_azone3_sh,tot_azone4_sh
58      integer tot_zone(1000)  !it was d.prec.
59c
60      integer iha,ihb !number of corr. alpha, beta holes
61      integer ipa,ipb !number of corr. alpha, beta particles
62c
63      integer mu,nu,rho,sigma
64      integer mu_lo,mu_hi
65      integer nu_lo,nu_hi
66      integer rho_lo,rho_hi
67      integer sigma_lo,sigma_hi
68      integer mu_range
69      integer nu_range
70      integer rho_range
71      integer sigma_range
72      integer mu1,nu1,rho1,sigma1
73      integer shift_mu,shift_nu
74      integer shift_rho,shift_sigma
75      integer work1,work2          ! Work array sizes
76      integer l_work1,k_work1      ! Work array 1
77      integer l_work2,k_work2      ! Work array 2
78      integer imu1,inu1,irho1,isigma1
79c
80      integer l_movecs_orb,k_movecs_orb
81c
82      integer l_integral,l_coeff
83      integer k_integral,k_coeff
84      integer l_aux,k_aux,size_aux
85      integer size_ic,size_icc,size_integral,size_coeff,max_na
86c
87      integer INDEX_PAIR,icol,irow
88c
89ccx      double precision tot_zone(1000)
90c
91      integer l_4af_offset,k_4af_offset,d_4af
92      integer sf_chunk,request
93      integer key_4af,offset_4af,size_4af
94      character*255 filename
95c
96      logical parallel
97c
98      INTEGER length
99      INTEGER next
100      INTEGER nprocs
101      INTEGER count
102      integer nnn, ilo,ihi,jlo,jhi
103      integer nxtask
104      external nxtask
105      logical nodezero
106      logical idiskl
107      logical compute_gggg
108c
109c
110      nodezero=(ga_nodeid().eq.0)
111c
112c
113      max_size_temp=imaxsize**4
114c
115      do ii=1,1000
116       tot_zone(ii)=0.0d0
117      enddo
118      if(atpart.gt.1000)
119     &  call errquit('tce_zones: atpart too big',1,MA_ERR)
120      sumx=0
121      do ii=1,atpart
122       tot_zone(ii)=sumx
123       sumx=sumx+nalength(ii)
124      enddo
125c
126c
127c this module is called only if intorb = .true.
128c N is the number of correlated orbitals
129        N = nmo(1) - nfc(1) - nfv(1)
130        iha = nocc(1)-nfc(1)
131        ihb = nocc(ipol)-nfc(ipol)
132        ipa = nmo(1)-nocc(1)-nfv(1)
133        ipb = nmo(ipol)-nocc(ipol)-nfv(ipol)
134c
135c     pre-compute size_loc_agaa & size_loc_aaaa
136
137      size_full_length=0
138      DO azone3 =1,atpart
139         size_full_length = size_full_length+nalength(azone3)
140      enddo
141c
142      DO azone2 = 1,atpart
143      size_loc_agaa=0
144      size_loc_aaaa=0
145      size_loc_agaa_f=0
146        DO azone1 = 1,atpart
147            nnn=nalength(azone1)*nalength(azone2)
148            DO g3b = 1,noa+nva
149                  size_loc_agaa_f= max(size_loc_agaa_f,
150     V                 nnn*
151     1                 size_full_length*int_mb(k_range_alpha+g3b-1))
152               DO azone3 = 1,atpart
153                  size_loc_agaa= max(size_loc_agaa,
154     V                 nnn*
155     1                 nalength(azone3)*int_mb(k_range_alpha+g3b-1))
156                  DO azone4=1,atpart
157
158                     size_loc_aaaa=max(size_loc_aaaa,
159     V                    nnn*
160     1                    nalength(azone3)*nalength(azone4))
161                  enddo
162               enddo
163            enddo
164         enddo
165c      enddo
166c     pre-compute size_loc_ggaa & size_loc_ggga
167c
168      size_loc_ggaa=0
169      size_loc_ggga=0
170c      DO azone2 = 1,atpart
171         DO g2b = 1,noa+nva
172            DO g3b = 1,noa+nva
173               DO g4b = g3b,noa+nva
174                  size_loc_ggga= max(size_loc_ggga,
175     N                 nalength(azone2)*int_mb(k_range_alpha+g2b-1)*
176     1 int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1))
177
178                  DO azone1=1,atpart ! azone1
179                     size_loc_ggaa=max(size_loc_ggaa,
180     N                    nalength(azone1)*nalength(azone2)*
181     1  int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1))
182
183                  enddo
184               enddo
185            enddo
186         enddo
187c      enddo
188c *** debug ***
189c        if(nodezero) then
190c         write(6,*) 'in step1'
191c         call util_flush(6)
192c        endif
193c *************
194c
195c     Offset for 4a file
196c
197      sf_chunk=(imaxsize)**4
198      call tce_4a_offset_ns_osu(l_4af_offset,k_4af_offset,size_4af,
199     &    azone2)
200
201      call tce_gacreate_sloc_osu(d_4af, size_4af,
202     S     size_loc_aaaa, 'd_4af')
203ccx      call ga_zero(d_4af)
204c *** debug ***
205c        if(nodezero) then
206c         write(6,*) 'in step2'
207c         call util_flush(6)
208c        endif
209c *************
210c
211c alpha orbitals only
212c
213      if (.not.ma_push_get(mt_dbl,nbf*(iha+ipa)
214     1  ,"sorted MO coeffs",
215     2  l_movecs_orb,k_movecs_orb))
216     3  call errquit('tce_mo2e_zone: MA problem 1',0,
217     2    BASIS_ERR)
218c      call dcopy(nbf*(iha+ipa),0.0d0, 0,dbl_mb(k_movecs_orb), 1)
219      call dcopy(nbf*iha,
220     c     dbl_mb(k_movecs_sorted),1,
221     C     dbl_mb(k_movecs_orb),1)
222      call dcopy(nbf*ipa,
223     c     dbl_mb(k_movecs_sorted+(iha+ihb)*nbf),1,
224     C     dbl_mb(k_movecs_orb+iha*nbf),1)
225c
226c
227      call int_mem_2e4c(work1,work2)
228      if (.not.ma_push_get(mt_dbl,work1,'work1',l_work1,k_work1))
229     1  call errquit('tce_ao2e: MA problem work1',0,MA_ERR)
230      if (.not.ma_push_get(mt_dbl,work2,'work2',l_work2,k_work2))
231     1  call errquit('tce_ao2e: MA problem work2',1,MA_ERR)
232c
233c
234c *** debug ***
235c        if(nodezero) then
236c         write(6,*) 'in step3'
237c         call util_flush(6)
238c        endif
239c *************
240c
241c 4af file formed here
242c
243c
244c       if(nodezero) then
245c          write(6,'(A,F20.2,A)') ' starting step 0 at ',
246c     C         util_wallsec(), ' secs '
247c          call util_flush(6)
248c       endif
249      nprocs = GA_NNODES()
250      count = 0
251      call  ga_distribution(d_4af,ga_nodeid(),ilo,ihi,jlo,jhi)
252      if(jlo.gt.0) then
253cold      next = NXTASK(nprocs, 1)
254      DO azone1 = 1,atpart      !nu
255c      DO azone2 = 1,atpart !mu
256      DO azone3 = 1,atpart      !sigma
257      DO azone4 = azone3,atpart !rho
258         size_4a = nalength(azone1)*nalength(azone2)*
259     1            nalength(azone3)*nalength(azone4)
260         key_4af=azone4 - 1 + atpart * (azone3 - 1 +
261     &          atpart * (azone2 - 1 + atpart * (azone1 - 1)))
262        call tce_hash(int_mb(k_4af_offset),key_4af,offset_4af)
263       if((ilo.le.offset_4af+size_4a).and.
264     A      (ihi.ge.offset_4af+1)) then
265cold      IF (next.eq.count) THEN
266c ---------------------------
267        if(.not.ma_push_get(mt_dbl,size_4a,'4a',l_4a,k_4a))
268     1     call errquit('tce_4af_zones1: MA problem',0,MA_ERR)
269        call dfill(size_4a, 0.0d0, dbl_mb(k_4a), 1)
270         shift_mu = 0
271         do mu    = a2length(azone2)+1,a2length(azone2+1)
272            if (.not.bas_cn2bfr(ao_bas_han,mu,mu_lo,mu_hi))
273     1      call errquit('tce_ao2e: basis fn range problem 1',0,
274     2      BASIS_ERR)
275            mu_range = mu_hi - mu_lo + 1
276         shift_nu = 0
277         do nu    = a2length(azone1)+1,a2length(azone1+1)
278            if (.not.bas_cn2bfr(ao_bas_han,nu,nu_lo,nu_hi))
279     1      call errquit('tce_ao2e: basis fn range problem 1',0,
280     2      BASIS_ERR)
281            nu_range = nu_hi - nu_lo + 1
282         shift_rho = 0
283         do rho   = a2length(azone4)+1,a2length(azone4+1)
284            if (.not.bas_cn2bfr(ao_bas_han,rho,rho_lo,rho_hi))
285     1      call errquit('tce_ao2e: basis fn range problem 1',0,
286     2      BASIS_ERR)
287            rho_range = rho_hi - rho_lo + 1
288         shift_sigma = 0
289         do sigma = a2length(azone3)+1,a2length(azone3+1)
290            if (.not.bas_cn2bfr(ao_bas_han,sigma,sigma_lo,sigma_hi))
291     1      call errquit('tce_ao2e: basis fn range problem 1',0,
292     2      BASIS_ERR)
293            sigma_range = sigma_hi - sigma_lo + 1
294            if (schwarz_shell(rho,sigma)*schwarz_shell(mu,nu)
295     1          .ge. tol2e) then
296            call int_2e4c(ao_bas_han,mu,nu,ao_bas_han,rho,sigma,
297     1           work2,dbl_mb(k_work2),work1,dbl_mb(k_work1))
298c
299            i=0
300            do mu1     = 1,mu_range
301               imu1=mu1+shift_mu
302               do nu1     = 1,nu_range
303                  inu1=nu1+shift_nu
304                  do rho1    = 1,rho_range
305                     irho1=rho1+shift_rho
306                     do sigma1  = 1,sigma_range
307                        i=i+1
308                        isigma1=sigma1+shift_sigma
309c     (isigma1,irho1|inu1, imu1)
310                        ipos1=(((imu1-1)*nalength(azone1)+inu1-1)*
311     1                       nalength(azone4)+irho1-1)*nalength(azone3)
312     2                       +isigma1
313                        dbl_mb(k_4a+ipos1-1)=dbl_mb(k_work1+i-1)
314                     enddo
315                  enddo
316               enddo
317            enddo
318            end if !schwarz  screening
319         shift_sigma = shift_sigma + sigma_range
320         enddo !sigma
321         shift_rho   = shift_rho + rho_range
322         enddo !rho
323         shift_nu    = shift_nu + nu_range
324         enddo !nu
325         shift_mu    = shift_mu + mu_range
326         enddo !mu
327c
328c fixing offsets and sf_writing
329        call ga_put(d_4af,offset_4af+1,offset_4af+size_4a,1,1,
330     1    dbl_mb(k_4a),1)
331c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "azone4 : ",azone4," azone3 : ",
332c     C  azone3,"azone2 : ",azone2," azone1 : ",azone1
333c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "azone4 : ",azone4," azone3 : ",
334c     C  azone3,"azone2 : ",azone2," azone1 : ",azone1
335c        do iter = 1,20
336c          write(6,'(F5.2,A,$)') dbl_mb(k_4a+iter),' , '
337c          call util_flush(6)
338c        enddo
339c closing l_4a file
340        if (.not.ma_pop_stack(l_4a))
341     1   call errquit('tcc_mo2e_4af2: l_4a',15,MA_ERR)
342c ---------------------------
343cold      next = NXTASK(nprocs, 1)
344      END IF
345      count = count + 1
346      ENDDO !azone4
347      ENDDO !azone3
348c      ENDDO !azone2
349      ENDDO !azone1
350      endif
351c
352cold      next = NXTASK(-nprocs, 1)
353      call ga_sync()
354c
355      call tce_offset_ggaa_osu(l_offset_ggaa,k_offset_ggaa,
356     &     size_ggaa,azone2)
357c     open offset l_ggga
358c      call tce_offset_ggga(l_offset_ggga,k_offset_ggga,size_ggga,azone2)
359c
360      max_na=0
361      do i=1,atpart
362       if(nalength(i).gt.max_na) max_na=nalength(i)
363      enddo
364      size_icc=tile_dim*max_na
365c
366       if (.not.ma_push_get(mt_dbl,size_icc,'l_coeff',
367     1  l_coeff,k_coeff))
368     1  call errquit('tce_4s: MA problem l_coeff',0,MA_ERR)
369c
370c     step1
371c
372c *** debug ***
373c        if(nodezero) then
374c         write(6,*) 'in step4'
375c         call util_flush(6)
376c        endif
377c *************
378c
379c     open ga d_agaa (size_agaa)
380c       if(nodezero) then
381c          write(6,'(A,F20.2,A)') ' starting step 1+2 at ',
382c     C         util_wallsec(), ' secs '
383c          call util_flush(6)
384c       endif
385      call tce_gacreate_sloc_osu(d_ggaa, size_ggaa,
386     S     size_loc_ggaa, 'd_ggaa')
387c
388      nprocs = GA_NNODES()
389      count = 0
390cold      next = NXTASK(nprocs, 1)
391      call  ga_distribution(d_ggaa,ga_nodeid(),ilo,ihi,jlo,jhi)
392c
393      if(jlo.gt.0) then
394       if (.not.ma_push_get(mt_dbl,size_loc_agaa_f,'loc_agaa_f',
395     1  l_loc_agaa_f,k_loc_agaa_f))
396     1  call errquit('step0:0',0,MA_ERR)
397       if (.not.ma_push_get(mt_dbl,size_loc_agaa,'loc_agaa',
398     1  l_loc_agaa,k_loc_agaa))
399     1  call errquit('step1:1',0,MA_ERR)
400        if (.not.ma_push_get(mt_dbl,size_loc_aaaa,'loc_aaaa',
401     1  l_loc_aaaa,k_loc_aaaa))
402     1  call errquit('step1:2',0,MA_ERR)
403       if (.not.ma_push_get(mt_dbl,
404     M   max(size_loc_aaaa,size_loc_agaa), 'auxaaaa',
405     1   l_aux,k_aux))
406     1   call errquit('step1:3',0,MA_ERR)
407c
408c *** debug ***
409c        if(nodezero) then
410c         write(6,*) 'in step5'
411c         call util_flush(6)
412c        endif
413c *************
414      nprocs = GA_NNODES()
415      count = 0
416cold      next = NXTASK(nprocs, 1)
417c do parallel
418c      DO azone2 = 1,atpart
419      DO azone1 = 1,atpart
420      DO g3b = 1,noa+nva
421c the above loops are loops fused for local memory enhancement
422       g4b=g3b
423       key_ggaa=g4b-1+(noa+nva)*(g3b-1+
424     &  (noa+nva)*(azone1-1+atpart*(azone2-1)))
425c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "Before g3b : ",g3b," g4b : ",
426c     C  g4b,"azone2 : ",azone2," azone1 : ",azone1
427c *** debug ***
428c        if(nodezero) then
429c         write(6,*) 'in step6'
430c         call util_flush(6)
431c        endif
432c *************
433       call tce_hash_n(int_mb(k_offset_ggaa),key_ggaa,offset_ggaa)
434c *** debug ***
435c        if(nodezero) then
436c         write(6,*) 'in step6a'
437c         call util_flush(6)
438c        endif
439c *************
440       size_loc_ggaa= nalength(azone1)*nalength(azone2)*
441     1 int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1)
442c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "After g3b : ",g3b," g4b : ",
443c     C  g4b,"azone2 : ",azone2," azone1 : ",azone1
444c the condition identifies, if this processor is responsible
445c for computing ggaa[azone2,azone1,g3b,*] tiles. If it is responsible for
446c combuting g4b=1, then it is responsible for computing all g4b
447       if((ilo.le.offset_ggaa+size_loc_ggaa).and.
448     A      (ihi.ge.offset_ggaa+1)) then
449c
450      call dfill(size_loc_agaa_f,0.0d0, dbl_mb(k_loc_agaa_f), 1)
451      k_offset_agaa_f = 0;
452      DO azone3 = 1,atpart
453       size_loc_agaa= nalength(azone1)*nalength(azone2)*
454     1 nalength(azone3)*int_mb(k_range_alpha+g3b-1)
455       call dfill(size_loc_agaa,0.0d0, dbl_mb(k_loc_agaa), 1)
456       DO azone4=1,atpart
457        tot_azone4_sh=tot_zone(azone4)
458        size_loc_aaaa=nalength(azone1)*nalength(azone2)*
459     1  nalength(azone3)*nalength(azone4)
460c --- C_[g3][azone4]
461        i = 0
462        do irho1 =  1,nalength(azone4)
463        do igi3    =  1,int_mb(k_range_alpha+g3b-1)
464        i = i+1
465        ipos1=(int_mb(k_offset_alpha+g3b-1)+igi3-1)*nbf+tot_azone4_sh
466     &       +irho1
467        dbl_mb(k_coeff+i-1)=dbl_mb(k_movecs_orb+ipos1-1)
468        enddo
469        enddo
470c ---
471        if(azone4.le.azone3) then ! azone4 <= azone3
472         key_aaaa=azone3-1+atpart*(azone4-1+
473     &          atpart*(azone2-1+atpart*(azone1-1)))
474         call tce_hash(int_mb(k_4af_offset),key_aaaa,offset_aaaa)
475         call ga_get(d_4af,offset_aaaa+1,offset_aaaa+size_loc_aaaa,1,1,
476     1    dbl_mb(k_loc_aaaa),1)
477c now C([g3][azone4]) ([azone4]<=[azone3]|[azone1]<=[azone2])
478        call  dgemm('N','N',
479     1   int_mb(k_range_alpha+g3b-1),  !m
480     1   nalength(azone1)*nalength(azone2)*nalength(azone3), !n
481     1   nalength(azone4), !k
482     1   1.0d0,
483     1   dbl_mb(k_coeff),                    !A
484     1   int_mb(k_range_alpha+g3b-1), !lda
485     1   dbl_mb(k_loc_aaaa),nalength(azone4),          ! B,ldb
486     1   1.0d0,dbl_mb(k_loc_agaa),
487     1   int_mb(k_range_alpha+g3b-1))
488        else    ! azone4 > azone3
489         key_aaaa=azone4-1+atpart*(azone3-1+
490     &          atpart*(azone2-1+atpart*(azone1-1)))
491         call tce_hash(int_mb(k_4af_offset),key_aaaa,offset_aaaa)
492         call ga_get(d_4af,offset_aaaa+1,offset_aaaa+size_loc_aaaa,1,1,
493     1    dbl_mb(k_loc_aaaa),1)
494         CALL TCE_SORT_4KG_(dbl_mb(k_loc_aaaa),dbl_mb(k_aux),
495     &   nalength(azone3),nalength(azone4),
496     &   nalength(azone1),nalength(azone2),
497     &   2,1,3,4,1.0d0)
498c C([g3][azone4]) ([azone4]>[azone3]|[azone1]<=[azone2])
499        call  dgemm('N','N',
500     1   int_mb(k_range_alpha+g3b-1),  !m
501     1   nalength(azone1)*nalength(azone2)*nalength(azone3), !n
502     1   nalength(azone4), !k
503     1   1.0d0,
504     1   dbl_mb(k_coeff),                    !A
505     1   int_mb(k_range_alpha+g3b-1), !lda
506     1   dbl_mb(k_aux),nalength(azone4),          ! B,ldb
507     1   1.0d0,dbl_mb(k_loc_agaa),
508     1   int_mb(k_range_alpha+g3b-1))
509        end if  ! azone4 <= azone3
510       ENDDO  !azone4
511C TRANSPOSITION HERE (g3 az3|az1 <=az2) => (az3 g3|az1 <=az2)
512       CALL TCE_SORT_4KG_(dbl_mb(k_loc_agaa),
513     &  dbl_mb(k_loc_agaa_f+k_offset_agaa_f),
514     &  int_mb(k_range_alpha+g3b-1),nalength(azone3),
515     &  nalength(azone1),nalength(azone2),
516     &  2,1,3,4,1.0d0)
517c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "g3b : ",g3b," azone3 : ",
518c     C  azone3,"azone2 : ",azone2," azone1 : ",azone1
519c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "g3b : ",g3b," azone3 : ",
520c     C  azone3,"azone2 : ",azone2," azone1 : ",azone1
521c        do iter = 1,20
522c          write(6,'(F5.2,A,$)') dbl_mb(k_aux+iter),' , '
523c          call util_flush(6)
524c        enddo
525c
526cold      next = NXTASK(nprocs, 1)
527      k_offset_agaa_f = k_offset_agaa_f + size_loc_agaa
528      count = count + 1
529      ENDDO !azone3
530c       if(nodezero) then
531c          write(6,'(A,F20.2,A)') ' starting step 2 at ',
532c     C         util_wallsec(), ' secs '
533c          call util_flush(6)
534c       endif
535
536      DO g4b = g3b,noa+nva
537cold      IF (next.eq.count) THEN
538       key_ggaa=g4b-1+(noa+nva)*(g3b-1+
539     &  (noa+nva)*(azone1-1+atpart*(azone2-1)))
540c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "g3b : ",g3b," g4b : ",
541c     C  g4b,"azone2 : ",azone2," azone1 : ",azone1
542c *** debug ***
543c        if(nodezero) then
544c         write(6,*) 'in step7'
545c         call util_flush(6)
546c        endif
547c *************
548       call tce_hash_n(int_mb(k_offset_ggaa),key_ggaa,offset_ggaa)
549c *** debug ***
550c        if(nodezero) then
551c         write(6,*) 'in step7a'
552c         call util_flush(6)
553c        endif
554c *************
555c       write(6,'(A,F20.2,A)') ' Done ',
556c     C      util_wallsec(), ' secs '
557c       call util_flush(6)
558       size_loc_ggaa= nalength(azone1)*nalength(azone2)*
559     1 int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1)
560       if (.not.ma_push_get(mt_dbl,size_loc_ggaa,'loc_ggaa',
561     1 l_loc_ggaa,k_loc_ggaa))
562     1 call errquit('step2:1',0,MA_ERR)
563      call dfill(size_loc_ggaa,0.0d0, dbl_mb(k_loc_ggaa), 1)
564      k_offset_agaa_f = 0
565c ----- Inner Loop for Computing ggaa ---------- c
566       DO azone3=1,atpart
567        tot_azone3_sh=tot_zone(azone3)
568        size_loc_agaa=nalength(azone1)*nalength(azone2)*
569     1  int_mb(k_range_alpha+g3b-1)*nalength(azone3)
570c        if (.not.ma_push_get(mt_dbl,size_loc_agaa,'loc_agaa',
571c     1  l_loc_agaa,k_loc_agaa))
572c     1 call errquit('step2:2',0,MA_ERR)
573c --- C_[g4][azone3]
574        i = 0
575        do isigma1 =  1,nalength(azone3)
576        do igi4    =  1,int_mb(k_range_alpha+g4b-1)
577        i = i+1
578        ipos1=(int_mb(k_offset_alpha+g4b-1)+igi4-1)*nbf+tot_azone3_sh
579     &       +isigma1
580        dbl_mb(k_coeff+i-1)=dbl_mb(k_movecs_orb+ipos1-1)
581        enddo
582        enddo
583c ---
584c  C([g4b][az3]) * ([az3][g3b]|[az1][az2])
585        call  dgemm('N','N',
586     1   int_mb(k_range_alpha+g4b-1),  !m
587     1   nalength(azone1)*nalength(azone2)*int_mb(k_range_alpha+g3b-1), !n
588     1   nalength(azone3), !k
589     1   1.0d0,
590     1   dbl_mb(k_coeff),                    !A
591     1   int_mb(k_range_alpha+g4b-1), !lda
592     1   dbl_mb(k_loc_agaa_f+k_offset_agaa_f),nalength(azone3),          ! B,ldb
593     1   1.0d0,dbl_mb(k_loc_ggaa),
594     1   int_mb(k_range_alpha+g4b-1))
595c
596       k_offset_agaa_f = k_offset_agaa_f+size_loc_agaa
597       ENDDO
598c
599       call ga_put(d_ggaa,offset_ggaa+1,offset_ggaa+size_loc_ggaa,1,1,
600     1  dbl_mb(k_loc_ggaa),1)
601
602       if (.not.ma_pop_stack(l_loc_ggaa))
603     1 call errquit('g4ind:MA1',15,MA_ERR)
604
605       ENDDO
606c
607c       if((.not.nodezero).and.(g3b.eq.1)) then
608c          write(6,'(A)') "a"
609c          write(6,'(A,I2,A,I2,A,I2,A,I2)') "g3b : ",g3b," g4b : ",
610c     C         g4b,"azone2 : ",azone2," azone1 : ",azone1
611cc
612c          do iter = 1,10
613c             write(6,'(F5.2,A,$)') dbl_mb(k_loc_ggaa+iter),' , '
614c             call util_flush(6)
615c          enddo
616c       end if
617cold      next = NXTASK(nprocs, 1)
618      END IF
619      ENDDO !g3b
620      ENDDO !azone1
621
622
623c      ENDDO !azone2
624        if (.not.ma_pop_stack(l_aux))
625     1  call errquit('l_aux:MA1',15,MA_ERR)
626        if (.not.ma_pop_stack(l_loc_aaaa))
627     1  call errquit('l_loc_aaaa:MA1',15,MA_ERR)
628       if (.not.ma_pop_stack(l_loc_agaa))
629     1 call errquit('l_loc_agaa:MA1',15,MA_ERR)
630       if (.not.ma_pop_stack(l_loc_agaa_f))
631     1 call errquit('l_loc_agaa_f:MA1',15,MA_ERR)
632cold      next = NXTASK(-nprocs, 1)
633
634      endif
635      call ga_sync()
636c     delete d_4af
637      call deletefile(d_4af)
638
639c
640c     step3
641c
642c     open ga d_ggga (size_ggga)
643c       if(nodezero) then
644c          write(6,'(A,F20.2,A)') ' starting step 3 at ',
645c     C         util_wallsec(), ' secs '
646c          call util_flush(6)
647c       endif
648      call  ga_distribution(d_v2,ga_nodeid(),ilo,ihi,jlo,jhi)
649      nprocs = GA_NNODES()
650      count = 0
651c      if(nodezero) then
652c         write(6,'(A,I5,A,I5,A,I5,A,I5)')
653c     C        "High Low :",ilo,", ",
654c     C        ihi,", ",jlo,", ",jhi
655c      endif
656      if(jlo.gt.0) then
657      DO g3b = 1,noa+nva
658      DO g4b = g3b,noa+nva
659      compute_gggg=.FALSE.
660c     Checking if these g3b and g4b are computed on this node
661      DO g1b = 1,noa+nva
662      DO g2b = g1b,noa+nva
663       IF (int_mb(k_spin_alpha+g3b-1)+int_mb(k_spin_alpha+g4b-1).eq.
664     &        int_mb(k_spin_alpha+g1b-1)
665     &        +int_mb(k_spin_alpha+g2b-1)) THEN
666          IF (ieor(int_mb(k_sym_alpha+g3b-1),
667     &         ieor(int_mb(k_sym_alpha+g4b-1)
668     &         ,ieor(int_mb(k_sym_alpha+g1b-1),
669     &         int_mb(k_sym_alpha+g2b-1)))) .eq.
670     &         irrep_v) THEN
671             IROW=INDEX_PAIR(g4b,g3b)
672             ICOL=INDEX_PAIR(g2b,g1b)
673             IF(IROW.GE.ICOL) THEN
674                size_loc_gggg= int_mb(k_range_alpha+g1b-1)*
675     1          int_mb(k_range_alpha+g2b-1)*
676     1          int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1)
677                key_gggg=g2b - 1 + (noa+nva) *
678     &               (g1b - 1 + (noa+nva) * (g4b-
679     &               1 + (noa+nva) * (g3b - 1)))
680                call tce_hash_v2(int_mb(k_v2_alpha_offset),
681     &               key_gggg,offset_gggg)
682                if((ilo.le.offset_gggg+size_loc_gggg).and.
683     A               (ihi.ge.offset_gggg+1)) then
684                   compute_gggg = .TRUE.
685c                   if(nodezero) then
686c                      write(6,'(I5,A,I5,A,I5,A,I5,A,I5,A,I5)')
687c     C                     ,g3b,", ",
688c     C                     g4b,", ",g1b,", ",g2b,") : ",offset_gggg,
689c     C                     ", ",offset_gggg+size_loc_gggg
690c                      endif
691c
692                endif
693                go to 10
694             endif
695          endif
696       endif
697      ENDDO
698      ENDDO
699 10   continue
700c The tile must be computed at this node
701      IF(compute_gggg) THEN
702      DO g2b = 1,noa+nva
703cold      IF (next.eq.count) THEN
704       size_loc_ggga= nalength(azone2)*int_mb(k_range_alpha+g2b-1)*
705     1 int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1)
706       if (.not.ma_push_get(mt_dbl,size_loc_ggga,'loc_ggga',
707     1 l_loc_ggga,k_loc_ggga))
708     1 call errquit('step3:1',0,MA_ERR)
709       call dfill(size_loc_ggga,0.0d0, dbl_mb(k_loc_ggga), 1)
710       DO azone1=1,atpart ! azone1
711        tot_azone1_sh=tot_zone(azone1)
712        size_loc_ggaa=nalength(azone1)*nalength(azone2)*
713     1  int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1)
714        if (.not.ma_push_get(mt_dbl,size_loc_ggaa,'loc_ggaa',
715     1  l_loc_ggaa,k_loc_ggaa))
716     1  call errquit('step3:2',0,MA_ERR)
717c --- C_[azone1][g2]
718        ii=0
719        do igi2 =  1,int_mb(k_range_alpha+g2b-1)
720        do inu1 =  1,nalength(azone1)
721         ii = ii+1
722         ipos1=(int_mb(k_offset_alpha+g2b-1)+igi2-1)*nbf+tot_azone1_sh
723     &           +inu1
724         dbl_mb(k_coeff+ii-1)=dbl_mb(k_movecs_orb+ipos1-1)
725        enddo
726        enddo
727c ---
728         size_loc_ggaa= nalength(azone1)*nalength(azone2)*
729     1   int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1)
730         key_ggaa=g4b-1+(noa+nva)*(g3b-1+
731     &   (noa+nva)*(azone1-1+atpart*(azone2-1)))
732         call tce_hash(int_mb(k_offset_ggaa),key_ggaa,offset_ggaa)
733         call ga_get(d_ggaa,offset_ggaa+1,offset_ggaa+size_loc_ggaa,1,1,
734     1   dbl_mb(k_loc_ggaa),1)
735         if (.not.ma_push_get(mt_dbl,size_loc_ggaa,'auxaaaa',
736     1   l_aux,k_aux))
737     1   call errquit('step3:3',0,MA_ERR)
738         CALL TCE_SORT_4KG_(dbl_mb(k_loc_ggaa),dbl_mb(k_aux),
739     &   int_mb(k_range_alpha+g4b-1),int_mb(k_range_alpha+g3b-1),
740     &   nalength(azone1),nalength(azone2),
741     &   1,2,4,3,1.0d0)
742c ([g4]=>[g3]|[azone2]>[azone1])*C_([azone1][g2])
743         call dgemm('N','N',
744     1   int_mb(k_range_alpha+g4b-1)*int_mb(k_range_alpha+g3b-1)*
745     1   nalength(azone2),                                        !m
746     1   int_mb(k_range_alpha+g2b-1),                             !n
747     3   nalength(azone1),                                        !k
748     4   1.0d0,dbl_mb(k_aux),
749     5   int_mb(k_range_alpha+g4b-1)*int_mb(k_range_alpha+g3b-1)*
750     5   nalength(azone2),
751     6   dbl_mb(k_coeff),nalength(azone1),
752     7   1.0d0,dbl_mb(k_loc_ggga),
753     8   int_mb(k_range_alpha+g4b-1)*int_mb(k_range_alpha+g3b-1)*
754     8   nalength(azone2))
755         if (.not.ma_pop_stack(l_aux))
756     1   call errquit('g4ind:MA1',15,MA_ERR)
757         if (.not.ma_pop_stack(l_loc_ggaa))
758     1   call errquit('g4ind:MA1',15,MA_ERR)
759       ENDDO ! azone1
760C TRANSPOSITION HERE (g4=>g3|az2 g2) => (g4=>g3|g2 az2)
761       if (.not.ma_push_get(mt_dbl,size_loc_ggga,'auxggga',
762     1  l_aux,k_aux))
763     1 call errquit('step3:4',0,MA_ERR)
764       CALL TCE_SORT_4KG_(dbl_mb(k_loc_ggga),dbl_mb(k_aux),
765     & int_mb(k_range_alpha+g4b-1),int_mb(k_range_alpha+g3b-1),
766     & nalength(azone2),int_mb(k_range_alpha+g2b-1),
767     & 1,2,4,3,1.0d0)
768c       if(.NOT.nodezero) then
769c          write(6,'(I5,A,I5,A,I5,A,I5)')
770c     C     ga_nodeid(), ": (",g3b,", ",
771c     C         g4b,", ",g2b
772c          write(6,'(A,F20.2,A)') ' starting step 4 at ',
773c     C         util_wallsec(), ' secs '
774c          call util_flush(6)
775c       endif
776      DO g1b = 1,g2b
777       IF (int_mb(k_spin_alpha+g3b-1)+int_mb(k_spin_alpha+g4b-1).eq.
778     & int_mb(k_spin_alpha+g1b-1)+int_mb(k_spin_alpha+g2b-1)) THEN
779       IF (ieor(int_mb(k_sym_alpha+g3b-1),ieor(int_mb(k_sym_alpha+g4b-1)
780     & ,ieor(int_mb(k_sym_alpha+g1b-1),int_mb(k_sym_alpha+g2b-1)))) .eq.
781     & irrep_v) THEN
782       IROW=INDEX_PAIR(g4b,g3b)
783       ICOL=INDEX_PAIR(g2b,g1b)
784       IF(IROW.GE.ICOL) THEN
785       size_loc_gggg= int_mb(k_range_alpha+g1b-1)*
786     1 int_mb(k_range_alpha+g2b-1)*
787     1 int_mb(k_range_alpha+g3b-1)*int_mb(k_range_alpha+g4b-1)
788       key_gggg=g2b - 1 + (noa+nva) *
789     &(g1b - 1 + (noa+nva) * (g4b-
790     & 1 + (noa+nva) * (g3b - 1)))
791       call tce_hash_v2(int_mb(k_v2_alpha_offset),key_gggg,offset_gggg)
792       if (.not.ma_push_get(mt_dbl,size_loc_gggg,'loc_gggg',
793     1 l_loc_gggg,k_loc_gggg))
794     1 call errquit('step4:1',0,MA_ERR)
795      call dfill(size_loc_gggg,0.0d0, dbl_mb(k_loc_gggg), 1)
796        tot_azone2_sh=tot_zone(azone2)
797c --- C_[azone2][g1]
798        ii=0
799        do igi1 =  1,int_mb(k_range_alpha+g1b-1)
800        do imu1 =  1,nalength(azone2)
801         ii = ii+1
802         ipos1=(int_mb(k_offset_alpha+g1b-1)+igi1-1)*nbf+tot_azone2_sh
803     &           +imu1
804         dbl_mb(k_coeff+ii-1)=dbl_mb(k_movecs_orb+ipos1-1)
805        enddo
806        enddo
807c ---
808        call dgemm('N','N',
809     1   int_mb(k_range_alpha+g4b-1)*int_mb(k_range_alpha+g3b-1)*
810     1   int_mb(k_range_alpha+g2b-1),                             !m
811     1   int_mb(k_range_alpha+g1b-1),                             !n
812     3   nalength(azone2),                                        !k
813     4   1.0d0,dbl_mb(k_aux),
814     5   int_mb(k_range_alpha+g4b-1)*int_mb(k_range_alpha+g3b-1)*
815     5   int_mb(k_range_alpha+g2b-1),
816     6   dbl_mb(k_coeff),nalength(azone2),
817     7   1.0d0,dbl_mb(k_loc_gggg),
818     8   int_mb(k_range_alpha+g4b-1)*int_mb(k_range_alpha+g3b-1)*
819     8   int_mb(k_range_alpha+g2b-1))
820c
821       call ga_acc(d_v2,offset_gggg+1,offset_gggg+size_loc_gggg,1,1,
822     1  dbl_mb(k_loc_gggg),1,1.0D0)
823c       write(6,'(A,I5,A,I5,A,I5,A,I5)') "g3b : ",g3b," g4b : ",
824c     C  g4b,"g2 : ",g2b," g1 : ",g1b
825c       if(.NOT.nodezero) then
826c          write(6,'(A,I5,A,I5,A,I5,A,I5,A,I5)') "Node ID : ",
827c     C     ga_nodeid(), "g3b : ",g3b," g4b : ",
828c     C         g4b,"g2 : ",g2b," g1 : ",g1b
829c          do iter = 1,20
830c             write(6,'(F5.2,A,$)') dbl_mb(k_loc_gggg+iter),' , '
831c             call util_flush(6)
832c          enddo
833c       endif
834       if (.not.ma_pop_stack(l_loc_gggg))
835     1 call errquit('g4ind:MA1',15,MA_ERR)
836cold      next = NXTASK(nprocs, 1)
837      END IF
838      END IF
839      END IF
840      ENDDO !g1b
841       if (.not.ma_pop_stack(l_aux))
842     1  call errquit('g4ind:MA1',15,MA_ERR)
843       if (.not.ma_pop_stack(l_loc_ggga))
844     1 call errquit('g4ind:MA1',15,MA_ERR)
845      ENDDO !g2b
846      ENDIF
847      ENDDO                     !g3b
848      ENDDO !g4b
849      endif
850      call ga_sync()
851c     delete l_coeff
852      if (.not.ma_pop_stack(l_coeff))
853     1  call errquit('tcc_off_4a: MA problem',15,MA_ERR)
854c     delete l_ggaa
855      if (.not.ma_pop_stack(l_offset_ggaa))
856     1  call errquit('ga4ind:ggaa',15,MA_ERR)
857c     delete l_agaa
858c      if (.not.ma_pop_stack(l_offset_agaa))
859c     1  call errquit('ga4ind:agaa',15,MA_ERR)
860c
861      if (.not.ma_pop_stack(l_work2))
862     1  call errquit('tcc_ao2e: MA problem',14,MA_ERR)
863c
864      if (.not.ma_pop_stack(l_work1))
865     1  call errquit('tcc_ao2e: MA problem',15,MA_ERR)
866c
867      if (.not.ma_pop_stack(l_movecs_orb))
868     1  call errquit('tcc_ao2e: MA problem',15,MA_ERR)
869c
870      if (.not.ma_pop_stack(l_4af_offset))
871     1  call errquit('ga4ind:4a',15,MA_ERR)
872
873      ENDDO
874c      if(nodezero) then
875c          write(6,'(A,F20.2,A)') ' done step 4 at ',
876c     C         util_wallsec(), ' secs '
877c          call util_flush(6)
878c       endif
879c
880c
881c
882      RETURN
883      END
884c
885c
886c
887c
888c
889c
890c
891c
892      SUBROUTINE tce_offset_agaa_osu(l_a_offset,k_a_offset,size,
893     & azone2_para)
894      IMPLICIT NONE
895#include "global.fh"
896#include "mafdecls.fh"
897#include "sym.fh"
898#include "errquit.fh"
899#include "tce.fh"
900#include "tce_main.fh"
901      INTEGER l_a_offset
902      INTEGER k_a_offset
903      INTEGER size
904      INTEGER length
905      INTEGER addr
906      INTEGER g3b,azone1,azone2,azone3,azone2_para
907c ([az3],[g3]|[az1]<=[az2])
908      length = 0
909      azone2 = azone2_para
910c      DO azone2 = 1,atpart
911      DO azone1 = 1,atpart
912      DO g3b = 1,noa+nva
913      DO azone3 = 1,atpart
914      length = length + 1
915      END DO
916      END DO
917c      END DO
918      END DO
919      IF (.not.MA_PUSH_GET(mt_int,2*length+1,'noname',l_a_offset,k_a_off
920     &set)) CALL ERRQUIT('tce_t2_offset',0,MA_ERR)
921      int_mb(k_a_offset) = length
922      addr = 0
923      size = 0
924c      DO azone2 = 1,atpart
925      DO azone1 = 1,atpart
926      DO g3b = 1,noa+nva
927      DO azone3 = 1,atpart
928      addr = addr + 1
929      int_mb(k_a_offset+addr) = azone3 - 1 + atpart * (g3b - 1 +
930     &  (noa+nva) * (azone1 - 1 + atpart * (azone2 - 1)))
931      int_mb(k_a_offset+length+addr) = size
932      size = size + nalength(azone2) * nalength(azone1) *
933     &  int_mb(k_range_alpha+g3b-1) * nalength(azone3)
934      END DO
935c      END DO
936      END DO
937      END DO
938      RETURN
939      END
940c
941c
942c
943      SUBROUTINE tce_offset_ggaa_osu(l_a_offset,k_a_offset,size,
944     & azone2_para)
945      IMPLICIT NONE
946#include "global.fh"
947#include "mafdecls.fh"
948#include "sym.fh"
949#include "errquit.fh"
950#include "tce.fh"
951#include "tce_main.fh"
952      INTEGER l_a_offset
953      INTEGER k_a_offset
954      INTEGER size
955      INTEGER length
956      INTEGER addr
957      INTEGER g3b,g4b,azone1,azone2,azone2_para
958      azone2 = azone2_para
959c ([g4]=>[g3]|[az1]<=[az2])
960      length = 0
961c      DO azone2 = 1,atpart
962      DO azone1 = 1,atpart
963      DO g3b = 1,noa+nva
964      DO g4b = g3b,noa+nva
965      length = length + 1
966      END DO
967      END DO
968c      END DO
969      END DO
970      IF (.not.MA_PUSH_GET(mt_int,2*length+1,'noname',l_a_offset,k_a_off
971     &set)) CALL ERRQUIT('tce_t2_offset',0,MA_ERR)
972      int_mb(k_a_offset) = length
973      addr = 0
974      size = 0
975c      DO azone2 = 1,atpart
976      DO azone1 = 1,atpart
977      DO g3b = 1,noa+nva
978      DO g4b = g3b,noa+nva
979      addr = addr + 1
980      int_mb(k_a_offset+addr) = g4b - 1 + (noa+nva) * (g3b - 1 +
981     &  (noa+nva) * (azone1 - 1 + atpart * (azone2 - 1)))
982      int_mb(k_a_offset+length+addr) = size
983      size = size + nalength(azone2) * nalength(azone1) *
984     &  int_mb(k_range_alpha+g3b-1) * int_mb(k_range_alpha+g4b-1)
985      END DO
986c     END DO
987      END DO
988      END DO
989      RETURN
990      END
991c
992c
993      SUBROUTINE tce_4a_offset_ns_osu(l_a_offset,k_a_offset,size,
994     & azone2_pa)
995C     $Id: tce_mo2e_zones_4a_disk_2s_new.F 19706 2010-10-29 17:52:31Z d3y133 $
996C     This is a Fortran77 program generated by Tensor Contraction Engine v.1.0
997C     Copyright (c) Battelle & Pacific Northwest National Laboratory (2002)
998      IMPLICIT NONE
999#include "global.fh"
1000#include "mafdecls.fh"
1001#include "sym.fh"
1002#include "errquit.fh"
1003#include "tce.fh"
1004#include "tce_main.fh"
1005      INTEGER l_a_offset
1006      INTEGER k_a_offset
1007      INTEGER size
1008      INTEGER length
1009      INTEGER addr
1010      INTEGER mu
1011      INTEGER nu
1012      INTEGER rho
1013      INTEGER sigma
1014      INTEGER azone1,azone2,azone3,azone4,azone2_pa
1015      length = 0
1016      azone2 = azone2_pa
1017      DO azone1 = 1,atpart      !nu
1018c      DO azone2 = 1,atpart !mu
1019      DO azone3 = 1,atpart      !sigma
1020      DO azone4 = azone3,atpart !rho
1021      length = length + 1
1022      END DO
1023      END DO
1024c     END DO
1025      END DO
1026      IF (.not.MA_PUSH_GET(mt_int,2*length+1,'noname',l_a_offset,k_a_off
1027     &set)) CALL ERRQUIT('tce_t2_offset',0,MA_ERR)
1028      int_mb(k_a_offset) = length
1029      addr = 0
1030      size = 0
1031      DO azone1 = 1,atpart      !nu
1032c      DO azone2 = 1,atpart !mu
1033      DO azone3 = 1,atpart      !sigma
1034      DO azone4 = azone3,atpart !rho
1035      addr = addr + 1
1036      int_mb(k_a_offset+addr) = azone4 - 1 + atpart * (azone3 - 1 +
1037     &  atpart * (azone2 - 1 + atpart * (azone1 - 1)))
1038      int_mb(k_a_offset+length+addr) = size
1039      size = size + nalength(azone1) * nalength(azone2) *
1040     &  nalength(azone3) * nalength(azone4)
1041      END DO
1042      END DO
1043c      END DO
1044      END DO
1045      RETURN
1046      END
1047
1048
1049
1050c
1051      SUBROUTINE tce_offset_ggga_osu(l_a_offset,k_a_offset,size,
1052     & azone2_para)
1053      IMPLICIT NONE
1054#include "global.fh"
1055#include "mafdecls.fh"
1056#include "sym.fh"
1057#include "errquit.fh"
1058#include "tce.fh"
1059#include "tce_main.fh"
1060      INTEGER l_a_offset
1061      INTEGER k_a_offset
1062      INTEGER size
1063      INTEGER length
1064      INTEGER addr
1065      INTEGER g3b,g4b,g2b,azone2, azone2_para
1066c ([g4]=>[g3]|[g2][az2])
1067      length = 0
1068      azone2 = azone2_para
1069c      DO azone2 = 1,atpart
1070      DO g2b = 1,noa+nva
1071      DO g3b = 1,noa+nva
1072      DO g4b = g3b,noa+nva
1073      length = length + 1
1074      END DO
1075      END DO
1076c      END DO
1077      END DO
1078      IF (.not.MA_PUSH_GET(mt_int,2*length+1,'noname',l_a_offset,k_a_off
1079     &set)) CALL ERRQUIT('tce_t2_offset',0,MA_ERR)
1080      int_mb(k_a_offset) = length
1081      addr = 0
1082      size = 0
1083c      DO azone2 = 1,atpart
1084      DO g2b = 1,noa+nva
1085      DO g3b = 1,noa+nva
1086      DO g4b = g3b,noa+nva
1087      addr = addr + 1
1088      int_mb(k_a_offset+addr) = g4b - 1 + (noa+nva) * (g3b - 1 +
1089     &  (noa+nva) * (g2b - 1 + (noa+nva) * (azone2 - 1)))
1090      int_mb(k_a_offset+length+addr) = size
1091      size = size + nalength(azone2) * int_mb(k_range_alpha+g2b-1) *
1092     &  int_mb(k_range_alpha+g3b-1) * int_mb(k_range_alpha+g4b-1)
1093      END DO
1094      END DO
1095c      END DO
1096      END DO
1097      RETURN
1098      END
1099c
1100      subroutine tce_gacreate_sloc_osu(g_a, size, size_loc, gname)
1101      implicit none
1102#include "global.fh"
1103#include "mafdecls.fh"
1104#include "errquit.fh"
1105#include "stdio.fh"
1106      integer g_a
1107      character*(*) gname
1108      integer size
1109      integer size_loc
1110c
1111      if (.not.ga_create(mt_dbl,size,1,gname,
1112     1     size_loc,1,g_a)) then
1113         write(LuOut,*) ' available GA memory ',
1114     1        ga_memory_avail(),' bytes'
1115         call errquit ('tce_gacreate: failed ga_create size/nproc bytes'
1116     S      ,   (size*ma_sizeof(mt_dbl,1,mt_byte))/ga_nnodes(),
1117     1        GA_ERR)
1118      endif
1119      return
1120      end
1121