Changeset 953 for trunk/LMDZ.VENUS/libf/phyvenus
- Timestamp:
- May 2, 2013, 10:33:18 AM (12 years ago)
- Location:
- trunk/LMDZ.VENUS/libf/phyvenus
- Files:
-
- 1 deleted
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.VENUS/libf/phyvenus/lw_venus_ve.F
r892 r953 17 17 C 18 18 c This routine uses the NER matrix 19 c (computed for a given cell and temp profile in radlwsw, 20 c from the initial matrixes computed in load_psi) 19 c (computed for a given cell and temp profile in radlwsw) 21 20 c to compute cooling rates and radiative fluxes. 22 21 c -
trunk/LMDZ.VENUS/libf/phyvenus/lwi.F
r892 r953 1 1 subroutine lwi(nl,netrad,deltapsi,deltap,temp,coolrate) 2 2 3 c !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!4 c §§§§!!! VERSION utilisable avec load_psi,5 c differente des versions *.1mat6 c !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!7 3 use dimphy 8 4 implicit none -
trunk/LMDZ.VENUS/libf/phyvenus/newstart.F
r927 r953 134 134 integer, dimension(4) :: start,counter 135 135 REAL phisinverse(iip1,jjp1) ! geopotentiel au sol avant inversion 136 logical topoflag,albedoflag 136 logical topoflag,albedoflag,razvitu 137 137 real albedo 138 138 … … 1040 1040 teta(iip1,:,:) = teta(1,:,:) 1041 1041 1042 ! RESETING U TO 0: may be done through run.def 1043 razvitu = . FALSE . 1044 CALL getin('razvitu',razvitu) 1045 1042 1046 c calcul des champ de vent; passage en vent covariant 1043 1047 write (*,*) 'uold ', uold (1,2,1) ! INFO … … 1056 1060 & rlonuold,rlatvold,rlonu,rlatv) 1057 1061 call scal_wind(us,vs,unat,vnat) 1062 ! Reseting u=0 1063 if (razvitu) then 1064 unat(:,:,:) = 0. 1065 endif 1058 1066 write (*,*) 'unat ', unat (1,2,1) ! INFO 1059 1067 do l=1,llm -
trunk/LMDZ.VENUS/libf/phyvenus/physiq.F
r892 r953 33 33 c lafin---input-L-variable logique indiquant le dernier passage 34 34 c rjour---input-R-numero du jour de l'experience 35 c gmtime--input-R- temps universel dans la journee (0 a RDAY s)35 c gmtime--input-R-fraction de la journee (0 a 1) 36 36 c pdtphys-input-R-pas d'integration pour la physique (seconde) 37 37 c paprs---input-R-pression pour chaque inter-couche (en Pa) … … 1087 1087 c==================================================================== 1088 1088 if (ballons.eq.1) then 1089 CALL ballon(30,pdtphys,rjourvrai,gmtime ,rlatd,rlond,1089 CALL ballon(30,pdtphys,rjourvrai,gmtime*RDAY,rlatd,rlond, 1090 1090 c C t,pplay,u,v,pphi) ! alt above surface (smoothed for GCM) 1091 1091 C t,pplay,u,v,zphi) ! alt above planet average radius … … 1119 1119 ENDDO 1120 1120 1121 CALL aaam_bud (27,klon,klev,rjourvrai,gmtime ,1121 CALL aaam_bud (27,klon,klev,rjourvrai,gmtime*RDAY, 1122 1122 C ra,rg,romega, 1123 1123 C rlatd,rlond,pphis, -
trunk/LMDZ.VENUS/libf/phyvenus/radlwsw.F
r892 r953 57 57 REAL swnet(klon,klev+1),lwnet(klon,klev+1) 58 58 c 59 INTEGER k, kk, i, j, nb_gr59 INTEGER k, kk, i, j, band 60 60 c 61 61 REAL PPB(klev+1) … … 71 71 cIM END 72 72 real,save,allocatable :: ksive(:,:,:,:) ! ksi matrixes in Vincent's file 73 real psimap(0:klev+1,0:klev+1,klon)74 real deltapsimap(0:klev+1,0:klev+1,klon) 73 real,save,allocatable :: ztop(:) ! in km 74 75 75 real psi(0:klev+1,0:klev+1) 76 76 real deltapsi(0:klev+1,0:klev+1) 77 77 real latdeg 78 real pt0(klon,0:klev+1) 79 real,save,allocatable :: ztop(:) ! in km 78 real pt0(0:klev+1) 79 real bplck(0:klev+1,nnuve) ! Planck luminances in table layers 80 real y(0:klev,nnuve) ! intermediaire Planck 81 real zdblay(0:klev+1,nnuve) ! gradient en temperature de planck 82 integer mat,mat0 83 real factp,factz,ksi 80 84 81 85 logical firstcall … … 83 87 save firstcall 84 88 85 c-------------------------------------------86 nb_gr = klon87 89 c------------------------------------------- 88 90 c Initialisations … … 116 118 117 119 endif ! firstcall 118 119 DO i = 1, klon 120 pt0(i,0) = tsol(i) 121 DO k = 1, klev 122 pt0(i,k) = t(i,k) 123 ENDDO 124 pt0(i,klev+1) = 0. 125 ENDDO !i 126 127 call load_psi(paprs(:,1),ztop,ksive,pt0,psimap,deltapsimap) 120 c------------------------------------------- 128 121 129 122 DO k = 1, klev 130 DO i = 1, klon123 DO i = 1, klon 131 124 heat(i,k)=0. 132 125 cool(i,k)=0. 126 ENDDO 133 127 ENDDO 134 ENDDO 135 c 128 136 129 c+++++++ BOUCLE SUR LA GRILLE +++++++++++++++++++++++++ 137 DO 99999 j = 1, nb_gr 138 130 DO 99999 j = 1, klon 131 132 c====================================================================== 133 c Initialisations 134 c --------------- 135 139 136 DO k = 1, klev 140 137 zheat(k) = 0.0 … … 154 151 zrmu0 = rmu0(j) 155 152 156 DO k = 1, klev+1153 DO k = 1, klev+1 157 154 PPB(k) = paprs(j,k)/1.e5 158 ENDDO 155 ENDDO 156 157 pt0(0) = tsol(j) 158 DO k = 1, klev 159 pt0(k) = t(j,k) 160 ENDDO 161 pt0(klev+1) = 0. 159 162 160 DO k = 0,klev+1161 DO i = 0,klev+1162 psi(i,k) = psimap(i,k,j)163 deltapsi(i,k) = deltapsimap(i,k,j)164 ENDDO165 ENDDO163 DO k = 0,klev+1 164 DO i = 0,klev+1 165 psi(i,k) = 0. ! positif quand nrj de i->k 166 deltapsi(i,k) = 0. 167 ENDDO 168 ENDDO 166 169 170 c====================================================================== 171 c Getting psi and deltapsi 172 c ------------------------ 173 174 c Planck function 175 c --------------- 176 do band=1,nnuve 177 do k=0,klev 178 c B(T,l) = al/(exp(bl/T)-1) 179 y(k,band) = exp(bl(band)/pt0(k))-1. 180 bplck(k,band) = al(band)/(y(k,band)) 181 zdblay(k,band)= al(band)*bl(band)*exp(bl(band)/pt0(k))/ 182 . ((pt0(k)*pt0(k))*(y(k,band)*y(k,band))) 183 enddo 184 bplck(klev+1,band) = 0.0 185 zdblay(klev+1,band)= 0.0 186 enddo 187 188 c finding the right matrixes 189 c -------------------------- 190 mat0 = 0 191 do mat=1,nbmat-nbztopve 192 if ( (psurfve(mat).ge.paprs(j,1)) 193 . .and.(psurfve(mat+nbztopve).lt.paprs(j,1)) 194 . .and.(ztopve(mat).lt.ztop(j)) 195 . .and.(ztopve(mat+1).ge.ztop(j)) ) then 196 mat0 = mat 197 c print*,'ig=',j,' mat0=',mat 198 factp = (paprs(j,1) -psurfve(mat)) 199 . /(psurfve(mat+nbztopve)-psurfve(mat)) 200 factz = (ztop(j) -ztopve(mat)) 201 . /(ztopve(mat+1)-ztopve(mat)) 202 exit 203 endif 204 enddo 205 if (mat0.eq.0) then 206 write(*,*) 'Finding the right matrix in radlwsw' 207 print*,'Probleme pour interpolation au point ig=',j 208 print*,'psurf = ',paprs(j,1),' ztop = ',ztop(j) 209 stop 210 endif 211 212 c interpolation of ksi and computation of psi,deltapsi 213 c ---------------------------------------------------- 214 do band=1,nnuve 215 do k=0,klev+1 216 do i=0,klev+1 217 ksi = ksive(i,k,band,mat0)*(1-factz)*(1-factp) 218 . +ksive(i,k,band,mat0+1)*factz *(1-factp) 219 . +ksive(i,k,band,mat0+nbztopve)*(1-factz)*factp 220 . +ksive(i,k,band,mat0+nbztopve+1)*factz *factp 221 psi(i,k) = psi(i,k) + 222 . ksi*(bplck(i,band)-bplck(k,band)) 223 deltapsi(i,k) = deltapsi(i,k) + ksi*zdblay(i,band) 224 enddo 225 enddo 226 enddo 227 167 228 c====================================================================== 168 229 c LW call
Note: See TracChangeset
for help on using the changeset viewer.