Changeset 1124


Ignore:
Timestamp:
Dec 10, 2013, 8:06:53 AM (11 years ago)
Author:
emillour
Message:

Mars GCM:

  • Improved 15um cooling routines, with also better handling of errors.

FGG

Location:
trunk/LMDZ.MARS
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/README

    r1119 r1124  
    19581958- Bug fix: Sun-Mars distance was not correctly taken into account in the
    19591959  solvarmod==1 (daily varying realistic EUV input) case.
     1960
     1961== 10/12/2013 == FGG
     1962- Improved 15um cooling routines, with also better handling of errors.
  • trunk/LMDZ.MARS/libf/phymars/nlte_aux.F

    r769 r1124  
     1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     2! Fast scheme for NLTE cooling rates at 15um by CO2 in a Martian GCM !
     3!                 Version dlvr11_03. 2012.                           !
     4! Software written and provided by IAA/CSIC, Granada, Spain,         !
     5! under ESA contract "Mars Climate Database and Physical Models"     !
     6! Person of contact: Miguel Angel Lopez Valverde  valverde@iaa.es    !
     7!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    18c**********************************************************************
    29
     
    198205c     ***********
    199206
    200       do 1, k=1,nbox_max
     207      do 1, k=1,nbox
    201208         temperatura = xtemp(k)
    202209         if (abs(xtemp(k)-thist(1)).le.0.01d0) then
     
    204211         elseif (abs(xtemp(k)-thist(nhist)).le.0.01d0) then
    205212            temperatura=thist(nhist)
     213         elseif (xtemp(k).lt.thist(1)) then
     214            temperatura=thist(1)
     215            write (*,*) ' WARNING intershphunt/ Too low atmosph Tk:'
     216            write (*,*) ' WARNING      k,xtemp = ', k,xtemp(k)
     217            write (*,*) ' Minimum Tk in histogram used : ', thist(1)
     218         elseif (xtemp(k).gt.thist(nhist)) then
     219            temperatura=thist(nhist)
     220            write (*,*) ' WARNING intershphunt/ Very high atmosph Tk:'
     221            write (*,*) ' WARNING      k,xtemp = ', k,xtemp(k)
     222            write (*,*) ' Max Tk in histogram used : ', thist(nhist)
    206223         endif
    207224         call huntdp ( thist,nhist, temperatura, i )
     
    209226            write (*,*) ' HUNT/ Limits input grid:',
    210227     @           thist(1),thist(nhist)
    211             write (*,*) ' HUNT/ location in new grid:', xtemp(k)
     228            write (*,*) ' HUNT/ location in grid:', xtemp(k)
    212229            stop ' INTERSHP/ Interpolation error. T out of Histogram.'
    213230         endif
     
    253270         elseif (abs(xtemp(k)-thist(nhist)).le.0.01d0) then
    254271            temperatura=thist(nhist)
     272         elseif (xtemp(k).lt.thist(1)) then
     273            temperatura=thist(1)
     274            write (*,*) ' WARNING interstrhunt/ Too low atmosph Tk:'
     275            write (*,*) ' WARNING     k,xtemp(k) = ', k,xtemp(k)
     276            write (*,*) ' Minimum Tk in histogram used : ', thist(1)
     277         elseif (xtemp(k).gt.thist(nhist)) then
     278            temperatura=thist(nhist)
     279            write (*,*) ' WARNING interstrhunt/ Very high atmosph Tk:'
     280            write (*,*) ' WARNING     k,xtemp(k) = ', k,xtemp(k)
     281            write (*,*) ' Max Tk in histogram used : ', thist(nhist)
    255282         endif
    256283         call huntdp ( thist,nhist, temperatura, i )
    257284         if ( i.eq.0 .or. i.eq.nhist ) then
    258             write(*,*)'HUNT/ Limits input grid:',thist(1),thist(nhist)
    259             write(*,*)'HUNT/ location in new grid:',xtemp(k)
     285            write(*,*)'HUNT/ Limits input grid:',
     286     $           thist(1),thist(nhist)
     287            write(*,*)'HUNT/ location in grid:',xtemp(k)
    260288            stop'INTERSTR/1/ Interpolation error. T out of Histogram.'
    261289         endif
     
    271299      elseif (abs(ts-thist(nhist)).le.0.01d0) then
    272300         temperatura=thist(nhist)
     301      elseif (ts.lt.thist(1)) then
     302         temperatura=thist(1)
     303         write (*,*) ' WARNING interstrhunt/ Too low atmosph Tk:'
     304         write (*,*) ' WARNING            ts = ', temperatura
     305         write (*,*) ' Minimum Tk in histogram used : ', thist(1)
     306      elseif (ts.gt.thist(nhist)) then
     307         temperatura=thist(nhist)
     308         write (*,*) ' WARNING interstrhunt/ Very high atmosph Tk:'
     309         write (*,*) ' WARNING            ts = ', temperatura
     310         write (*,*) ' Max Tk in histogram used : ', thist(nhist)
    273311      endif
    274312      call huntdp ( thist,nhist, temperatura, i )
     
    276314         write (*,*) ' HUNT/ Limits input grid:',
    277315     @        thist(1),thist(nhist)
    278          write (*,*) ' HUNT/ location in new grid:', ts
     316         write (*,*) ' HUNT/ location in grid:', ts
    279317         stop ' INTERSTR/2/ Interpolat error. T out of Histogram.'
    280318      endif
     
    23132351
    23142352c     local variables
    2315       integer   j, r, mm
    2316       real*8          xx
     2353      integer         j, r
    23172354
    23182355c     ***************
     
    23502387
    23512388c     local variables
    2352       integer   j, r, mm
     2389      integer         j, r
    23532390      real*8          xx
    23542391
     
    23692406      endif
    23702407
    2371       do j=1,mm_stored(ihist)
     2408      do j=mm_stored(ihist),1,-1
    23722409         read(3,*) thist_stored(ihist,j)
    23732410         do r=1,nbox_stored(ihist)
  • trunk/LMDZ.MARS/libf/phymars/nlte_calc.F

    r757 r1124  
     1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     2! Fast scheme for NLTE cooling rates at 15um by CO2 in a Martian GCM !
     3!                 Version dlvr11_03. 2012.                           !
     4! Software written and provided by IAA/CSIC, Granada, Spain,         !
     5! under ESA contract "Mars Climate Database and Physical Models"     !
     6! Person of contact: Miguel Angel Lopez Valverde  valverde@iaa.es    !
     7!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    18c**********************************************************************
    29c     
     
    1926
    2027c***********************************************************************
    21       subroutine MZESC110 (nl_cts_real, nzy_cts_real)
     28      subroutine MZESC110 (ig,nl_cts_real, nzy_cts_real,ierr,varerr)
    2229c***********************************************************************
    2330
     
    3037c     arguments
    3138      integer     nl_cts_real, nzy_cts_real ! i
     39      integer     ig
    3240
    3341c     old arguments
     
    3644
    3745c     local variables and constants
    38       integer   i, in, ir, iaquiHIST , iaquiZ
    39       integer         ib, isot
    40       real*8            argumento
    41       real*8            tauinf(nl_cts)
    42       real*8            con(nzy_cts), coninf
    43       real*8            c1, c2 , ccc
    44       real*8            t1, t2
    45       real*8            p1, p2
    46       real*8            mr1, mr2
    47       real*8            st1, st2
    48       real*8            c1box(nbox_max), c2box(nbox_max)
    49       real*8            ff      ! to avoid too small numbers
    50       real*8            st, beta, ts
     46      integer   i, iaquiHIST , iaquiZ
     47      integer   isot
     48      real*8    argumento
     49      real*8    tauinf(nl_cts)
     50      real*8    con(nzy_cts), coninf
     51      real*8    c1, c2 , ccc
     52      real*8    t1, t2
     53      real*8    p1, p2
     54      real*8    mr1, mr2
     55      real*8    st1, st2
     56      real*8    c1box(nbox_max), c2box(nbox_max)
     57      real*8    ff      ! to avoid too small numbers
     58      real*8    st, beta, ts
    5159      real*8    tyd(nzy_cts)
    52       real*8            correc
    53       real*8            deltanudbl, deltazdbl
    54       real*8          yy
     60      real*8    correc
     61      real*8    deltanudbl, deltazdbl
     62      real*8    yy
    5563
    5664c     external function
    57       external        we_clean
    58       real*8          we_clean
    59 
    60 c***********************************************************************
    61 
    62 c     
    63       ib = 1
     65      external  we_clean
     66      real*8    we_clean
     67
     68c***********************************************************************
     69      ierr = 0
     70      varerr = 0.d0
     71c     
    6472      beta = 1.8d5
    6573      ibcode1 = '1'
     
    7785         mr_cts(i) = dble(co2y_cts(i)/nty_cts(i))
    7886      end do
    79       coninf = dble( con(nzy_cts_real) /
    80      @     log( con(nzy_cts_real-1) / con(nzy_cts_real) ) )
    81                                 ! Correccion pequeña para la FB, nos la saltamos
    82                                 !call mztf_correccion_cts ( coninf, con, ib )
    83 
     87      if ( con(nzy_cts_real) .le. 0.0d0 ) then
     88         ierr = 33
     89         varerr = con(nzy_cts_real)
     90         return
     91      elseif ( con(nzy_cts_real-1) .le. con(nzy_cts_real) ) then
     92         write (*,*) ' WARNING in MZESC110 '
     93         write (*,*) '    [CO2] growing with altitude at TOA.'
     94         write (*,*) '    [CO2] @ TOA = ', con(nzy_cts_real)
     95         coninf = dble( con(nzy_cts_real) )
     96      else
     97         coninf = dble( con(nzy_cts_real) /
     98     @        log( con(nzy_cts_real-1) / con(nzy_cts_real) ) )
     99      endif
    84100ccc   
    85101      call gethist_03 ( 1 )
     
    187203         argumento = eqw / deltanudbl
    188204         tauinf(i) = dexp( - argumento )
    189 
    190205         if (i.eq.nl_cts_real) then
    191206            taustar11_cts(i) = 0.0d0
     
    223238
    224239c     local variables and constants
    225       integer   i, in, ir, iaquiHIST , iaquiZ
     240      integer         i, in, ir, iaquiHIST , iaquiZ
    226241      integer         ib, isot
    227       real*8            tau(nl,nl), argumento
    228       real*8            tauinf(nl)
    229       real*8            con(nzy), coninf
    230       real*8            c1, c2
    231       real*8            t1, t2
    232       real*8            p1, p2
    233       real*8            mr1, mr2
    234       real*8            st1, st2
    235       real*8            c1box(nbox_max), c2box(nbox_max)
    236       real*8            ff      ! to avoid too small numbers
    237       real*8            tvtbs(nzy)
    238       real*8            st, beta, ts
    239       real*8    zld(nl), zyd(nzy), deltazdbl
    240       real*8            correc
    241       real*8            deltanudbl
     242      real*8          tau(nl,nl), argumento
     243      real*8          tauinf(nl)
     244      real*8          con(nzy), coninf
     245      real*8          c1, c2
     246      real*8          t1, t2
     247      real*8          p1, p2
     248      real*8          mr1, mr2
     249      real*8          st1, st2
     250      real*8          c1box(nbox_max), c2box(nbox_max)
     251      real*8          ff      ! to avoid too small numbers
     252      real*8          tvtbs(nzy)
     253      real*8          st, beta, ts
     254      real*8          zld(nl), zyd(nzy), deltazdbl
     255      real*8          correc
     256      real*8          deltanudbl
    242257      real*8          maxtau, yy
    243258
     
    248263c***********************************************************************
    249264
     265      ierr = 0
     266      varerr = 0.d0
    250267c     
    251268      ib = 1
     
    271288         mr(i) = dble(co2y(i)/nty(i))
    272289      end do
    273       coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     290      if ( con(nzy) .le. 0.0d0 ) then
     291         ierr = 43
     292         varerr = con(nzy)
     293         return
     294      elseif ( con(nzy-1) .le. con(nzy) ) then
     295         write (*,*) ' WARNING in MZTUD110 '
     296         write (*,*) '    [CO2] grows with height at CurtisMatrix top.'
     297         write (*,*) '    [CO2] @ top = ', con(nzy)
     298         coninf = dble( con(nzy) )
     299      else
     300         coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     301      endif
    274302      call mztf_correccion ( coninf, con, ib )
    275303
     
    294322            end do
    295323            call interstrhunt (iaquiHIST, st2,t2,ka,ta)
    296                                 ! Check interpolation errors :
     324            ! Check interpolation errors :
    297325            if (c2.le.0.0d0) then
    298                ierr=15
     326               ierr=45
    299327               varerr=c2
    300328               return
    301329            elseif (p2.le.0.0d0) then
    302                ierr=16
     330               ierr=46
    303331               varerr=p2
    304332               return
    305333            elseif (mr2.le.0.0d0) then
    306                ierr=17
     334               ierr=47
    307335               varerr=mr2
    308336               return
    309337            elseif (t2.le.0.0d0) then
    310                ierr=18
     338               ierr=48
    311339               varerr=t2
    312340               return
    313341            elseif (st2.le.0.0d0) then
    314                ierr=19
     342               ierr=49
    315343               varerr=st2
    316344               return
     
    337365               c1box(kr) = c1 * ka(kr) * deltazdbl
    338366            end do
     367            ! Check interpolation errors :
     368            if (c1.le.0.0d0) then
     369               ierr=75
     370               varerr=c1
     371               return
     372            elseif (p1.le.0.0d0) then
     373               ierr=76
     374               varerr=p1
     375               return
     376            elseif (mr1.le.0.0d0) then
     377               ierr=77
     378               varerr=mr1
     379               return
     380            elseif (t1.le.0.0d0) then
     381               ierr=78
     382               varerr=t1
     383               return
     384            elseif (st1.le.0.0d0) then
     385               ierr=79
     386               varerr=st1
     387               return
     388            endif
     389            !
    339390            c1 = c1 * st1 * deltazdbl
    340391            aa = aa + ( p1*mr1*(c1*ff) + p2*mr2*(c2*ff)) / 2.d0
     
    549600      end do
    550601      if (maxtau .gt. 1.0d0) then
    551          ierr = 13
     602         ierr = 42
    552603         varerr = maxtau
    553604         return
     
    673724                                ! local variables
    674725
    675       real*8 cax1(nl,nl)
    676       real*8 v1(nl), cm_factor, vc_factor
    677       real nuaux1, nuaux2, nuaux3
    678       real*8 faux2,faux3, daux2,daux3
     726      real*8  cax1(nl,nl)
     727      real*8  v1(nl), cm_factor, vc_factor
     728      real    nuaux1, nuaux2, nuaux3
     729      real*8  faux2,faux3, daux2,daux3
     730      real*8  varerr
    679731
    680732      integer i,j,ik,ib
     733      integer ierr     
    681734
    682735************************************************************************
     
    700753         cax1(1:nl,1:nl)=0.d0
    701754!         call zerom (cax1,nl)
    702          call MZTUD121 ( cax1,v1, ib )
     755         call MZTUD121 ( cax1,v1, ib, ierr, varerr )
     756         if (ierr .gt. 0) call ERRORS (ierr,varerr)
    703757
    704758         do i=1,nl
     
    734788
    735789c***********************************************************************
    736       subroutine MZTUD121 ( cf,vc, ib )
     790      subroutine MZTUD121 ( cf,vc, ib, ierr, varerr )
    737791c***********************************************************************
    738792
     
    745799
    746800c     arguments
    747       real*8    cf(nl,nl)       ! o.
    748       real*8            vc(nl)  ! o
    749       integer           ib      ! i
     801      real*8          cf(nl,nl) ! o
     802      real*8          vc(nl)    ! o
     803      integer         ib        ! i
     804      integer         ierr      ! o
     805      real*8          varerr    ! o
    750806
    751807
    752808c     local variables and constants
    753       integer   i, in, ir, iaquiHIST, iaquiZ
    754       integer   nmu
    755       parameter         (nmu = 8)
    756       real*8            tau(nl,nl), argumento, deltazdbl
    757       real*8            tauinf(nl)
    758       real*8            con(nzy), coninf
    759       real*8            c1, c2
    760       real*8            t1, t2
    761       real*8            p1, p2
    762       real*8            mr1, mr2
    763       real*8            st1, st2
    764       real*8            c1box(70), c2box(70)
    765       real*8            ff      ! to avoid too small numbers
    766       real*8            tvtbs(nzy)
    767       real*8            st, beta, ts
    768 
    769       real*8    zld(nl), zyd(nzy)
    770       real*8            correc
    771       real*8    deltanudbl
    772       integer         isot
     809      integer         i, in, ir, iaquiHIST, iaquiZ
     810      integer         isot
     811      real*8          tau(nl,nl), argumento, deltazdbl
     812      real*8          tauinf(nl)
     813      real*8          con(nzy), coninf
     814      real*8          c1, c2
     815      real*8          t1, t2
     816      real*8          p1, p2
     817      real*8          mr1, mr2
     818      real*8          st1, st2
     819      real*8          c1box(nbox_max), c2box(nbox_max)
     820      real*8          ff      ! to avoid too small numbers
     821      real*8          tvtbs(nzy)
     822      real*8          st, beta, ts
     823      real*8          zld(nl), zyd(nzy)
     824      real*8          correc
     825      real*8          deltanudbl
    773826      real*8          yy
    774827
     
    781834 101  format(i1)
    782835c***********************************************************************
     836
     837      ierr = 0
     838      varerr = 0.d0
    783839
    784840c     some values
     
    809865      end do
    810866
    811       coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     867      if ( con(nzy) .le. 0.0d0 ) then
     868         ierr = 83
     869         varerr = con(nzy)
     870         return
     871      elseif ( con(nzy-1) .le. con(nzy) ) then
     872         write (*,*) ' WARNING in MZTUD121 '
     873         write (*,*) '    [CO2] grows with height at CurtisMatrix top.'
     874         write (*,*) '    [CO2] @ top = ', con(nzy)
     875         coninf = dble( con(nzy) )
     876      else
     877         coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     878      endif
    812879      call mztf_correccion ( coninf, con, ib )
    813880
     
    852919               c1box(kr) = c1 * ka(kr) * deltazdbl
    853920            end do
     921            ! Check interpolation errors :
     922            if (c1.le.0.0d0) then
     923               ierr=85
     924               varerr=c1
     925               return
     926            elseif (p1.le.0.0d0) then
     927               ierr=86
     928               varerr=p1
     929               return
     930            elseif (mr1.le.0.0d0) then
     931               ierr=87
     932               varerr=mr1
     933               return
     934            elseif (t1.le.0.0d0) then
     935               ierr=88
     936               varerr=t1
     937               return
     938            elseif (st1.le.0.0d0) then
     939               ierr=89
     940               varerr=st1
     941               return
     942            endif
     943            !
    854944            c1 = c1 * st1 * deltazdbl
    855945            aa = aa + ( p1*mr1*(c1*ff) + p2*mr2*(c2*ff)) / 2.d0
     
    11791269
    11801270c     local variables
    1181       integer   i
     1271      integer         i,ierr
    11821272      real*8          factor0200, factor0220, factor1000
    11831273      real*8          aux_0200(nl), aux2_0200(nl)
    11841274      real*8          aux_0220(nl), aux2_0220(nl)
    11851275      real*8          aux_1000(nl), aux2_1000(nl)
     1276      real*8          varerr
    11861277
    11871278c***********************************************************************
     
    11931284      call zero2v(aux_1000,aux2_1000, nl)
    11941285
    1195       call MZESC121sub (aux_0200,aux2_0200, 2 )
    1196       call MZESC121sub (aux_0220,aux2_0220, 3 )
    1197       call MZESC121sub (aux_1000,aux2_1000, 4 )
     1286      call MZESC121sub (aux_0200,aux2_0200, 2 , ierr, varerr)
     1287      if (ierr .gt. 0) call ERRORS (ierr,varerr)
     1288      call MZESC121sub (aux_0220,aux2_0220, 3 , ierr, varerr)
     1289      if (ierr .gt. 0) call ERRORS (ierr,varerr)
     1290      call MZESC121sub (aux_1000,aux2_1000, 4 , ierr, varerr)
     1291      if (ierr .gt. 0) call ERRORS (ierr,varerr)
    11981292
    11991293      factor0220 = 1.d0
     
    12181312c***********************************************************************
    12191313
    1220       subroutine MZESC121sub (taustar,tauinf, ib )
     1314      subroutine MZESC121sub (taustar,tauinf, ib, ierr, varerr )
    12211315
    12221316c***********************************************************************
     
    12301324
    12311325c     arguments
    1232       real*8            taustar(nl) ! o
    1233       real*8            tauinf(nl) ! o
    1234       integer           ib      ! i
     1326      real*8          taustar(nl) ! o
     1327      real*8          tauinf(nl)  ! o
     1328      integer         ib          ! i
     1329      integer         ierr        ! o
     1330      real*8          varerr      ! o
    12351331
    12361332
    12371333c     local variables and constants
    1238       integer   i, iaquiHIST, iaquiZ, isot
    1239       real*8            con(nzy), coninf
    1240       real*8            c1, c2, ccc
    1241       real*8            t1, t2
    1242       real*8            p1, p2
    1243       real*8            mr1, mr2
    1244       real*8            st1, st2
    1245       real*8            c1box(70), c2box(70)
    1246       real*8            ff      ! to avoid too small numbers
    1247       real*8            tvtbs(nzy)
    1248       real*8            st, beta, ts
    1249       real*8    zld(nl), zyd(nzy)
    1250       real*8            correc
    1251       real*8            deltanudbl, deltazdbl
     1334      integer         i, iaquiHIST, iaquiZ, isot
     1335      real*8          con(nzy), coninf
     1336      real*8          c1, c2, ccc
     1337      real*8          t1, t2
     1338      real*8          p1, p2
     1339      real*8          mr1, mr2
     1340      real*8          st1, st2
     1341      real*8          c1box(70), c2box(70)
     1342      real*8          ff      ! to avoid too small numbers
     1343      real*8          tvtbs(nzy)
     1344      real*8          st, beta, ts
     1345      real*8          zld(nl), zyd(nzy)
     1346      real*8          correc
     1347      real*8          deltanudbl, deltazdbl
    12521348      real*8          yy
    12531349
    12541350c     external function
    12551351      external        we_clean
    1256       real*8  we_clean
     1352      real*8          we_clean
    12571353
    12581354c     formats
     
    12611357c***********************************************************************
    12621358
     1359      ierr = 0
     1360      varerr = 0.d0
    12631361c     
    12641362      beta = 1.8d5
     
    12851383         mr(i) = dble(co2y(i)/nty(i))
    12861384      end do
    1287 
    1288       coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     1385      if ( con(nzy) .le. 0.0d0 ) then
     1386         ierr = 63
     1387         varerr = con(nzy)
     1388         return
     1389      elseif ( con(nzy-1) .le. con(nzy) ) then
     1390         write (*,*) ' WARNING in MZESC121sub '
     1391         write (*,*) '    [CO2] grows with height at CurtisMatrix top.'
     1392         write (*,*) '    [CO2] @ top = ', con(nzy)
     1393         coninf = dble( con(nzy) )
     1394      else
     1395         coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     1396      endif
    12891397      call mztf_correccion ( coninf, con, ib )
    12901398
     
    13091417            end do
    13101418            call interstrhunt (iaquiHIST, st2,t2,ka,ta)
     1419            ! Check interpolation errors :
     1420            if (c2.le.0.0d0) then
     1421               ierr=65
     1422               varerr=c2
     1423               return
     1424            elseif (p2.le.0.0d0) then
     1425               ierr=66
     1426               varerr=p2
     1427               return
     1428            elseif (mr2.le.0.0d0) then
     1429               ierr=67
     1430               varerr=mr2
     1431               return
     1432            elseif (t2.le.0.0d0) then
     1433               ierr=68
     1434               varerr=t2
     1435               return
     1436            elseif (st2.le.0.0d0) then
     1437               ierr=69
     1438               varerr=st2
     1439               return
     1440            endif
     1441            !
    13111442            aa = p2 * coninf * mr2 * (st2 * ff)
    13121443            cc = coninf * st2
     
    14011532
    14021533
    1403 !     arguments
    14041534      integer ierr
    14051535      real*8 varerr
     
    14481578
    14491579c***********************************************************************
    1450 c     mztf.f
    1451 c***********************************************************************
    14521580
    14531581      subroutine MZTVC121sub  ( vc, ib,  ierr, varerr )
     
    14631591
    14641592c     arguments
    1465       real*8            vc(nl)  ! o
    1466       integer           ib      ! i
    1467       integer ierr              ! o
    1468       real*8  varerr            ! o
     1593      real*8        vc(nl)  ! o
     1594      integer       ib      ! i
     1595      integer       ierr    ! o
     1596      real*8        varerr  ! o
    14691597
    14701598c     local variables and constants
    1471       integer   i, in, ir, iaquiHIST , iaquiZ, isot
    1472       integer   nmu
    1473       parameter         (nmu = 8)
    1474       real*8            tau(nl,nl), argumento
    1475       real*8            con(nzy), coninf
    1476       real*8            c1, c2
    1477       real*8            t1, t2
    1478       real*8            p1, p2
    1479       real*8            mr1, mr2
    1480       real*8            st1, st2
    1481       real*8            c1box(70), c2box(70)
    1482       real*8            ff      ! to avoid too small numbers
    1483       real*8            tvtbs(nzy)
    1484       real*8            st, beta, ts
    1485       real*8    zld(nl), zyd(nzy), deltazdbl
    1486       real*8            correc
    1487       real*8            deltanudbl, pideltanu,pi
    1488       real*8          yy
    1489       real*8          minvc, maxtau
     1599      integer       i, in, ir, iaquiHIST , iaquiZ, isot
     1600      real*8        tau(nl,nl), argumento
     1601      real*8        con(nzy), coninf
     1602      real*8        c1, c2
     1603      real*8        t1, t2
     1604      real*8        p1, p2
     1605      real*8        mr1, mr2
     1606      real*8        st1, st2
     1607      real*8        c1box(70), c2box(70)
     1608      real*8        ff      ! to avoid too small numbers
     1609      real*8        tvtbs(nzy)
     1610      real*8        st, beta, ts
     1611      real*8        zld(nl), zyd(nzy), deltazdbl
     1612      real*8        correc
     1613      real*8        deltanudbl, pideltanu,pi
     1614      real*8        yy
     1615      real*8        minvc, maxtau
    14901616
    14911617c     external function
    1492       external        we_clean
    1493       real*8          we_clean
     1618      external      we_clean
     1619      real*8        we_clean
    14941620
    14951621c     formats
     
    14971623
    14981624c***********************************************************************
    1499 
     1625     
     1626      ierr = 0
     1627      varerr = 0.d0
    15001628c     
    15011629      pi=3.141592
     
    15271655      end do
    15281656
    1529       coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     1657      if ( con(nzy) .le. 0.0d0 ) then
     1658         ierr = 53
     1659         varerr = con(nzy)
     1660         return
     1661      elseif ( con(nzy-1) .le. con(nzy) ) then
     1662         write (*,*) ' WARNING in MZTVC121sub '
     1663         write (*,*) '    [CO2] grows with height at CurtisMatrix top.'
     1664         write (*,*) '    [CO2] @ top = ', con(nzy)
     1665         coninf = dble( con(nzy) )
     1666      else
     1667         coninf = dble( con(nzy) / log( con(nzy-1) / con(nzy) ) )
     1668      endif
    15301669      call mztf_correccion ( coninf, con, ib )
    15311670
     
    15571696                                ! Check interpolation errors :
    15581697      if (c1.le.0.0d0) then
    1559          ierr=15
     1698         ierr=55
    15601699         varerr=c1
    15611700         return
    15621701      elseif (p1.le.0.0d0) then
    1563          ierr=16
     1702         ierr=56
    15641703         varerr=p1
    15651704         return
    15661705      elseif (mr1.le.0.0d0) then
    1567          ierr=17
     1706         ierr=57
    15681707         varerr=mr1
    15691708         return
    15701709      elseif (t1.le.0.0d0) then
    1571          ierr=18
     1710         ierr=58
    15721711         varerr=t1
    15731712         return
    15741713      elseif (st1.le.0.0d0) then
    1575          ierr=19
     1714         ierr=59
    15761715         varerr=st1
    15771716         return
     
    16561795      end do
    16571796      if (maxtau .gt. 1.0d0) then
    1658          ierr = 13
     1797         ierr = 52
    16591798         varerr = maxtau
    16601799         return
    16611800      else if (minvc .lt. 0.0d0) then
    1662          ierr = 14
     1801         ierr = 51
    16631802         varerr = minvc
    16641803         return
  • trunk/LMDZ.MARS/libf/phymars/nlte_tcool.F

    r1047 r1124  
     1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     2! Fast scheme for NLTE cooling rates at 15um by CO2 in a Martian GCM !
     3!                 Version dlvr11_03. 2012.                           !
     4! Software written and provided by IAA/CSIC, Granada, Spain,         !
     5! under ESA contract "Mars Climate Database and Physical Models"     !
     6! Person of contact: Miguel Angel Lopez Valverde  valverde@iaa.es    !
     7!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    18c**********************************************************************
    29c     
     
    613c     -NLTEdlvr11_CZALU_03
    714c     -NLTEdlvr11_FB626CTS_03
     15c     -NLTEdlvr11_ERRORS
    816c
    917c     
     
    1826     $     p_gcm, t_gcm, z_gcm,
    1927     $     co2vmr_gcm, n2vmr_gcm, covmr_gcm, o3pvmr_gcm,
    20      $     q15umco2_gcm )
     28     $     q15umco2_gcm , ierr, varerr)
    2129
    2230c***********************************************************************
     
    8088
    8189                                ! From GCM's grid to NLTE's grid
    82          call NLTEdlvr11_ZGRID_02 (n_gcm,
     90         call NLTEdlvr11_ZGRID (n_gcm,
    8391     $        p_ig, t_ig, z_ig,
    8492     $        co2_ig, n2_ig, co_ig, o3p_ig,
     
    9199
    92100                                ! Tstar para NLTE-CTS
    93          call MZESC110 ( nl_cts_real, nzy_cts_real )
    94 
    95                                 ! 626FB C.M.
     101         call MZESC110 ( ig,nl_cts_real, nzy_cts_real,ierr,varerr )
     102         if (ierr .gt. 0) call ERRORS (ierr,varerr)
     103
     104         ! 626FB C.M.
    96105         call leetvt
    97106         c110(1:nl,1:nl)=0.d0
     
    99108         call zero2v (vc110,taustar11, nl)
    100109         call MZTUD110 ( ierr, varerr )
    101          if (ierr .gt. 0) goto 900
     110         if (ierr .gt. 0) call ERRORS (ierr,varerr)
    102111
    103112         input_cza = 0
    104          call NLTEdlvr11_CZALU
     113         call NLTEdlvr11_CZALU(ierr,varerr)
     114         if (ierr .gt. 0) call ERRORS (ierr,varerr)
    105115
    106116         input_cza = 1
    107          call NLTEdlvr11_CZALU
     117         call NLTEdlvr11_CZALU(ierr,varerr)
     118         if (ierr .gt. 0) call ERRORS (ierr,varerr)
    108119
    109120                                !  call NLTEdlvr11_FB626CTS
     
    137148                                !   hacemos la media entre hrTotal y hr110CTS :
    138149         i=nl-1
    139          q15umco2_nltot(i) = 0.5*( q15umco2_nltot(i) + hr110CTS(1) )
     150         q15umco2_nltot(i) = 0.5d0*( q15umco2_nltot(i) + hr110CTS(1) )
    140151         do i=2,nl_cts_real
    141152            indice = (nl-2) + i
     
    188199c     end subroutine
    189200      return
    190 
    191 c     Error messages
    192  900  write (*,*) ' ERROR in MZTUD (banda 110).    ierr=',ierr
    193       write (*,*) ' VAR available : ', varerr
    194       return
    195 
    196  901  write (*,*) ' ERROR in MZTVC for vc210.    ierr=',ierr
    197       write (*,*) ' VAR available : ', varerr
    198       return
    199 
    200  902  write (*,*) ' ERROR in MZTVC for vc310.    ierr=',ierr
    201       write (*,*) ' VAR available : ', varerr
    202       return
    203 
    204  903  write (*,*) ' ERROR in MZTVC for vc410.    ierr=',ierr
    205       write (*,*) ' VAR available : ', varerr
    206       return
    207 
    208  904  write (*,*) ' ERROR in mzescape_fb    ierr=',ierr
    209       write (*,*) ' VAR available : ', varerr
    210       return
    211      
    212  930  write (*,*) ' ERROR in mztvc3iso. Temp is NaN'
    213       write (*,*) ' ierr , VAR =', ierr, varerr
    214       if (ierr.eq.30) write (*,*) ' During computation of VC210.'
    215       if (ierr.eq.31) write (*,*) ' During computation of VC310.'
    216       if (ierr.eq.32) write (*,*) ' During computation of VC410.'
    217       return
    218 
    219 c     end subroutine
    220201      end
    221202
     
    223204c***********************************************************************
    224205
    225       subroutine NLTEdlvr11_ZGRID_02 (n_gcm,
     206      subroutine NLTEdlvr11_ZGRID (n_gcm,
    226207     $     p_gcm, t_gcm, z_gcm, co2vmr_gcm, n2vmr_gcm,
    227208     $     covmr_gcm, o3pvmr_gcm, mmean_gcm,cpnew_gcm,
     
    244225      real cpnew_gcm(n_gcm)     ! I
    245226      integer   nl_cts_real, nzy_cts_real ! O
     227      real zaux_gcm(n_gcm)
    246228
    247229c     local variables
     
    250232      real  zmin, zmax
    251233      real  mmean_nlte(n_gcm),cpnew_nlte(n_gcm)
     234      real  gg,masa,radio,kboltzman
    252235
    253236c     functions
     
    264247!     Primero, construimos escala z_gcm
    265248
    266 !     z_gcm(1) = zmin_gcm             ! [km]
    267 
    268 !     do iz = 2, n_gcm
    269 !     meanm = ( co2vmr_gcm(iz)*44. + o3pvmr_gcm(iz)*16.
    270 !     @               + n2vmr_gcm(iz)*28. + covmr_gcm(iz)*28. )
    271 !     meanm = meanm / n_avog
    272 !     distancia = ( radio + z_gcm(iz-1) )*1.e5
    273 !     gz = gg * masa / ( distancia * distancia )
    274 !     Hkm = 0.5*( t_gcm(iz)+t_gcm(iz-1) ) / ( meanm * gz )
    275 !     Hkm = kboltzman * Hkm *1e-5                           ! [km]
    276 !     z_gcm(iz) = z_gcm(iz-1) - Hkm * log( p_gcm(iz)/p_gcm(iz-1) )
    277 !     enddo
     249!      zaux_gcm(1) = z_gcm(1)             ! [km]
     250!      gg=6.67259e-8
     251!      masa=6.4163e26
     252!      radio=3390.
     253!      kboltzman=1.381e-16
     254
     255!      do iz = 2, n_gcm
     256!         distancia = ( radio + zaux_gcm(iz-1) )*1.e5
     257!         gz = gg * masa / ( distancia * distancia )
     258!         Hkm = 0.5*( t_gcm(iz)+t_gcm(iz-1) ) /
     259!     $        ( mmean_gcm(iz)/6.023e23 * gz )
     260!         Hkm = kboltzman * Hkm *1e-5 ! [km]
     261!         zaux_gcm(iz) = zaux_gcm(iz-1) -
     262!     $        Hkm * log( p_gcm(iz)/p_gcm(iz-1) )
     263!      enddo
     264     
    278265
    279266!     Segundo, definimos los límites de los 2 modelos de NLTE.
     
    352339
    353340      do i=1,nl
    354          if (t(i) .gt. 400.0) then
     341         if (t(i) .gt. 500.0) then
    355342            write (*,*) '!!!! WARNING    Temp higher than Histogram.'
    356343            write (*,*) ' Histogram will be extrapolated. '
     
    496483
    497484
    498       subroutine NLTEdlvr11_CZALU
     485      subroutine NLTEdlvr11_CZALU(ierr,varerr)
    499486
    500487c***********************************************************************
     
    506493      include 'nlte_paramdef.h'
    507494      include 'nlte_commons.h'
     495
     496
     497c     Arguments
     498
     499      integer ierr
     500      real*8 varerr
    508501
    509502
     
    578571      real*8 pl11, pl12, pl21, pl31, pl41
    579572
     573      real*8 minvt11, minvt21, minvt31, minvt41
    580574
    581575c     local constants and indexes
     
    596590!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!start program
    597591
     592      ierr = 0
     593      varerr = 0.d0
    598594
    599595      call zero4v( aa11, aa21, aa31, aa41, nl)
     
    691687            c121(1:nl,1:nl)=0.d0
    692688            call MZESC121
    693             call MZTVC121
     689            call MZTVC121( ierr,varerr )
     690            if (ierr .gt. 0) call ERRORS (ierr,varerr)
    694691
    695692         endif
     
    814811         end if
    815812
    816 
    817                                 ! For ITT=13,15
    818 
     813 
     814         ! For ITT=13,15
    819815         a21_einst(i) = a2_010_000 * 1.8d0 / 4.d0 * taustar21(i)
    820816         a31_einst(i) = a3_010_000 * 1.8d0 / 4.d0 * taustar31(i)
     
    825821         l41 = l41 + a41_einst(i)
    826822
    827                                 ! For ITT=13
     823         ! For ITT=13
    828824         if (input_cza.ge.1 .and. itt_cza.eq.13) then
    829825            a12_einst(i) = a1_020_010/3.d0 * 1.8d0/4.d0 * taustar12(i)
     
    832828
    833829
    834                                 !
     830         ! Checking for collisional severe errors
     831         if (l11 .le. 0.0d0) then
     832            ierr = 21
     833            varerr = l11
     834            return
     835         elseif (l21 .le. 0.0d0) then
     836            ierr = 22
     837            varerr = l21
     838            return
     839         elseif (l31 .le. 0.0d0) then
     840            ierr = 23
     841            varerr = l31
     842            return
     843         elseif (l41 .le. 0.0d0) then
     844            ierr = 24
     845            varerr = l41
     846            return
     847         endif
     848         if (input_cza.ge.1) then
     849            if (l12 .lt. 0.0d0) then
     850               ierr = 25
     851               varerr = l12
     852               return
     853            endif
     854         endif
     855         !   
    835856
    836857         a11(i) = gamma*nu11**3.d0 * 1.d0/2.d0 * (p11) /
     
    11541175      if (input_cza.lt.1) then
    11551176
     1177         minvt11 = 1.d6
     1178         minvt21 = 1.d6
     1179         minvt31 = 1.d6
     1180         minvt41 = 1.d6
    11561181         do i=1,nl
    11571182            pl11 = el11(i)/( gamma * nu11**3.0d0  * 1.d0/2.d0 /n10(i) )
     
    11661191            hr310(i) = sl310(i) -hplanck*vlight*nu31 *a31_einst(i)*pl31
    11671192            hr410(i) = sl410(i) -hplanck*vlight*nu41 *a41_einst(i)*pl41
     1193
     1194            minvt11 = min( minvt11,vt11(i) )
     1195            minvt21 = min( minvt21,vt21(i) )
     1196            minvt31 = min( minvt31,vt31(i) )
     1197            minvt41 = min( minvt41,vt41(i) )
    11681198         enddo
     1199
     1200         ! Checking for errors in Tvibs
     1201         if (minvt11 .le. 0.d0) then
     1202            ierr = 26
     1203            varerr = minvt11
     1204            return
     1205         elseif (minvt21 .le. 0.d0) then
     1206            ierr = 27
     1207            varerr = minvt21
     1208            return
     1209         elseif (minvt31 .le. 0.d0) then
     1210            ierr = 28
     1211            varerr = minvt31
     1212            return
     1213         elseif (minvt41 .le. 0.d0) then
     1214            ierr = 29
     1215            varerr = minvt41
     1216            return
     1217         endif
    11691218
    11701219         v626t1(1:nl)=vt11(1:nl)
     
    13291378      return                                 
    13301379      end       
     1380
     1381
     1382c     *** Old NLTEdlvr11_ERRORS ***
     1383c     
     1384c***********************************************************************
     1385
     1386
     1387
     1388      subroutine ERRORS (ierr,varerr)
     1389
     1390c***********************************************************************
     1391
     1392      implicit none
     1393
     1394c Arguments
     1395      integer ierr
     1396      real*8 varerr
     1397     
     1398c***************
     1399
     1400      if (ierr .eq. 15) then
     1401         write (*,*) ' ERROR in MZESC110.    ierr=',ierr
     1402         write (*,*) '                VAR available=', varerr
     1403         write (*,*) '   c2 < 0 after INTZHUNT_CTS'
     1404         
     1405      elseif (ierr .eq. 16) then
     1406         write (*,*) ' ERROR in MZESC110.    ierr=',ierr
     1407         write (*,*) '                VAR available=', varerr
     1408         write (*,*) '   p2 < 0 after INTZHUNT_CTS'
     1409         
     1410      elseif (ierr .eq. 17) then
     1411         write (*,*) ' ERROR in MZESC110.    ierr=',ierr
     1412         write (*,*) '                VAR available=', varerr
     1413         write (*,*) '   mr2 < 0 after INTZHUNT_CTS'
     1414         
     1415      elseif (ierr .eq. 18) then
     1416         write (*,*) ' ERROR in MZESC110.    ierr=',ierr
     1417         write (*,*) '                VAR available=', varerr
     1418         write (*,*) '   t2 < 0 after INTZHUNT_CTS'
     1419         
     1420      elseif (ierr .eq. 19) then
     1421         write (*,*) ' ERROR in MZESC110.    ierr=',ierr
     1422         write (*,*) '                VAR available=', varerr
     1423         write (*,*) '   st2 < 0 after INTZHUNT_CTS'
     1424         
     1425      elseif (ierr .eq. 33) then
     1426         write (*,*) ' ERROR in MZESC110.    ierr=',ierr
     1427         write (*,*) '                VAR available=', varerr
     1428         write (*,*) '   [CO2] < 0 at TOA.'
     1429         
     1430      elseif (ierr .eq. 42) then
     1431         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1432         write (*,*) '                VAR available=', varerr
     1433         write (*,*) '   Atmospheric transmittance too large. '
     1434         
     1435      elseif (ierr .eq. 43) then
     1436         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1437         write (*,*) '                VAR available=', varerr
     1438         write (*,*) '   [CO2] < 0 at  CurtisMatrix top.'
     1439         
     1440      elseif (ierr .eq. 45) then
     1441         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1442         write (*,*) '                VAR available=', varerr
     1443         write (*,*) '   c2 < 0 after INTZHUNT'
     1444         
     1445      elseif (ierr .eq. 46) then
     1446         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1447         write (*,*) '                VAR available=', varerr
     1448         write (*,*) '   p2 < 0 after INTZHUNT'
     1449         
     1450      elseif (ierr .eq. 47) then
     1451         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1452         write (*,*) '                VAR available=', varerr
     1453         write (*,*) '   mr2 < 0 after INTZHUNT'
     1454         
     1455      elseif (ierr .eq. 48) then
     1456         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1457         write (*,*) '                VAR available=', varerr
     1458         write (*,*) '   t2 < 0 after INTZHUNT'
     1459         
     1460      elseif (ierr .eq. 49) then
     1461         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1462         write (*,*) '                VAR available=', varerr
     1463         write (*,*) '   st2 < 0 after INTZHUNT'
     1464         
     1465      elseif (ierr .eq. 75) then
     1466         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1467         write (*,*) '                VAR available=', varerr
     1468         write (*,*) '   c1 < 0 after INTZHUNT'
     1469         
     1470      elseif (ierr .eq. 76) then
     1471         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1472         write (*,*) '                VAR available=', varerr
     1473         write (*,*) '   p1 < 0 after INTZHUNT'
     1474         
     1475      elseif (ierr .eq. 77) then
     1476         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1477         write (*,*) '                VAR available=', varerr
     1478         write (*,*) '   mr1 < 0 after INTZHUNT'
     1479         
     1480      elseif (ierr .eq. 78) then
     1481         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1482         write (*,*) '                VAR available=', varerr
     1483         write (*,*) '   t1 < 0 after INTZHUNT'
     1484         
     1485      elseif (ierr .eq. 79) then
     1486         write (*,*) ' ERROR in MZTUD110.    ierr=',ierr
     1487         write (*,*) '                VAR available=', varerr
     1488         write (*,*) '   st1 < 0 after INTZHUNT'
     1489         
     1490      elseif (ierr .eq. 83) then
     1491         write (*,*) ' ERROR in MZTUD121.    ierr=',ierr
     1492         write (*,*) '                VAR available=', varerr
     1493         write (*,*) '   [CO2] < 0 at  CurtisMatrix top.'
     1494         
     1495      elseif (ierr .eq. 85) then
     1496         write (*,*) ' ERROR in MZTUD121.    ierr=',ierr
     1497         write (*,*) '                VAR available=', varerr
     1498         write (*,*) '   c1 < 0 after INTZHUNT'
     1499         
     1500      elseif (ierr .eq. 86) then
     1501         write (*,*) ' ERROR in MZTUD121.    ierr=',ierr
     1502         write (*,*) '                VAR available=', varerr
     1503         write (*,*) '   p1 < 0 after INTZHUNT'
     1504         
     1505      elseif (ierr .eq. 87) then
     1506         write (*,*) ' ERROR in MZTUD121.    ierr=',ierr
     1507         write (*,*) '                VAR available=', varerr
     1508         write (*,*) '   mr1 < 0 after INTZHUNT'
     1509         
     1510      elseif (ierr .eq. 88) then
     1511         write (*,*) ' ERROR in MZTUD121.    ierr=',ierr
     1512         write (*,*) '                VAR available=', varerr
     1513         write (*,*) '   t1 < 0 after INTZHUNT'
     1514         
     1515      elseif (ierr .eq. 89) then
     1516         write (*,*) ' ERROR in MZTUD121.    ierr=',ierr
     1517         write (*,*) '                VAR available=', varerr
     1518         write (*,*) '   st1 < 0 after INTZHUNT'
     1519         
     1520      elseif (ierr .eq. 51) then
     1521         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1522         write (*,*) '                VAR available=', varerr
     1523         write (*,*) '   Ground transmittance vector VC < 0 '
     1524         
     1525      elseif (ierr .eq. 52) then
     1526         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1527         write (*,*) '                VAR available=', varerr
     1528         write (*,*) '   Atmospheric transmittance too large. '
     1529         
     1530      elseif (ierr .eq. 53) then
     1531         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1532         write (*,*) '                VAR available=', varerr
     1533         write (*,*) '   [CO2] < 0 at  CurtisMatrix top.'
     1534         
     1535      elseif (ierr .eq. 55) then
     1536         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1537         write (*,*) '                VAR available=', varerr
     1538         write (*,*) '   c2 < 0 after INTZHUNT'
     1539         
     1540      elseif (ierr .eq. 56) then
     1541         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1542         write (*,*) '                VAR available=', varerr
     1543         write (*,*) '   p2 < 0 after INTZHUNT'
     1544         
     1545      elseif (ierr .eq. 57) then
     1546         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1547         write (*,*) '                VAR available=', varerr
     1548         write (*,*) '   mr2 < 0 after INTZHUNT'
     1549         
     1550      elseif (ierr .eq. 58) then
     1551         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1552         write (*,*) '                VAR available=', varerr
     1553         write (*,*) '   t2 < 0 after INTZHUNT'
     1554         
     1555      elseif (ierr .eq. 59) then
     1556         write (*,*) ' ERROR in MZTVC121.    ierr=',ierr
     1557         write (*,*) '                VAR available=', varerr
     1558         write (*,*) '   st2 < 0 after INTZHUNT'
     1559         
     1560      elseif (ierr .eq. 63) then
     1561         write (*,*) ' ERROR in MZESC121sub.    ierr=',ierr
     1562         write (*,*) '                VAR available=', varerr
     1563         write (*,*) '   [CO2] < 0 at  CurtisMatrix top.'
     1564         
     1565      elseif (ierr .eq. 65) then
     1566         write (*,*) ' ERROR in MZESC121sub.    ierr=',ierr
     1567         write (*,*) '                VAR available=', varerr
     1568         write (*,*) '   c2 < 0 after INTZHUNT'
     1569         
     1570      elseif (ierr .eq. 66) then
     1571         write (*,*) ' ERROR in MZESC121sub.    ierr=',ierr
     1572         write (*,*) '                VAR available=', varerr
     1573         write (*,*) '   p2 < 0 after INTZHUNT'
     1574         
     1575      elseif (ierr .eq. 67) then
     1576         write (*,*) ' ERROR in MZESC121sub.    ierr=',ierr
     1577         write (*,*) '                VAR available=', varerr
     1578         write (*,*) '   mr2 < 0 after INTZHUNT'
     1579         
     1580      elseif (ierr .eq. 68) then
     1581         write (*,*) ' ERROR in MZESC121sub.    ierr=',ierr
     1582         write (*,*) '                VAR available=', varerr
     1583         write (*,*) '   t2 < 0 after INTZHUNT'
     1584         
     1585      elseif (ierr .eq. 69) then
     1586         write (*,*) ' ERROR in MZESC121sub.    ierr=',ierr
     1587         write (*,*) '                VAR available=', varerr
     1588         write (*,*) '   st2 < 0 after INTZHUNT'
     1589         
     1590      elseif (ierr .eq. 21) then
     1591         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1592         write (*,*) '                VAR available=', varerr
     1593         write (*,*) '   l11 < 0 '
     1594         
     1595      elseif (ierr .eq. 22) then
     1596         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1597         write (*,*) '                VAR available=', varerr
     1598         write (*,*) '   l21 < 0 '
     1599         
     1600      elseif (ierr .eq. 23) then
     1601         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1602         write (*,*) '                VAR available=', varerr
     1603         write (*,*) '   l31 < 0 '
     1604         
     1605      elseif (ierr .eq. 24) then
     1606         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1607         write (*,*) '                VAR available=', varerr
     1608         write (*,*) '   l41 < 0 '
     1609         
     1610      elseif (ierr .eq. 25) then
     1611         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1612         write (*,*) '                VAR available=', varerr
     1613         write (*,*) '   l12 < 0 '
     1614         
     1615      elseif (ierr .eq. 26) then
     1616         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1617         write (*,*) '                VAR available=', varerr
     1618         write (*,*) '   Negative vibr.temp   xvt11 < 0 '
     1619         
     1620      elseif (ierr .eq. 27) then
     1621         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1622         write (*,*) '                VAR available=', varerr
     1623         write (*,*) '   Negative vibr.temp   xvt21 < 0 '
     1624         
     1625      elseif (ierr .eq. 28) then
     1626         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1627         write (*,*) '                VAR available=', varerr
     1628         write (*,*) '   Negative vibr.temp   xvt31 < 0 '
     1629         
     1630      elseif (ierr .eq. 29) then
     1631         write (*,*) ' ERROR in CZA.    ierr=',ierr
     1632         write (*,*) '                VAR available=', varerr
     1633         write (*,*) '   Negative vibr.temp   xvt41 < 0 '
     1634         
     1635
     1636      endif
     1637
     1638
     1639      stop ' Stopped in NLTE scheme due to severe error.'
     1640      end
  • trunk/LMDZ.MARS/libf/phymars/physiq.F

    r1114 r1124  
    364364      real ovmr_gcm(ngrid,nlayer)
    365365      real covmr_gcm(ngrid,nlayer)
    366 
     366      integer ierr_nlte
     367      real*8  varerr
    367368
    368369c Variables for PBL
     
    647648                 CALL nlte_tcool(ngrid,nlayer,zplay*9.869e-6,
    648649     $                pt,zzlay,co2vmr_gcm, n2vmr_gcm, covmr_gcm,
    649      $                ovmr_gcm,  zdtnlte )
     650     $                ovmr_gcm,  zdtnlte,ierr_nlte,varerr )
     651                 if(ierr_nlte.gt.0) then
     652                    write(*,*)
     653     $                'WARNING: nlte_tcool output with error message',
     654     $                'ierr_nlte=',ierr_nlte,'varerr=',varerr
     655                    write(*,*)'I will continue anyway'
     656                 endif
    650657
    651658                 zdtnlte(1:ngrid,1:nlayer)=
     
    876883
    877884               zdtdif(ig,l)=zdhdif(ig,l)*zpopsk(ig,l) ! for diagnostic only
    878 
    879885            ENDDO
    880886         ENDDO
Note: See TracChangeset for help on using the changeset viewer.