Changeset 2346 for LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd
- Timestamp:
- Aug 21, 2015, 5:13:46 PM (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd/iniphysiq.F90
r2343 r2346 3 3 4 4 5 SUBROUTINE iniphysiq(ii,jj,nlayer,punjours, pdayref,ptimestep, 6 rlatu,rl onv,aire,cu,cv,&5 SUBROUTINE iniphysiq(ii,jj,nlayer,punjours, pdayref,ptimestep, & 6 rlatu,rlatv,rlonu,rlonv,aire,cu,cv, & 7 7 prad,pg,pr,pcpp,iflag_phys) 8 8 USE dimphy, ONLY: klev ! number of atmospheric levels … … 35 35 USE phystokenc_mod, ONLY: init_phystokenc 36 36 USE phyaqua_mod, ONLY: iniaqua 37 USE regular_lonlat_mod, ONLY : init_regular_lonlat, & 38 east, west, north, south, & 39 north_east, north_west, & 40 south_west, south_east 37 41 IMPLICIT NONE 38 42 … … 44 48 include "dimensions.h" 45 49 include "comvert.h" 50 include "comconst.h" 46 51 include "iniprint.h" 47 52 include "temps.h" … … 57 62 INTEGER, INTENT (IN) :: jj ! number of atompsheric columns along latitudes 58 63 REAL, INTENT (IN) :: rlatu(jj+1) ! latitudes of the physics grid 64 REAL, INTENT (IN) :: rlatv(jj) ! latitude boundaries of the physics grid 59 65 REAL, INTENT (IN) :: rlonv(ii+1) ! longitudes of the physics grid 66 REAL, INTENT (IN) :: rlonu(ii+1) ! longitude boundaries of the physics grid 60 67 REAL, INTENT (IN) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2) 61 68 REAL, INTENT (IN) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u) … … 71 78 REAL :: total_area_phy, total_area_dyn 72 79 80 ! boundaries, on global grid 81 REAL,ALLOCATABLE :: boundslon_reg(:,:) 82 REAL,ALLOCATABLE :: boundslat_reg(:,:) 73 83 74 84 ! global array, on full physics grid: … … 87 97 !call init_phys_lmdz(ii,jj+1,llm,1,(/(jj-1)*ii+2/)) 88 98 99 ! init regular global longitude-latitude grid points and boundaries 100 ALLOCATE(boundslon_reg(ii,2)) 101 ALLOCATE(boundslat_reg(jj+1,2)) 102 103 DO i=1,ii 104 boundslon_reg(i,east)=rlonu(i) 105 boundslon_reg(i,west)=rlonu(i+1) 106 ENDDO 107 108 boundslat_reg(1,north)= PI/2 109 boundslat_reg(1,south)= rlatv(1) 110 DO j=2,jj 111 boundslat_reg(j,north)=rlatv(j-1) 112 boundslat_reg(j,south)=rlatv(j) 113 ENDDO 114 boundslat_reg(jj+1,north)= rlatv(jj) 115 boundslat_reg(jj+1,south)= -PI/2 116 117 ! Write values in module regular_lonlat_mod 118 CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, & 119 boundslon_reg, boundslat_reg) 120 89 121 ! Generate global arrays on full physics grid 90 122 ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
Note: See TracChangeset
for help on using the changeset viewer.