TABLE OF CONTENTS


ABINIT/m_dtfil [ Modules ]

[ Top ] [ Modules ]

NAME

  m_dtfil

FUNCTION

   object and procedures dealing with input/output filenames

COPYRIGHT

  Copyright (C) 2008-2018 ABINIT group (XG, MT)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

PARENTS

CHILDREN

SOURCE

21 #if defined HAVE_CONFIG_H
22 #include "config.h"
23 #endif
24 
25 #include "abi_common.h"
26 
27 module m_dtfil
28 
29  use defs_basis
30  use defs_abitypes
31  use m_abicore
32  use m_errors
33  use m_xmpi
34  use m_build_info
35 
36  use m_clib,         only : clib_rename
37  use m_fstrings,     only : int2char4
38  use m_io_tools,     only : open_file
39  use m_libpaw_tools, only : libpaw_log_flag_set
40 
41  implicit none
42 
43  private

ABINIT/status [ Functions ]

[ Top ] [ Functions ]

NAME

 status

FUNCTION

 Routine for description of the status of the calculation
 Eventually open the status file, write different information,
 and close the file. The output rate and shift are governed by istat

INPUTS

  counter=value of the loop counter at that level
  file (optional argument)=name of the status file
  istat=gives the rate or shift of output. The status file will be opened
     and written only once every "istatr" calls.
     This variable is saved at the first call (just after the first
     call to invars0. The shift "istatshft" is saved at the second call.
     In subsequent calls, istat has no meaning.
  level=number of the level of the calling subroutine (see the description later)
  routine=string of 14 characters indicating the status inside the level

OUTPUT

  (only writing)

NOTES

 Warning : The string "routine" can have any size but
 it is truncated to a size of 14.
 because of the behaviour of some compilers, the string
 "routine" should always have 14 characters in the calling subroutine

PARENTS

      abinit,dfpt_accrho,dfpt_looppert,dfpt_rhofermi,dfpt_scfcv,dfpt_vtowfk
      dfpt_wfkfermi,dfptnl_loop,dfptnl_mv,dfptnl_resp,driver,gstate,gstateimg
      m_ab7_invars_f90,mover,nonlinear,respfn,scfcv,vtorhotf

CHILDREN

      timab

SOURCE

1495 subroutine status(counter,filstat,istat,level,routine)
1496 
1497  use m_time,       only : timab
1498 
1499 !This section has been created automatically by the script Abilint (TD).
1500 !Do not modify the following lines by hand.
1501 #undef ABI_FUNC
1502 #define ABI_FUNC 'status'
1503 !End of the abilint section
1504 
1505  implicit none
1506 
1507 !Arguments ------------------------------------
1508 !scalars
1509  integer,intent(in) :: counter,istat,level
1510  character(len=*),intent(in) :: routine
1511  character(len=*),intent(in) :: filstat
1512 
1513 !Local variables-------------------------------
1514 !scalars
1515  integer,parameter :: mcounter=2,mlevel=120
1516  integer,save :: output_rate=1,shift_rate=1,statnu=0
1517  integer :: ilevel,temp_unit
1518  character(len=12) :: headwr
1519  character(len=500) :: message
1520 !arrays
1521  integer,save :: active(mlevel),actual_counter(mlevel,mcounter)
1522  integer,save :: ncounter(mlevel)
1523  integer,save :: list_level(29)=&
1524 &  (/1,2,3,100,101,102,110,111,112,113,114,10,11,12,13,14,15,16,17,18,20,30,31,40,41,50,51,52,90/)
1525  real(dp) :: tsec(2)
1526  character(len=20),save :: nm_levels(mlevel),nm_routine(mlevel)
1527  character(len=12),save :: nm_counter(mlevel,mcounter)
1528 
1529 !***********************************************************************
1530 
1531  if (.not.do_write_status .or. output_rate==0) return
1532 
1533  call timab(73,1,tsec)
1534 
1535 !Note : all processors have their own file, so no special
1536 !attention must be paid to the parallel case.
1537 !Initialisation
1538  if(statnu==0)then
1539    nm_routine(:)='                  '
1540    active(:)=0
1541    actual_counter(:,:)=0
1542 
1543 !  List of names for each level
1544 !  Numbers from 1 to 9 are for abinit and driver
1545 !  Numbers from 100 to 120 are for optdriver=0 routines (GS)
1546 !  Numbers between 10 and 19 are for optdriver=1 routines (RF)
1547 !  Numbers between 20 and 29 are for optdriver=2 routines (suscep)
1548 !  Numbers between 30 and 39 are for optdriver=3 routines (screening)
1549 !  Numbers between 40 and 49 are for optdriver=4 routines (sigma)
1550 !  Numbers between 50 and 59 are for optdriver=5 routines (nonlinear)
1551 !  When you add a level number, or modify one, do not forget to change list_level
1552 
1553    nm_levels(1)   ='abinit              '
1554    ncounter(1)=0
1555    nm_counter(1,1)='            '
1556 
1557    nm_levels(2)   ='driver              '
1558    ncounter(2)=1
1559    nm_counter(2,1)='jdtset     ='
1560 
1561    nm_levels(3)   ='ab7_invars_load     '
1562    ncounter(3)=0
1563    nm_counter(3,1)='            '
1564 
1565 !  Optdriver=0
1566    nm_levels(100)   ='gstateimg           '
1567    ncounter(100)=2
1568    nm_counter(100,1)='idynimage  ='
1569    nm_counter(100,2)='itimimage  ='
1570 
1571    nm_levels(101)   ='gstate              '
1572    ncounter(101)=1
1573    nm_counter(101,1)='itime      ='
1574 
1575    nm_levels(102)   ='mover               '
1576    ncounter(102)=2
1577    nm_counter(102,1)='icycle     ='
1578    nm_counter(102,2)='itime      ='
1579 
1580    nm_levels(110)   ='scfcv               '
1581    ncounter(110)=1
1582    nm_counter(110,1)='istep      ='
1583 
1584    nm_levels(111)   ='vtorho(tf)          '
1585    ncounter(111)=2
1586    nm_counter(111,1)='isppol     ='
1587    nm_counter(111,2)='ikpt       ='
1588 
1589    nm_levels(112)   ='vtowfk              '
1590    ncounter(112)=2
1591    nm_counter(112,1)='inonsc     ='
1592    nm_counter(112,2)='iband      ='
1593 
1594    nm_levels(113)   ='cgwf                '
1595    ncounter(113)=1
1596    nm_counter(113,1)='iline      ='
1597 
1598    nm_levels(114)   ='getghc              '
1599    ncounter(114)=0
1600    nm_counter(114,1)='            '
1601 
1602 
1603 !  Optdriver=1
1604    nm_levels(10)   ='respfn              '
1605    ncounter(10)=0
1606    nm_counter(10,1)='            '
1607 
1608    nm_levels(11)   ='dfpt_looppert       '
1609    ncounter(11)=1
1610    nm_counter(11,1)='respcase   ='
1611 
1612    nm_levels(12)   ='dfpt_scfcv          '
1613    ncounter(12)=1
1614    nm_counter(12,1)='istep      ='
1615 
1616    nm_levels(13)   ='dfpt_vtorho         '
1617    ncounter(13)=2
1618    nm_counter(13,1)='isppol     ='
1619    nm_counter(13,2)='ikpt       ='
1620 
1621    nm_levels(14)   ='dfpt_vtowfk         '
1622    ncounter(14)=2
1623    nm_counter(14,1)='inonsc     ='
1624    nm_counter(14,2)='iband      ='
1625 
1626    nm_levels(15)   ='dfpt_cgwf           '
1627    ncounter(15)=1
1628    nm_counter(15,1)='iline      ='
1629 
1630    nm_levels(16)   ='getgh1c             '
1631    ncounter(16)=0
1632    nm_counter(16,1)='            '
1633 
1634    nm_levels(17)   ='dfpt_rhofermi       '
1635    ncounter(17)=2
1636    nm_counter(17,1)='isppol     ='
1637    nm_counter(17,2)='ikpt       ='
1638 
1639    nm_levels(18)   ='dfpt_wfkfermi       '
1640    ncounter(18)=2
1641    nm_counter(18,1)='inonsc     ='
1642    nm_counter(18,2)='iband      ='
1643 
1644 
1645 !  Optdriver=2
1646    nm_levels(20)   ='suscep              '
1647    ncounter(20)=0
1648    nm_counter(20,1)='            '
1649 
1650 
1651 !  Optdriver=3
1652    nm_levels(30)   ='screening           '
1653    ncounter(30)=1
1654    nm_counter(30,1)='iqpt       ='
1655 
1656 !  Optdriver=4
1657    nm_levels(40)   ='sigma               '
1658    ncounter(40)=1
1659    nm_counter(40,1)='ikpt_gw    ='
1660 
1661 !  Optdriver=5
1662    nm_levels(50)   ='nonlinear           '
1663    ncounter(50)=0
1664    nm_counter(50,1)='            '
1665 
1666    nm_levels(51)   ='dfptnl_loop         '
1667    ncounter(51)=2
1668    nm_counter(51,1)='pert1case  ='
1669    nm_counter(51,2)='pert3case  ='
1670 
1671    nm_levels(52)   ='mv_/dfptnl_resp     '
1672    ncounter(52)=2
1673    nm_counter(52,2)='ikpt       ='
1674 
1675 !  Optdriver=9
1676    nm_levels(90)   ='bethe_salpether     '
1677    ncounter(90)=0
1678    nm_counter(90,1)='            '
1679 
1680    if(istat<0)then
1681      write(message, '(a,i7,a,a,a,a,a)' )&
1682 &     'the value of the input variable istatr is',istat,' ,',ch10,&
1683 &     'while it must be a positive number.',ch10,&
1684 &     'Action : change istatr in your input file.'
1685      MSG_ERROR(message)
1686    end if
1687    output_rate=istat
1688  end if
1689 
1690 !The true input variable "shift_rate" is only available at the second call
1691  if(statnu==1)then
1692    if(istat<0 .or. istat>output_rate)then
1693      write(message, '(a,i7,3a,i7,2a)' )&
1694 &     'the value of the input variable istatshft is',istat,' ,',ch10,&
1695 &     'while it must be a positive number smaller or equal to istatr=',output_rate,ch10,&
1696 &     'Action: change istatshft in your input file.'
1697      MSG_ERROR(message)
1698    end if
1699    shift_rate=istat
1700    if(shift_rate==output_rate)shift_rate=0
1701 
1702 !  At this second call, also feed information that the abinit routine called ab7_invars_load
1703    write(unit=nm_routine(1),fmt='(a20)') 'call ab7_invars_load'
1704    active(1)=1
1705  end if
1706 
1707 !Check the value of level
1708  if( minval(abs(list_level(:)-level)) /= 0)then
1709    write(message, '(a,i5,3a)' )&
1710 &   '  The value of level in the calling routine is',level,' ,',ch10,&
1711 &   '  which is not an allowed value.'
1712    MSG_BUG(message)
1713  end if
1714 
1715 !Assign the info about the actual routine
1716  write(unit=nm_routine(level),fmt='(a20)') trim(adjustl(routine))
1717  if(trim(adjustl(nm_routine(level)))=='exit')then
1718 !  The value 2 will be changed to 0 at the end of the routine.
1719    active(level)=2
1720  else if(trim(adjustl(nm_routine(level)))=='')then
1721    active(level)=0
1722  else
1723    active(level)=1
1724  end if
1725 
1726 !Assign the info about the actual counter
1727  if(counter>=0)then
1728    if(ncounter(level)==1)then
1729      actual_counter(level,1)=counter
1730    else if(ncounter(level)==2)then
1731      actual_counter(level,2)=counter/100
1732 !    The counter number 1 does not allow more than 99 passes
1733      actual_counter(level,1)=counter-(counter/100)*100
1734    end if
1735  end if
1736 
1737 !============================================================
1738 
1739 !After treatment of present information, output of the status
1740  statnu=statnu+1
1741 
1742 !DEBUG
1743 ! write(std_out,*)' status : statnu, output_rate, shift_rate=',statnu,output_rate, shift_rate
1744 ! write(std_out,*)'level,routine=',level,routine
1745 ! write(std_out,*)'active(level)=',active(level)
1746 ! write(std_out,*)'counter,actual_counter(level,1:2)=',counter,actual_counter(level,1:2)
1747 ! write(std_out,*)'List of active levels :'
1748 ! do ilevel=1,mlevel
1749 !   if(active(ilevel)/=0)write(std_out,*)' Active level number=',ilevel
1750 ! end do
1751 !ENDDEBUG
1752 
1753  if(statnu > 2 )then
1754    if( mod(statnu,output_rate)==shift_rate )then
1755 
1756      if (open_file(filstat,message,newunit=temp_unit,form='formatted',status='unknown') /= 0) then
1757        MSG_ERROR(message)
1758      end if
1759 
1760      rewind temp_unit
1761      write(temp_unit,*)
1762 
1763      headwr='(a,i4,a,i6 )'
1764      if(statnu>=100000)   headwr='(a,i4,a,i9 )'
1765      if(output_rate>=1000)headwr='(a,i6,a,i6 )'
1766      if(statnu>=100000 .and. output_rate>=1000)   headwr='(a,i6,a,i9 )'
1767      if(statnu>=100000000)headwr='(a,i6,a,i12)'
1768      write(temp_unit,headwr)' Status file, with repetition rate',output_rate,', status number',statnu
1769      write(temp_unit,*)
1770 
1771 !    Treat every level one after the other
1772      do ilevel=1,mlevel
1773 !      This level must be activated in order to have a corresponding output
1774        if(active(ilevel)==1 .or. active(ilevel)==2)then
1775 
1776          write(temp_unit,'(4a)')'  Level ',nm_levels(ilevel),' : ',adjustl(nm_routine(ilevel))
1777 
1778 !        Is there a counter for this level ?
1779          if(ncounter(ilevel)>=1)then
1780 
1781            if(actual_counter(ilevel,1)>0)then
1782              write(temp_unit,'(a,a,i5)')'  ',nm_counter(ilevel,1),actual_counter(ilevel,1)
1783            end if
1784            if(ncounter(ilevel)==2)then
1785              if(actual_counter(ilevel,2)>0)then
1786                write(temp_unit,'(a,a,i5)')'  ',nm_counter(ilevel,2),actual_counter(ilevel,2)
1787              end if
1788            end if
1789 
1790          end if
1791        end if ! End of the check on activation of the level
1792      end do ! End of the loop on the levels
1793 
1794      close (temp_unit)
1795    end if ! End of the repetition rate check
1796  end if ! statnu > 2
1797 
1798  if (active(level)==2)then
1799    active(level)=0
1800    nm_routine(level)='              '
1801  end if
1802 
1803  call timab(73,2,tsec)
1804 
1805 end subroutine status

m_dtfil/dtfil_init [ Functions ]

[ Top ] [ m_dtfil ] [ Functions ]

NAME

 dtfil_init

FUNCTION

 Initialize most of the dtfil structured variable
 (what is left should be initialized inside the itimimage,
 iimage and itime loops).

INPUTS

 dtset=<type datasets_type>contain all input variables for the current dataset
 filnam(5)=character strings giving file names
 filstat=character strings giving name of status file
 idtset=number of the dataset
 image_index= -optional argument- index of image to be used when appending
             "_IMGxxx" string to file names. To be used only when an algorithm
             using images of the cell is activated
 jdtset_(0:ndtset)=actual index of the datasets
 mpi_enreg=information about MPI parallelization
 ndtset=number of datasets

OUTPUT

 dtfil=<type datafiles_type>infos about file names, file unit numbers
  (part of which were initialized previously)

NOTES

 The array filnam is used for the name of input and output files,
 and roots for generic input, output or temporary files.
 Pseudopotential file names are set in pspini and pspatm,
 using another name. The name filstat will be needed beyond gstate to check
 the appearance of the "exit" flag, to make a hasty exit, as well as
 in order to output the status of the computation.

PARENTS

      driver,gstateimg

CHILDREN

      appdig,int2char4,mkfilename

SOURCE

101 subroutine dtfil_init(dtfil,dtset,filnam,filstat,idtset,jdtset_,mpi_enreg,ndtset,&
102 &                      image_index) ! optional argument
103 
104 
105 !This section has been created automatically by the script Abilint (TD).
106 !Do not modify the following lines by hand.
107 #undef ABI_FUNC
108 #define ABI_FUNC 'dtfil_init'
109 !End of the abilint section
110 
111  implicit none
112 
113 !Arguments ------------------------------------
114 !scalars
115  integer, intent(in) :: idtset,ndtset
116  integer, optional, intent(in) :: image_index
117  character(len=fnlen),intent(in) :: filstat
118  type(MPI_type),intent(in) :: mpi_enreg
119  type(datafiles_type),intent(inout) :: dtfil !vz_i
120 !arrays
121  integer :: jdtset_(0:ndtset)
122  character(len=fnlen),intent(in) :: filnam(5)
123  type(dataset_type),intent(in) :: dtset
124 
125 !Local variables-------------------------------
126 !scalars
127 ! Define input and output unit numbers (do not forget, unit 5 and 6 are standard input and output)
128 ! Also, unit number 21, 22 and 23 are used in dfpt_nstdy, for the 3 dot wavefunctions.
129 ! Unit 50,51,52 and 53 are used in dfpt_looppert (for ipert=natom+2, ipert=natom+10 and ipert=natom+11).
130 ! Others unit numbers will be used in the case of the variational and 2n+1 expressions.
131 ! In defs_basis, one defines :
132 !   std_in=5, ab_in=5, std_out=6, ab_out=7, tmp_unit=9, tmp_unit2=10
133  integer,parameter :: unchi0=42,unddb=16,unddk=50,undkdk=54,undkde=55,unkg1=19,unkg=17,unkgq=18
134  integer,parameter :: unpaw=26,unpaw1=27,unpawq=28,unpos=30
135  integer,parameter :: unwff1=1,unwff2=2,unwff3=8,unwffgs=3,unwfkq=4,unwft1=11
136  integer,parameter :: unwft2=12,unwft3=15,unwftgs=13,unwftkq=14,unylm=24,unylm1=25
137  integer,parameter :: unkss=40,unscr=41,unqps=43
138  integer :: ii,iimage,ireadden,ireadkden,ireadwf,ixx,jdtset,will_read
139  character(len=10) :: appen,tag
140  character(len=9) :: stringvar
141  character(len=15) :: stringfile
142  character(len=500) :: message
143  character(len=fnlen) :: filsus,filddbsin,fildens1in,fildensin,filpawdensin,filkdensin,filqps,filscr,fil_efmas
144  character(len=fnlen) :: fnamewff1,fnamewffddk,fnamewffdelfd,fnamewffdkdk,fnamewffdkde,fnamewffk,fnamewffq
145  character(len=fnlen) :: filbseig,filfft,filhaydock,fil_bsreso,fil_bscoup
146  character(len=fnlen) :: filwfkfine
147  character(len=fnlen) :: filnam_ds(5)
148  character(len=fnlen) :: tmpfil(14)
149  integer :: idtmpfil(14)
150 
151 !******************************************************************
152 
153  DBG_ENTER("COLL")
154 
155  iimage=0;if (present(image_index)) iimage=image_index
156 
157  dtfil%unchi0 =unchi0
158  dtfil%unddb  =unddb
159  dtfil%unddk  =unddk
160  dtfil%undkde =undkde
161  dtfil%undkdk =undkdk
162  dtfil%unkg   =unkg
163  dtfil%unkgq  =unkgq
164  dtfil%unkg1  =unkg1
165  dtfil%unkss  =unkss
166  dtfil%unqps  =unqps
167  dtfil%unscr  =unscr
168  dtfil%unwff1 =unwff1
169  dtfil%unwff2 =unwff2
170  dtfil%unwff3 =unwff3
171  dtfil%unwffgs=unwffgs
172  dtfil%unwffkq=unwfkq
173  dtfil%unwft1 =unwft1
174  dtfil%unwft2 =unwft2
175  dtfil%unwft3 =unwft3
176  dtfil%unwftgs=unwftgs
177  dtfil%unwftkq=unwftkq
178  dtfil%unylm  =unylm
179  dtfil%unylm1 =unylm1
180  dtfil%unpaw  =unpaw
181  dtfil%unpaw1 =unpaw1
182  dtfil%unpawq =unpawq
183  dtfil%unpos  =unpos
184  filnam_ds(1:5)=filnam(1:5)
185  jdtset=dtset%jdtset
186 
187 !If multi dataset mode, special treatment of filenames 3 and 4 (density and
188 !wavefunctions input and output, as well as other output files)
189  if(ndtset>0)then
190    call appdig(jdtset,'',appen)
191    filnam_ds(3)=trim(filnam(3))//'_DS'//trim(appen)
192    filnam_ds(4)=trim(filnam(4))//'_DS'//trim(appen)
193  end if
194 
195 !If multi image mode (nimage>1), special treatment of filenames 4 and 5
196  if(iimage>0)then
197    call appdig(iimage,'',appen)
198    filnam_ds(4)=trim(filnam_ds(4))//'_IMG'//trim(appen)
199    filnam_ds(5)=trim(filnam_ds(5))//'_IMG'//trim(appen)
200  end if
201 
202 !According to getwfk and irdwfk, build _WFK file name, referred as fnamewffk
203  if (iimage>0.and.dtfil%getwfk_from_image/=0) then
204    if (dtfil%getwfk_from_image==-1) then
205      call appdig(iimage,'',appen)
206    else
207      call appdig(dtfil%getwfk_from_image,'',appen)
208    end if
209    stringfile='_IMG'//trim(appen)//'_WFK'
210  else
211    stringfile='_WFK'
212  end if
213  stringvar='wfk'
214  call mkfilename(filnam,fnamewffk,dtset%getwfk,idtset,dtset%irdwfk,jdtset_,&
215 & ndtset,stringfile,stringvar,will_read)
216 
217  if(dtset%optdriver/=RUNL_RESPFN)ireadwf=will_read
218  if(ndtset/=0 .and. dtset%optdriver==RUNL_RESPFN .and. will_read==0)then
219    write(message, '(5a,i3,3a,i3,a,i3,3a)' )&
220 &   'At least one of the input variables irdwfk and getwfk ',ch10,&
221 &   'must refer to a valid _WFK file, in the response function',ch10,&
222 &   'case, while for idtset=',idtset,',',ch10,&
223 &   'they are irdwfk=',dtset%irdwfk,', and getwfk=',dtset%getwfk,'.',ch10,&
224 &   'Action: correct irdwfk or getwfk in your input file.'
225    MSG_ERROR(message)
226  end if
227 
228 !Treatment of the other get wavefunction variable, if response function case or nonlinear case
229  if ( ANY(dtset%optdriver == (/RUNL_RESPFN, RUNL_NONLINEAR, RUNL_EPH/)) ) then
230 
231 !  According to getwfq and irdwfq, build _WFQ file name, referred as fnamewffq
232    stringfile='_WFQ' ; stringvar='wfq'
233    call mkfilename(filnam,fnamewffq,dtset%getwfq,idtset,dtset%irdwfq,jdtset_,&
234 &   ndtset,stringfile,stringvar,will_read)
235 !  If fnamewffq is not initialized thanks to getwfq or irdwfq, use fnamewffk
236    if(will_read==0)fnamewffq=fnamewffk
237 
238 !  According to get1wf and ird1wf, build _1WF file name, referred as fnamewff1
239    stringfile='_1WF' ; stringvar='1wf'
240    call mkfilename(filnam,fnamewff1,dtset%get1wf,idtset,dtset%ird1wf,jdtset_,&
241 &   ndtset,stringfile,stringvar,will_read)
242    ireadwf=will_read
243 
244 !  According to getddk and irdddk, build _1WF file name, referred as fnamewffddk
245    stringfile='_1WF' ; stringvar='ddk'
246    call mkfilename(filnam,fnamewffddk,dtset%getddk,idtset,dtset%irdddk,jdtset_,&
247 &   ndtset,stringfile,stringvar,will_read)
248 
249 !  According to getdelfd, build _1WF file name, referred as fnamewffdelfd
250    stringfile='_1WF' ; stringvar='delfd'
251    call mkfilename(filnam,fnamewffdelfd,dtset%getdelfd,idtset,0,jdtset_,&
252 &   ndtset,stringfile,stringvar,will_read)
253 
254 !  According to getdkdk, build _1WF file name, referred as fnamewffdkdk
255    stringfile='_1WF' ; stringvar='dkdk'
256    call mkfilename(filnam,fnamewffdkdk,dtset%getdkdk,idtset,0,jdtset_,&
257 &   ndtset,stringfile,stringvar,will_read)
258 
259 !  According to getdkde, build _1WF file name, referred as fnamewffdkde
260    stringfile='_1WF' ; stringvar='dkde'
261    call mkfilename(filnam,fnamewffdkde,dtset%getdkde,idtset,0,jdtset_,&
262 &   ndtset,stringfile,stringvar,will_read)
263 
264  end if
265 
266 !-------------------------------------------------------------------------------------------
267 !Build name of files from dtfil%filnam_ds(3)
268 
269 !SP :According to getddb, build _DDB file name, referred as filddbsin
270  stringfile='_DDB'
271  stringvar='ddb'
272  call mkfilename(filnam,filddbsin,dtset%getddb,idtset,dtset%irdddb,jdtset_,&
273 & ndtset,stringfile,stringvar,will_read)
274 
275 !According to getden, build _DEN file name, referred as fildensin
276 !A default is available if getden is 0
277  if (iimage>0.and.dtfil%getden_from_image/=0) then
278    if (dtfil%getden_from_image==-1) then
279      call appdig(iimage,'',appen)
280    else
281      call appdig(dtfil%getden_from_image,'',appen)
282    end if
283    stringfile='_IMG'//trim(appen)//'_DEN'
284  else
285    stringfile='_DEN'
286  end if
287  stringvar='den'
288  call mkfilename(filnam,fildensin,dtset%getden,idtset,dtset%irdden,jdtset_,&
289 & ndtset,stringfile,stringvar,will_read)
290 
291  if(will_read==0)fildensin=trim(filnam_ds(3))//'_DEN'
292  ireadden=will_read
293 
294  if ((dtset%optdriver==RUNL_GWLS.or.dtset%optdriver==RUNL_GSTATE) &
295 & .and.dtset%iscf<0) ireadden=1
296 !if (optdriver==RUNL_GSTATE.and.ireadwf/=0) ireadden=0
297 
298 !According to getpawden, build _PAWDEN file name, referred as filpawdensin
299 !A default is available if getden is 0
300  if (iimage>0.and.dtfil%getpawden_from_image/=0) then
301    if (dtfil%getpawden_from_image==-1) then
302      call appdig(iimage,'',appen)
303    else
304      call appdig(dtfil%getpawden_from_image,'',appen)
305    end if
306    stringfile='_IMG'//trim(appen)//'_PAWDEN'
307  else
308    stringfile='_PAWDEN'
309  end if
310  stringvar='pawden'
311  call mkfilename(filnam,filpawdensin,dtset%getpawden,idtset,dtset%irdden,jdtset_,&
312 & ndtset,stringfile,stringvar,will_read)
313  if(will_read==0)filpawdensin=trim(filnam_ds(3))//'_PAWDEN'
314 
315 !According to getden and usekden, build _KDEN file name, referred as filkdensin
316 !A default is available if getden is 0
317  if(dtset%usekden==1)then
318    if (iimage>0.and.dtfil%getden_from_image/=0) then
319      if (dtfil%getden_from_image==-1) then
320        call appdig(iimage,'',appen)
321      else
322        call appdig(dtfil%getden_from_image,'',appen)
323      end if
324      stringfile='_IMG'//trim(appen)//'_KDEN'
325    else
326      stringfile='_KDEN'
327    end if
328    stringvar='kden'
329    call mkfilename(filnam,filkdensin,dtset%getden,idtset,dtset%irdden,jdtset_,&
330 &   ndtset,stringfile,stringvar,will_read)
331    if(will_read==0)filkdensin=trim(filnam_ds(3))//'_KDEN'
332    ireadkden=will_read
333    if ((dtset%optdriver==RUNL_GSTATE.or.dtset%optdriver==RUNL_GWLS).and.dtset%iscf<0) ireadkden=1
334  else
335    ireadkden=0
336  end if
337 
338 !According to get1den, build _DEN file name, referred as fildens1in
339 !A default is available if get1den is 0
340  stringfile='_DEN' ; stringvar='1den'
341  call mkfilename(filnam,fildens1in,dtset%get1den,idtset,dtset%ird1den,jdtset_,&
342 & ndtset,stringfile,stringvar,will_read)
343  if(will_read==0)fildens1in=trim(filnam_ds(3))//'_DEN'
344 
345 !According to getefmas and irdefmas, build _EFMAS file name, referred as fil_efmas
346 !A default is available if getefmas is 0
347  stringfile='_EFMAS.nc' ; stringvar='efmas'
348  call mkfilename(filnam,fil_efmas,dtset%getefmas,idtset,dtset%irdefmas,jdtset_,&
349 & ndtset,stringfile,stringvar,will_read)
350  if(will_read==0)fil_efmas=trim(filnam_ds(3))//'_EFMAS.nc'
351 
352 !According to getscr and irdscr, build _SCR file name, referred as filscr
353 !A default is available if getscr is 0
354  stringfile='_SCR' ; stringvar='scr'
355  call mkfilename(filnam,filscr,dtset%getscr,idtset,dtset%irdscr,jdtset_,&
356 & ndtset,stringfile,stringvar,will_read)
357  if(will_read==0)filscr=trim(filnam_ds(3))//'_SCR'
358 
359 !According to getsuscep and irdsuscep, build _SUS file name, referred as filsus
360 !A default is available if getsuscep is 0
361  stringfile='_SUS' ; stringvar='sus'
362  call mkfilename(filnam,filsus,dtset%getsuscep,idtset,dtset%irdsuscep,jdtset_,&
363 & ndtset,stringfile,stringvar,will_read)
364  if(will_read==0)filsus=TRIM(filnam_ds(3))//'_SUS'
365 
366 !According to getqps and irdqps, build _QPS file name, referred as filqps
367 !A default is available if getqps is 0
368  stringfile='_QPS' ; stringvar='qps'
369  call mkfilename(filnam,filqps,dtset%getqps,idtset,dtset%irdqps,jdtset_,&
370 & ndtset,stringfile,stringvar,will_read)
371  if(will_read==0)filqps=trim(filnam_ds(3))//'_QPS'
372 
373 !According to getbseig and irdbseig, build _BSEIG file name, referred as filbseig
374 !A default is available if getbseig is 0
375  stringfile='_BSEIG' ; stringvar='bseig'
376  call mkfilename(filnam,filbseig,dtset%getbseig,idtset,dtset%irdbseig,jdtset_,ndtset,stringfile,stringvar,will_read)
377  if(will_read==0)filbseig=trim(filnam_ds(3))//'_BSEIG'
378 
379 !According to gethaydock and irdhaydock, build _HAYD file name, referred as filhaydock.
380 !A default is available if gethaydock is 0
381  stringfile='_HAYDR_SAVE' ; stringvar='haydock'
382  call mkfilename(filnam,filhaydock,dtset%gethaydock,idtset,dtset%irdhaydock,jdtset_,ndtset,stringfile,stringvar,will_read)
383  if(will_read==0)filhaydock=trim(filnam_ds(3))//'_HAYDR_SAVE'
384 
385 !According to getbsr and irdbsr, build _BSR file name, referred as fil_bsreso
386 !A default is available if getbsr is 0
387  stringfile='_BSR' ; stringvar='bsreso'
388  call mkfilename(filnam,fil_bsreso,dtset%getbsreso,idtset,dtset%irdbsreso,jdtset_,ndtset,stringfile,stringvar,will_read)
389  if(will_read==0) fil_bsreso=trim(filnam_ds(3))//'_BSR'
390 
391 !According to getbsc and irdbsc, build _BSC file name, referred as fil_bscoup
392 !A default is available if getbsc is 0
393  stringfile='_BSC' ; stringvar='bscoup'
394  call mkfilename(filnam,fil_bscoup,dtset%getbscoup,idtset,dtset%irdbscoup,jdtset_,ndtset,stringfile,stringvar,will_read)
395  if(will_read==0)fil_bscoup=trim(filnam_ds(3))//'_BSC'
396 
397 !According to getwfkfine and irdwfkfine, build _WFK file name, referred as filwfkfine
398 !A default is avaible if getwfkfine is 0
399  stringfile='_WFK' ; stringvar='wfkfine'
400  call mkfilename(filnam,filwfkfine,dtset%getwfkfine,idtset,dtset%irdwfkfine,jdtset_,&
401 & ndtset,stringfile,stringvar,will_read)
402  if(will_read==0)filwfkfine=trim(filnam_ds(3))//'_WFK'
403 
404  dtfil%ireadden      =ireadden
405  dtfil%ireadkden     =ireadkden
406  dtfil%ireadwf       =ireadwf
407  dtfil%filnam_ds(1:5)=filnam_ds(1:5)
408 
409  dtfil%fnameabi_bsham_reso=fil_bsreso
410  dtfil%fnameabi_bsham_coup=fil_bscoup
411  dtfil%fnameabi_bseig=filbseig
412  dtfil%fnameabi_haydock=filhaydock
413  dtfil%fnameabi_sus  =filsus
414  dtfil%fnameabi_qps  =filqps
415  dtfil%fnameabi_scr  =filscr
416  dtfil%fnameabi_efmas=fil_efmas
417  dtfil%filddbsin     =filddbsin
418  dtfil%fildensin     =fildensin
419  dtfil%fildens1in    =fildens1in
420  dtfil%filkdensin    =filkdensin
421  dtfil%filpawdensin  =filpawdensin
422  dtfil%fnameabi_wfkfine = filwfkfine
423  dtfil%filstat       =filstat
424  dtfil%fnamewffk     =fnamewffk
425  dtfil%fnamewffq     =fnamewffq
426  dtfil%fnamewffddk   =fnamewffddk
427  dtfil%fnamewffdelfd =fnamewffdelfd
428  dtfil%fnamewffdkdk  =fnamewffdkdk
429  dtfil%fnamewffdkde  =fnamewffdkde
430  dtfil%fnamewff1     =fnamewff1
431 
432  dtfil%fnameabi_hes=trim(dtfil%filnam_ds(3))//'_HES'
433  dtfil%fnameabi_phfrq=trim(dtfil%filnam_ds(3))//'_PHFRQ'
434  dtfil%fnameabi_phvec=trim(dtfil%filnam_ds(3))//'_PHVEC'
435 
436 !-------------------------------------------------------------------------------------------
437 !Build name of files from dtfil%filnam_ds(4)
438 
439  dtfil%fnameabo_ddb=trim(dtfil%filnam_ds(4))//'_DDB'
440  dtfil%fnameabo_den=trim(dtfil%filnam_ds(4))//'_DEN'
441  dtfil%fnameabo_dos=trim(dtfil%filnam_ds(4))//'_DOS'
442  dtfil%fnameabo_eelf=trim(dtfil%filnam_ds(4))//'_EELF'
443  dtfil%fnameabo_eig=trim(dtfil%filnam_ds(4))//'_EIG'
444  dtfil%fnameabo_eigi2d=trim(dtfil%filnam_ds(4))//'_EIGI2D'
445  dtfil%fnameabo_eigr2d=trim(dtfil%filnam_ds(4))//'_EIGR2D'
446  dtfil%fnameabo_em1=trim(dtfil%filnam_ds(4))//'_EM1'
447  dtfil%fnameabo_em1_lf=trim(dtfil%filnam_ds(4))//'_EM1_LF'
448  dtfil%fnameabo_em1_nlf=trim(dtfil%filnam_ds(4))//'_EM1_NLF'
449  dtfil%fnameabo_fan=trim(dtfil%filnam_ds(4))//'_FAN'
450  dtfil%fnameabo_gkk=trim(dtfil%filnam_ds(4))//'_GKK'
451  dtfil%fnameabo_gw=trim(dtfil%filnam_ds(4))//'_GW' ! TODO change name
452  dtfil%fnameabo_gwdiag=trim(dtfil%filnam_ds(4))//'_GWDIAG'
453  dtfil%fnameabo_kss=trim(dtfil%filnam_ds(4))//'_KSS'
454  dtfil%fnameabo_moldyn=trim(dtfil%filnam_ds(4))//'_MOLDYN'
455  dtfil%fnameabo_pot=trim(dtfil%filnam_ds(4))//'_POT'
456  dtfil%fnameabo_qps=trim(dtfil%filnam_ds(4))//'_QPS'
457  dtfil%fnameabo_qp_den=trim(dtfil%filnam_ds(4))//'_QP_DEN'
458  dtfil%fnameabo_qp_pawden=trim(dtfil%filnam_ds(4))//'_QP_PAWDEN'
459  dtfil%fnameabo_qp_dos=trim(dtfil%filnam_ds(4))//'_QP_DOS'
460  dtfil%fnameabo_qp_eig=trim(dtfil%filnam_ds(4))//'_QP_DB.nc' ! TODO change name
461  dtfil%fnameabo_rpa=trim(dtfil%filnam_ds(4))//'_RPA'
462  dtfil%fnameabo_scr=trim(dtfil%filnam_ds(4))//'_SCR'
463  dtfil%fnameabo_sgm=trim(dtfil%filnam_ds(4))//'_SGM'
464  dtfil%fnameabo_sgr=trim(dtfil%filnam_ds(4))//'_SGR'
465  dtfil%fnameabo_sig=trim(dtfil%filnam_ds(4))//'_SIG'
466  dtfil%fnameabo_spcur=trim(dtfil%filnam_ds(4))//'_SPCUR'
467  dtfil%fnameabo_sus=trim(dtfil%filnam_ds(4))//'_SUS'
468  dtfil%fnameabo_vha=trim(dtfil%filnam_ds(4))//'_VHA'
469  dtfil%fnameabo_vpsp=trim(dtfil%filnam_ds(4))//'_VPSP'
470  dtfil%fnameabo_vso=trim(dtfil%filnam_ds(4))//'_VSO'
471  dtfil%fnameabo_vxc=trim(dtfil%filnam_ds(4))//'_VXC'
472  dtfil%fnameabo_wan=trim(dtfil%filnam_ds(4))//'_WAN'
473  dtfil%fnameabo_wfk=trim(dtfil%filnam_ds(4))//'_WFK'
474  dtfil%fnameabo_wfq=trim(dtfil%filnam_ds(4))//'_WFQ'
475  dtfil%fnameabo_w90=trim(dtfil%filnam_ds(4))//'_w90'
476  dtfil%fnameabo_1wf=trim(dtfil%filnam_ds(4))//'_1WF'
477  dtfil%fnameabo_nlcc_derivs=trim(dtfil%filnam_ds(4))//'_nlcc_derivs_'
478  dtfil%fnameabo_pspdata=trim(dtfil%filnam_ds(4))//'_pspdata_'
479 
480 !-------------------------------------------------------------------------------------------
481 !Build name of files from dtfil%filnam_ds(5)
482 
483  dtfil%fnametmp_eig=trim(dtfil%filnam_ds(5))//'_EIG'
484  dtfil%fnametmp_1wf1_eig=trim(dtfil%filnam_ds(5))//'_1WF1_EIG' ! This appendix should be changed !
485  dtfil%fnametmp_kgs=trim(dtfil%filnam_ds(5))//'_KGS'
486  dtfil%fnametmp_sustr=trim(dtfil%filnam_ds(5))//'_SUSTR'
487  dtfil%fnametmp_tdexcit=trim(dtfil%filnam_ds(5))//'_TDEXCIT'
488  dtfil%fnametmp_tdwf=trim(dtfil%filnam_ds(5))//'_TDWF'
489  dtfil%fnametmp_cg=trim(dtfil%filnam_ds(5))//'_cg'
490  dtfil%fnametmp_cprj=trim(dtfil%filnam_ds(5))//'_cprj'
491 
492 !'_WF1' -> dtfil%unwft1
493 !'_WF2' -> dtfil%unwft2
494 !'_KG' ->  dtfil%unkg
495 !'_DUM' -> tmp_unit (real dummy name)
496 !'_YLM' -> dtfil%unylm
497 !'_PAW' -> dtfil%unpaw
498 
499  tmpfil(1)=trim(dtfil%filnam_ds(5))//'_WF1'  ! tmpfil(1)
500  tmpfil(2)=trim(dtfil%filnam_ds(5))//'_WF2'  ! tmpfil(2)
501 
502  tmpfil(3)=trim(dtfil%filnam_ds(5))//'_KG'   ! tmpfil(3)
503  tmpfil(4)=trim(dtfil%filnam_ds(5))//'_DUM'  ! tmpfil(4)
504  tmpfil(5)=' '  ! to avoid Valgrind complain.
505  tmpfil(6)=trim(dtfil%filnam_ds(5))//'_YLM'  ! tmpfil(6)
506  tmpfil(7)=trim(dtfil%filnam_ds(5))//'_PAW'  ! tmpfil(7)
507 
508  if(xmpi_paral==1)then ! parallel case : the index of the processor must be appended
509    call int2char4(mpi_enreg%me,tag)
510    ABI_CHECK((tag(1:1)/='#'),'Bug: string length too short!')
511    ixx=1
512    if (xmpi_mpiio == 1 .and. dtset%iomode == IO_MODE_MPI ) ixx=3
513    do ii=ixx,7
514      tmpfil(ii)=trim(tmpfil(ii))//'_P-'//trim(tag)
515    end do
516  end if
517 
518  dtfil%fnametmp_wf1=trim(tmpfil(1))
519  dtfil%fnametmp_wf2=trim(tmpfil(2))
520 
521  dtfil%fnametmp_kg =trim(tmpfil(3))
522  dtfil%fnametmp_dum=trim(tmpfil(4))
523  dtfil%fnametmp_ylm=trim(tmpfil(6))
524  dtfil%fnametmp_paw=trim(tmpfil(7))
525 
526 !Create names for the temporary files based on dtfil%filnam_ds(5)
527 !by appending adequate string.
528 !'_1WF1' -> dtfil%unwft1
529 !'_1WF2' -> dtfil%unwft2
530 !'_KG'   -> dtfil%unkg
531 !'_KGQ'  -> dtfil%unkgq (not used for the time being)
532 !'_KG1'  -> dtfil%unkg1
533 !'_DUM'  -> tmp_unit (real dummy name)
534 !'_WFGS' -> dtfil%unwftgs
535 !'_WFKQ' -> dtfil%unwftkq
536 !'_YLM'  -> dtfil%unylm
537 !'_YLM1' -> dtfil%unylm1
538 !'_PAW'  -> dtfil%unpaw
539 !'_PAW1' -> dtfil%unpaw1
540 !'_PAWQ' -> dtfil%unpawq
541  tmpfil(1) =trim(dtfil%filnam_ds(5))//'_1WF1'
542  tmpfil(2) =trim(dtfil%filnam_ds(5))//'_1WF2'
543  tmpfil(3) =trim(dtfil%filnam_ds(5))//'_KG'
544  tmpfil(4) =trim(dtfil%filnam_ds(5))//'_KGQ'
545  tmpfil(5) =trim(dtfil%filnam_ds(5))//'_KG1'
546  tmpfil(6) =trim(dtfil%filnam_ds(5))//'_DUM'
547  tmpfil(7) =trim(dtfil%filnam_ds(5))//'_WFGS'
548  tmpfil(8) =trim(dtfil%filnam_ds(5))//'_WFKQ'
549  tmpfil(9) =' ' ! for Valgrind, to avoid uninitialized
550  tmpfil(10)=trim(dtfil%filnam_ds(5))//'_YLM'
551  tmpfil(11)=trim(dtfil%filnam_ds(5))//'_YLM1'
552  tmpfil(12)=trim(dtfil%filnam_ds(5))//'_PAW'
553  tmpfil(13)=trim(dtfil%filnam_ds(5))//'_PAW1'
554  tmpfil(14)=trim(dtfil%filnam_ds(5))//'_PAWQ'
555 
556  if(xmpi_paral==1) then
557    idtmpfil(:)=0
558    do ii=1,14
559      idtmpfil(ii)=ii
560    end do
561    if (xmpi_mpiio==1.and.dtset%iomode==IO_MODE_MPI)then
562      idtmpfil(1)=0              !_1wf1
563      idtmpfil(2)=0              ! s1wf2
564      idtmpfil(7)=0              !  WFGS
565      idtmpfil(8)=0              !  WFKQ
566    end if
567    call int2char4(mpi_enreg%me,tag)
568    ABI_CHECK((tag(1:1)/='#'),'Bug: string length too short!')
569    do ii=1,14
570      if(idtmpfil(ii) /= 0) tmpfil(ii)=trim(tmpfil(ii))//'_P-'//trim(tag)
571    end do
572  end if
573 
574  dtfil%fnametmp_1wf1=trim(tmpfil(1))
575  dtfil%fnametmp_1wf2=trim(tmpfil(2))
576  dtfil%fnametmp_kg  =trim(tmpfil(3))
577  dtfil%fnametmp_kgq =trim(tmpfil(4))
578  dtfil%fnametmp_kg1 =trim(tmpfil(5))
579  dtfil%fnametmp_dum =trim(tmpfil(6))
580  dtfil%fnametmp_wfgs=trim(tmpfil(7))
581  dtfil%fnametmp_wfkq=trim(tmpfil(8))
582  dtfil%fnametmp_ylm =trim(tmpfil(10))
583  dtfil%fnametmp_ylm1=trim(tmpfil(11))
584  dtfil%fnametmp_paw =trim(tmpfil(12))
585  dtfil%fnametmp_paw1=trim(tmpfil(13))
586  dtfil%fnametmp_pawq=trim(tmpfil(14))
587 
588 !Prepare the name of the _FFT file
589  filfft=trim(dtfil%filnam_ds(5))//'_FFT'
590 !There is a definite problem in the treatment of // by CPP ...
591  if(xmpi_paral==1 .or. mpi_enreg%paral_kgb==1)then
592    call int2char4(mpi_enreg%me,tag)
593    ABI_CHECK((tag(1:1)/='#'),'Bug: string length too short!')
594    filfft=trim(filfft)//'_P-'//trim(tag)
595  end if
596  dtfil%fnametmp_fft=filfft
597 
598 !These keywords are only used in algorithms using images of the cell
599  if (iimage==0) then
600    dtfil%getwfk_from_image   =0
601    dtfil%getden_from_image   =0
602    dtfil%getpawden_from_image=0
603  end if
604 
605  DBG_EXIT("COLL")
606 
607 end subroutine dtfil_init

m_dtfil/dtfil_init_img [ Functions ]

[ Top ] [ m_dtfil ] [ Functions ]

NAME

 dtfil_init_img

FUNCTION

 Initialize few scalars in the dtfil structured variable
 when an alogrithm using image of the cell is selected.
 (initialize index of images from which read files)

INPUTS

  dtset=<type datasets_type>=input variables for the current dataset
  dtsets(0:ndtset_alloc)=<type datasets_type>=input variables for all datasets
  idtset=number of the dataset
  jdtset(0:ndtset)=actual index of the datasets
  ndtset=number of datasets
  ndtset_alloc=number of datasets, corrected for allocation of at least one data set

OUTPUT

SIDE EFFECTS

 dtfil=<type datafiles_type>= only getxxx_from_image flags are modified

PARENTS

      driver

CHILDREN

SOURCE

846 subroutine dtfil_init_img(dtfil,dtset,dtsets,idtset,jdtset,ndtset,ndtset_alloc)
847 
848 
849 !This section has been created automatically by the script Abilint (TD).
850 !Do not modify the following lines by hand.
851 #undef ABI_FUNC
852 #define ABI_FUNC 'dtfil_init_img'
853 !End of the abilint section
854 
855  implicit none
856 
857 !Arguments ------------------------------------
858 !scalars
859  integer, intent(in) :: idtset,ndtset,ndtset_alloc
860  type(datafiles_type),intent(out) :: dtfil
861  type(dataset_type),intent(in) :: dtset
862 !arrays
863  integer :: jdtset(0:ndtset)
864  type(dataset_type),intent(in) :: dtsets(0:ndtset_alloc)
865 
866 !Local variables -------------------------
867 !scalars
868  integer :: iget
869 !arrays
870 
871 ! *********************************************************************
872 
873  DBG_ENTER("COLL")
874 
875 !Default values
876  dtfil%getwfk_from_image   =0 ! Get standard WFK from previous dataset
877  dtfil%getden_from_image   =0 ! Get standard DEN from previous dataset
878  dtfil%getpawden_from_image=0 ! Get standard PAWDEN from previous dataset
879 
880  if (dtset%optdriver==RUNL_GSTATE.and.dtset%nimage>1) then
881 
882 !  Define getwfk_from_image
883    if (dtset%getwfk/=0.or.dtset%irdwfk/=0) then
884      iget=-1
885      if(dtset%getwfk<0) iget=jdtset(idtset+dtset%getwfk)
886      if(dtset%getwfk>0) iget=dtset%getwfk
887      if(dtset%irdwfk>0) iget=0
888      if (iget>=0) then
889        if (iget==0.or.dtsets(iget)%nimage==dtset%nimage) then
890          dtfil%getwfk_from_image=-1     ! Get WFK from the same image of previous dataset
891        else if (dtsets(iget)%nimage>1) then
892          dtfil%getwfk_from_image=1      ! Get WFK from the first image of previous dataset
893        end if
894      end if
895    end if
896 
897 !  Define getden_from_image
898    if (dtset%getden/=0.or.dtset%irdden/=0) then
899      iget=-1
900      if(dtset%getden<0) iget=jdtset(idtset+dtset%getden)
901      if(dtset%getden>0) iget=dtset%getden
902      if(dtset%irdden>0) iget=0
903      if (iget>=0) then
904        if (iget==0.or.dtsets(iget)%nimage==dtset%nimage) then
905          dtfil%getden_from_image=-1     ! Get DEN from the same image of previous dataset
906        else if (dtsets(iget)%nimage>1) then
907          dtfil%getden_from_image=1      ! Get DEN from the first image of previous dataset
908        end if
909      end if
910    end if
911 
912 !  Define getpawden_from_image
913    if (dtset%getpawden/=0.or.dtset%irdpawden/=0) then
914      iget=-1
915      if(dtset%getpawden<0) iget=jdtset(idtset+dtset%getpawden)
916      if(dtset%getpawden>0) iget=dtset%getpawden
917      if(dtset%irdpawden>0) iget=0
918      if (iget>=0) then
919        if (iget==0.or.dtsets(iget)%nimage==dtset%nimage) then
920          dtfil%getpawden_from_image=-1     ! Get PAWDEN from the same image of previous dataset
921        else if (dtsets(iget)%nimage>1) then
922          dtfil%getpawden_from_image=1      ! Get PAWDEN from the first image of previous dataset
923        end if
924      end if
925    end if
926  end if
927 
928  DBG_EXIT("COLL")
929 
930 end subroutine dtfil_init_img

m_dtfil/dtfil_init_time [ Functions ]

[ Top ] [ m_dtfil ] [ Functions ]

NAME

 dtfil_init_time

FUNCTION

 Inside the itimimage, iimage and itime loops (this is only needed for optdriver=0),
 initialize the remaining parts of dtfil.

INPUTS

 iapp=indicates the eventual suffix to be appended to the generic output root
         if 0 : no suffix to be appended (called directly from gstate)
         if positive : append "_TIM//iapp" (called from move or brdmin)
         if -1 : append "_TIM0" (called from brdmin)
         if -2, -3, -4, -5: append "_TIMA", ... ,"_TIMD", (called from move)

OUTPUT

SIDE EFFECTS

 dtfil=<type datafiles_type>infos about file names, file unit numbers
  (part of which were initialized previously)

PARENTS

      gstate,mover

CHILDREN

      fappnd

SOURCE

640 subroutine dtfil_init_time(dtfil,iapp)
641 
642 
643 !This section has been created automatically by the script Abilint (TD).
644 !Do not modify the following lines by hand.
645 #undef ABI_FUNC
646 #define ABI_FUNC 'dtfil_init_time'
647 !End of the abilint section
648 
649  implicit none
650 
651 !Arguments ------------------------------------
652 !scalars
653  integer, intent(in) :: iapp
654  type(datafiles_type),intent(inout) :: dtfil
655 
656 !Local variables-------------------------------
657 !scalars
658  character(len=fnlen) :: filapp,filprot
659 
660 !******************************************************************
661 
662  DBG_ENTER("COLL")
663 
664 !--------------------------------------------------------
665 !Names based on dtfil%filnam_ds(4)+iapp
666 
667 !Prepare the name of the auxiliary files DOS, EIG...
668  call fappnd(filapp,dtfil%filnam_ds(4),iapp)
669  dtfil%fnameabo_app=trim(filapp)
670  dtfil%fnameabo_app_atmden_core=trim(filapp)//'_ATMDEN_CORE'
671  dtfil%fnameabo_app_atmden_val=trim(filapp)//'_ATMDEN_VAL'
672  dtfil%fnameabo_app_atmden_full=trim(filapp)//'_ATMDEN_FULL'
673  dtfil%fnameabo_app_n_tilde=trim(filapp)//'_N_TILDE'
674  dtfil%fnameabo_app_n_one=trim(filapp)//'_N_ONE'
675  dtfil%fnameabo_app_nt_one=trim(filapp)//'_NT_ONE'
676  dtfil%fnameabo_app_bxsf=trim(filapp)//'_BXSF'
677  dtfil%fnameabo_app_cif=trim(filapp)//'.cif'
678  dtfil%fnameabo_app_den=trim(filapp)//'_DEN'
679  dtfil%fnameabo_app_dos=trim(filapp)//'_DOS'
680  dtfil%fnameabo_app_eig=trim(filapp)//'_EIG'
681  dtfil%fnameabo_app_elf=trim(filapp)//'_ELF'
682  dtfil%fnameabo_app_elf_down=trim(filapp)//'_ELF_DOWN'
683  dtfil%fnameabo_app_elf_up=trim(filapp)//'_ELF_UP'
684  dtfil%fnameabo_app_fatbands=trim(filapp)//'_FATBANDS'
685  dtfil%fnameabo_app_gden1=trim(filapp)//'_GDEN1'
686  dtfil%fnameabo_app_gden2=trim(filapp)//'_GDEN2'
687  dtfil%fnameabo_app_gden3=trim(filapp)//'_GDEN3'
688  dtfil%fnameabo_app_geo=trim(filapp)//'_GEO'
689  dtfil%fnameabo_app_kden=trim(filapp)//'_KDEN'
690  dtfil%fnameabo_app_lden=trim(filapp)//'_LDEN'
691  dtfil%fnameabo_app_nesting=trim(filapp)//'_NEST'
692  dtfil%fnameabo_app_opt=trim(filapp)//'_OPT'
693  dtfil%fnameabo_app_opt2=trim(filapp)//'_OPT2'
694  dtfil%fnameabo_app_pawden=trim(filapp)//'_PAWDEN'
695  dtfil%fnameabo_app_pot=trim(filapp)//'_POT'
696  dtfil%fnameabo_app_stm=trim(filapp)//'_STM'
697  dtfil%fnameabo_app_vclmb=trim(filapp)//'_VCLMB'
698  dtfil%fnameabo_app_vha=trim(filapp)//'_VHA'
699  dtfil%fnameabo_app_vhxc=trim(filapp)//'_VHXC'
700  dtfil%fnameabo_app_vpsp=trim(filapp)//'_VPSP'
701  dtfil%fnameabo_app_vxc=trim(filapp)//'_VXC'
702  dtfil%fnameabo_app_wfk=trim(filapp)//'_WFK'
703  dtfil%fnameabo_app_vha_1dm=trim(filapp)//'_VHA_1DM'
704  dtfil%fnameabo_app_vclmb_1dm=trim(filapp)//'_VCLMB_1DM'
705  dtfil%fnameabo_app_1dm=trim(filapp)//'_1DM'
706 
707 !--------------------------------------------------------
708 !Names based on dtfil%filnam_ds(5)+iapp
709 
710 !Prepare the name of the auxiliary files for protection
711  call fappnd(filprot,dtfil%filnam_ds(5),iapp)
712  dtfil%fnametmp_app_den=trim(filprot)//'_DEN'
713  dtfil%fnametmp_app_kden=trim(filprot)//'_KDEN'
714 
715  DBG_EXIT("COLL")
716 
717 end subroutine dtfil_init_time

m_dtfil/fappnd [ Functions ]

[ Top ] [ m_dtfil ] [ Functions ]

NAME

 fappnd

FUNCTION

 Create the modified root name to be used for output of density, potential,
 and geometry files. See the description of the iapp input variable.

INPUTS

 filnam= generic output root name
 iapp=indicates the eventual suffix to be appended to the generic output root
      (the suffixe depends on the presence of the suff (optional) argument.
        if 0 : no suffix to be appended (called directly from gstate)
        if positive : append "_SUF//iapp" (called from move or brdmin)
        if -1 : append "_SUF0" (called from brdmin)
        if -2, -3, -4, -5: append "_SUFA", ... ,"_SUFD", (called from move)
      SUF can be TIM (default) or IMG
 suff= --optional argument--indicates the suffixe to be appended:
       SUF=TIM (default) or SUF=IMG or ...

OUTPUT

 filapp= filename with appended string

PARENTS

      dtfil_init_time

CHILDREN

SOURCE

751 subroutine fappnd(filapp,filnam,iapp,&
752 &                 suff) ! optional argument
753 
754 
755 !This section has been created automatically by the script Abilint (TD).
756 !Do not modify the following lines by hand.
757 #undef ABI_FUNC
758 #define ABI_FUNC 'fappnd'
759 !End of the abilint section
760 
761  implicit none
762 
763 !Arguments ------------------------------------
764 !scalars
765  integer,intent(in) :: iapp
766  character(len=fnlen),intent(in) :: filnam
767  character(len=fnlen),intent(out) :: filapp
768  character(len=3),optional,intent(in) :: suff
769 
770 !Local variables-------------------------------
771 !scalars
772  integer :: ndig
773  character(len=3) :: suffixe
774  character(len=8) :: nchar
775  character(len=500) :: msg
776 
777 ! *************************************************************************
778 
779  if(iapp==0)then
780    filapp=trim(filnam)
781  else
782    suffixe="TIM"
783    if (present(suff)) suffixe=trim(suff(1:3))
784    if(iapp>0)then
785 !    Create character string for filename. Determine the number of digits in iapp.
786      ndig=int(log10(dble(iapp)+0.5_dp))+1
787 !    Make integer format field of exact size (internal write)
788 !    for assumed nchar string of 8 characters
789      write(nchar, '(i8)' ) iapp
790      if (ndig>8) then
791        write(msg,'(5a,i0,2a,i0,2a)')&
792 &       'Requested file name extension has more than the allowed 8 digits.',ch10,&
793 &       'Action: resubmit the job with smaller value for ntime.',ch10,&
794 &       'Value computed here was ndig=',ndig,ch10,&
795 &       'iapp= ',iapp,' filnam=',TRIM(filnam)
796        MSG_ERROR(msg)
797      end if
798 !    Concatenate into character string, picking off exact number of digits
799 !    The potential or density label will be appended in ioarr
800      filapp=trim(filnam)//'_'//suffixe(1:3)//nchar(9-ndig:8)
801    else if(iapp==-1)then
802      filapp=trim(filnam)//'_'//suffixe(1:3)//'0'
803    else if(iapp==-2)then
804      filapp=trim(filnam)//'_'//suffixe(1:3)//'A'
805    else if(iapp==-3)then
806      filapp=trim(filnam)//'_'//suffixe(1:3)//'B'
807    else if(iapp==-4)then
808      filapp=trim(filnam)//'_'//suffixe(1:3)//'C'
809    else if(iapp==-5)then
810      filapp=trim(filnam)//'_'//suffixe(1:3)//'D'
811    end if
812  end if
813 
814 end subroutine fappnd

m_dtfil/iofn1 [ Functions ]

[ Top ] [ m_dtfil ] [ Functions ]

NAME

 iofn1

FUNCTION

 Begin by eventual redefinition of unit std_in and std_out
 Then, print greetings for interactive user.
 Next, Read filenames from unit std_in, AND check that new
 output file does not already exist.

INPUTS

  comm=MPI communicator.

OUTPUT

  character(len=fnlen) :: filnam(5)=character strings giving file names
  character(len=fnlen) :: filstat=character strings giving name of status file

NOTES

 If it does exist, isfile will create a new name to avoid overwriting the output file.
 Also create name of status file

 File names refer to following files, in order:
  (1) Formatted input file  (std_in)
  (2) Formatted output file (std_out)
  (3) Root name for generic input files (wavefunctions, potential, density ...)
  (4) Root name for generic output files (wavefunctions, potential, density,
                                          DOS, hessian ...)
  (5) Root name for generic temporary files (wftmp1,wftmp2,kgunit,status ...)

PARENTS

      abinit

CHILDREN

      abi_log_status_state,int2char4,isfile,libpaw_log_flag_set,xmpi_barrier
      xmpi_bcast

SOURCE

1243 subroutine iofn1(filnam,filstat,comm)
1244 
1245 
1246 !This section has been created automatically by the script Abilint (TD).
1247 !Do not modify the following lines by hand.
1248 #undef ABI_FUNC
1249 #define ABI_FUNC 'iofn1'
1250 !End of the abilint section
1251 
1252  implicit none
1253 
1254 !Arguments ------------------------------------
1255  integer,intent(in) :: comm
1256  character(len=fnlen), intent(out) :: filstat
1257  character(len=fnlen), intent(out) :: filnam(5)
1258 
1259 !Local variables-------------------------------
1260  character(len=1) :: blank
1261  integer,parameter :: master=0
1262  integer :: me,ios,nproc,ierr
1263  logical :: ex
1264  character(len=fnlen) :: fillog,tmpfil
1265  character(len=10) :: tag
1266  character(len=500) :: message,errmsg
1267 
1268 !*************************************************************************
1269 
1270  ! NOTE: In this routine it's very important to perform tests
1271  ! on possible IO errors (err=10, iomsg) because we are initializing the IO stuff
1272  ! It there's some problem with the hardware or some misconfiguration,
1273  ! it's very likely that the code will crash here and we should try to give useful error messages.
1274 
1275  blank = ' '; tmpfil = ''
1276 
1277 !Determine who I am in comm
1278  me = xmpi_comm_rank(comm)
1279  nproc = xmpi_comm_size(comm)
1280 
1281 !Define values of do_write_log and do_write_status parameters
1282 !if a _NOLOG file exists no LOG file and no STATUS file are created for each cpu core
1283 !if a _LOG file exists, a LOG file and a STATUS file are created for each cpu core
1284 !if the #_of_cpu_core>NPROC_NO_EXTRA_LOG OR presence of ABI_MAIN_LOG_FILE, LOG file is only created for master proc
1285 !if the #_of_cpu_core>NPROC_NO_EXTRA_STATUS OR presence of ABI_MAIN_LOG_FILE, STATUS file is only created for master proc
1286  inquire(file=ABI_NO_LOG_FILE,iostat=ios,exist=ex)
1287  if (ios/=0) ex=.false.
1288  if (ex) then
1289    do_write_log=.false. ; do_write_status=.false.
1290    call abi_log_status_state(new_do_write_log=.false.,new_do_write_status=.false.)
1291    call libpaw_log_flag_set(.false.)
1292  else
1293    inquire(file=ABI_ENFORCE_LOG_FILE,iostat=ios,exist=ex)
1294    if (ios/=0) ex=.false.
1295    if (ex) then
1296      do_write_log=.true. ; do_write_status=.true.
1297      call abi_log_status_state(new_do_write_log=.true.,new_do_write_status=.true.)
1298      call libpaw_log_flag_set(.true.)
1299    else
1300      inquire(file=ABI_MAIN_LOG_FILE,iostat=ios,exist=ex)
1301      if (ios/=0) ex=.false.
1302      if (ex.and.me/=0) then
1303        do_write_log=.false. ; do_write_status=.false.
1304        call abi_log_status_state(new_do_write_log=.false.,new_do_write_status=.false.)
1305        call libpaw_log_flag_set(.false.)
1306      else
1307        if (me/=0) then
1308          do_write_log= (nproc<NPROC_NO_EXTRA_LOG)
1309          call abi_log_status_state(new_do_write_log=(nproc<NPROC_NO_EXTRA_LOG))
1310          call libpaw_log_flag_set((nproc<NPROC_NO_EXTRA_LOG))
1311          do_write_status= (nproc<NPROC_NO_EXTRA_STATUS)
1312          call abi_log_status_state(new_do_write_status=(nproc<NPROC_NO_EXTRA_STATUS))
1313        end if
1314      end if
1315    end if
1316  end if
1317 
1318  if (me==master) then
1319 
1320 !  Eventually redefine standard input and standard output
1321 
1322    if (do_write_log) then
1323 #if defined READ_FROM_FILE
1324 !    Take care of the output file
1325      tmpfil(1:fnlen)=blank
1326      tmpfil(1:3)='log'
1327      call isfile(tmpfil,'new')
1328      close(std_out, err=10, iomsg=errmsg)
1329      if (open_file(tmpfil,message,unit=std_out,form='formatted',status='new',action="write") /= 0) then
1330        MSG_ERROR(message)
1331      end if
1332 #endif
1333    else
1334 !    Redirect standard output to null
1335      close(std_out, err=10, iomsg=errmsg)
1336      if (open_file(NULL_FILE,message,unit=std_out,action="write") /= 0) then
1337        MSG_ERROR(message)
1338      end if
1339    end if
1340 
1341 #if defined READ_FROM_FILE
1342 !  Now take care of the "files" file
1343    tmpfil(1:fnlen)=blank
1344    tmpfil(1:9)='ab.files'
1345    write(message, '(4a)' )&
1346 &   'Because of cpp option READ_FROM_FILE,',ch10,&
1347 &   'read file "ab.files" instead of standard input ' ,ch10
1348    MSG_COMMENT(message)
1349    call isfile(tmpfil,'old')
1350    close(std_in, err=10, iomsg=errmsg)
1351    if (open_file(tmpfil,message,unit=std_in,form='formatted',status='old',action="read") /= 0) then
1352      MSG_ERROR(message)
1353    end if
1354 #endif
1355 
1356 !  Print greetings for interactive user
1357    write(std_out,*,err=10,iomsg=errmsg)' ABINIT ',trim(abinit_version)
1358    write(std_out,*,err=10,iomsg=errmsg)' '
1359 
1360 !  Read name of input file (std_in):
1361    write(std_out,*,err=10,iomsg=errmsg)' Give name for formatted input file: '
1362    read(std_in, '(a)',err=10,iomsg=errmsg ) filnam(1)
1363    write(std_out, '(a)',err=10,iomsg=errmsg ) trim(filnam(1))
1364    write(std_out,*)' Give name for formatted output file:'
1365    read (std_in, '(a)',err=10,iomsg=errmsg ) filnam(2)
1366    write (std_out, '(a)',err=10,iomsg=errmsg ) trim(filnam(2))
1367    write(std_out,*)' Give root name for generic input files:'
1368    read (std_in, '(a)',err=10,iomsg=errmsg ) filnam(3)
1369    write (std_out, '(a)',err=10,iomsg=errmsg ) trim(filnam(3))
1370    write(std_out,*, err=10, iomsg=errmsg )' Give root name for generic output files:'
1371    read (std_in, '(a)', err=10, iomsg=errmsg ) filnam(4)
1372    write (std_out, '(a)', err=10, iomsg=errmsg ) trim(filnam(4))
1373    write(std_out,*, err=10, iomsg=errmsg)' Give root name for generic temporary files:'
1374    read (std_in, '(a)', err=10, iomsg=errmsg ) filnam(5)
1375    write (std_out, '(a)', err=10, iomsg=errmsg ) trim(filnam(5))
1376 
1377 !  Check that old input file exists
1378    call isfile(filnam(1),'old')
1379 
1380 !  Check that new output file does NOT exist
1381    call isfile(filnam(2),'new')
1382 
1383 !  Check that root name for generic input and output differ
1384    if ( trim(filnam(3))==trim(filnam(4)) ) then
1385      write(message, '(a,a,a)' )&
1386 &     'Root name for generic input and output files must differ ',ch10,&
1387 &     'Action: correct your "file" file.'
1388      MSG_ERROR(message)
1389    end if
1390 
1391 !  Check that root names are at least 20 characters less than fnlen
1392    if ( len_trim(filnam(3)) >= (fnlen-20) ) then
1393      write(message, '(a,a,a,a,a,i4,a,i4,a,a)' )&
1394 &     'Root name for generic input files is too long. ',ch10,&
1395 &     'It must be 20 characters less than the maximal allowed ',ch10,&
1396 &     'length of names, that is ',fnlen,', while it is ',len_trim(filnam(3)),ch10,&
1397 &     'Action : correct your "file" file.'
1398      MSG_ERROR(message)
1399    end if
1400    if ( len_trim(filnam(4)) >= (fnlen-20) ) then
1401      write(message, '(a,a,a,a,a,i4,a,i4,a,a)' )&
1402 &     'Root name for generic output files is too long. ',ch10,&
1403 &     'It must be 20 characters less than the maximal allowed ',ch10,&
1404 &     'length of names, that is ',fnlen,', while it is ',len_trim(filnam(4)),ch10,&
1405 &     'Action: correct your "file" file.'
1406      MSG_ERROR(message)
1407    end if
1408    if ( len_trim(filnam(5)) >= (fnlen-20) ) then
1409      write(message, '(a,a,a,a,a,i4,a,i4,a,a)' )&
1410 &     'Root name for generic temporary files is too long. ',ch10,&
1411 &     'It must be 20 characters less than the maximal allowed ',ch10,&
1412 &     'length of names, that is ',fnlen,', while it is ',len_trim(filnam(5)),ch10,&
1413 &     'Action: correct your "file" file.'
1414      MSG_ERROR(message)
1415    end if
1416 
1417  end if ! master only
1418 
1419 !Communicate filenames to all processors
1420  call xmpi_bcast(filnam,master,comm,ierr)
1421 
1422 !Check
1423 !Create a name for the status file, based on filnam(5)
1424  filstat=trim(filnam(5))//'_STATUS'
1425 
1426 !Redefine the log unit if not the master
1427  if (me/=master) then
1428    call int2char4(me,tag)
1429    ABI_CHECK((tag(1:1)/='#'),'Bug: string length too short!')
1430    filstat=trim(filstat)//'_P-'//trim(tag)
1431    if (do_write_log) then
1432      fillog=trim(filnam(5))//'_LOG_'//trim(tag)
1433      close(std_out, err=10, iomsg=errmsg)
1434      if (open_file(fillog,message,unit=std_out,status='unknown',action="write") /= 0) then
1435        MSG_ERROR(message)
1436      end if
1437    else
1438      close(std_out, err=10, iomsg=errmsg)
1439      if (open_file(NULL_FILE,message,unit=std_out,action="write") /= 0) then
1440        MSG_ERROR(message)
1441      end if
1442    end if
1443  end if
1444 
1445  call xmpi_barrier(comm)
1446  return
1447 
1448 ! Handle possibe IO errors
1449  10 continue
1450  MSG_ERROR(errmsg)
1451 
1452 end subroutine iofn1

m_dtfil/isfile [ Functions ]

[ Top ] [ m_dtfil ] [ Functions ]

NAME

 isfile

FUNCTION

 Inquire Status of FILE
 Checks that for status =
 'old': file already exists
 'new': file does not exist; if file exists,
 filnam is modified to filnam.A or filnam.B,....

INPUTS

 filnam=character string to specify filename
 status='old' or 'new'

OUTPUT

 stops processing if old file does not exist; changes name
 and returns new name in redefined filnam if new file already exists.

PARENTS

      anaddb,iofn1,m_effective_potential,m_polynomial_coeff,m_vcoul
      multibinit,ujdet

CHILDREN

      clib_rename,int2char4

SOURCE

1096 subroutine isfile(filnam, status)
1097 
1098 
1099 !This section has been created automatically by the script Abilint (TD).
1100 !Do not modify the following lines by hand.
1101 #undef ABI_FUNC
1102 #define ABI_FUNC 'isfile'
1103 !End of the abilint section
1104 
1105  implicit none
1106 
1107 !Arguments ------------------------------------
1108 !scalars
1109  character(len=3),intent(in) :: status
1110  character(len=fnlen),intent(inout) :: filnam
1111 
1112 !Local variables-------------------------------
1113 !scalars
1114  logical :: ex,found
1115  integer :: ii,ios, ioserr
1116  character(len=500) :: message
1117  character(len=fnlen) :: filnam_tmp
1118  character(len=fnlen) :: trialnam
1119 
1120 ! *************************************************************************
1121 
1122  filnam_tmp=filnam
1123 
1124  if (status=='old') then !  Check that old file exists
1125    inquire(file=filnam,iostat=ios,exist=ex)
1126 
1127    if (ios/=0) then
1128      write(message,'(4a,i0,2a)')&
1129 &     'Checks for existence of file  ',trim(filnam),ch10,&
1130 &     'but INQUIRE statement returns error code',ios,ch10,&
1131 &     'Action: identify which problem appears with this file.'
1132      MSG_ERROR(message)
1133    else if (.not.ex) then
1134      write(message, '(5a)' )&
1135 &     'Checks for existence of file  ',trim(filnam),ch10,&
1136 &     'but INQUIRE finds file does not exist.',&
1137 &     'Action: check file name and re-run.'
1138      MSG_ERROR(message)
1139    end if
1140 
1141  else if (status=='new') then
1142 
1143    ! Check that new output file does NOT exist
1144    ioserr = 0
1145    trialnam = filnam
1146    ii = 0
1147    inquire(file=trim(trialnam),iostat=ios,exist=ex)
1148    if ( ios /= 0 ) then
1149      write(message,'(4a)') 'Something is wrong with permissions for ', &
1150 &     'reading/writing on this filesystem.',ch10,&
1151 &     'Action : Check permissions.'
1152      MSG_ERROR(message)
1153    end if
1154 
1155    if ( ex .eqv. .true. ) then
1156      write(message,'(4a)')'Output file ',trim(trialnam),ch10,' already exists.'
1157      MSG_COMMENT(message)
1158      found=.false.
1159 
1160      ii=1
1161      do while ( (found .eqv. .false.) .and. (ii < 10000) )
1162        call int2char4(ii,message)
1163        trialnam=trim(trim(filnam_tmp)//message)
1164        inquire(file=trim(trialnam),iostat=ios,exist=ex)
1165        if ( (ex .eqv. .false.) .and. (ios == 0)) then
1166          found  = .true.
1167        end if
1168        if ( ios /= 0 )  ioserr=ioserr+1
1169        if ( ioserr > 10 ) then
1170 !        There is a problem => stop
1171          write(message, '(2a,i0,2a)' )&
1172 &         'Check for permissions of reading/writing files on the filesystem', &
1173 &         '10 INQUIRE statements returned an error code like ',ios,ch10,&
1174 &         'Action: Check permissions'
1175          MSG_ERROR(message)
1176        end if
1177        ii=ii+1
1178      end do
1179      if ( found .eqv. .true. ) then
1180        write(message,'(4a)') 'Renaming old ',trim(filnam),' to ',trim(trialnam)
1181        MSG_COMMENT(message)
1182        call clib_rename(filnam,trialnam,ioserr)
1183        if ( ioserr /= 0 ) then
1184          write(message,'(4a)') 'Failed to rename file ', trim(filnam),' to ',trim(trialnam)
1185          MSG_ERROR(message)
1186        end if
1187      else
1188        write(message,'(3a)')&
1189 &       'Have used all names of the form filenameXXXX, X in [0-9]',ch10,&
1190 &       'Action: clean up your directory and start over.'
1191        MSG_ERROR(message)
1192      end if
1193    end if
1194 
1195    ! if ii > 0 we iterated so rename abi_out to abi_outXXXX
1196    ! and just write to abi_out
1197  else ! status not recognized
1198    write(message,'(3a)')'  Input status= ',status,' not recognized.'
1199    MSG_BUG(message)
1200  end if
1201 
1202 end subroutine isfile

m_dtfil/mkfilename [ Functions ]

[ Top ] [ m_dtfil ] [ Functions ]

NAME

 mkfilename

FUNCTION

 From the root (input or output) file names, produce a real file name.

INPUTS

 character(len=fnlen):: filnam(5)=the root file names
  (only filnam(3) and filnam(4) are really needed)
 get=input 'get variable', if 1, must get the file from another dataset
 idtset=number of the dataset
 ird=input 'iread variable', if 1, must get the file from the input root
 jdtset_(0:ndtset)=actual index of the dataset
 ndtset=number of datasets
 stringfil character(len=*)=the string of characters to be appended e.g. '_WFK' or '_DEN'
 stringvar tcharacter(len=*)=the string of characters to be appended
   that defines the 'get' or 'ird' variables, e.g. 'wfk' or 'ddk'

OUTPUT

 character(len=fnlen):: filnam_out=the new file name
 will_read=1 if the file must be read ; 0 otherwise (ird and get were zero)

PARENTS

      dtfil_init,finddistrproc

CHILDREN

      appdig,wrtout

SOURCE

 965 subroutine mkfilename(filnam,filnam_out,get,idtset,ird,jdtset_,ndtset,stringfil,stringvar,will_read)
 966 
 967 
 968 !This section has been created automatically by the script Abilint (TD).
 969 !Do not modify the following lines by hand.
 970 #undef ABI_FUNC
 971 #define ABI_FUNC 'mkfilename'
 972 !End of the abilint section
 973 
 974  implicit none
 975 
 976 !Arguments ------------------------------------
 977 !scalars
 978  integer,intent(in) :: get,idtset,ird,ndtset
 979  integer,intent(out) :: will_read
 980  character(len=*),intent(in) :: stringfil
 981  character(len=*),intent(in) :: stringvar
 982  character(len=fnlen),intent(out) :: filnam_out
 983 !arrays
 984  integer,intent(in) :: jdtset_(0:ndtset)
 985  character(len=fnlen),intent(in) :: filnam(5)
 986 
 987 !Local variables-------------------------------
 988 !scalars
 989  integer :: jdtset,jget
 990  character(len=4) :: appen
 991  character(len=500) :: message
 992  character(len=fnlen) :: filnam_appen
 993 
 994 ! *************************************************************************
 995 
 996 !Here, defaults if no get variable
 997  will_read=ird
 998 
 999  filnam_appen=trim(filnam(3))
1000  if(ndtset>0)then
1001    jdtset=jdtset_(idtset)
1002    call appdig(jdtset,'',appen)
1003    filnam_appen=trim(filnam_appen)//'_DS'//appen
1004  end if
1005  filnam_out=trim(filnam_appen)//trim(stringfil)
1006 
1007 !Treatment of the multi-dataset case  (get is not relevant otherwise)
1008  if(ndtset/=0)then
1009 
1010    if(ndtset==1.and.get<0.and.(jdtset_(1)+get>0))then
1011      write(message, '(7a,i3,a,i3,5a)' )&
1012 &     'You cannot use a negative value of get',trim(stringvar),' with only 1 dataset!',ch10, &
1013 &     ' If you want to refer to a previously computed dataset,',ch10, &
1014 &     ' you should give the absolute index of it (i.e. ', &
1015 &     jdtset_(idtset)+get,' instead of ',get,').',ch10, &
1016 &     'Action: correct get',trim(stringvar),' in your input file.'
1017      MSG_ERROR(message)
1018    end if
1019 
1020    if(idtset+get<0)then
1021      write(message, '(a,a,a,a,a,i3,a,a,a,i3,a,a,a,a)' )&
1022 &     'The sum of idtset and get',trim(stringvar),' cannot be negative,',ch10,&
1023 &     'while they are idtset=',idtset,', and get',trim(stringvar),'=',get,ch10,&
1024 &     'Action: correct get',trim(stringvar),' in your input file.'
1025      MSG_ERROR(message)
1026    end if
1027 
1028    if(get>0 .or. (get<0 .and. idtset+get>0) )then
1029 
1030      if(ird/=0 .and. get/=0)then
1031        write(message, '(a,a,a,a,a,a,a,a,a,a,a,i3,a,i3,a,a,a,a,a,a,a)' )&
1032 &       'The input variables ird',trim(stringvar),' and get',trim(stringvar),' cannot be',ch10,&
1033 &       'simultaneously non-zero, while for idtset=',idtset,',',ch10,&
1034 &       'they are ',ird,', and ',get,'.',ch10,&
1035 &       'Action: correct ird',trim(stringvar),' or get',trim(stringvar),' in your input file.'
1036        MSG_ERROR(message)
1037      end if
1038 
1039      will_read=1
1040 
1041 !    Compute the dataset from which to take the file, and the corresponding index
1042      if(get<0 .and. idtset+get>0) jget=jdtset_(idtset+get)
1043      if(get>0) jget=get
1044      call appdig(jget,'',appen)
1045 
1046 !    Note use of output filename (filnam(4))
1047      filnam_out=trim(filnam(4))//'_DS'//trim(appen)//trim(stringfil)
1048 
1049      if(jdtset>=100)then
1050        write(message, '(a,a,a,a,a,i5,a,a)' )&
1051 &       ' mkfilename : get',trim(stringvar) ,'/=0, take file ',trim(stringfil),&
1052 &       ' from output of DATASET ',jget,'.',ch10
1053      else
1054        write(message, '(a,a,a,a,a,i3,a,a)' )&
1055 &       ' mkfilename : get',trim(stringvar) ,'/=0, take file ',trim(stringfil),&
1056 &       ' from output of DATASET ',jget,'.',ch10
1057      end if
1058      call wrtout(ab_out,message,'COLL')
1059      call wrtout(std_out,message,'COLL')
1060    end if ! conditions on get and idtset
1061 
1062  end if ! ndtset/=0
1063 
1064 end subroutine mkfilename