source: trunk/LMDZ.MARS/libf/phymars/getslopes.F90 @ 3118

Last change on this file since 3118 was 2399, checked in by emillour, 4 years ago

Mars GCM:
More code cleanup: use "call abort_physic()" instead of "stop" or "call abort"
EM

File size: 3.2 KB
RevLine 
[1047]1subroutine getslopes(ngrid,geopot)
[900]2   
[1543]3use geometry_mod, only: longitude, latitude ! in radians
[1047]4use slope_mod, only: theta_sl, psi_sl
[1528]5use comcstfi_h, only: g, rad, pi
6use mod_phys_lmdz_para, only: is_parallel
7use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat
[900]8implicit none
9
10
11! This routine computes slope inclination and orientation for the GCM (callslope=.true. in callphys.def)
12! It works fine with a non-regular grid for zoomed simulations.
13! slope inclination angle (deg)  0 == horizontal, 90 == vertical
14! slope orientation angle (deg)  0 == Northward,  90 == Eastward, 180 == Southward, 270 == Westward
15! TN 04/1013
16
[1047]17integer,intent(in) :: ngrid ! nnumber of atmospheric columns
18real,intent(in) :: geopot(ngrid)     ! geopotential on phy grid
[1528]19real topogrid(nbp_lon,nbp_lat) ! topography on lat/lon grid with poles and only one -180/180 point
20real latigrid(nbp_lon,nbp_lat),longgrid(nbp_lon,nbp_lat) ! meshgrid of latitude and longitude values (radians)
[900]21real theta_val ! slope inclination
22real psi_val   ! slope orientation
[1528]23real gradx(nbp_lon,nbp_lat) ! x: latitude-wise topography gradient,  increasing northward
24real grady(nbp_lon,nbp_lat) ! y: longitude-wise topography gradient, increasing westward
[900]25integer i,j,ig0
[998]26integer id2,idm1 ! a trick to compile testphys1d with debug option
[900]27
[1528]28if (is_parallel) then
29  ! This routine only works in serial mode so stop now.
30  write(*,*) "getslopes Error: this routine is not designed to run in parallel"
[2399]31  call abort_physic("getslopes",'cannot be run in parallel',1)
[1528]32endif
33
[998]34id2  = 2
[1528]35idm1 = nbp_lon-1
[900]36
37! rearrange topography on a 2d array
[1528]38do j=2,nbp_lat-1
39   ig0= 1+(j-2)*nbp_lon
40   do i=1,nbp_lon
[900]41      topogrid(i,j)=geopot(ig0+i)/g
[1541]42      latigrid(i,j)=latitude(ig0+i)
43      longgrid(i,j)=longitude(ig0+i)
[900]44   enddo
45enddo
46!poles :
47topogrid(:,1) = geopot(1)/g
[1541]48latigrid(:,1) = latitude(1)
49longgrid(:,1) = longitude(1)
[1528]50topogrid(:,nbp_lat) = geopot(ngrid)/g
[1541]51latigrid(:,nbp_lat) = latitude(ngrid)
52longgrid(:,nbp_lat) = longitude(ngrid)
[900]53
54
55
56! compute topography gradient
57! topogrid and rad are both in meters
[1528]58do j=2,nbp_lat-1
59   do i=1,nbp_lon
[900]60     gradx(i,j) = (topogrid(i,j+1) - topogrid(i,j-1)) / (latigrid(i,j+1)-latigrid(i,j-1))
61     gradx(i,j) = gradx(i,j) / rad
62   enddo
[1528]63   grady(1,j) = (topogrid(id2,j) - topogrid(nbp_lon,j)) / (2*pi+longgrid(id2,j)-longgrid(nbp_lon,j))
[900]64   grady(1,j) = grady(1,j) / rad
[1528]65   grady(nbp_lon,j) = (topogrid(1,j) - topogrid(idm1,j)) / (2*pi+longgrid(1,j)-longgrid(idm1,j))
66   grady(nbp_lon,j) = grady(nbp_lon,j) / rad
67   do i=2,nbp_lon-1
[900]68     grady(i,j) = (topogrid(i+1,j) - topogrid(i-1,j)) / (longgrid(i+1,j)-longgrid(i-1,j))
69     grady(i,j) = grady(i,j) / rad
70   enddo
71enddo
72! poles :
73gradx(:,1) = 0.
74grady(:,1) = 0.
[1528]75gradx(:,nbp_lat) = 0.
76grady(:,nbp_lat) = 0.
[900]77
78
79
80! compute slope inclination and orientation :
81theta_sl(:) = 0.
82psi_sl(:)   = 0.
[1528]83do j=2,nbp_lat-1
84   do i=1,nbp_lon
[900]85   
[1528]86     ig0= 1+(j-2)*nbp_lon
[900]87   
88     theta_val=atan(sqrt( (gradx(i,j))**2 + (grady(i,j))**2 ))
89     
90     psi_val=0.
91     if (gradx(i,j) .ne. 0.) psi_val= -pi/2. - atan(grady(i,j)/gradx(i,j))
92     if (gradx(i,j) .ge. 0.) psi_val= psi_val - pi
93     psi_val = 3*pi/2. - psi_val
94     psi_val = psi_val*180./pi
95     psi_val = MODULO(psi_val,360.)
96     
97     theta_sl(ig0+i) = theta_val
98     psi_sl(ig0+i)   = psi_val
99         
100   enddo
101enddo
102
103
104
105end subroutine getslopes
Note: See TracBrowser for help on using the repository browser.