[1017] | 1 | module sponge_mod |
---|
| 2 | |
---|
| 3 | implicit none |
---|
| 4 | |
---|
| 5 | ! sponge parameters (set/read via conf_gcm.F) |
---|
| 6 | logical,save :: callsponge ! do we use a sponge on upper layers |
---|
| 7 | integer,save :: mode_sponge ! sponge mode |
---|
| 8 | integer,save :: nsponge ! number of sponge layers |
---|
| 9 | real,save :: tetasponge ! sponge time scale (s) at topmost layer |
---|
| 10 | |
---|
| 11 | |
---|
| 12 | contains |
---|
| 13 | |
---|
| 14 | subroutine sponge(ucov,vcov,h,ps,dt,mode) |
---|
| 15 | |
---|
| 16 | ! Sponge routine: Quench ucov, vcov and potential temperature near the |
---|
| 17 | ! top of the model |
---|
| 18 | ! Depending on 'mode' relaxation of variables is towards: |
---|
| 19 | ! mode = 0 : h -> h_mean , ucov -> 0 , vcov -> 0 |
---|
| 20 | ! mode = 1 : h -> h_mean , ucov -> ucov_mean , vcov -> 0 |
---|
| 21 | ! mode >= 2 : h -> h_mean , ucov -> ucov_mean , vcov -> vcov_mean |
---|
| 22 | ! Number of layer over which sponge is applied is 'nsponge' (read from def file) |
---|
| 23 | ! Time scale for quenching at top level is given by 'tetasponge' (read from |
---|
| 24 | ! def file) and doubles as level indexes decrease. |
---|
| 25 | ! Quenching is modeled as: A(t)=Am+A0exp(-lambda*t) |
---|
| 26 | ! where Am is the zonal average of the field (or zero), and lambda the inverse |
---|
| 27 | ! of the characteristic quenching/relaxation time scale |
---|
| 28 | ! Thus, assuming Am to be time-independent, field at time t+dt is given by: |
---|
| 29 | ! A(t+dt)=A(t)-(A(t)-Am)*(1-exp(-lambda*dt)) |
---|
| 30 | |
---|
[1422] | 31 | USE comvert_mod, ONLY: ap,bp,preff,scaleheight |
---|
[1017] | 32 | |
---|
| 33 | implicit none |
---|
| 34 | #include "dimensions.h" |
---|
| 35 | #include "paramet.h" |
---|
| 36 | #include "comdissip.h" |
---|
| 37 | #include "comgeom2.h" |
---|
| 38 | #include "iniprint.h" |
---|
| 39 | |
---|
| 40 | ! Arguments: |
---|
| 41 | !------------ |
---|
| 42 | real,intent(inout) :: ucov(iip1,jjp1,llm) ! covariant zonal wind |
---|
| 43 | real,intent(inout) :: vcov(iip1,jjm,llm) ! covariant meridional wind |
---|
| 44 | real,intent(inout) :: h(iip1,jjp1,llm) ! potential temperature |
---|
| 45 | ! real,intent(in) :: pext(iip1,jjp1) ! extensive pressure |
---|
| 46 | real,intent(in) :: ps(iip1,jjp1) ! surface pressure (Pa) |
---|
| 47 | real,intent(in) :: dt ! time step |
---|
| 48 | integer,intent(in) :: mode ! sponge mode |
---|
| 49 | |
---|
| 50 | ! Local: |
---|
| 51 | ! ------ |
---|
| 52 | |
---|
| 53 | real,save :: sig_s(llm) !sigma au milieu des couches |
---|
| 54 | REAL vm,um,hm,ptot(jjp1) |
---|
| 55 | real,save :: cst(llm) |
---|
| 56 | real :: pext(iip1,jjp1) ! extensive pressure |
---|
| 57 | |
---|
| 58 | INTEGER l,i,j |
---|
| 59 | integer,save :: l0 ! layer down to which sponge is applied |
---|
| 60 | |
---|
| 61 | real ssum |
---|
| 62 | |
---|
| 63 | real zkm |
---|
| 64 | logical,save :: firstcall=.true. |
---|
| 65 | |
---|
| 66 | |
---|
| 67 | |
---|
| 68 | if (firstcall) then |
---|
| 69 | |
---|
| 70 | ! build approximative sigma levels at midlayer |
---|
| 71 | do l=1,llm |
---|
| 72 | sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2. |
---|
| 73 | enddo |
---|
| 74 | |
---|
| 75 | l0=llm-nsponge+1 |
---|
| 76 | |
---|
| 77 | write(lunout,*) |
---|
| 78 | write(lunout,*)'sponge mode',mode |
---|
| 79 | write(lunout,*)'nsponge tetasponge ',nsponge,tetasponge |
---|
| 80 | write(lunout,*)'Coeffs for the sponge layer' |
---|
| 81 | write(lunout,*)'Z (km) tau cst' |
---|
| 82 | do l=llm,l0,-1 |
---|
| 83 | ! double time scale with every level, starting from the top |
---|
| 84 | cst(l)=1.-exp(-dt/(tetasponge*2**(llm-l))) |
---|
| 85 | enddo |
---|
| 86 | |
---|
| 87 | do l=l0,llm |
---|
| 88 | zkm=-scaleheight*log(sig_s(l)) |
---|
| 89 | print*,zkm,tetasponge*2**(llm-l),cst(l) |
---|
| 90 | enddo |
---|
| 91 | PRINT* |
---|
| 92 | |
---|
| 93 | firstcall=.false. |
---|
| 94 | endif ! of if (firstcall) |
---|
| 95 | |
---|
| 96 | !----------------------------------------------------------------------- |
---|
| 97 | ! compute sponge relaxation: |
---|
| 98 | ! ------------------------- |
---|
| 99 | |
---|
| 100 | pext(1:iip1,1:jjp1)=ps(1:iip1,1:jjp1)*aire(1:iip1,1:jjp1) |
---|
| 101 | |
---|
| 102 | do j=1,jjp1 |
---|
| 103 | ptot(j)=sum(pext(1:iim,j)) |
---|
| 104 | enddo |
---|
| 105 | |
---|
| 106 | ! potential temperature |
---|
| 107 | do l=l0,llm |
---|
| 108 | do j=1,jjp1 |
---|
| 109 | ! compute zonal average |
---|
| 110 | hm=0. |
---|
| 111 | do i=1,iim |
---|
| 112 | hm=hm+h(i,j,l)*pext(i,j) |
---|
| 113 | enddo |
---|
| 114 | hm=hm/ptot(j) |
---|
| 115 | ! update h() |
---|
| 116 | do i=1,iim |
---|
| 117 | h(i,j,l)=h(i,j,l)-cst(l)*(h(i,j,l)-hm) |
---|
| 118 | enddo |
---|
| 119 | h(iip1,j,l)=h(1,j,l) |
---|
| 120 | enddo |
---|
| 121 | enddo |
---|
| 122 | |
---|
| 123 | ! zonal wind |
---|
| 124 | do l=l0,llm |
---|
| 125 | do j=2,jjm |
---|
| 126 | um=0. |
---|
| 127 | if(mode.ge.1) then ! compute zonal average |
---|
| 128 | do i=1,iim |
---|
| 129 | um=um+0.5*ucov(i,j,l)*(pext(i,j)+pext(i+1,j))/cu(i,j) |
---|
| 130 | enddo |
---|
| 131 | um=um/ptot(j) |
---|
| 132 | endif |
---|
| 133 | ! update ucov() |
---|
| 134 | do i=1,iim |
---|
| 135 | ucov(i,j,l)=ucov(i,j,l)-cst(l)*(ucov(i,j,l)-um*cu(i,j)) |
---|
| 136 | enddo |
---|
| 137 | ucov(iip1,j,l)=ucov(1,j,l) |
---|
| 138 | enddo |
---|
| 139 | enddo |
---|
| 140 | |
---|
| 141 | ! meridional wind |
---|
| 142 | do l=l0,llm |
---|
| 143 | do j=1,jjm |
---|
| 144 | vm=0. |
---|
| 145 | if(mode.ge.2) then ! compute zonal average |
---|
| 146 | do i=1,iim |
---|
| 147 | vm=vm+vcov(i,j,l)*(pext(i,j)+pext(i,j+1))/cv(i,j) |
---|
| 148 | enddo |
---|
| 149 | vm=vm/(ptot(j)+ptot(j+1)) |
---|
| 150 | endif |
---|
| 151 | ! update vcov |
---|
| 152 | do i=1,iim |
---|
| 153 | vcov(i,j,l)=vcov(i,j,l)-cst(l)*(vcov(i,j,l)-vm*cv(i,j)) |
---|
| 154 | enddo |
---|
| 155 | vcov(iip1,j,l)=vcov(1,j,l) |
---|
| 156 | enddo |
---|
| 157 | enddo |
---|
| 158 | |
---|
| 159 | end subroutine sponge |
---|
| 160 | |
---|
| 161 | end module sponge_mod |
---|
| 162 | |
---|