source: trunk/LMDZ.MARS/libf/aeronomars/inichim_newstart.F90 @ 1416

Last change on this file since 1416 was 1381, checked in by emillour, 10 years ago

Mars GCM:

EM

File size: 20.2 KB
RevLine 
[1047]1      subroutine inichim_newstart(ngrid, nq, pq, qsurf, ps, &
2                                  flagh2o, flagthermo)
[38]3
[1036]4      use tracer_mod
[38]5      implicit none
6
[618]7!=======================================================================
8!
9!  Purpose:
10!  --------
11!
12!  Initialization of the chemistry for use in the building of a new start file
13!  used by program newstart.F
14!  also used by program testphys1d.F
15!
16!  Authors:
17!  -------
18!  Initial version 11/2002 by Sebastien Lebonnois
19!  Revised 07/2003 by Monica Angelats-i-Coll to use input files
20!  Modified 10/2008 Identify tracers by their names Ehouarn Millour
21!  Modified 11/2011 Addition of methane Franck Lefevre
22!  Rewritten 04/2012 Franck Lefevre
23!
24!  Arguments:
25!  ----------
26!
[1036]27!  pq(iip1,jjp1,llm,nq)  Advected fields, ie chemical species here
[1047]28!  qsurf(ngrid,nq)     Amount of tracer on the surface (kg/m2)
[618]29!  ps(iip1,jjp1)           Surface pressure (Pa)
30!  flagh2o                 flag for initialisation of h2o (1: yes / 0: no)
31!  flagthermo              flag for initialisation of thermosphere only (1: yes / 0: no)
32!
33!=======================================================================
[38]34
35#include "dimensions.h"
36#include "paramet.h"
[618]37#include "comvert.h"
[38]38#include "callkeys.h"
39#include "datafile.h"
40
[618]41! inputs :
[38]42
[1047]43      integer,intent(in) :: ngrid         ! number of atmospheric columns in the physics
[1036]44      integer,intent(in) :: nq                    ! number of tracers
[618]45      real,intent(in) :: ps(iip1,jjp1)            ! surface pressure in the gcm (Pa)   
46      integer,intent(in) :: flagh2o               ! flag for h2o initialisation
47      integer,intent(in) :: flagthermo            ! flag for thermosphere initialisation only
[38]48
[618]49! outputs :
[38]50
[1036]51      real,intent(out) :: pq(iip1,jjp1,llm,nq)  ! advected fields, ie chemical species
[1047]52      real,intent(out) :: qsurf(ngrid,nq)     ! surface values (kg/m2) of tracers
[38]53
[618]54! local :
[38]55
[618]56      integer :: iq, i, j, l, n, nbqchem
57      integer :: count, ierr, dummy
58      real    :: mmean(iip1,jjp1,llm)             ! mean molecular mass (g)
59      real    :: pgcm                             ! pressure at each layer in the gcm (Pa)
[38]60
[618]61      integer, parameter         :: nalt = 252    ! number of altitudes in the initialization files
[655]62      integer                    :: nspe          ! number of species in the initialization files
63      integer, allocatable       :: niq(:)        ! local index of species in initialization files
[618]64      real, dimension(nalt)      :: tinit, zzfile ! temperature in initialization files
65      real, dimension(nalt)      :: pinit         ! pressure in initialization files
66      real, dimension(nalt)      :: densinit      ! total number density in initialization files
[655]67      real, allocatable          :: vmrinit(:,:)  ! mixing ratios in initialization files
68      real, allocatable          :: vmrint(:)     ! mixing ratio interpolated onto the gcm vertical grid
[618]69      real                       :: vmr
[38]70
[618]71      character(len=20)          :: txt           ! to store some text
[655]72      logical                    :: flagnitro     ! checks if N species present
[38]73
[618]74! 1. identify tracers by their names: (and set corresponding values of mmol)
[38]75
[618]76! 1.1 initialize tracer indexes to zero:
[1036]77      nqmx=nq ! initialize value of nqmx
78     
[618]79      do iq = 1,nqmx
80        igcm_dustbin(iq) = 0
81      end do
[38]82
[618]83      igcm_dust_mass   = 0
84      igcm_dust_number = 0
85      igcm_ccn_mass    = 0
86      igcm_ccn_number  = 0
87      igcm_h2o_vap     = 0
88      igcm_h2o_ice     = 0
89      igcm_co2         = 0
90      igcm_co          = 0
91      igcm_o           = 0
92      igcm_o1d         = 0
93      igcm_o2          = 0
94      igcm_o3          = 0
95      igcm_h           = 0
96      igcm_h2          = 0
97      igcm_oh          = 0
98      igcm_ho2         = 0
99      igcm_h2o2        = 0
100      igcm_ch4         = 0
101      igcm_n2          = 0
102      igcm_ar          = 0
103      igcm_ar_n2       = 0
104      igcm_co2plus     = 0
105      igcm_oplus       = 0
106      igcm_o2plus      = 0
107      igcm_coplus      = 0
108      igcm_cplus       = 0
109      igcm_nplus       = 0
110      igcm_noplus      = 0
111      igcm_n2plus      = 0
112      igcm_hplus       = 0
[655]113      igcm_hco2plus    = 0
[618]114      igcm_elec        = 0
[38]115
[618]116! 1.2 find dust tracers
[38]117
[618]118      count = 0
[38]119
[618]120      if (dustbin > 0) then
121         do iq = 1,nqmx
122            txt = " "
123            write(txt,'(a4,i2.2)') 'dust', count + 1
124            if (noms(iq) == txt) then
125               count = count + 1
126               igcm_dustbin(count) = iq
127               mmol(iq) = 100.
128            end if
129         end do !do iq=1,nqmx
130      end if ! of if (dustbin.gt.0)
[38]131
132      if (doubleq) then
[618]133         do iq = 1,nqmx
134            if (noms(iq) == "dust_mass") then
135               igcm_dust_mass = iq
136               count = count + 1
137            end if
138            if (noms(iq) == "dust_number") then
139               igcm_dust_number = iq
140               count = count + 1
141            end if
142         end do
143      end if ! of if (doubleq)
144
145      if (scavenging) then
146         do iq = 1,nqmx
147            if (noms(iq) == "ccn_mass") then
148               igcm_ccn_mass = iq
149               count = count + 1
150            end if
151            if (noms(iq) == "ccn_number") then
152               igcm_ccn_number = iq
153               count = count + 1
154            end if
155         end do
156      end if ! of if (scavenging)
157
158      if (submicron) then
159         do iq=1,nqmx
160            if (noms(iq) == "dust_submicron") then
161               igcm_dust_submicron = iq
162               mmol(iq) = 100.
163               count = count + 1
164            end if
165         end do
166      end if ! of if (submicron)
167
168! 1.3 find chemistry and water tracers
169
170      nbqchem = 0
171      do iq = 1,nqmx
172         if (noms(iq) == "co2") then
173            igcm_co2 = iq
174            mmol(igcm_co2) = 44.
175            count = count + 1
176            nbqchem = nbqchem + 1
177        end if
178        if (noms(iq) == "co") then
179           igcm_co = iq
180           mmol(igcm_co) = 28.
181           count = count + 1
182           nbqchem = nbqchem + 1
183        end if
184        if (noms(iq) == "o") then
185           igcm_o = iq
186           mmol(igcm_o) = 16.
187           count = count + 1
188           nbqchem = nbqchem + 1
189        end if
190        if (noms(iq) == "o1d") then
191           igcm_o1d = iq
192           mmol(igcm_o1d) = 16.
193           count = count + 1
194           nbqchem = nbqchem + 1
195        end if
196        if (noms(iq) == "o2") then
197           igcm_o2 = iq
198           mmol(igcm_o2) = 32.
199           count = count + 1
200           nbqchem = nbqchem + 1
201        end if
202        if (noms(iq) == "o3") then
203           igcm_o3 = iq
204           mmol(igcm_o3) = 48.
205           count = count + 1
206           nbqchem = nbqchem + 1
207        end if
208        if (noms(iq) == "h") then
209           igcm_h = iq
210           mmol(igcm_h) = 1.
211           count = count + 1
212           nbqchem = nbqchem + 1
213        end if
214        if (noms(iq) == "h2") then
215           igcm_h2 = iq
216           mmol(igcm_h2) = 2.
217           count = count + 1
218           nbqchem = nbqchem + 1
219        end if
220        if (noms(iq) == "oh") then
221           igcm_oh = iq
222           mmol(igcm_oh) = 17.
223           count = count + 1
224           nbqchem = nbqchem + 1
225        end if
226        if (noms(iq) == "ho2") then
227           igcm_ho2 = iq
228           mmol(igcm_ho2) = 33.
229           count = count + 1
230           nbqchem = nbqchem + 1
231        end if
232        if (noms(iq) == "h2o2") then
233           igcm_h2o2 = iq
234           mmol(igcm_h2o2) = 34.
235           count = count + 1
236           nbqchem = nbqchem + 1
237        end if
238        if (noms(iq) == "ch4") then
239           igcm_ch4 = iq
240           mmol(igcm_ch4) = 16.
241           count = count + 1
242           nbqchem = nbqchem + 1
243        end if
244        if (noms(iq) == "n2") then
245           igcm_n2 = iq
246           mmol(igcm_n2) = 28.
247           count = count + 1
248           nbqchem = nbqchem + 1
249        end if
250        if (noms(iq) == "n") then
251           igcm_n = iq
252           mmol(igcm_n) = 14.
253           count = count + 1
254           nbqchem = nbqchem + 1
255        end if
256        if (noms(iq) == "n2d") then
257           igcm_n2d = iq
258           mmol(igcm_n2d) = 14.
259           count = count + 1
260           nbqchem = nbqchem + 1
261        end if
262        if (noms(iq) == "no") then
263           igcm_no = iq
264           mmol(igcm_no) = 30.
265           count = count + 1
266           nbqchem = nbqchem + 1
267        end if
268        if (noms(iq) == "no2") then
269           igcm_no2 = iq
270           mmol(igcm_no2) = 46.
271           count = count + 1
272           nbqchem = nbqchem + 1
273        end if
274        if (noms(iq) == "ar") then
275           igcm_ar = iq
276           mmol(igcm_ar) = 40.
277           count = count + 1
278           nbqchem = nbqchem + 1
279        end if
280        if (noms(iq) == "h2o_vap") then
281           igcm_h2o_vap = iq
282           mmol(igcm_h2o_vap) = 18.
283           count = count + 1
284           nbqchem = nbqchem + 1
285        end if
286        if (noms(iq) == "h2o_ice") then
287           igcm_h2o_ice = iq
288           mmol(igcm_h2o_ice) = 18.
289           count = count + 1
290           nbqchem = nbqchem + 1
291        end if
292
293! 1.4 find ions
294
295        if (noms(iq) == "co2plus") then
296           igcm_co2plus = iq
297           mmol(igcm_co2plus) = 44.
298           count = count + 1
299           nbqchem = nbqchem + 1
300        end if
301        if (noms(iq) == "oplus") then
302           igcm_oplus = iq
303           mmol(igcm_oplus) = 16.
304           count = count + 1
305           nbqchem = nbqchem + 1
306        end if
307        if (noms(iq) == "o2plus") then
308           igcm_o2plus = iq
309           mmol(igcm_o2plus) = 32.
310           count = count + 1
311           nbqchem = nbqchem + 1
312        end if
313        if (noms(iq) == "coplus") then
314           igcm_coplus = iq
315           mmol(igcm_coplus) = 28.
316           count = count + 1
317           nbqchem = nbqchem + 1
318        end if
319        if (noms(iq) == "cplus") then
320           igcm_cplus = iq
321           mmol(igcm_cplus) = 12.
322           count = count + 1
323           nbqchem = nbqchem + 1
324        end if
325        if (noms(iq) == "nplus") then
326           igcm_nplus = iq
327           mmol(igcm_nplus) = 14.
328           count = count + 1
329           nbqchem = nbqchem + 1
330        end if
331        if (noms(iq) == "noplus") then
332           igcm_noplus = iq
333           mmol(igcm_noplus) = 30.
334           count = count + 1
335           nbqchem = nbqchem + 1
336        end if
337        if (noms(iq) == "n2plus") then
338           igcm_n2plus = iq
339           mmol(igcm_n2plus) = 28.
340           count = count + 1
341           nbqchem = nbqchem + 1
342        end if
343        if (noms(iq) == "hplus") then
344           igcm_hplus = iq
345           mmol(igcm_hplus) = 1.
346           count = count + 1
347           nbqchem = nbqchem + 1
348        end if
[655]349        if (noms(iq) == "hco2plus") then
350           igcm_hco2plus = iq
351           mmol(igcm_hco2plus) = 45.
352           count = count + 1
353           nbqchem = nbqchem + 1
354        end if
[618]355        if (noms(iq) == "elec") then
356           igcm_elec = iq
357           mmol(igcm_elec) = 1./1822.89
358           count = count + 1
359        end if
360
361! 1.5 find idealized non-condensible tracer
362
363        if (noms(iq) == "Ar_N2") then
364           igcm_ar_n2 = iq
365           mmol(igcm_ar_n2) = 30.
366           count = count + 1
367        end if
368
369      end do ! of do iq=1,nqmx
[38]370     
[618]371! 1.6 check that we identified all tracers:
372
373      if (count /= nqmx) then
374         write(*,*) "inichim_newstart: found only ",count," tracers"
375         write(*,*) "                  expected ",nqmx
376         do iq = 1,count
377            write(*,*) '      ', iq, ' ', trim(noms(iq))
378         end do
379         stop
[38]380      else
[618]381         write(*,*) "inichim_newstart: found all expected tracers"
382         do iq = 1,nqmx
383            write(*,*) '      ', iq, ' ', trim(noms(iq))
384         end do
385      end if
[38]386
[655]387! 1.7 check if nitrogen species are present:
388
389      if(igcm_no == 0) then
390         !check that no N species is in traceur.def
391         if(igcm_n /= 0 .or. igcm_no2 /= 0 .or. igcm_n2d /= 0) then
392            write(*,*)'inichim_newstart error:'
393            write(*,*)'N, NO2 and/or N2D are in traceur.def, but not NO'
394            write(*,*)'stop'
395            stop
396         endif
397         flagnitro = .false.
398         nspe = 14
399      else
400         !check that all N species are in traceur.def
401         if(igcm_n == 0 .or. igcm_no2 == 0 .or. igcm_n2d == 0) then
402            write(*,*)'inichim_newstart error:'
403            write(*,*)'if NO is in traceur.def, N, NO2 and N2D must also be'
404            write(*,*)'stop'
405            stop
406         endif
407         flagnitro = .true.
408         nspe = 18
409      endif
410
411! 1.8 allocate arrays
412
413      allocate(niq(nspe))
414      allocate(vmrinit(nalt,nspe))
415      allocate(vmrint(nspe))
416
[618]417! 2. load in chemistry data for initialization:
[38]418
[618]419! order of major species in initialization file:
420!
421!    1: co2
422!    2: ar
423!    3: n2 
424!    4: o2 
425!    5: co 
426!    6: o   
427!    7: h2
428!
429! order of minor species in initialization file:
430!
431!    1: h 
432!    2: oh
433!    3: ho2
434!    4: h2o
435!    5: h2o2
436!    6: o1d
437!    7: o3
[655]438!
439! order of nitrogen species in initialization file:
440!
441!    1: n
442!    2: no
443!    3: no2
444!    4: n2d
[38]445
[618]446! major species:
[38]447
[618]448      niq(1) = igcm_co2
449      niq(2) = igcm_ar
450      niq(3) = igcm_n2
451      niq(4) = igcm_o2
452      niq(5) = igcm_co
453      niq(6) = igcm_o
454      niq(7) = igcm_h2
[38]455
[618]456! minor species:
[38]457
[618]458      niq(8)  = igcm_h
459      niq(9)  = igcm_oh
460      niq(10) = igcm_ho2
461      niq(11) = igcm_h2o_vap
462      niq(12) = igcm_h2o2
463      niq(13) = igcm_o1d
464      niq(14) = igcm_o3
[38]465
[655]466! nitrogen species:
467
468      if (flagnitro) then
469         niq(15) = igcm_n
470         niq(16) = igcm_no
471         niq(17) = igcm_no2
472         niq(18) = igcm_n2d         
473      end if
474
[618]475! 2.1 open initialization files
[38]476
[655]477      open(210, iostat=ierr,file=trim(datafile)//'/atmosfera_LMD_may.dat')
[618]478      if (ierr /= 0) then
479         write(*,*)'Error : cannot open file atmosfera_LMD_may.dat '
[655]480         write(*,*)'(in aeronomars/inichim_newstart.F)'
[618]481         write(*,*)'It should be in :', trim(datafile),'/'
482         write(*,*)'1) You can change this path in callphys.def with'
483         write(*,*)'   datadir=/path/to/datafiles/'
484         write(*,*)'2) If necessary atmosfera_LMD_may.dat (and others)'
485         write(*,*)'   can be obtained online on:'
[1381]486         write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir'
[618]487         stop
488      end if
[655]489      open(220, iostat=ierr,file=trim(datafile)//'/atmosfera_LMD_min.dat')
[618]490      if (ierr /= 0) then
491         write(*,*)'Error : cannot open file atmosfera_LMD_min.dat '
[655]492         write(*,*)'(in aeronomars/inichim_newstart.F)'
[618]493         write(*,*)'It should be in :', trim(datafile),'/'
494         write(*,*)'1) You can change this path in callphys.def with'
495         write(*,*)'   datadir=/path/to/datafiles/'
496         write(*,*)'2) If necessary atmosfera_LMD_min.dat (and others)'
497         write(*,*)'   can be obtained online on:'
[1381]498         write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir'
[618]499         stop
500      end if
[655]501      if(flagnitro) then
502         open(230, iostat=ierr,file=trim(datafile)//'/atmosfera_LMD_nitr.dat')
503         if (ierr.ne.0) then
504            write(*,*)'Error : cannot open file atmosfera_LMD_nitr.dat '
505            write(*,*)'(in aeronomars/inichim_newstart.F)'
506            write(*,*)'It should be in :', datafile
507            write(*,*)'1) You can change this directory address in '
508            write(*,*)'   file phymars/datafile.h'
509            write(*,*)'2) If necessary atmosfera_LMD_nitr.dat (and others)'
510            write(*,*)'   can be obtained online on:'
[1381]511            write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir'
[655]512            STOP
513         endif
514      endif   ! Of if(flagnitro)
[38]515
[618]516! 2.2 read initialization files
[38]517
[618]518! major species
[38]519
[618]520      read(210,*)
521      do l = 1,nalt
522         read(210,*) dummy, tinit(l), pinit(l), densinit(l), &
523                     (vmrinit(l,n), n = 1,7)
524         pinit(l) = pinit(l)*100.              ! conversion in Pa
525         pinit(l) = log(pinit(l))              ! for the vertical interpolation
526      end do
527      close(210)
[38]528
[618]529! minor species
[38]530
[618]531      read(220,*)
532      do l = 1,nalt
533         read(220,*) dummy, (vmrinit(l,n), n = 8,14)
534      end do
535      close(220)
[38]536
[655]537! nitrogen species
538
539      if (flagnitro) then
540         read(230,*)
541         do l = 1,nalt
542            read(230,*) dummy, (vmrinit(l,n), n = 15,18)
543         end do
544         close(230)
545      end if
546     
[618]547! 3. initialization of tracers
548
549      do i = 1,iip1
550         do j = 1,jjp1
551            do l = 1,llm
552
553               pgcm = aps(l) + bps(l)*ps(i,j)  ! gcm pressure
554               pgcm = log(pgcm)                ! for the vertical interpolation
555               mmean(i,j,l) = 0.
556
557! 3.1 vertical interpolation
558
559               do n = 1,nspe
560                  call intrplf(pgcm,vmr,pinit,vmrinit(:,n),nalt)
561                  vmrint(n) = vmr
562                  iq = niq(n)
563                  mmean(i,j,l) = mmean(i,j,l) + vmrint(n)*mmol(iq)
564               end do
565
566! 3.2 attribute mixing ratio: - all layers or only thermosphere
567!                             - with our without h2o
568
[655]569               if (flagthermo == 0 .or. (flagthermo == 1 .and. exp(pgcm) < 0.1)) then
[618]570                  do n = 1,nspe
571                     iq = niq(n)
572                     if (iq /= igcm_h2o_vap .or. flagh2o == 1) then
573                        pq(i,j,l,iq) = vmrint(n)*mmol(iq)/mmean(i,j,l)
574                     end if
575                  end do
576               end if
577
578            end do
579         end do
580      end do
581
[655]582! set surface values of chemistry tracers to zero
583
[618]584      if (flagthermo == 0) then
[655]585         ! NB: no problem for "surface water vapour" tracer which is always 0
586         do n = 1,nspe
587            iq = niq(n)
[1047]588            qsurf(1:ngrid,iq) = 0.
[655]589         end do
590      end if
[618]591
592! 3.3 initialization of tracers not contained in the initialization files
593
594! methane : 10 ppbv
595
596      if (igcm_ch4 /= 0) then
597         vmr = 10.e-9       
598         do i = 1,iip1
599            do j = 1,jjp1
600               do l = 1,llm
601                  pq(i,j,l,igcm_ch4) = vmr*mmol(igcm_ch4)/mmean(i,j,l)
602               end do
603            end do
604         end do
605         ! set surface value to zero
[1047]606         qsurf(1:ngrid,igcm_ch4) = 0.
[618]607      end if
608
[655]609! ions: 0
610
611      if (igcm_co2plus /= 0) then
612         !check that all required ions are in traceur.def
613         if (igcm_o2plus == 0 .or. igcm_oplus == 0 .or. igcm_coplus == 0          &
614              .or. igcm_cplus == 0 .or. igcm_nplus == 0 .or. igcm_noplus == 0    &
615              .or. igcm_n2plus == 0 .or. igcm_hplus == 0 .or. igcm_hco2plus == 0 &
616              .or. igcm_elec == 0) then
617            write(*,*)'inichim_newstart error:'
618            write(*,*)'if co2plus is in traceur.def, all other ions must also be'
619            write(*,*)'o2plus, oplus, coplus, cplus, nplus, noplus, n2plus'
620            write(*,*)'hplus, hco2plus and elec'
621            write(*,*)'stop'
622            stop
623         end if
624
625         do i = 1,iip1
626            do j = 1,jjp1
627               do l = 1,llm
628                  ! all ions to 0     
629                  pq(i,j,l,igcm_co2plus)  = 0.
630                  pq(i,j,l,igcm_o2plus)   = 0.
631                  pq(i,j,l,igcm_oplus)    = 0.
632                  pq(i,j,l,igcm_coplus)   = 0.
633                  pq(i,j,l,igcm_cplus)    = 0.
634                  pq(i,j,l,igcm_nplus)    = 0.
635                  pq(i,j,l,igcm_noplus)   = 0.
636                  pq(i,j,l,igcm_n2plus)   = 0.
637                  pq(i,j,l,igcm_hplus)    = 0.
638                  pq(i,j,l,igcm_hco2plus) = 0.
639                  pq(i,j,l,igcm_elec)     = 0.
640               end do
641            end do
642         end do
643
644         ! surface value to 0
645
[1047]646         qsurf(1:ngrid,igcm_co2plus)  = 0.
647         qsurf(1:ngrid,igcm_o2plus)   = 0.
648         qsurf(1:ngrid,igcm_oplus)    = 0.
649         qsurf(1:ngrid,igcm_coplus)   = 0.
650         qsurf(1:ngrid,igcm_cplus)    = 0.
651         qsurf(1:ngrid,igcm_nplus)    = 0.
652         qsurf(1:ngrid,igcm_noplus)   = 0.
653         qsurf(1:ngrid,igcm_n2plus)   = 0.
654         qsurf(1:ngrid,igcm_hplus)    = 0.
655         qsurf(1:ngrid,igcm_hco2plus) = 0.
656         qsurf(1:ngrid,igcm_elec)     = 0.
[655]657
658      else
659
660         if (igcm_o2plus /= 0 .or. igcm_oplus /= 0 .or. igcm_coplus /= 0          &
661              .or. igcm_cplus /= 0 .or. igcm_nplus /= 0 .or. igcm_noplus /= 0    &
662              .or. igcm_n2plus /= 0 .or. igcm_hplus /= 0 .or. igcm_hco2plus /= 0 &
663              .or. igcm_elec /= 0) then
664            write(*,*)'inichim_newstart error:'
665            write(*,*)'some ions are in traceur.def, but not co2plus'
666            write(*,*)'stop'
667            stop
668         end if
669      end if    ! of if(igcm_co2 /= 0)
670     
671      ! deallocations
672
673      deallocate(niq)
674      deallocate(vmrinit)
675      deallocate(vmrint)
676
[38]677      end
Note: See TracBrowser for help on using the repository browser.