Changeset 575 for trunk/LMDZ.MARS/libf/dyn3d
- Timestamp:
- Mar 12, 2012, 2:10:08 PM (13 years ago)
- Location:
- trunk/LMDZ.MARS/libf/dyn3d
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/dyn3d/defrun_new.F
r38 r575 457 457 458 458 WRITE(tapeout,*) "" 459 WRITE(tapeout,*) "Sponge: mode0(u=v=0), mode1(u=umoy,v=0), ", 460 & "mode2(u=umoy,v=vmoy)" 459 WRITE(tapeout,*) "Sponge: number of layers over which", 460 & " sponge extends" 461 nsponge=3 ! default value 462 call getin("nsponge",nsponge) 463 WRITE(tapeout,*)" nsponge = ",nsponge 464 465 WRITE(tapeout,*)"" 466 WRITE(tapeout,*)"Sponge mode: (forcing is towards ..." 467 WRITE(tapeout,*)" over upper nsponge layers)" 468 WRITE(tapeout,*)" 0: (h=hmean,u=v=0)" 469 WRITE(tapeout,*)" 1: (h=hmean,u=umean,v=0)" 470 WRITE(tapeout,*)" 2: (h=hmean,u=umean,v=vmean)" 461 471 mode_sponge=2 ! default value 462 472 call getin("mode_sponge",mode_sponge) … … 464 474 465 475 WRITE(tapeout,*) "" 466 WRITE(tapeout,*) "Sponge: hauteur de sponge (km)" 467 hsponge=90.0 ! default value 468 call getin("hsponge",hsponge) 469 WRITE(tapeout,*)" hsponge = ",hsponge 470 471 WRITE(tapeout,*) "" 472 WRITE(tapeout,*) "Sponge: tetasponge (secondes)" 473 tetasponge=50000.0 476 WRITE(tapeout,*) "Sponge: characteristice time scale tetasponge" 477 WRITE(tapeout,*) "(seconds) at topmost layer (time scale then " 478 WRITE(tapeout,*) " doubles with decreasing layer index)." 479 tetasponge=30000.0 474 480 call getin("tetasponge",tetasponge) 475 481 WRITE(tapeout,*)" tetasponge = ",tetasponge -
trunk/LMDZ.MARS/libf/dyn3d/sponge.F
r38 r575 1 1 subroutine sponge(ucov,vcov,h,pext,dt,mode) 2 3 ! Sponge routine: Quench ucov, vcov and potential temperature near the 4 ! top of the model 5 ! Depending on 'mode' relaxation of variables is towards: 6 ! mode = 0 : h -> h_mean , ucov -> 0 , vcov -> 0 7 ! mode = 1 : h -> h_mean , ucov -> ucov_mean , vcov -> 0 8 ! mode >= 2 : h -> h_mean , ucov -> ucov_mean , vcov -> vcov_mean 9 ! Number of layer over which sponge is applied is 'nsponge' (read from def file) 10 ! Time scale for quenching at top level is given by 'tetasponge' (read from 11 ! def file) and doubles as level indexes decrease. 12 2 13 implicit none 3 14 #include "dimensions.h" … … 7 18 #include "comgeom2.h" 8 19 #include "sponge.h" 9 real ucov(iip1,jjp1,llm),vcov(iip1,jjm,llm) 10 real h(iip1,jjp1,llm),pext(iip1,jjp1) 11 real dt 12 real sig_s(llm) !sigma au milieu des couches 13 save sig_s 20 21 ! Arguments: 22 !------------ 23 real,intent(inout) :: ucov(iip1,jjp1,llm) ! covariant zonal wind 24 real,intent(inout) :: vcov(iip1,jjm,llm) ! covariant meridional wind 25 real,intent(inout) :: h(iip1,jjp1,llm) ! potential temperature 26 real,intent(in) :: pext(iip1,jjp1) ! extensive pressure 27 real,intent(in) :: dt ! time step 28 integer,intent(in) :: mode ! sponge mode 14 29 15 30 c Local: 16 31 c ------ 17 32 33 real,save :: sig_s(llm) !sigma au milieu des couches 18 34 REAL vm,um,hm,ptot(jjp1) 19 real cst(llm) 20 integer mode 35 real,save :: cst(llm) 21 36 22 INTEGER l,i,j,l0 37 INTEGER l,i,j 38 integer,save :: l0 ! layer down to which sponge is applied 23 39 24 40 real ssum 25 41 26 42 real echelle,zkm 27 logical firstcall28 save firstcall,cst,l0 29 data firstcall/.true./ 43 logical,save :: firstcall=.true. 44 45 30 46 31 47 if (firstcall) then 32 do l=1,llm 48 49 ! build approximative sigma levels at midlayer 50 do l=1,llm 33 51 sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2. 34 enddo52 enddo 35 53 36 IF (mode.eq.3) then 37 l0=llm-2 38 echelle=10. 54 l0=llm-nsponge+1 39 55 40 PRINT* 41 print*,'sponge mode',mode 42 print*,'hsponge',hsponge 43 print*,'tetasponge n intervient pas' 44 print*,'Constantes de dissipations fixees comme les anglais' 45 print*,'Coeffs pour la couche en eponge' 46 print*,'Z (km) tau' 47 cst(llm)=dt/(0.5*88775) 48 cst(llm-1)=dt/(88775) 49 cst(l0)=dt/(2*88775) 50 do l=l0,llm 51 zkm=-echelle*log(sig_s(l)) 52 print*,zkm,dt/cst(l),cst(l) 53 enddo 54 firstcall=.false. 55 PRINT* 56 ELSE 57 l0=1 58 echelle=10. 59 60 PRINT* 61 print*,'sponge mode',mode 62 print*,'hsponge tetasponge ',hsponge,tetasponge 63 print*,'Coeffs pour la couche en eponge' 64 print*,'Z (km) tau' 65 do l=1,llm 66 zkm=-echelle*log(sig_s(l)) 67 cst(l)=.5*(1.+tanh((zkm-hsponge)/echelle)) 68 cst(l)= max(tetasponge*1.e-15,cst(l)) 69 print*,zkm,1./cst(l)*tetasponge,cst(l)*dt/tetasponge 70 cst(l)=cst(l)*dt/tetasponge 71 enddo 72 firstcall=.false. 73 PRINT* 74 ENDIF 75 endif 56 PRINT* 57 print*,'sponge mode',mode 58 print*,'nsponge tetasponge ',nsponge,tetasponge 59 print*,'Coeffs for the sponge layer' 60 print*,'Z (km) tau cst' 61 do l=llm,l0,-1 62 ! double time scale with every level, starting from the top 63 cst(l)=dt/(tetasponge*2**(llm-l)) 64 enddo 65 66 echelle=10. 67 do l=l0,llm 68 zkm=-echelle*log(sig_s(l)) 69 print*,zkm,dt/cst(l),cst(l) 70 enddo 71 PRINT* 72 73 firstcall=.false. 74 endif ! of if (firstcall) 76 75 77 76 c----------------------------------------------------------------------- … … 83 82 enddo 84 83 85 c temperature potentielle84 c potential temperature 86 85 do l=l0,llm 87 86 do j=1,jjp1 … … 98 97 enddo 99 98 100 c vent zonal99 c zonal wind 101 100 do l=l0,llm 102 101 do j=2,jjm … … 116 115 enddo 117 116 118 c vent meridien117 c meridional wind 119 118 do l=l0,llm 120 119 do j=1,jjm … … 134 133 enddo 135 134 136 RETURN137 135 end -
trunk/LMDZ.MARS/libf/dyn3d/sponge.h
r38 r575 2 2 c INCLUDE 'sponge.h' 3 3 4 COMMON/com_sponge/callsponge,mode_sponge,hsponge,tetasponge 4 COMMON/com_sponge_l/callsponge 5 common/com_sponge_i/mode_sponge,nsponge 6 common/com_sponge_r/tetasponge 5 7 6 LOGICAL callsponge 7 INTEGER mode_sponge 8 REAL hsponge,tetasponge 8 LOGICAL callsponge ! do we use a sponge on upper layers 9 INTEGER mode_sponge ! sponge mode 10 INTEGER nsponge ! number of sponge layers 11 REAL tetasponge ! sponge time scale (s) at topmost layer 9 12 c-----------------------------------------------------------------------
Note: See TracChangeset
for help on using the changeset viewer.