- Timestamp:
- Mar 28, 2016, 5:27:51 PM (9 years ago)
- File:
-
- 1 moved
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/dynlonlat_phylonlat/phymars/iniphysiq_mod.F90
r1520 r1523 1 MODULE iniphysiq_mod 2 3 CONTAINS 4 1 5 subroutine iniphysiq(ii,jj,nlayer,punjours, pdayref,ptimestep, & 2 rlatu,rl onv,aire,cu,cv,&6 rlatu,rlatv,rlonu,rlonv,aire,cu,cv, & 3 7 prad,pg,pr,pcpp,iflag_phys) 4 8 … … 19 23 use infotrac, only : nqtot ! number of advected tracers 20 24 use comgeomfi_h, only: ini_fillgeom 25 use regular_lonlat_mod, only: init_regular_lonlat, & 26 east, west, north, south, & 27 north_east, north_west, & 28 south_west, south_east 21 29 22 30 implicit none … … 33 41 integer,intent(in) :: ii ! number of atmospheric coulumns along longitudes 34 42 integer,intent(in) :: jj ! number of atompsheric columns along latitudes 35 real,intent(in) :: rlatu(jj+1) ! latitudes of the dynamics U grid 36 real,intent(in) :: rlonv(ii+1) ! longitudes of the dynamics V grid 43 real,intent(in) :: rlatu(jj+1) ! latitudes of the physics grid 44 real,intent(in) :: rlatv(jj) ! latitude boundaries of the physics grid 45 real,intent(in) :: rlonv(ii+1) ! longitudes of the physics grid 46 real,intent(in) :: rlonu(ii+1) ! longitude boundaries of the physics grid 37 47 real,intent(in) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2) 38 48 real,intent(in) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u) … … 47 57 character(len=80) :: abort_message 48 58 real :: total_area_phy, total_area_dyn 59 real :: pi 49 60 61 ! boundaries, on global grid 62 real,allocatable :: boundslon_reg(:,:) 63 real,allocatable :: boundslat_reg(:,:) 50 64 51 65 ! global array, on full physics grid: … … 55 69 real,allocatable :: cvfi(:) 56 70 real,allocatable :: airefi(:) 71 72 pi=2.*asin(1.0) 57 73 58 74 IF (nlayer.NE.klev) THEN … … 73 89 ! CALL abort_gcm (modname,abort_message,1) 74 90 !ENDIF 91 92 ! init regular global longitude-latitude grid points and boundaries 93 ALLOCATE(boundslon_reg(ii,2)) 94 ALLOCATE(boundslat_reg(jj+1,2)) 95 96 DO i=1,ii 97 boundslon_reg(i,east)=rlonu(i) 98 boundslon_reg(i,west)=rlonu(i+1) 99 ENDDO 100 101 boundslat_reg(1,north)= PI/2 102 boundslat_reg(1,south)= rlatv(1) 103 DO j=2,jj 104 boundslat_reg(j,north)=rlatv(j-1) 105 boundslat_reg(j,south)=rlatv(j) 106 ENDDO 107 boundslat_reg(jj+1,north)= rlatv(jj) 108 boundslat_reg(jj+1,south)= -PI/2 109 110 ! Write values in module regular_lonlat_mod 111 CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, & 112 boundslon_reg, boundslat_reg) 75 113 76 114 ! Generate global arrays on full physics grid … … 141 179 142 180 end subroutine iniphysiq 181 182 183 END MODULE iniphysiq_mod
Note: See TracChangeset
for help on using the changeset viewer.