TABLE OF CONTENTS


ABINIT/invars0 [ Functions ]

[ Top ] [ Functions ]

NAME

 invars0

FUNCTION

 Initialisation phase : prepare the main input subroutine call by
 reading most of the NO MULTI variables, as well as natom, nimage, and ntypat,
 needed for allocating some input arrays in abinit, and also useri
 and userr. The variable usewvl is also read here for later reading
 of input path for the atomic orbital file (if required).

COPYRIGHT

 Copyright (C) 1999-2018 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

  lenstr=actual length of string
  ndtset= number of datasets to be read; if 0, no multi-dataset mode
  ndtset_alloc=number of datasets, corrected for allocation of at least
               one data set.
  string*(*)=string of characters containing all input variables and data

OUTPUT

  dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables,
   some of which are initialized here :
   cpus,jdtset,natom,nimage,npsp,ntypat,useri*,userr*
  istatr=repetition rate for status file
  istatshft=shift of the repetition rate for status file
  msym=maximal value of input msym for all the datasets
  mxnatom=maximal value of input natom for all the datasets
  mxnimage=maximal value of input nimage for all the datasets
  mxntypat=maximal value of input ntypat for all the datasets
  npsp=number of pseudopotentials

PARENTS

      m_ab7_invars_f90

CHILDREN

      get_ndevice,intagm

SOURCE

 47 #if defined HAVE_CONFIG_H
 48 #include "config.h"
 49 #endif
 50 
 51 #include "abi_common.h"
 52 
 53 subroutine invars0(dtsets,istatr,istatshft,lenstr,&
 54 & msym,mxnatom,mxnimage,mxntypat,ndtset,ndtset_alloc,npsp,papiopt,timopt,string)
 55 
 56  use defs_basis
 57  use defs_abitypes
 58  use m_errors
 59  use m_profiling_abi
 60 #if defined HAVE_GPU_CUDA
 61  use m_initcuda, only : Get_ndevice
 62 #endif
 63 
 64  use m_parser,  only : intagm
 65 
 66 !This section has been created automatically by the script Abilint (TD).
 67 !Do not modify the following lines by hand.
 68 #undef ABI_FUNC
 69 #define ABI_FUNC 'invars0'
 70 !End of the abilint section
 71 
 72  implicit none
 73 
 74 !Arguments ------------------------------------
 75 !scalars
 76  integer,intent(in) :: lenstr,ndtset,ndtset_alloc
 77  integer,intent(out) :: istatr,istatshft,msym,mxnatom,mxnimage,mxntypat,npsp,papiopt
 78  integer,intent(inout) :: timopt
 79  character(len=*),intent(in) :: string
 80 !arrays
 81  type(dataset_type),intent(inout) :: dtsets(0:ndtset_alloc) !vz_i
 82 
 83 !Local variables-------------------------------
 84 !scalars
 85  integer :: i1,i2,idtset,ii,jdtset,marr,tjdtset,tread,treadh,treadm
 86  integer :: treads,use_gpu_cuda
 87  real(dp) :: cpus
 88  character(len=500) :: message
 89 !arrays
 90  integer,allocatable :: intarr(:)
 91  real(dp),allocatable :: dprarr(:)
 92 
 93 !******************************************************************
 94 
 95 !Set ii to avoid warning of uninitialised variable
 96  ii = 0
 97 
 98  marr=max(ndtset_alloc,2)
 99  ABI_ALLOCATE(dprarr,(marr))
