[3727] | 1 | module soil_mod |
---|
| 2 | |
---|
| 3 | implicit none |
---|
| 4 | |
---|
| 5 | contains |
---|
| 6 | |
---|
[38] | 7 | subroutine soil(ngrid,nsoil,firstcall, |
---|
| 8 | & therm_i, |
---|
| 9 | & timestep,tsurf,tsoil, |
---|
| 10 | & capcal,fluxgrd) |
---|
[1224] | 11 | |
---|
| 12 | use comsoil_h, only: layer, mlayer, volcapa, |
---|
| 13 | & mthermdiff, thermdiff, coefq, |
---|
[2887] | 14 | & coefd, alph, beta, mu,flux_geo |
---|
[1047] | 15 | use surfdat_h, only: watercaptag, inert_h2o_ice |
---|
[2900] | 16 | use comslope_mod, ONLY: nslope |
---|
[3726] | 17 | use callkeys_mod, only: surfaceice_tifeedback, poreice_tifeedback |
---|
[38] | 18 | implicit none |
---|
| 19 | |
---|
| 20 | !----------------------------------------------------------------------- |
---|
| 21 | ! Author: Ehouarn Millour |
---|
| 22 | ! |
---|
| 23 | ! Purpose: Compute soil temperature using an implict 1st order scheme |
---|
| 24 | ! |
---|
| 25 | ! Note: depths of layers and mid-layers, soil thermal inertia and |
---|
[1047] | 26 | ! heat capacity are commons in comsoil_h |
---|
[38] | 27 | !----------------------------------------------------------------------- |
---|
| 28 | |
---|
| 29 | c----------------------------------------------------------------------- |
---|
| 30 | ! arguments |
---|
| 31 | ! --------- |
---|
| 32 | ! inputs: |
---|
[3727] | 33 | integer,intent(in) :: ngrid ! number of (horizontal) grid-points |
---|
| 34 | integer,intent(in) :: nsoil ! number of soil layers |
---|
| 35 | logical,intent(in) :: firstcall ! identifier for initialization call |
---|
| 36 | real,intent(in) :: therm_i(ngrid,nsoil,nslope) ! thermal inertia |
---|
| 37 | real,intent(in) :: timestep ! time step |
---|
| 38 | real,intent(in) :: tsurf(ngrid,nslope) ! surface temperature |
---|
[38] | 39 | ! outputs: |
---|
[3727] | 40 | real,intent(out) :: tsoil(ngrid,nsoil,nslope) ! soil (mid-layer) temperature |
---|
| 41 | real,intent(out) :: capcal(ngrid,nslope) ! surface specific heat |
---|
| 42 | real,intent(out) :: fluxgrd(ngrid,nslope) ! surface diffusive heat flux |
---|
[38] | 43 | |
---|
| 44 | ! local variables: |
---|
[2900] | 45 | integer ig,ik,islope |
---|
[38] | 46 | |
---|
| 47 | ! 0. Initialisations and preprocessing step |
---|
[3230] | 48 | if(firstcall.or.surfaceice_tifeedback.or.poreice_tifeedback) then |
---|
[38] | 49 | ! note: firstcall is set to .true. or .false. by the caller |
---|
| 50 | ! and not changed by soil.F |
---|
| 51 | ! 0.1 Build mthermdiff(:), the mid-layer thermal diffusivities |
---|
| 52 | do ig=1,ngrid |
---|
[2900] | 53 | do islope = 1,nslope |
---|
[283] | 54 | if (watercaptag(ig)) then |
---|
| 55 | do ik=0,nsoil-1 |
---|
| 56 | ! If we have permanent ice, we use the water ice thermal inertia from ground to last layer. |
---|
[2900] | 57 | mthermdiff(ig,ik,islope)=inert_h2o_ice* |
---|
| 58 | & inert_h2o_ice/volcapa |
---|
[283] | 59 | enddo |
---|
| 60 | else |
---|
| 61 | do ik=0,nsoil-1 |
---|
[2900] | 62 | mthermdiff(ig,ik,islope)=therm_i(ig,ik+1,islope)* |
---|
| 63 | & therm_i(ig,ik+1,islope)/volcapa |
---|
[283] | 64 | enddo |
---|
| 65 | endif |
---|
[38] | 66 | enddo |
---|
[2900] | 67 | enddo |
---|
[38] | 68 | |
---|
[285] | 69 | #ifdef MESOSCALE |
---|
| 70 | do ig=1,ngrid |
---|
[2900] | 71 | do islope = 1,nslope |
---|
| 72 | if ( therm_i(ig,1,islope) .ge. inert_h2o_ice ) then |
---|
| 73 | print *, "limit max TI ", therm_i(ig,1,islope), inert_h2o_ice |
---|
[285] | 74 | do ik=0,nsoil-1 |
---|
[2900] | 75 | mthermdiff(ig,ik,islope)=inert_h2o_ice* |
---|
| 76 | & inert_h2o_ice/volcapa |
---|
[285] | 77 | enddo |
---|
| 78 | endif |
---|
[2900] | 79 | enddo |
---|
[285] | 80 | enddo |
---|
| 81 | #endif |
---|
| 82 | |
---|
[38] | 83 | ! 0.2 Build thermdiff(:), the "interlayer" thermal diffusivities |
---|
| 84 | do ig=1,ngrid |
---|
[2900] | 85 | do islope = 1,nslope |
---|
[38] | 86 | do ik=1,nsoil-1 |
---|
[2900] | 87 | thermdiff(ig,ik,islope)=((layer(ik)-mlayer(ik-1)) |
---|
| 88 | & *mthermdiff(ig,ik,islope) |
---|
| 89 | & +(mlayer(ik)-layer(ik)) |
---|
| 90 | & *mthermdiff(ig,ik-1,islope)) |
---|
[38] | 91 | & /(mlayer(ik)-mlayer(ik-1)) |
---|
| 92 | ! write(*,*),'soil: ik: ',ik,' thermdiff:',thermdiff(ig,ik) |
---|
| 93 | enddo |
---|
[2900] | 94 | enddo |
---|
[38] | 95 | enddo |
---|
| 96 | |
---|
| 97 | ! 0.3 Build coefficients mu, q_{k+1/2}, d_k, alpha_k and capcal |
---|
| 98 | ! mu |
---|
| 99 | mu=mlayer(0)/(mlayer(1)-mlayer(0)) |
---|
| 100 | |
---|
| 101 | ! q_{1/2} |
---|
| 102 | coefq(0)=volcapa*layer(1)/timestep |
---|
| 103 | ! q_{k+1/2} |
---|
| 104 | do ik=1,nsoil-1 |
---|
| 105 | coefq(ik)=volcapa*(layer(ik+1)-layer(ik)) |
---|
| 106 | & /timestep |
---|
| 107 | enddo |
---|
| 108 | |
---|
| 109 | do ig=1,ngrid |
---|
[2900] | 110 | do islope = 1,nslope |
---|
[38] | 111 | ! d_k |
---|
| 112 | do ik=1,nsoil-1 |
---|
[2900] | 113 | coefd(ig,ik,islope)=thermdiff(ig,ik,islope) |
---|
| 114 | & /(mlayer(ik)-mlayer(ik-1)) |
---|
[38] | 115 | enddo |
---|
| 116 | |
---|
| 117 | ! alph_{N-1} |
---|
[2900] | 118 | alph(ig,nsoil-1,islope)=coefd(ig,nsoil-1,islope)/ |
---|
| 119 | & (coefq(nsoil-1)+coefd(ig,nsoil-1,islope)) |
---|
[38] | 120 | ! alph_k |
---|
| 121 | do ik=nsoil-2,1,-1 |
---|
[2900] | 122 | alph(ig,ik,islope)=coefd(ig,ik,islope)/ |
---|
| 123 | & (coefq(ik)+coefd(ig,ik+1,islope)* |
---|
| 124 | & (1.-alph(ig,ik+1,islope))+coefd(ig,ik,islope)) |
---|
[38] | 125 | enddo |
---|
| 126 | |
---|
| 127 | ! capcal |
---|
| 128 | ! Cstar |
---|
[2900] | 129 | capcal(ig,islope)=volcapa*layer(1)+ |
---|
| 130 | & (thermdiff(ig,1,islope)/(mlayer(1)-mlayer(0)))* |
---|
| 131 | & (timestep*(1.-alph(ig,1,islope))) |
---|
[38] | 132 | ! Cs |
---|
[2900] | 133 | capcal(ig,islope)=capcal(ig,islope)/ |
---|
| 134 | & (1.+mu*(1.0-alph(ig,1,islope))* |
---|
| 135 | & thermdiff(ig,1,islope)/mthermdiff(ig,0,islope)) |
---|
[38] | 136 | ! write(*,*)'soil: ig=',ig,' capcal(ig)=',capcal(ig) |
---|
[2900] | 137 | enddo ! islope |
---|
[38] | 138 | enddo ! of do ig=1,ngrid |
---|
| 139 | |
---|
[833] | 140 | endif ! of if (firstcall.or.tifeedback) |
---|
[38] | 141 | |
---|
| 142 | ! 1. Compute soil temperatures |
---|
[833] | 143 | IF (.not.firstcall) THEN |
---|
[38] | 144 | ! First layer: |
---|
[2900] | 145 | do islope = 1,nslope |
---|
[38] | 146 | do ig=1,ngrid |
---|
[2900] | 147 | tsoil(ig,1,islope)=(tsurf(ig,islope)+mu*beta(ig,1,islope)* |
---|
| 148 | & thermdiff(ig,1,islope)/mthermdiff(ig,0,islope))/ |
---|
| 149 | & (1.+mu*(1.0-alph(ig,1,islope))* |
---|
| 150 | & thermdiff(ig,1,islope)/mthermdiff(ig,0,islope)) |
---|
[38] | 151 | enddo |
---|
| 152 | ! Other layers: |
---|
| 153 | do ik=1,nsoil-1 |
---|
| 154 | do ig=1,ngrid |
---|
[2900] | 155 | tsoil(ig,ik+1,islope)=alph(ig,ik,islope)* |
---|
| 156 | & tsoil(ig,ik,islope)+beta(ig,ik,islope) |
---|
[38] | 157 | enddo |
---|
| 158 | enddo |
---|
[2900] | 159 | enddo ! islope |
---|
[833] | 160 | ENDIF! of if (.not.firstcall) |
---|
[38] | 161 | |
---|
| 162 | ! 2. Compute beta coefficients (preprocessing for next time step) |
---|
| 163 | ! Bottom layer, beta_{N-1} |
---|
[2900] | 164 | do islope = 1,nslope |
---|
[38] | 165 | do ig=1,ngrid |
---|
[2900] | 166 | beta(ig,nsoil-1,islope)=coefq(nsoil-1)*tsoil(ig,nsoil,islope) |
---|
| 167 | & /(coefq(nsoil-1)+coefd(ig,nsoil-1,islope)) |
---|
[2919] | 168 | & +flux_geo(ig,islope)/ |
---|
| 169 | & (coefq(nsoil-1) + coefd(ig,nsoil-1,islope)) |
---|
[38] | 170 | enddo |
---|
[2887] | 171 | |
---|
[38] | 172 | ! Other layers |
---|
| 173 | do ik=nsoil-2,1,-1 |
---|
| 174 | do ig=1,ngrid |
---|
[2900] | 175 | beta(ig,ik,islope)=(coefq(ik)*tsoil(ig,ik+1,islope)+ |
---|
| 176 | & coefd(ig,ik+1,islope)*beta(ig,ik+1,islope))/ |
---|
| 177 | & (coefq(ik)+coefd(ig,ik+1,islope)* |
---|
| 178 | & (1.0-alph(ig,ik+1,islope)) |
---|
| 179 | & +coefd(ig,ik,islope)) |
---|
[38] | 180 | enddo |
---|
| 181 | enddo |
---|
| 182 | |
---|
| 183 | ! 3. Compute surface diffusive flux & calorific capacity |
---|
| 184 | do ig=1,ngrid |
---|
| 185 | ! Cstar |
---|
| 186 | ! capcal(ig)=volcapa(ig,1)*layer(ig,1)+ |
---|
| 187 | ! & (thermdiff(ig,1)/(mlayer(ig,1)-mlayer(ig,0)))* |
---|
| 188 | ! & (timestep*(1.-alph(ig,1))) |
---|
| 189 | ! Fstar |
---|
[2900] | 190 | fluxgrd(ig,islope)=(thermdiff(ig,1,islope)/ |
---|
| 191 | & (mlayer(1)-mlayer(0)))* |
---|
| 192 | & (beta(ig,1,islope)+(alph(ig,1,islope)-1.0) |
---|
| 193 | & *tsoil(ig,1,islope)) |
---|
[38] | 194 | |
---|
| 195 | ! mu=mlayer(ig,0)/(mlayer(ig,1)-mlayer(ig,0)) |
---|
| 196 | ! capcal(ig)=capcal(ig)/(1.+mu*(1.0-alph(ig,1))* |
---|
| 197 | ! & thermdiff(ig,1)/mthermdiff(ig,0)) |
---|
| 198 | ! Fs |
---|
[2900] | 199 | fluxgrd(ig,islope)=fluxgrd(ig,islope)+(capcal(ig,islope) |
---|
| 200 | & /timestep)* |
---|
| 201 | & (tsoil(ig,1,islope)* |
---|
| 202 | & (1.+mu*(1.0-alph(ig,1,islope))* |
---|
| 203 | & thermdiff(ig,1,islope)/mthermdiff(ig,0,islope)) |
---|
| 204 | & -tsurf(ig,islope)-mu*beta(ig,1,islope)* |
---|
| 205 | & thermdiff(ig,1,islope)/mthermdiff(ig,0,islope)) |
---|
[38] | 206 | enddo |
---|
[2900] | 207 | enddo ! islope |
---|
[38] | 208 | |
---|
[3727] | 209 | end subroutine soil |
---|
| 210 | |
---|
| 211 | end module soil_mod |
---|
| 212 | |
---|