1 | subroutine soil_pem(ngrid,nsoil,firstcall, & |
---|
2 | therm_i, & |
---|
3 | timestep,tsurf,tsoil,alph_PEM,beta_PEM) |
---|
4 | |
---|
5 | |
---|
6 | use comsoil_h_PEM, only: layer_PEM, mlayer_PEM, & |
---|
7 | mthermdiff_PEM, thermdiff_PEM, coefq_PEM, & |
---|
8 | coefd_PEM, mu_PEM,fluxgeo |
---|
9 | use comsoil_h,only: volcapa |
---|
10 | implicit none |
---|
11 | |
---|
12 | !----------------------------------------------------------------------- |
---|
13 | ! Author: LL |
---|
14 | ! Purpose: Compute soil temperature using an implict 1st order scheme |
---|
15 | ! |
---|
16 | ! Note: depths of layers and mid-layers, soil thermal inertia and |
---|
17 | ! heat capacity are commons in comsoil_PEM.h |
---|
18 | ! A convergence loop is added until equilibrium |
---|
19 | !----------------------------------------------------------------------- |
---|
20 | |
---|
21 | #include "dimensions.h" |
---|
22 | !#include "dimphys.h" |
---|
23 | |
---|
24 | !#include"comsoil.h" |
---|
25 | |
---|
26 | |
---|
27 | !----------------------------------------------------------------------- |
---|
28 | ! arguments |
---|
29 | ! --------- |
---|
30 | ! inputs: |
---|
31 | integer,intent(in) :: ngrid ! number of (horizontal) grid-points |
---|
32 | integer,intent(in) :: nsoil ! number of soil layers |
---|
33 | logical,intent(in) :: firstcall ! identifier for initialization call |
---|
34 | real,intent(in) :: therm_i(ngrid,nsoil) ! thermal inertia |
---|
35 | real,intent(in) :: timestep ! time step |
---|
36 | real,intent(in) :: tsurf(ngrid) ! surface temperature |
---|
37 | |
---|
38 | ! outputs: |
---|
39 | real,intent(inout) :: tsoil(ngrid,nsoil) ! soil (mid-layer) temperature |
---|
40 | real,intent(inout) :: alph_PEM(ngrid,nsoil-1) |
---|
41 | real,intent(inout) :: beta_PEM(ngrid,nsoil-1) |
---|
42 | |
---|
43 | |
---|
44 | |
---|
45 | |
---|
46 | |
---|
47 | ! local variables: |
---|
48 | integer ig,ik |
---|
49 | |
---|
50 | |
---|
51 | |
---|
52 | |
---|
53 | ! 0. Initialisations and preprocessing step |
---|
54 | if (firstcall) then |
---|
55 | |
---|
56 | ! 0.1 Build mthermdiff_PEM(:), the mid-layer thermal diffusivities |
---|
57 | do ig=1,ngrid |
---|
58 | do ik=0,nsoil-1 |
---|
59 | mthermdiff_PEM(ig,ik)=therm_i(ig,ik+1)*therm_i(ig,ik+1)/volcapa |
---|
60 | |
---|
61 | enddo |
---|
62 | enddo |
---|
63 | |
---|
64 | |
---|
65 | ! 0.2 Build thermdiff(:), the "interlayer" thermal diffusivities |
---|
66 | do ig=1,ngrid |
---|
67 | do ik=1,nsoil-1 |
---|
68 | thermdiff_PEM(ig,ik)=((layer_PEM(ik)-mlayer_PEM(ik-1))*mthermdiff_PEM(ig,ik) & |
---|
69 | +(mlayer_PEM(ik)-layer_PEM(ik))*mthermdiff_PEM(ig,ik-1)) & |
---|
70 | /(mlayer_PEM(ik)-mlayer_PEM(ik-1)) |
---|
71 | |
---|
72 | enddo |
---|
73 | enddo |
---|
74 | |
---|
75 | ! 0.3 Build coefficients mu_PEM, q_{k+1/2}, d_k, alph_PEMa_k and capcal |
---|
76 | ! mu_PEM |
---|
77 | mu_PEM=mlayer_PEM(0)/(mlayer_PEM(1)-mlayer_PEM(0)) |
---|
78 | |
---|
79 | ! q_{1/2} |
---|
80 | coefq_PEM(0)=volcapa*layer_PEM(1)/timestep |
---|
81 | ! q_{k+1/2} |
---|
82 | do ik=1,nsoil-1 |
---|
83 | coefq_PEM(ik)=volcapa*(layer_PEM(ik+1)-layer_PEM(ik)) & |
---|
84 | /timestep |
---|
85 | enddo |
---|
86 | |
---|
87 | do ig=1,ngrid |
---|
88 | ! d_k |
---|
89 | do ik=1,nsoil-1 |
---|
90 | coefd_PEM(ig,ik)=thermdiff_PEM(ig,ik)/(mlayer_PEM(ik)-mlayer_PEM(ik-1)) |
---|
91 | enddo |
---|
92 | |
---|
93 | ! alph_PEM_{N-1} |
---|
94 | alph_PEM(ig,nsoil-1)=coefd_PEM(ig,nsoil-1)/ & |
---|
95 | (coefq_PEM(nsoil-1)+coefd_PEM(ig,nsoil-1)) |
---|
96 | ! alph_PEM_k |
---|
97 | do ik=nsoil-2,1,-1 |
---|
98 | alph_PEM(ig,ik)=coefd_PEM(ig,ik)/(coefq_PEM(ik)+coefd_PEM(ig,ik+1)* & |
---|
99 | (1.-alph_PEM(ig,ik+1))+coefd_PEM(ig,ik)) |
---|
100 | enddo |
---|
101 | |
---|
102 | |
---|
103 | enddo ! of do ig=1,ngrid |
---|
104 | |
---|
105 | |
---|
106 | |
---|
107 | endif ! of if (firstcall) |
---|
108 | |
---|
109 | |
---|
110 | |
---|
111 | IF (.not.firstcall) THEN |
---|
112 | ! 2. Compute soil temperatures |
---|
113 | ! First layer: |
---|
114 | |
---|
115 | |
---|
116 | |
---|
117 | |
---|
118 | do ig=1,ngrid |
---|
119 | tsoil(ig,1)=(tsurf(ig)+mu_PEM*beta_PEM(ig,1)* & |
---|
120 | thermdiff_PEM(ig,1)/mthermdiff_PEM(ig,0))/ & |
---|
121 | (1.+mu_PEM*(1.0-alph_PEM(ig,1))*& |
---|
122 | thermdiff_PEM(ig,1)/mthermdiff_PEM(ig,0)) |
---|
123 | |
---|
124 | |
---|
125 | |
---|
126 | |
---|
127 | |
---|
128 | ! Other layers: |
---|
129 | do ik=1,nsoil-1 |
---|
130 | |
---|
131 | tsoil(ig,ik+1)=alph_PEM(ig,ik)*tsoil(ig,ik)+beta_PEM(ig,ik) |
---|
132 | |
---|
133 | |
---|
134 | enddo |
---|
135 | |
---|
136 | enddo |
---|
137 | ENDIF |
---|
138 | |
---|
139 | |
---|
140 | |
---|
141 | |
---|
142 | ! 2. Compute beta_PEM coefficients (preprocessing for next time step) |
---|
143 | ! Bottom layer, beta_PEM_{N-1} |
---|
144 | do ig=1,ngrid |
---|
145 | beta_PEM(ig,nsoil-1)=coefq_PEM(nsoil-1)*tsoil(ig,nsoil) & |
---|
146 | /(coefq_PEM(nsoil-1)+coefd_PEM(ig,nsoil-1)) & |
---|
147 | + fluxgeo/(coefq_PEM(nsoil-1)+coefd_PEM(ig,nsoil-1)) |
---|
148 | enddo |
---|
149 | ! Other layers |
---|
150 | do ik=nsoil-2,1,-1 |
---|
151 | do ig=1,ngrid |
---|
152 | beta_PEM(ig,ik)=(coefq_PEM(ik)*tsoil(ig,ik+1)+ & |
---|
153 | coefd_PEM(ig,ik+1)*beta_PEM(ig,ik+1))/ & |
---|
154 | (coefq_PEM(ik)+coefd_PEM(ig,ik+1)*(1.0-alph_PEM(ig,ik+1)) & |
---|
155 | +coefd_PEM(ig,ik)) |
---|
156 | enddo |
---|
157 | enddo |
---|
158 | |
---|
159 | |
---|
160 | |
---|
161 | end |
---|
162 | |
---|