100  ABI_ALLOCATE(intarr,(marr))
101 
102 !Set up jdtset
103  if(ndtset/=0)then
104 
105 !  Default values
106    dtsets(0)%jdtset = -1 ! unused value
107    dtsets(1:ndtset_alloc)%jdtset=(/ (ii,ii=1,ndtset_alloc) /)
108 
109 !  Read explicitely the jdtset array
110    call intagm(dprarr,intarr,0,marr,ndtset,string(1:lenstr),'jdtset',tjdtset,'INT')
111    if(tjdtset==1) dtsets(1:ndtset)%jdtset=intarr(1:ndtset)
112 
113 !  Read the udtset array
114    call intagm(dprarr,intarr,0,marr,2,string(1:lenstr),'udtset',tread,'INT')
115 
116 !  jdtset and udtset cannot be defined together
117    if(tjdtset==1 .and. tread==1)then
118      write(message, '(3a)' )&
119 &     'jdtset and udtset cannot be defined both in the input file.',ch10,&
120 &     'Action: remove one of them from your input file.'
121      MSG_ERROR(message)
122    end if
123 
124 !  Check values of udtset
125    if(tread==1)then
126      if(intarr(1)<1 .or. intarr(1)>999)then
127        write(message, '(a,i0,3a)' )&
128 &       'udtset(1) must be between 1 and 999, but it is ',intarr(1),'.',ch10,&
129 &       'Action: change the value of udtset(1) in your input file.'
130        MSG_ERROR(message)
131      end if
132      if(intarr(2)<1 .or. intarr(2)>9)then
133        write(message, '(a,i0,3a)' )&
134 &       'udtset(2) must be between 1 and 9, but it is ',intarr(2),'.',ch10,&
135 &       'Action: change the value of udtset(2) in your input file.'
136        MSG_ERROR(message)
137      end if
138      if(intarr(1)*intarr(2) /= ndtset)then
139        write(message, '(3a,i0,3a,i0,a,i0,3a,i0,3a)' )&
140 &       'udtset(1)*udtset(2) must be equal to ndtset,',ch10,&
141 &       'but it is observed that udtset(1) = ',intarr(1),',',ch10,&
142 &       'and udtset(2) = ',intarr(2),' so that their product is ',intarr(1)*intarr(2),',',ch10,&
143 &       'while ndtset is ',ndtset,'.',ch10,&
144 &       'Action: change udtset or ndtset in your input file.'
145        MSG_ERROR(message)
146      end if
147      idtset=0
148      do i1=1,intarr(1)
149        do i2=1,intarr(2)
150          idtset=idtset+1
151          dtsets(idtset)%jdtset=i1*10+i2
152        end do
153      end do
154    end if
155 
156 !  Final check on the jdtset values
157    do idtset=1,ndtset
158      if(dtsets(idtset)%jdtset<1 .or. dtsets(idtset)%jdtset>9999)then
159        write(message, '(3a,i0,a,i0,a,a)' )&
160 &       'The components of jdtset must be between 1 and 9999.',ch10,&
161 &       'However, the input value of the component ',idtset,' of jdtset is ',dtsets(idtset)%jdtset,ch10,&
162 &       'Action : correct jdtset in your input file.'
163        MSG_ERROR(message)
164      end if
165    end do
166 
167  else
168    dtsets(1)%jdtset=0
169  end if
170 
171  papiopt = 0
172  call intagm(dprarr,intarr,0,1,1,string(1:lenstr),'papiopt',tread,'INT')
173  if(tread==1) papiopt=intarr(1)
174 
175 !Read timopt and pass it to timab
176  call intagm(dprarr,intarr,0,1,1,string(1:lenstr),'timopt',tread,'INT')
177  if(tread==1) timopt=intarr(1)
178 
179  istatr=0
180  dtsets(0)%istatr=istatr
181  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'istatr',tread,'INT')
182  if(tread==1) istatr=intarr(1)
183  dtsets(1:)%istatr=istatr
184 
185  istatshft=1
186  dtsets(0)%istatshft=istatshft
187  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'istatshft',tread,'INT')
188  if(tread==1) istatshft=intarr(1)
189  dtsets(1:)%istatshft=istatshft
190 
191  cpus=zero
192  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'cpus ',treads,'DPR')
193  if(treads==1) cpus=dprarr(1)
194  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'cpum ',treadm,'DPR')
195  if(treadm==1) cpus=dprarr(1)*60.0_dp
196  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'cpuh ',treadh,'DPR')
197 
198  if(treadh==1) cpus=dprarr(1)*3600.0_dp
199  if(treads+treadm+treadh>1)then
200    write(message, '(5a)' )&
201 &   'More than one input variable is used to defined the CPU time limit.',ch10,&
202 &   'This is not allowed.',ch10,&
203 &   'Action: in the input file, suppress either cpus, cpum or cpuh.'
204    MSG_ERROR(message)
205  end if
206  dtsets(:)%cpus=cpus
207 
208 !Default for natom, nimage, ntypat, useri and userr
209  dtsets(:)%natom=1
210  dtsets(:)%nimage=1
211  dtsets(:)%ntypat=1 ; dtsets(0)%ntypat=0    ! Will always echo ntypat
212  dtsets(:)%macro_uj=0
213  dtsets(:)%maxnsym=384
214  dtsets(:)%useria=0
215  dtsets(:)%userib=0
216  dtsets(:)%useric=0
217  dtsets(:)%userid=0
218  dtsets(:)%userie=0
219  dtsets(:)%userra=zero
220  dtsets(:)%userrb=zero
221  dtsets(:)%userrc=zero
222  dtsets(:)%userrd=zero
223  dtsets(:)%userre=zero
224  dtsets(:)%usewvl = 0
225  dtsets(:)%plowan_compute=0
226 
227 !Loop on datasets, to find natom and mxnatom, as well as useri and userr
228  do idtset=1,ndtset_alloc
229    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
230 
231 !  Read natom from string
232    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natom',tread,'INT')
233 !  Might initialize natom from XYZ file
234    if(tread==0)then
235      call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'_natom',tread,'INT')
236    end if
237 
238    if(tread==1)then
239      dtsets(idtset)%natom=intarr(1)
240    else
241      write(message, '(a,i0,2a)' )&
242 &     'Input natom must be defined, but was absent for dataset ',jdtset,ch10,&
243 &     'Action: check the input file.'
244      MSG_ERROR(message)
245    end if
246 !  Check that natom is greater than 0
247    if (dtsets(idtset)%natom<=0) then
248      write(message, '(a,i0,2a,i0,3a)' )&
249 &     'Input natom must be > 0, but was ',dtsets(idtset)%natom,ch10,&
250 &     'for dataset ',jdtset,'. This is not allowed.',ch10,&
251 &     'Action: check the input file.'
252      MSG_ERROR(message)
253    end if
254 
255    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nimage',tread,'INT')
256    if(tread==1) dtsets(idtset)%nimage=intarr(1)
257 
258 !  Check that nimage is greater than 0
259    if (dtsets(idtset)%nimage<=0) then
260      write(message, '(a,i0,4a)' )&
261 &     'nimage must be > 0, but was ',dtsets(idtset)%nimage,ch10,&
262 &     'This is not allowed.',ch10,&
263 &     'Action: check the input file.'
264      MSG_ERROR(message)
265    end if
266 
267    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ntypat',tread,'INT')
268    if(tread==1)dtsets(idtset)%ntypat=intarr(1)
269 !  Check that ntypat is greater than 0
270    if (dtsets(idtset)%ntypat<=0) then
271      write(message, '(a,i0,2a,i0,3a)' )&
272 &     'Input ntypat must be > 0, but was ',dtsets(idtset)%ntypat,ch10,&
273 &     'for dataset ',jdtset,'. This is not allowed.',ch10,&
274 &     'Action: check the input file.'
275      MSG_ERROR(message)
276    end if
277 
278 !  Read msym from string
279    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'maxnsym',tread,'INT')
280    if(tread==1)dtsets(idtset)%maxnsym=intarr(1)
281 !  Check that maxnsym is greater than 1
282    if (dtsets(idtset)%maxnsym<1) then
283      write(message, '(a,i0,2a,i0,3a)' )&
284 &     'Input maxnsym must be > 1, but was ',dtsets(idtset)%maxnsym,ch10,&
285 &     'for dataset ',jdtset,'. This is not allowed.',ch10,&
286 &     'Action: check the input file.'
287      MSG_ERROR(message)
288    end if
289 
290 ! Read plowan_compute
291    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'plowan_compute',tread,'INT')
292    if(tread==1) dtsets(idtset)%plowan_compute=intarr(1)
293 
294    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'useria',tread,'INT')
295    if(tread==1) dtsets(idtset)%useria=intarr(1)
296    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userib',tread,'INT')
297    if(tread==1) dtsets(idtset)%userib=intarr(1)
298    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'useric',tread,'INT')
299    if(tread==1) dtsets(idtset)%useric=intarr(1)
300    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userid',tread,'INT')
301    if(tread==1) dtsets(idtset)%userid=intarr(1)
302    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userie',tread,'INT')
303    if(tread==1) dtsets(idtset)%userie=intarr(1)
304 
305    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userra',tread,'DPR')
306    if(tread==1) dtsets(idtset)%userra=dprarr(1)
307    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userrb',tread,'DPR')
308    if(tread==1) dtsets(idtset)%userrb=dprarr(1)
309    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userrc',tread,'DPR')
310    if(tread==1) dtsets(idtset)%userrc=dprarr(1)
311    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userrd',tread,'DPR')
312    if(tread==1) dtsets(idtset)%userrd=dprarr(1)
313    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userre',tread,'DPR')
314    if(tread==1) dtsets(idtset)%userre=dprarr(1)
315 
316    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'usewvl',tread,'INT')
317    if(tread==1) dtsets(idtset)%usewvl=intarr(1)
318 
319  end do
320 
321 !mxnatom =maxval(dtsets(1:ndtset_alloc)%natom)
322 !mxntypat =maxval(dtsets(1:ndtset_alloc)%ntypat)
323 !msym =maxval(dtsets(1:ndtset_alloc)%maxnsym)
324 !There is a bug in the HP compiler, the following should execute properly
325  mxnatom=dtsets(1)%natom ; mxnimage=dtsets(1)%nimage
326  mxntypat=dtsets(1)%ntypat ; msym=dtsets(1)%maxnsym
327  if(ndtset_alloc>1)then
328    do idtset=2,ndtset_alloc
329      mxnatom =max(dtsets(idtset)%natom,mxnatom)
330      mxnimage=max(dtsets(idtset)%nimage,mxnimage)
331      mxntypat=max(dtsets(idtset)%ntypat,mxntypat)
332      msym    =max(dtsets(idtset)%maxnsym,msym)
333    end do
334  end if
335 
336  if(mxnimage>1)then
337    do idtset=2,ndtset_alloc
338      if(mxnatom/=dtsets(idtset)%natom)then
339        write(message,'(5a,i0,a,i0,3a,i0,a)')&
340 &       'When there exist one dataset with more than one image,',ch10,&
341 &       'the number of atoms in each dataset must be the same.',ch10,&
342 &       'However, it has been found that for dataset= ',idtset,ch10,&
343 &       'natom= ',dtsets(idtset)%natom,' differs from the maximum number',ch10,&
344 &       'of atoms, mxnatom= ',mxnatom,&
345 &       'Action: check the input variables natom for different datasets.'
346        MSG_ERROR(message)
347      end if
348    end do
349  end if
350 
351 !Set up npsp
352  npsp=mxntypat   ! Default value
353  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'npsp',tread,'INT')
354  if(tread==1)then
355    npsp=intarr(1)
356  else
357    if(ndtset_alloc>1)then
358      do idtset=1,ndtset_alloc
359        if(dtsets(idtset)%ntypat/=mxntypat)then
360          write(message, '(5a,i0,a,i0,2a,i0,2a)' )&
361 &         'When npsp is not defined, the input variable ntypat must be',ch10,&
362 &         'the same for all datasets. However, it has been found that for',ch10,&
363 &         'jdtset: ',dtsets(idtset)%jdtset,', ntypat= ',dtsets(idtset)%ntypat,ch10,&
364 &         'differs from the maximum value of ntypat= ',mxntypat,ch10,&
365 &         'Action: check the input variables npsp and ntypat.'
366          MSG_ERROR(message)
367        end if
368      end do
369    end if
370  end if
371  dtsets(0)%npsp = mxntypat   ! Default value
372  dtsets(1:ndtset_alloc)%npsp = npsp
373 
374 !KGB parallelism information (needed at this stage)
375  dtsets(:)%paral_kgb=0
376  do idtset=1,ndtset_alloc
377    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
378    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'paral_kgb',tread,'INT')
379    if(tread==1)dtsets(idtset)%paral_kgb=intarr(1)
380 
381    if (dtsets(idtset)%paral_kgb<0 .or. dtsets(idtset)%paral_kgb>1) then
382      write(message,'(a,i0,2a,i0,3a)')&
383 &     'Input paral_kgb must be 0 or 1, but was ',dtsets(idtset)%paral_kgb,ch10,&
384 &     'for dataset',jdtset,'. This is not allowed.',ch10,&
385 &     'Action: check the input file.'
386      MSG_ERROR(message)
387    end if
388  end do
389 
390 !GPU information
391  use_gpu_cuda=0
392  dtsets(:)%use_gpu_cuda=0
393 #if defined HAVE_GPU_CUDA && defined HAVE_GPU_CUDA_DP
394  call Get_ndevice(ii)
395  if (ii>0) then
396    do i1=1,ndtset_alloc
397      dtsets(i1)%use_gpu_cuda=-1
398    end do
399  end if
400 #endif
401  do idtset=1,ndtset_alloc
402    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
403    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'use_gpu_cuda',tread,'INT')
404    if(tread==1)dtsets(idtset)%use_gpu_cuda=intarr(1)
405    if (dtsets(idtset)%use_gpu_cuda==1) use_gpu_cuda=1
406  end do
407  if (use_gpu_cuda==1) then
408 #if defined HAVE_GPU_CUDA && defined HAVE_GPU_CUDA_DP
409    if (ii<=0) then
410      write(message,'(3a)')&
411 &     'Input variables use_gpu_cuda is on',ch10,&
412 &     'but no available GPU device has been detected !',ch10,&
413 &     'Action: change the input variable use_gpu_cuda.'
414      MSG_ERROR(message)
415    end if
416 #else
417    write(message,'(7a)')&
418 &   'Input variables use_gpu_cuda is on but abinit hasn''t been built',ch10,&
419 &   'with (double precision) gpu mode enabled !',ch10,&
420 &   'Action: change the input variable use_gpu_cuda',ch10,&
421 &   '        or re-compile ABINIT with double-precision Cuda enabled.'
422    MSG_ERROR(message)
423 #endif
424  end if
425 
426  ABI_DEALLOCATE(dprarr)
427  ABI_DEALLOCATE(intarr)
428 
429 !We allocate the internal array, depending on the computed values.
430 !WARNING : do not forget to deallocate these arrays in the routine dtset_free
431 !(should make a separate subroutine for allocating/deallocating these records)
432  do idtset=0,ndtset_alloc
433    ABI_ALLOCATE(dtsets(idtset)%acell_orig,(3,mxnimage))
434    ABI_ALLOCATE(dtsets(idtset)%algalch,(mxntypat))
435    ABI_ALLOCATE(dtsets(idtset)%amu_orig,(mxntypat,mxnimage))
436    ABI_ALLOCATE(dtsets(idtset)%corecs,(mxntypat))
437    ABI_ALLOCATE(dtsets(idtset)%densty,(mxntypat,4))
438    ABI_ALLOCATE(dtsets(idtset)%dynimage,(mxnimage))
439    ABI_ALLOCATE(dtsets(idtset)%iatfix,(3,mxnatom))
440    ABI_ALLOCATE(dtsets(idtset)%f4of2_sla,(mxntypat))
441    ABI_ALLOCATE(dtsets(idtset)%f6of2_sla,(mxntypat))
442    ABI_ALLOCATE(dtsets(idtset)%jpawu,(mxntypat,mxnimage))
443    ABI_ALLOCATE(dtsets(idtset)%kberry,(3,20))
444    ABI_ALLOCATE(dtsets(idtset)%lexexch,(mxntypat))
445    ABI_ALLOCATE(dtsets(idtset)%ldaminushalf,(mxntypat))
446    ABI_ALLOCATE(dtsets(idtset)%lpawu,(mxntypat))
447    ABI_ALLOCATE(dtsets(idtset)%mixalch_orig,(npsp,mxntypat,mxnimage))
448    ABI_ALLOCATE(dtsets(idtset)%nucdipmom,(3,mxnatom))
449    ABI_ALLOCATE(dtsets(idtset)%pimass,(mxntypat))
450    ABI_ALLOCATE(dtsets(idtset)%ptcharge,(mxntypat))
451    ABI_ALLOCATE(dtsets(idtset)%prtatlist,(mxnatom))
452    ABI_ALLOCATE(dtsets(idtset)%quadmom,(mxntypat))
453    ABI_ALLOCATE(dtsets(idtset)%ratsph,(mxntypat))
454    ABI_ALLOCATE(dtsets(idtset)%rprim_orig,(3,3,mxnimage))
455    ABI_ALLOCATE(dtsets(idtset)%rprimd_orig,(3,3,mxnimage))
456    ABI_ALLOCATE(dtsets(idtset)%so_psp,(npsp))
457    ABI_ALLOCATE(dtsets(idtset)%spinat,(3,mxnatom))
458    ABI_ALLOCATE(dtsets(idtset)%shiftk,(3,210))
459    ABI_ALLOCATE(dtsets(idtset)%typat,(mxnatom))
460    ABI_ALLOCATE(dtsets(idtset)%upawu,(mxntypat,mxnimage))
461 !   if (dtsets(idtset)%plowan_compute>0) then
462    ABI_ALLOCATE(dtsets(idtset)%plowan_iatom,(mxnatom))
463    ABI_ALLOCATE(dtsets(idtset)%plowan_it,(100*3))
464    ABI_ALLOCATE(dtsets(idtset)%plowan_nbl,(mxnatom))
465    ABI_ALLOCATE(dtsets(idtset)%plowan_lcalc,(12*mxnatom))
466    ABI_ALLOCATE(dtsets(idtset)%plowan_projcalc,(12*mxnatom))
467 !   endif
468    ABI_ALLOCATE(dtsets(idtset)%vel_orig,(3,mxnatom,mxnimage))
469    ABI_ALLOCATE(dtsets(idtset)%vel_cell_orig,(3,3,mxnimage))
470    ABI_ALLOCATE(dtsets(idtset)%xred_orig,(3,mxnatom,mxnimage))
471    ABI_ALLOCATE(dtsets(idtset)%ziontypat,(mxntypat))
472    ABI_ALLOCATE(dtsets(idtset)%znucl,(npsp))
473  end do
474 
475 !DEBUG
476 !write(std_out,*)' invars0 : nimage, mxnimage = ',dtsets(:)%nimage, mxnimage
477 !write(std_out,*)' invars0 : natom = ',dtsets(:)%natom
478 !write(std_out,*)' invars0 : mxnatom = ',mxnatom
479 !ENDDEBUG
480 
481 end subroutine invars0