Ignore:
Timestamp:
Nov 8, 2011, 3:14:39 PM (13 years ago)
Author:
acolaitis
Message:

Added possibility for Zrecast to work if there is no temperature, but potential temperature.

Temperature is recomputed from potential temperature and pressure field, using constant r/cp = 0.256793.
R can be taken constant for this conversion from teta to temp because this is what is done in the model in the output section of physiq.F

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/util/zrecast.F90

    r281 r360  
    101101real,dimension(:,:,:,:),allocatable :: press ! GCM atmospheric pressure
    102102real,dimension(:,:,:,:),allocatable :: temp ! GCM atmospheric temperature
     103real,dimension(:,:,:,:),allocatable :: teta ! GCM atmospheric potential temperature
    103104real,dimension(:,:,:,:),allocatable :: rho ! GCM atmospheric density
    104105real,dimension(:,:,:,:),allocatable :: za_gcm ! GCM above areoid levels (m)
     
    773774endif
    774775
     776allocate(teta(lonlength,latlength,altlength,timelength),stat=ierr)
     777if (ierr.ne.0) then
     778  write(*,*) "Error: Failed to allocate teta(lonlength,latlength,altlength,timelength)"
     779  write(*,*) "       lonlength=",lonlength," latlength=",latlength
     780  write(*,*) "       altlength=",altlength," timelength=",timelength
     781  stop
     782endif
     783
     784
    775785ierr=NF_INQ_VARID(infid,"temp",tmpvarid)
    776786if (ierr.ne.NF_NOERR) then
     
    779789  ierr=NF_INQ_VARID(infid,"t",tmpvarid)
    780790  if (ierr.ne.NF_NOERR) then
    781     stop "Error: Failed to get t ID"
    782   else
    783     ierr=NF_GET_VAR_REAL(infid,tmpvarid,temp)
    784     if (ierr.ne.NF_NOERR) then
    785       stop "Error: Failed reading atmospheric temperature"
    786     endif
     791     ierr=NF_INQ_VARID(infid,"teta",tmpvarid)
     792     if (ierr.ne.NF_NOERR) then
     793       stop "Error: Failed to get t or teta ID"
     794     endif
     795       ierr=NF_GET_VAR_REAL(infid,tmpvarid,teta)
     796       if (ierr.ne.NF_NOERR) then
     797         stop "Error: Failed reading atmospheric temperature"
     798       endif
     799     
     800  do itim=1,timelength
     801    do ilev=1,altlength
     802      do ilat=1,latlength
     803        do ilon=1,lonlength
     804          temp(ilon,ilat,ilev,itim)=teta(ilon,ilat,ilev,itim)*(press(ilon,ilat,ilev,itim)/ps(ilon,ilat,itim))**(.256793)
     805        enddo
     806      enddo
     807    enddo
     808  enddo
     809
    787810  endif
    788811else
Note: See TracChangeset for help on using the changeset viewer.