1 | module sponge_mod |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | ! sponge parameters (set/read via defrun_new.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,pext,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*t)) |
---|
30 | |
---|
31 | |
---|
32 | implicit none |
---|
33 | #include "dimensions.h" |
---|
34 | #include "paramet.h" |
---|
35 | #include "comdissip.h" |
---|
36 | #include "comvert.h" |
---|
37 | #include "comgeom2.h" |
---|
38 | |
---|
39 | ! Arguments: |
---|
40 | !------------ |
---|
41 | real,intent(inout) :: ucov(iip1,jjp1,llm) ! covariant zonal wind |
---|
42 | real,intent(inout) :: vcov(iip1,jjm,llm) ! covariant meridional wind |
---|
43 | real,intent(inout) :: h(iip1,jjp1,llm) ! potential temperature |
---|
44 | real,intent(in) :: pext(iip1,jjp1) ! extensive pressure |
---|
45 | real,intent(in) :: dt ! time step |
---|
46 | integer,intent(in) :: mode ! sponge mode |
---|
47 | |
---|
48 | ! Local: |
---|
49 | ! ------ |
---|
50 | |
---|
51 | real,save :: sig_s(llm) !sigma au milieu des couches |
---|
52 | REAL vm,um,hm,ptot(jjp1) |
---|
53 | real,save :: cst(llm) |
---|
54 | |
---|
55 | INTEGER l,i,j |
---|
56 | integer,save :: l0 ! layer down to which sponge is applied |
---|
57 | |
---|
58 | real ssum |
---|
59 | |
---|
60 | real echelle,zkm |
---|
61 | logical,save :: firstcall=.true. |
---|
62 | |
---|
63 | |
---|
64 | |
---|
65 | if (firstcall) then |
---|
66 | |
---|
67 | ! build approximative sigma levels at midlayer |
---|
68 | do l=1,llm |
---|
69 | sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2. |
---|
70 | enddo |
---|
71 | |
---|
72 | l0=llm-nsponge+1 |
---|
73 | |
---|
74 | PRINT* |
---|
75 | print*,'sponge mode',mode |
---|
76 | print*,'nsponge tetasponge ',nsponge,tetasponge |
---|
77 | print*,'Coeffs for the sponge layer' |
---|
78 | print*,'Z (km) tau cst' |
---|
79 | do l=llm,l0,-1 |
---|
80 | ! double time scale with every level, starting from the top |
---|
81 | cst(l)=1.-exp(-dt/(tetasponge*2**(llm-l))) |
---|
82 | enddo |
---|
83 | |
---|
84 | echelle=10. |
---|
85 | do l=l0,llm |
---|
86 | zkm=-echelle*log(sig_s(l)) |
---|
87 | print*,zkm,tetasponge*2**(llm-l),cst(l) |
---|
88 | enddo |
---|
89 | PRINT* |
---|
90 | |
---|
91 | firstcall=.false. |
---|
92 | endif ! of if (firstcall) |
---|
93 | |
---|
94 | !----------------------------------------------------------------------- |
---|
95 | ! compute sponge relaxation: |
---|
96 | ! ------------------------- |
---|
97 | |
---|
98 | do j=1,jjp1 |
---|
99 | ptot(j)=ssum(iim,pext(1,j),1) |
---|
100 | enddo |
---|
101 | |
---|
102 | ! potential temperature |
---|
103 | do l=l0,llm |
---|
104 | do j=1,jjp1 |
---|
105 | ! compute zonal average |
---|
106 | hm=0. |
---|
107 | do i=1,iim |
---|
108 | hm=hm+h(i,j,l)*pext(i,j) |
---|
109 | enddo |
---|
110 | hm=hm/ptot(j) |
---|
111 | ! update h() |
---|
112 | do i=1,iim |
---|
113 | h(i,j,l)=h(i,j,l)-cst(l)*(h(i,j,l)-hm) |
---|
114 | enddo |
---|
115 | h(iip1,j,l)=h(1,j,l) |
---|
116 | enddo |
---|
117 | enddo |
---|
118 | |
---|
119 | ! zonal wind |
---|
120 | do l=l0,llm |
---|
121 | do j=2,jjm |
---|
122 | um=0. |
---|
123 | if(mode.ge.1) then ! compute zonal average |
---|
124 | do i=1,iim |
---|
125 | um=um+0.5*ucov(i,j,l)*(pext(i,j)+pext(i+1,j))/cu(i,j) |
---|
126 | enddo |
---|
127 | um=um/ptot(j) |
---|
128 | endif |
---|
129 | ! update ucov() |
---|
130 | do i=1,iim |
---|
131 | ucov(i,j,l)=ucov(i,j,l)-cst(l)*(ucov(i,j,l)-um*cu(i,j)) |
---|
132 | enddo |
---|
133 | ucov(iip1,j,l)=ucov(1,j,l) |
---|
134 | enddo |
---|
135 | enddo |
---|
136 | |
---|
137 | ! meridional wind |
---|
138 | do l=l0,llm |
---|
139 | do j=1,jjm |
---|
140 | vm=0. |
---|
141 | if(mode.ge.2) then ! compute zonal average |
---|
142 | do i=1,iim |
---|
143 | vm=vm+vcov(i,j,l)*(pext(i,j)+pext(i,j+1))/cv(i,j) |
---|
144 | enddo |
---|
145 | vm=vm/(ptot(j)+ptot(j+1)) |
---|
146 | endif |
---|
147 | ! update vcov |
---|
148 | do i=1,iim |
---|
149 | vcov(i,j,l)=vcov(i,j,l)-cst(l)*(vcov(i,j,l)-vm*cv(i,j)) |
---|
150 | enddo |
---|
151 | vcov(iip1,j,l)=vcov(1,j,l) |
---|
152 | enddo |
---|
153 | enddo |
---|
154 | |
---|
155 | end subroutine sponge |
---|
156 | |
---|
157 | end module sponge_mod |
---|