Changeset 826 for trunk/LMDZ.TITAN/Tools


Ignore:
Timestamp:
Oct 30, 2012, 12:18:36 PM (12 years ago)
Author:
slebonnois
Message:

SL: amelioration angmom.F90 (Venus, Titan Tools) + correction bug GW non-oro (Venus)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.TITAN/Tools/angmom.F90

    r819 r826  
    1616!                                       or if simple friction
    1717! tdyn  -- 1D -- Dynamics torque IF dudyn is present
     18! tajs  -- 1D -- Torque from convective adjustment IF dudyn is present
    1819! tgwo  -- 1D -- Orographic GW torque IF dugwo is present
    1920! tgwno -- 1D -- Non-Orographic GW torque IF dugwno is present
     
    2324! - surface pressure and surface geopotential
    2425! - zonal wind
    25 ! Optional: dudyn, duvdf, dugwo, dugwno (acceleration terms from physiq param)
     26! Optional: dudyn, duvdf, duajs, dugwo, dugwno (acceleration terms from physiq param)
    2627
    2728implicit none
     
    6768real,dimension(:,:,:,:),allocatable :: duvdf  ! Friction in BL
    6869real,dimension(:,:,:,:),allocatable :: dudyn  ! Dynamics
     70real,dimension(:,:,:,:),allocatable :: duajs  ! Convective adjustment
    6971real,dimension(:,:,:,:),allocatable :: dugwo  ! Orographic Gravity Waves
    7072real,dimension(:,:,:,:),allocatable :: dugwno ! Non-Orographic Gravity Waves
     
    8486real,dimension(:),allocatable :: tmou ! mountain torque (kg m2/s2)
    8587real,dimension(:),allocatable :: tdyn ! dynamics torque (kg m2/s2)
     88real,dimension(:),allocatable :: tajs ! convective adjustment torque (kg m2/s2)
    8689real,dimension(:),allocatable :: tbls ! friction torque (kg m2/s2)
    8790real,dimension(:),allocatable :: tgwo ! oro GW torque (kg m2/s2)
     
    99102integer i,j,ilon,ilat,ilev,itim ! for loops
    100103integer idlsurf ! for option ideal surface
    101 logical :: flag_duvdf,flag_dudyn,flag_dugwo,flag_dugwno,lmdflag,dzflag
     104logical :: flag_duvdf,flag_dudyn,flag_duajs,flag_dugwo,flag_dugwno,lmdflag,dzflag
    102105
    103106real :: deltalat,deltalon ! lat and lon intervals in radians
     
    236239allocate(phis(lonlength,latlength))
    237240
    238 if(lmdflag) then
    239 text="psol"
    240 else
    241241text="PS"
    242 endif
    243242call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2)
    244243if (ierr1.ne.NF_NOERR) then
    245   write(*,*) "  looking for ps instead... "
    246   text="ps"
     244  write(*,*) "  looking for psol instead... "
     245  text="psol"
    247246  call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2)
    248   if (ierr1.ne.NF_NOERR) stop "Error: Failed to get ps ID"
     247  if (ierr1.ne.NF_NOERR) stop "Error: Failed to get psol ID"
    249248endif
    250249if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface pressure"
     
    374373else   ! lmdflag
    375374  print*,"dugwo and dugwno not in CAM simulations"
     375  flag_dugwo = .false.
     376  flag_dugwno = .false.
     377endif  ! lmdflag
     378
     379!===============================================================================
     380! 2.1.65 Accelerations from convective adjustment
     381!===============================================================================
     382allocate(duajs(lonlength,latlength,altlength,timelength))
     383
     384if(lmdflag) then
     385
     386text="duajs"
     387call get_var4d(infid,lonlength,latlength,altlength,timelength,text,duajs,ierr1,ierr2)
     388if (ierr1.ne.NF_NOERR) then
     389  write(*,*) "Failed to get duajs ID"
     390  flag_duajs = .false.
     391else
     392  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading duajs"
     393  flag_duajs = .true.
     394endif
     395
     396else   ! lmdflag
     397  print*,"duajs not in CAM simulations"
     398  flag_duajs = .false.
    376399endif  ! lmdflag
    377400
     
    518541
    519542!===============================================================================
    520 ! 2.2.5 Mountain, friction, dynamics and GW torques
     543! 2.2.5 Mountain, friction, convective adjustment, dynamics and GW torques
    521544!===============================================================================
    522545allocate(tmou(timelength))
    523546if (flag_dudyn)  allocate(tdyn(timelength))
     547if (flag_duajs)  allocate(tajs(timelength))
    524548if (flag_duvdf)  allocate(tbls(timelength))
    525549if (flag_dugwo)  allocate(tgwo(timelength))
     
    559583enddo
    560584  if (tdyn(itim).eq.0.) tdyn(itim) = miss_val
     585endif
     586
     587if (flag_duajs) then
     588tajs(itim) = 0.
     589do ilon=1,lonlength
     590 do ilat=1,latlength
     591  do ilev=1,altlength
     592    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
     593      tajs(itim) = tajs(itim) + duajs(ilon,ilat,ilev,itim) &
     594       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))      &
     595       * dmass(ilon,ilat,ilev,itim) &
     596       / hadley
     597    endif
     598  enddo
     599 enddo
     600enddo
     601  if (tajs(itim).eq.0.) tajs(itim) = miss_val
    561602endif
    562603
     
    682723!===============================================================================
    683724
     725! Check variables to output
     726
     727do itim=1,timelength
     728  if (flag_dudyn .and.( tdyn(itim).eq.miss_val)) flag_dudyn =.false.
     729  if (flag_duajs .and.( tajs(itim).eq.miss_val)) flag_duajs =.false.
     730  if (flag_duvdf .and.( tvdf(itim).eq.miss_val)) flag_duvdf =.false.
     731  if (flag_dugwo .and.( tgwo(itim).eq.miss_val)) flag_dugwo =.false.
     732  if (flag_dugwno.and.(tgwno(itim).eq.miss_val)) flag_dugwno=.false.
     733enddo ! timelength
     734
     735if(dzflag) then
     736do itim=1,timelength
     737  if (tdyndz(itim).eq.miss_val) dzflag=.false.
     738  if (tdisdz(itim).eq.miss_val) dzflag=.false.
     739  if (tspgdz(itim).eq.miss_val) dzflag=.false.
     740  if (tphydz(itim).eq.miss_val) dzflag=.false.
     741enddo ! timelength
     742endif ! dzflag
     743
     744
    684745! 1D Variables
    685746
     
    702763                "tdyn      ", "Dynamics torque     ","E18kgm2s-2",miss_val,&
    703764                 tdyn )
     765endif
     766
     767if (flag_duajs) then
     768call write_var1d(outfid,datashape1d,timelength,&
     769                "tajs      ", "Dynamics torque     ","E18kgm2s-2",miss_val,&
     770                 tajs )
    704771endif
    705772
Note: See TracChangeset for help on using the changeset viewer.