TABLE OF CONTENTS


ABINIT/memana [ Functions ]

[ Top ] [ Functions ]

NAME

 memana

FUNCTION

 Analysis of the memory and disk space needed for the job,
 thanks to the data computed in the calling routine : for each
 array, the number of blocks of size mpw or nfft bytes, and the
 additional memory occupation;
 the list of arrays that are used for each chain.

 According to the value of the option variable,
 the routine will eventually try to allocate this amount of memory,
 and if it fails, estimate the maximum value nfft compatible with
 the available memory.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (XG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  cadd(marrays)= count of bytes needed in addition of cmpw, cfftc and cfft.
  cfft(marrays) =for each array, count of blocks of size nfft bytes (coarse grid, if PAW)
  cfftf(marrays)=for each array, count of blocks of size nfft bytes (fine grid, if PAW)
  chain(marrays,nchain)=logical variable, that informs whether an array
    belongs to a given chain.
  cmpw(marrays)=for each array, count of blocks of size mpw bytes.
  dttyp(marrays)=datatype of the array : 4 for integers, 8 for real(dp)
  iout=unit number for output of formatted data.
  iprcel=govern the choice of preconditioner for the SCF cycle
  iscf=governs the choice of SCF algorithm, or non-SCF calculation.
  marrays=maximal number of arrays (or group of arrays) to be monitored.
  mbcg=number of MB needed for the cg array.
  mbdiskpd=number of MB needed to store a density or potential file on disk
  mbdiskwf=number of MB needed to store a wavefunction file on disk
  mbf_fftgr=number of MB needed for the f_fftgr array.
  mbgylm=number of MB needed for the pawfgrtab%gylm array (paw only)
  mffmem =governs the number of FFT arrays which are fit in core memory
  mpw   =maximum number of planewaves in basis sphere (large number)
  natom =number of atoms in unit cell
  nchain=number of chains to be used in the estimation of memory.
  nfft =(effective) number of FFT grid points (for one processor) (coarse grid, if PAW)
  nfftf=(effective) number of FFT grid points (for one processor) (fine grid, if PAW)
  occopt=option for occupation numbers. If 3<=occopt<=8, varying occupation
  option : if 0 , no test of available memory
           if 1 , the routine tries to allocate the estimated memory, for testing
                    purposes, and if a failure occurs, the routine stops.
           if 2 , like 1, but before stopping, the routine will provide
                    an estimation of the available memory.
  prtvol=control print volume

OUTPUT

  (only writing)

PARENTS

      memorf,memory

CHILDREN

      wrtout

SOURCE

 67 #if defined HAVE_CONFIG_H
 68 #include "config.h"
 69 #endif
 70 
 71 #include "abi_common.h"
 72 
 73 subroutine memana(cadd,cfft,cfftf,chain,cmpw,dttyp,iout,iprcel,iscf,&
 74 & marrays,mbcg,mbdiskpd,mbdiskwf,mbf_fftgr,mbgylm,mffmem,&
 75 & mpw,natom,nchain,nfft,nfftf,occopt,option,prtvol)
 76 
 77  use defs_basis
 78  use m_profiling_abi
 79  use m_errors
 80  use m_xmpi
 81 
 82 !This section has been created automatically by the script Abilint (TD).
 83 !Do not modify the following lines by hand.
 84 #undef ABI_FUNC
 85 #define ABI_FUNC 'memana'
 86  use interfaces_14_hidewrite
 87 !End of the abilint section
 88 
 89  implicit none
 90 
 91 !Arguments ------------------------------------
 92 !scalars
 93  integer,intent(in) :: iout,iprcel,iscf,marrays,mffmem,mpw,natom,nchain
 94  integer,intent(in) :: nfft,nfftf,occopt,option,prtvol
 95  real(dp),intent(in) :: mbcg,mbdiskpd,mbdiskwf,mbf_fftgr,mbgylm
 96 !arrays
 97  integer,intent(in) :: dttyp(marrays)
 98  logical,intent(in) :: chain(marrays,nchain)
 99  real(dp),intent(in) :: cadd(marrays),cfft(marrays),cfftf(marrays),cmpw(marrays)
100 
101 !Local variables-------------------------------
102 !scalars
103  integer :: biggest,ichain,ier,ier1,ier2,ier3,ier4,ier5,ier6,ier7,ier8,ii
104 !integer :: jj,kk
105  integer :: mu,nmbytes,nquarter_mbytes,quit
106  real(dp) :: mbbigarr,mbbiggest
107  character(len=500) :: message
108 !arrays
109  real(dp),allocatable :: bigarray(:,:),bigarray1(:,:),bigarray2(:,:)
110  real(dp),allocatable :: bigarray3(:,:),bigarray4(:,:),bigarray5(:,:)
111  real(dp),allocatable :: bigarray6(:,:),bigarray7(:,:),bigarray8(:,:)
112  real(dp),allocatable :: cdpadd(:),cdpfft(:),cdpfftf(:),cdpmpw(:)
113  real(dp),allocatable :: cintfft(:),cintfftf(:),cintmpw(:),cintadd(:)
114  real(dp),allocatable :: mbdpadd(:),mbdpfft(:),mbdpfftf(:)
115  real(dp),allocatable :: mbdpmpw(:),mbintadd(:),mbintfft(:),mbintfftf(:)
116  real(dp),allocatable :: mbintmpw(:),mbother(:),mbtot(:)
117 
118 ! **************************************************************************
119 
120 !write(std_out,*)' memana : nchain=',nchain
121 
122  ABI_ALLOCATE(cdpfftf,(nchain))
123  ABI_ALLOCATE(cdpfft,(nchain))
124  ABI_ALLOCATE(cdpmpw,(nchain))
125  ABI_ALLOCATE(cintfftf,(nchain))
126  ABI_ALLOCATE(cintfft,(nchain))
127  ABI_ALLOCATE(cintmpw,(nchain))
128  ABI_ALLOCATE(cdpadd,(nchain))
129  ABI_ALLOCATE(cintadd,(nchain))
130  ABI_ALLOCATE(mbdpadd,(nchain))
131  ABI_ALLOCATE(mbdpfftf,(nchain))
132  ABI_ALLOCATE(mbdpfft,(nchain))
133  ABI_ALLOCATE(mbdpmpw,(nchain))
134  ABI_ALLOCATE(mbintadd,(nchain))
135  ABI_ALLOCATE(mbintfftf,(nchain))
136  ABI_ALLOCATE(mbintfft,(nchain))
137  ABI_ALLOCATE(mbintmpw,(nchain))
138  ABI_ALLOCATE(mbother,(nchain))
139  ABI_ALLOCATE(mbtot,(nchain))
140 
141  biggest=0
142  mbbiggest=0.0_dp
143 
144 !For each chain, compute the number of bytes
145  do ichain=1,nchain
146 
147 !  First, the number of integer or real(dp), fft, mpw or add blocks
148    cdpmpw(ichain) =sum(cmpw(:),MASK=(dttyp(:)==8).and.chain(:,ichain))
149    cintmpw(ichain)=sum(cmpw(:),MASK=(dttyp(:)==4).and.chain(:,ichain))
150    cdpfftf(ichain) =sum(cfftf(:),MASK=(dttyp(:)==8).and.chain(:,ichain))
151    cintfftf(ichain)=sum(cfftf(:),MASK=(dttyp(:)==4).and.chain(:,ichain))
152    cdpfft(ichain) =sum(cfft(:),MASK=(dttyp(:)==8).and.chain(:,ichain))
153    cintfft(ichain)=sum(cfft(:),MASK=(dttyp(:)==4).and.chain(:,ichain))
154    cdpadd(ichain) =sum(cadd(:),MASK=(dttyp(:)==8).and.chain(:,ichain))
155    cintadd(ichain)=sum(cadd(:),MASK=(dttyp(:)==4).and.chain(:,ichain))
156 
157 !  Compute the corresponding number of Mbytes
158    mbdpmpw(ichain) =8*cdpmpw(ichain) *dble(mpw) /1024._dp**2
159    mbintmpw(ichain)=4*cintmpw(ichain)*dble(mpw) /1024._dp**2
160    mbdpfftf(ichain) =8*cdpfftf(ichain) *dble(nfftf)/1024._dp**2
161    mbintfftf(ichain)=4*cintfftf(ichain)*dble(nfftf)/1024._dp**2
162    mbdpfft(ichain) =8*cdpfft(ichain) *dble(nfft)/1024._dp**2
163    mbintfft(ichain)=4*cintfft(ichain)*dble(nfft)/1024._dp**2
164    mbdpadd(ichain) =8*cdpadd(ichain)              /1024._dp**2
165    mbintadd(ichain)=4*cintadd(ichain)             /1024._dp**2
166    mbother(ichain) =dble(231+6*natom)/1024._dp
167    if(3<=occopt .and. occopt<=8)mbother(ichain)=dble(991+natom)/1024._dp
168 
169 !  Compute the total number of Mbytes
170    mbtot(ichain)=mbdpmpw(ichain)+mbintmpw(ichain)&
171 &   +mbdpfftf(ichain)+mbintfftf(ichain)&
172 &   +mbdpfft(ichain)+mbintfft(ichain)&
173 &   +mbdpadd(ichain)+mbintadd(ichain)+mbother(ichain)
174 
175 !  Select the biggest chain
176    if(mbtot(ichain)>mbbiggest)then
177      mbbiggest=mbtot(ichain)
178      biggest=ichain
179    end if
180  end do
181 !When iprcel<20, the biggest chains cannot be number 8 or 9 ...
182  if(modulo(iprcel,100)<20 .and. (biggest==8 .or. biggest==9))then
183    write(message,'(a,a,a,a,i3,a,a,a)') ch10,&
184 &   ' memana: BUG -',ch10,&
185 &   '  The biggest chain is number',biggest,' while iprcel==20.',ch10,&
186 &   '  This is not allowed.'
187    call wrtout(std_out,message,'COLL')
188  end if
189 
190  write(message, '(a,f11.3,a)' ) &
191 & 'P This job should need less than                 ',&
192 & mbbiggest+tol10,' Mbytes of memory. '
193  call wrtout(std_out,message,'COLL')
194  call wrtout(iout,message,'COLL')
195 
196  if(prtvol>=10)then
197    if(biggest==1)write(message,'(a)')'P Max. in main chain + fourwf.f '
198    if(biggest==2)write(message,'(a)')'P Max. in main chain + nonlop.f + opernl.f '
199    if(biggest==3)write(message,'(a)')'P Max. in XC chain '
200    if(biggest==4)write(message,'(a)')'P Max. in mkrho chain '
201    if(biggest==5)write(message,'(a)')'P Max. in fourdp chain '
202    if(biggest==6)write(message,'(a)')'P Max. in parallel k-point chain '
203    if(biggest==7)write(message,'(a)')'P Max. in newvtr chain '
204    if(biggest==8)write(message,'(a)')'P Max. in suscep chain '
205    if(biggest==9)write(message,'(a)')'P Max. in dielmt chain '
206    if(biggest==10)write(message,'(a)')'P Max. in tddft chain '
207    call wrtout(iout,message,'COLL')
208 
209    write(message, '(a,i13,a,f11.3,a)' )&
210 &   'P',nint(cintmpw(biggest)),' blocks of mpw  integer numbers, for',&
211 &   mbintmpw(biggest)+tol10,' Mbytes. '
212    call wrtout(iout,message,'COLL')
213    write(message, '(a,i13,a,f11.3,a)' )&
214 &   'P',nint(cdpmpw(biggest)),' blocks of mpw  real(dp)  numbers, for',&
215 &   mbdpmpw(biggest)+tol10,' Mbytes. '
216    call wrtout(iout,message,'COLL')
217    if (nfft==nfftf) then
218      if(mbintfft(biggest)+mbintfftf(biggest)>0.001)then
219        write(message, '(a,i13,a,f11.3,a)' )&
220 &       'P',nint(cintfft(biggest)+cintfftf(biggest)),' blocks of nfft integer numbers, for',&
221 &       mbintfft(biggest)+mbintfftf(biggest)+tol10,' Mbytes. '
222        call wrtout(iout,message,'COLL')
223      end if
224      write(message, '(a,i13,a,f11.3,a)' )&
225 &     'P',nint(cdpfft(biggest)+cdpfftf(biggest)),' blocks of nfft real(dp)  numbers, for',&
226 &     mbdpfft(biggest)+mbdpfftf(biggest)+tol10,' Mbytes. '
227      call wrtout(iout,message,'COLL')
228    else
229      if(mbintfftf(biggest)>0.001)then
230        write(message, '(a,i13,a,f11.3,a)' )&
231 &       'P',nint(cintfftf(biggest)),' blocks of nfft (fine grid) integer numbers, for',&
232 &       mbintfftf(biggest)+tol10,' Mbytes. '
233        call wrtout(iout,message,'COLL')
234      end if
235      write(message, '(a,i13,a,f11.3,a)' )&
236 &     'P',nint(cdpfftf(biggest)),' blocks of nfft (fine grid) real(dp)  numbers, for',&
237 &     mbdpfftf(biggest)+tol10,' Mbytes. '
238      call wrtout(iout,message,'COLL')
239      if(mbintfft(biggest)>0.001)then
240        write(message, '(a,i13,a,f11.3,a)' )&
241 &       'P',nint(cintfft(biggest)),' blocks of nfft (coarse grid) integer numbers, for',&
242 &       mbintfft(biggest)+tol10,' Mbytes. '
243        call wrtout(iout,message,'COLL')
244      end if
245      write(message, '(a,i13,a,f11.3,a)' )&
246 &     'P',nint(cdpfft(biggest)),' blocks of nfft (coarse grid) real(dp)  numbers, for',&
247 &     mbdpfft(biggest)+tol10,' Mbytes. '
248      call wrtout(iout,message,'COLL')
249    end if
250    if(mbintadd(biggest)>0.001)then
251      write(message, '(a,13x,a,f11.3,a)' )&
252 &     'P',' Additional     integer numbers, for',mbintadd(biggest)+tol10,' Mbytes. '
253      call wrtout(iout,message,'COLL')
254    end if
255    write(message, '(a,13x,a,f11.3,a)' )&
256 &   'P',' Additional     real(dp)  numbers, for',mbdpadd(biggest)+tol10,' Mbytes. '
257    call wrtout(iout,message,'COLL')
258    write(message, '(a,13x,a,f11.3,a)' )&
259 &   'P',' With residue estimated to be       ',mbother(biggest)+tol10,' Mbytes. '
260    call wrtout(iout,message,'COLL')
261    write(message, '(a)' )'P'
262    call wrtout(iout,message,'COLL')
263    write(message, '(a)' )&
264 &   'P Comparison of the memory needs of different chains'
265    call wrtout(iout,message,'COLL')
266 
267    write(message, '(a,f11.3,a)' )&
268 &   'P Main chain + fourwf.f           ',mbtot(1)+tol10,' Mbytes. '
269    call wrtout(iout,message,'COLL')
270    write(message, '(a,f11.3,a)' )&
271 &   'P Main chain + nonlop.f + opernl.f',mbtot(2)+tol10,' Mbytes. '
272    call wrtout(iout,message,'COLL')
273 
274 !  The next chains are not defined in the RF case.
275    if(nchain>2)then
276      write(message, '(a,f11.3,a)' )&
277 &     'P XC chain                        ',mbtot(3)+tol10,' Mbytes. '
278      call wrtout(iout,message,'COLL')
279      write(message, '(a,f11.3,a)' )&
280 &     'P mkrho chain                     ',mbtot(4)+tol10,' Mbytes. '
281      call wrtout(iout,message,'COLL')
282      write(message, '(a,f11.3,a)' )&
283 &     'P fourdp chain                    ',mbtot(5)+tol10,' Mbytes. '
284      call wrtout(iout,message,'COLL')
285      if(xmpi_paral==1)then
286        write(message, '(a,f11.3,a)' )&
287 &       '- parallel k-point chain          ',mbtot(6)+tol10,' Mbytes. '
288        call wrtout(iout,message,'COLL')
289      end if
290      write(message, '(a,f11.3,a)' )&
291 &     'P newvtr chain                    ',mbtot(7)+tol10,' Mbytes. '
292      call wrtout(iout,message,'COLL')
293      if(modulo(iprcel,100)>=20.and.modulo(iprcel,100)<70)then
294        write(message, '(a,f11.3,a)' )&
295 &       'P suscep chain                    ',mbtot(8)+tol10,' Mbytes. '
296        call wrtout(iout,message,'COLL')
297        write(message, '(a,f11.3,a)' )&
298 &       'P dielmt chain                    ',mbtot(9)+tol10,' Mbytes. '
299        call wrtout(iout,message,'COLL')
300      end if
301      if(iscf==-1)then
302        write(message, '(a,f11.3,a)' )&
303 &       'P tddft  chain                    ',mbtot(10)+tol10,' Mbytes. '
304      end if
305    end if ! nchain>2
306 
307  end if
308 
309 !--------------------------------------------------------------------
310 
311  write(message, '(a)' ) &
312 & '  Rough estimation (10% accuracy) of disk space for files :'
313  call wrtout(iout,message,'COLL')
314  call wrtout(std_out,message,'COLL')
315 
316  write(message, '(a,f11.3,a,a,f11.3,a)' ) &
317 & '_ WF disk file :',mbdiskwf+tol10,' Mbytes ;',&
318 & ' DEN or POT disk file :',mbdiskpd+tol10,' Mbytes.'
319  call wrtout(iout,message,'COLL')
320  call wrtout(std_out,message,'COLL')
321 
322  if(mffmem==0 .and. iscf>0)then
323    if(iscf==1)then
324      write(message, '(a,a,a)' )&
325 &     '  mffmem==0, iscf==1 => use of 1 FFT temporary disk file,',ch10,&
326 &     '                       5 times bigger than a DEN file.'
327    else if(iscf==2.or.iscf==12)then
328      write(message, '(a,a,a)' )&
329 &     '  mffmem==0, iscf==2 => use of 1 FFT temporary disk file,',ch10,&
330 &     '                       3 times bigger than a DEN file.'
331    else if(iscf==3.or.iscf==13)then
332      write(message, '(a,a,a)' )&
333 &     '  mffmem==0, iscf==3 => use of 1 FFT temporary disk file,',ch10,&
334 &     '                       4 times bigger than a DEN file.'
335    else if(iscf==4.or.iscf==14)then
336      write(message, '(a,a,a)' )&
337 &     '  mffmem==0, iscf==4 => use of 1 FFT temporary disk file,',ch10,&
338 &     '                       6 times bigger than a DEN file.'
339    else if(iscf==5)then
340      write(message, '(a,a,a)' )&
341 &     '  mffmem==0, iscf==5 => use of 1 FFT temporary disk file,',ch10,&
342 &     '                       10 times bigger than a DEN file.'
343    else if(iscf==6)then
344      write(message, '(a,a,a)' )&
345 &     '  mffmem==0, iscf==6 => use of 1 FFT temporary disk file,',ch10,&
346 &     '                       10 times bigger than a DEN file.'
347    else if(iscf==7.or.iscf==17)then
348      write(message, '(a,a,a)' )&
349 &     '  mffmem==0, iscf==7 => use of 1 FFT temporary disk file,',ch10,&
350 &     '                       (2+2*npulayit) times bigger than a DEN file.'
351    end if
352    call wrtout(iout,message,'COLL')
353    call wrtout(std_out,message,'COLL')
354  end if
355 
356 !Temporary message - estimation of PAW specific data has to be done...
357 !Have to add the usepaw argument to use this.
358 !if (usepaw==1) then
359 !write(message,'(5a)') '  WARNING: You are using PAW formalism;',ch10,&
360 !&       '           Above estimations do not take PAW',ch10,&
361 !&       '           specific data into account !'
362 !call wrtout(iout,message,'COLL')
363 !call wrtout(std_out,message,'COLL')
364 !end if
365 
366  write(message,'(80a,a)') ('=',mu=1,80),ch10
367  call wrtout(iout,message,'COLL')
368  call wrtout(std_out,message,'COLL')
369 
370 !--------------------------------------------------------------------
371 !Here, each processor must test its memory, so use
372 !the PERS mode for error messages, followed by synchronisation
373 
374  mbbigarr=max(mbf_fftgr,mbcg,mbgylm)
375  if(mbbigarr==mbcg) then
376    write(message, '(a,f12.4,a)' ) &
377 &   ' Biggest array : cg(disk), with',mbcg+tol10,' MBytes.'
378  else if (mbbigarr==mbf_fftgr) then
379    write(message, '(a,f12.4,a)' ) &
380 &   ' Biggest array : f_fftgr(disk), with',mbf_fftgr+tol10,' MBytes.'
381  else if (mbbigarr==mbgylm)then
382    write(message, '(a,f12.4,a)' ) &
383 &   ' Biggest array : pawfgrtab%gylm(gr), with',mbgylm+tol10,' MBytes.'
384  end if
385  call wrtout(std_out,message,'COLL')
386 
387 !if (mpi_enreg%my_nimage>1) then
388 !write(message, '(a,f12.4,a)' ) &
389 !&   ' These estimations take the distribution over replicas (images) of the cell into account.'
390 !call wrtout(std_out,message,'COLL')
391 !end if
392 
393  quit=0
394 
395  if(option>=1)then
396 
397 !  Test the ability to allocate the biggest array
398    nquarter_mbytes=4.0_dp*mbbigarr+1.0_dp
399    ABI_STAT_ALLOCATE(bigarray,(32*1024,nquarter_mbytes), ier)
400    if(ier/=0)then
401      write(message,'(a,f11.3,a,a,a,a,a,a,a)')&
402 &     'Test failed to allocate an array of',mbbigarr,' Mbytes',ch10,&
403 &     'It is not worth to continue ',ch10,&
404 &     'Action: modify input variable to fit the available memory,',ch10,&
405 &     'increase limit on maximal array size or set mem_test to 0 to disable this test.'
406      call wrtout(std_out,message,'PERS')
407      if(option==1)then
408        MSG_ERROR_CLASS(message, "MemanaError")
409      else
410        MSG_WARNING(message)
411        quit=1
412      end if
413    end if
414    if(quit==0)then
415      write(message,'(a,f11.3,a)')&
416 &     ' memana : allocated an array of',mbbigarr+tol10,' Mbytes, for testing purposes. '
417      call wrtout(std_out,message,'COLL')
418    end if
419    if(allocated(bigarray)) then
420      ABI_DEALLOCATE(bigarray)
421    end if
422 
423 !  Test the ability to allocate the needed total memory : use 8 segments,
424 !  hoping that the maximal segment size is not so much smaller than the
425 !  total memory
426    nquarter_mbytes=0.5_dp*mbbiggest+1.0_dp
427    ABI_STAT_ALLOCATE(bigarray1,(32*1024,nquarter_mbytes), ier1)
428    ABI_STAT_ALLOCATE(bigarray2,(32*1024,nquarter_mbytes), ier2)
429    ABI_STAT_ALLOCATE(bigarray3,(32*1024,nquarter_mbytes), ier3)
430    ABI_STAT_ALLOCATE(bigarray4,(32*1024,nquarter_mbytes), ier4)
431    ABI_STAT_ALLOCATE(bigarray5,(32*1024,nquarter_mbytes), ier5)
432    ABI_STAT_ALLOCATE(bigarray6,(32*1024,nquarter_mbytes), ier6)
433    ABI_STAT_ALLOCATE(bigarray7,(32*1024,nquarter_mbytes), ier7)
434    ABI_STAT_ALLOCATE(bigarray8,(32*1024,nquarter_mbytes), ier8)
435 
436    if(ier1/=0 .or. ier2/=0 .or. ier3/=0 .or. ier4/=0 .or.&
437 &   ier5/=0 .or. ier6/=0 .or. ier7/=0 .or. ier8/=0) then
438      write(message,'(a,f11.3,a,a,a,a,a,a,a)')&
439 &     'Test failed to allocate ',mbbiggest,' Mbytes',ch10,&
440 &     'It is not worth to continue ',ch10,&
441 &     'Action: modify input variables or submission parameters to fit the available memory,',ch10,&
442 &     'increase limit on available memory or set mem_test to 0 to disable this test.'
443      if(option==1)then
444        MSG_ERROR_CLASS(message, "MemanaError")
445      else
446        MSG_WARNING(message)
447        quit=1
448      end if
449    end if
450 
451    if(quit==0)then
452      write(message,'(a,f11.3,a,a,a)')&
453 &     ' memana: allocated ',mbbiggest,'Mbytes, for testing purposes. ',ch10,&
454 &     ' The job will continue.'
455      call wrtout(std_out,message,'COLL')
456    end if
457    if(allocated(bigarray1)) then
458      ABI_DEALLOCATE(bigarray1)
459    end if
460    if(allocated(bigarray2)) then
461      ABI_DEALLOCATE(bigarray2)
462    end if
463    if(allocated(bigarray3)) then
464      ABI_DEALLOCATE(bigarray3)
465    end if
466    if(allocated(bigarray4)) then
467      ABI_DEALLOCATE(bigarray4)
468    end if
469    if(allocated(bigarray5)) then
470      ABI_DEALLOCATE(bigarray5)
471    end if
472    if(allocated(bigarray6)) then
473      ABI_DEALLOCATE(bigarray6)
474    end if
475    if(allocated(bigarray7)) then
476      ABI_DEALLOCATE(bigarray7)
477    end if
478    if(allocated(bigarray8)) then
479      ABI_DEALLOCATE(bigarray8)
480    end if
481 
482  end if
483 
484 !--------------------------------------------------------------------
485 
486  if(option==2 .and. quit==1 )then
487 
488 !  Estimation of the available memory
489 !  
490 !  A quarter of Mbyte is 256*1024/8 real(dp) numbers,
491 !  that is 32*1024 dp numbers.
492 !  One begins with the allocation of 4 Mbytes. If successful,
493 !  one increases that number, until the allocation is not successfull
494 !  any more. Unfortunately, on a P6 with the pghpf compiler, the
495 !  allocate instruction generate a core dump, instead of returning
496 !  an error code, so that this part of code has been made optional.
497 
498    nquarter_mbytes=16
499    nmbytes=nquarter_mbytes/4.0_dp
500 
501 !  With an increase ratio of 1.25_dp (see below), ii=5 leads to 9 MB,
502 !  ii=10 leads to 28 MB, ii=15 leads to 85 MB, ii=18 leads to 165 MB,
503 !  ii=30 is over 2 GB
504    do ii=1,30
505      ABI_STAT_ALLOCATE(bigarray,(32*1024,nquarter_mbytes), ier)
506      if(ier/=0)then
507        write(message,'(a,i0,a)')' memana : failed to allocate ',nmbytes,' Mbytes'
508        call wrtout(std_out,message,'PERS')
509        exit
510      end if
511      write(message,'(a,i0,a)')' memana : succeeded to allocate ',nmbytes,' Mbytes'
512      call wrtout(std_out,message,'PERS')
513 !    Here really test the space
514 !    do kk=1,nquarter_mbytes
515 !    do jj=1,32*1024,37
516 !    bigarray(jj,kk)=0.0_dp
517 !    end do
518 !    write(std_out,*)' memana : wrote ',kk,' quarter of mbytes'
519 !    end do
520      ABI_DEALLOCATE(bigarray)
521      nquarter_mbytes=dble(nquarter_mbytes)*1.25_dp
522      nmbytes=nquarter_mbytes/4.0_dp
523    end do
524    if(allocated(bigarray)) then
525      ABI_DEALLOCATE(bigarray)
526    end if
527 
528    MSG_ERROR_CLASS("in memana with option==2 .and. quit==1", "MemanaError")
529  end if !  End the test of the available memory
530 
531 !--------------------------------------------------------------------
532 
533  ABI_DEALLOCATE(cdpfftf)
534  ABI_DEALLOCATE(cdpfft)
535  ABI_DEALLOCATE(cdpmpw)
536  ABI_DEALLOCATE(cintfftf)
537  ABI_DEALLOCATE(cintfft)
538  ABI_DEALLOCATE(cintmpw)
539  ABI_DEALLOCATE(cdpadd)
540  ABI_DEALLOCATE(cintadd)
541  ABI_DEALLOCATE(mbdpadd)
542  ABI_DEALLOCATE(mbdpfftf)
543  ABI_DEALLOCATE(mbdpfft)
544  ABI_DEALLOCATE(mbdpmpw)
545  ABI_DEALLOCATE(mbintadd)
546  ABI_DEALLOCATE(mbintfftf)
547  ABI_DEALLOCATE(mbintfft)
548  ABI_DEALLOCATE(mbintmpw)
549  ABI_DEALLOCATE(mbother)
550  ABI_DEALLOCATE(mbtot)
551 
552 end subroutine memana