Changeset 1707 for LMDZ5/branches/testing/libf/dyn3d
- Timestamp:
- Jan 11, 2013, 10:19:19 AM (12 years ago)
- Location:
- LMDZ5/branches/testing
- Files:
-
- 3 deleted
- 11 edited
- 1 copied
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/branches/testing
- Property svn:mergeinfo changed
/LMDZ5/trunk merged: 1670-1692,1694-1703,1705-1706
- Property svn:mergeinfo changed
-
LMDZ5/branches/testing/libf/dyn3d/calfis.F
r1669 r1707 507 507 . debut_split, !! firstcall 508 508 . lafin_split, !! lastcall 509 . float(day_ini), !! pday <-- day_ini (dans temps.h)509 . jD_cur, !! pday. see leapfrog 510 510 . jH_cur_split, !! ptime "fraction of day" 511 511 . zdt_split, !! ptimestep -
LMDZ5/branches/testing/libf/dyn3d/comconst.h
r1505 r1707 21 21 REAL dtdiss ! (s) time step for the dissipation 22 22 REAL rad ! (m) radius of the planet 23 REAL r ! Gas constant R=8.31 J.K-1.mol-1 24 REAL cpp ! Cp 23 REAL r ! Reduced Gas constant r=R/mu 24 ! with R=8.31.. J.K-1.mol-1, mu: mol mass of atmosphere (kg/mol) 25 REAL cpp ! Specific heat Cp (J.kg-1.K-1) 25 26 REAL kappa ! kappa=R/Cp 26 27 REAL cotot -
LMDZ5/branches/testing/libf/dyn3d/comdissnew.h
r1319 r1707 12 12 13 13 COMMON/comdissnew/ lstardis,nitergdiv,nitergrot,niterh,tetagdiv, & 14 & tetagrot,tetatemp,coefdis 14 & tetagrot,tetatemp,coefdis, vert_prof_dissip 15 15 16 16 LOGICAL lstardis 17 17 INTEGER nitergdiv, nitergrot, niterh 18 19 integer vert_prof_dissip ! vertical profile of horizontal dissipation 20 ! Allowed values: 21 ! 0: rational fraction, function of pressure 22 ! 1: tanh of altitude 23 18 24 REAL tetagdiv, tetagrot, tetatemp, coefdis 19 25 -
LMDZ5/branches/testing/libf/dyn3d/conf_gcm.F
r1665 r1707 14 14 #endif 15 15 USE infotrac, ONLY : type_trac 16 use assert_m, only: assert 17 16 18 IMPLICIT NONE 17 19 c----------------------------------------------------------------------- … … 93 95 CALL getin('lunout', lunout) 94 96 IF (lunout /= 5 .and. lunout /= 6) THEN 95 OPEN(lunout,FILE='lmdz.out') 97 OPEN(UNIT=lunout,FILE='lmdz.out',ACTION='write', 98 & STATUS='unknown',FORM='formatted') 96 99 ENDIF 97 100 … … 173 176 174 177 !Config Key = nsplit_phys 175 !Config Desc = nombre de pas par jour176 !Config Def = 1177 !Config Help = nombre de pas par jour (multiple de iperiod) (178 !Config ici pour dt = 1 min )179 178 nsplit_phys = 1 180 179 CALL getin('nsplit_phys',nsplit_phys) … … 625 624 CALL getin('ok_dyn_ave',ok_dyn_ave) 626 625 627 628 626 write(lunout,*)' #########################################' 629 627 write(lunout,*)' Configuration des parametres du gcm: ' … … 635 633 write(lunout,*)' day_step = ', day_step 636 634 write(lunout,*)' iperiod = ', iperiod 635 write(lunout,*)' nsplit_phys = ', nsplit_phys 637 636 write(lunout,*)' iconser = ', iconser 638 637 write(lunout,*)' iecri = ', iecri … … 805 804 !Config Desc = sortie des transports zonaux dans la dynamique 806 805 !Config Def = n 807 !Config Help = 806 !Config Help = Permet de mettre en route le calcul des transports 808 807 !Config 809 ok_dynzon = .FALSE.810 808 ok_dynzon = .FALSE. 809 CALL getin('ok_dynzon',ok_dynzon) 811 810 812 811 !Config Key = ok_dyn_ins … … 838 837 write(lunout,*)'STOP !!!' 839 838 write(lunout,*)'use_filtre_fft n est pas implemente dans dyn3d' 840 STOP 839 STOP 1 841 840 ENDIF 842 841 … … 848 847 ok_strato=.FALSE. 849 848 CALL getin('ok_strato',ok_strato) 849 850 vert_prof_dissip = merge(1, 0, ok_strato .and. llm==39) 851 CALL getin('vert_prof_dissip', vert_prof_dissip) 852 call assert(vert_prof_dissip == 0 .or. vert_prof_dissip == 1, 853 $ "bad value for vert_prof_dissip") 850 854 851 855 !Config Key = ok_gradsfile -
LMDZ5/branches/testing/libf/dyn3d/fxhyp.F
r1403 r1707 68 68 xzoom = xzoomdeg * pi/180. 69 69 c 70 if (iim==1) then 71 72 rlonm025(1)=-pi/2. 73 rlonv(1)=0. 74 rlonu(1)=pi 75 rlonp025(1)=pi/2. 76 rlonm025(2)=rlonm025(1)+depi 77 rlonv(2)=rlonv(1)+depi 78 rlonu(2)=rlonu(1)+depi 79 rlonp025(2)=rlonp025(1)+depi 80 81 xprimm025(:)=1. 82 xprimv(:)=1. 83 xprimu(:)=1. 84 xprimp025(:)=1. 85 champmin=depi 86 champmax=depi 87 return 88 89 endif 90 70 91 decalx = .75 71 92 IF( grossism.EQ.1..AND.scal180 ) THEN … … 286 307 287 308 309 288 310 IF(ik.EQ.1.and.grossism.EQ.1.) THEN 289 311 xvrai(1) = xvrai(iip1)-depi 290 312 xxprim(1) = xxprim(iip1) 291 313 ENDIF 314 292 315 DO i = 1 , iim 293 316 xlon(i) = xvrai(i) -
LMDZ5/branches/testing/libf/dyn3d/gcm.F
r1665 r1707 405 405 406 406 CALL inidissip( lstardis, nitergdiv, nitergrot, niterh , 407 * tetagdiv, tetagrot , tetatemp 407 * tetagdiv, tetagrot , tetatemp, vert_prof_dissip) 408 408 409 409 c----------------------------------------------------------------------- … … 433 433 ! Physics: 434 434 #ifdef CPP_PHYS 435 CALL iniphysiq(ngridmx,llm,daysec,day_ini,dtphys/nsplit_phys , 436 , latfi,lonfi,airefi,zcufi,zcvfi,rad,g,r,cpp ) 435 CALL iniphysiq(ngridmx,llm,daysec,day_ini,dtphys/nsplit_phys, 436 & latfi,lonfi,airefi,zcufi,zcvfi,rad,g,r,cpp, 437 & iflag_phys) 437 438 #endif 438 439 call_iniphys=.false. … … 457 458 301 FORMAT('1'/,15x,'run du ', i2,'/',i2,'/',i4) 458 459 302 FORMAT('1'/,15x,' au ', i2,'/',i2,'/',i4) 459 #endif460 461 #ifdef CPP_PHYS462 ! Create start file (startphy.nc) and boundary conditions (limit.nc)463 ! for the Earth verstion464 if (iflag_phys>=100) then465 call iniaqua(ngridmx,latfi,lonfi,iflag_phys)466 endif467 460 #endif 468 461 -
LMDZ5/branches/testing/libf/dyn3d/groupe.F
r524 r1707 38 38 integer i,j,l 39 39 40 logical firstcall 41 save firstcall 40 logical firstcall,groupe_ok 41 save firstcall,groupe_ok 42 42 43 43 data firstcall/.true./ 44 data groupe_ok/.true./ 45 46 if (iim==1) then 47 groupe_ok=.false. 48 endif 44 49 45 50 if (firstcall) then 46 if(mod(iim,2**ngroup).ne.0) stop'probleme du nombre ede point' 51 if (groupe_ok) then 52 if(mod(iim,2**ngroup).ne.0) stop'probleme du nombre de point' 53 endif 47 54 firstcall=.false. 48 55 endif 56 49 57 50 58 c Champs 1D … … 52 60 call convflu(pbaru,pbarv,llm,zconvm) 53 61 54 c55 62 call scopy(ijp1llm,zconvm,1,zconvmm,1) 56 63 call scopy(ijmllm,pbarv,1,pbarvm,1) 57 64 58 c 65 if (groupe_ok) then 59 66 call groupeun(jjp1,llm,zconvmm) 60 67 call groupeun(jjm,llm,pbarvm) 61 68 62 69 c Champs 3D 63 64 70 do l=1,llm 65 71 do j=2,jjm … … 74 80 enddo 75 81 enddo 82 83 else 84 pbarum(:,:,:)=pbaru(:,:,:) 85 pbarvm(:,:,:)=pbarv(:,:,:) 86 endif 76 87 77 88 c integration de la convergence de masse de haut en bas ...... -
LMDZ5/branches/testing/libf/dyn3d/inidissip.F90
r1665 r1707 3 3 ! 4 4 SUBROUTINE inidissip ( lstardis,nitergdiv,nitergrot,niterh , & 5 tetagdiv,tetagrot,tetatemp 5 tetagdiv,tetagrot,tetatemp, vert_prof_dissip) 6 6 !======================================================================= 7 7 ! initialisation de la dissipation horizontale … … 25 25 INTEGER,INTENT(in) :: nitergdiv,nitergrot,niterh 26 26 REAL,INTENT(in) :: tetagdiv,tetagrot,tetatemp 27 28 integer, INTENT(in):: vert_prof_dissip 29 ! Vertical profile of horizontal dissipation 30 ! Allowed values: 31 ! 0: rational fraction, function of pressure 32 ! 1: tanh of altitude 27 33 28 34 ! Local variables: … … 167 173 ! -------------------------------------------------- 168 174 169 if ( ok_strato .and. llm==39) then175 if (vert_prof_dissip == 1) then 170 176 do l=1,llm 171 177 pseudoz=8.*log(preff/presnivs(l)) -
LMDZ5/branches/testing/libf/dyn3d/leapfrog.F
r1669 r1707 383 383 jD_cur = jD_ref + day_ini - day_ref + & 384 384 & itau/day_step 385 386 IF (planet_type .eq."generic") THEN 387 ! AS: we make jD_cur to be pday 388 jD_cur = int(day_ini + itau/day_step) 389 ENDIF 390 385 391 jH_cur = jH_ref + start_time + & 386 392 & mod(itau,day_step)/float(day_step) -
LMDZ5/branches/testing/libf/dyn3d/paramet.h
r792 r1707 17 17 INTEGER jcfil,jcfllm 18 18 19 PARAMETER( iip1= iim+1 -1/iim,iip2=iim+2,iip3=iim+3&19 PARAMETER( iip1= iim+1,iip2=iim+2,iip3=iim+3 & 20 20 & ,jjp1=jjm+1-1/jjm) 21 21 PARAMETER( llmp1 = llm+1, llmp2 = llm+2, llmm1 = llm-1 )
Note: See TracChangeset
for help on using the changeset viewer.