1 | SUBROUTINE computeice_table(timelen,ngrid,nslope,nsoil_GCM,nsoil_PEM,tsoil,tsurf,q_co2,q_h2o,ps,ice_table) |
---|
2 | USE comsoil_h, only: inertiedat, volcapa |
---|
3 | USE vertical_layers_mod, ONLY: ap,bp |
---|
4 | USE comsoil_h_PEM, only: layer_PEM |
---|
5 | implicit none |
---|
6 | |
---|
7 | integer,intent(in) :: timelen,ngrid,nslope,nsoil_PEM,nsoil_GCM |
---|
8 | real,intent(in) :: tsoil(ngrid,nsoil_PEM,nslope,timelen) ! soil temperature, timeseries [K] |
---|
9 | real,intent(in) :: tsurf(ngrid,nslope,timelen) ! surface temperature [K] |
---|
10 | real,intent(in) :: q_co2(ngrid,timelen) ! MMR tracer co2 [kg/kg] |
---|
11 | real,intent(in) :: q_h2o(ngrid,timelen) ! MMR tracer h2o [kg/kg] |
---|
12 | real,intent(in) :: ps(ngrid,timelen) ! surface pressure [Pa] |
---|
13 | real,intent(out) :: ice_table(ngrid,nslope) ! ice table [m] |
---|
14 | |
---|
15 | real :: m_h2o = 18.01528E-3 |
---|
16 | real :: m_co2 = 44.01E-3 |
---|
17 | real :: m_noco2 = 33.37E-3 |
---|
18 | real :: A,B,z1,z2 |
---|
19 | real :: alpha = -6143.7 |
---|
20 | real :: beta = 29.9074 |
---|
21 | |
---|
22 | integer ig, islope,isoil,it |
---|
23 | real,allocatable :: mass_mean(:,:) ! mean mass above the surface |
---|
24 | real,allocatable :: zplev_mean(:,:) ! pressure above the surface |
---|
25 | real,allocatable :: pvapor(:,:) ! partial pressure above the surface |
---|
26 | |
---|
27 | real,allocatable :: rhovapor(:,:,:) |
---|
28 | real,allocatable :: rhovapor_avg(:,:) ! mean vapor_density at the surface yearly averaged |
---|
29 | |
---|
30 | real :: psv_surf |
---|
31 | real,allocatable :: rho_soil(:,:,:,:) ! water vapor in the soil |
---|
32 | real,allocatable :: rho_soil_avg(:,:,:) ! water vapor in the soil yearly averaged |
---|
33 | |
---|
34 | real,allocatable :: diff_rho(:,:,:) ! difference of vapor content |
---|
35 | |
---|
36 | allocate(rhovapor(ngrid,nslope,timelen)) |
---|
37 | allocate(rhovapor_avg(ngrid,nslope)) |
---|
38 | allocate(pvapor(ngrid,timelen)) |
---|
39 | |
---|
40 | allocate(mass_mean(ngrid,timelen)) |
---|
41 | allocate(zplev_mean(ngrid,timelen)) |
---|
42 | |
---|
43 | ! 0. Some initializations |
---|
44 | |
---|
45 | A =(1/m_co2 - 1/m_noco2) |
---|
46 | B=1/m_noco2 |
---|
47 | ! 1. Compute rho surface yearly averaged |
---|
48 | |
---|
49 | ! 1.1 Compute the partial pressure of vapor |
---|
50 | !a. the molecular mass into the column |
---|
51 | do ig = 1,ngrid |
---|
52 | mass_mean(ig,:) = 1/(A*q_co2(ig,:) +B) |
---|
53 | enddo |
---|
54 | |
---|
55 | ! b. pressure level |
---|
56 | do it = 1,timelen |
---|
57 | do ig = 1,ngrid |
---|
58 | zplev_mean(ig,it) = ap(1) + bp(1)*ps(ig,it) |
---|
59 | enddo |
---|
60 | enddo |
---|
61 | |
---|
62 | ! c. Vapor pressure |
---|
63 | pvapor(:,:) = mass_mean(:,:)/m_h2o*q_h2o(:,:)*zplev_mean(:,:) |
---|
64 | |
---|
65 | deallocate(mass_mean) |
---|
66 | deallocate(zplev_mean) |
---|
67 | |
---|
68 | ! 1.2 Check if there is frost at the surface and then compute the density |
---|
69 | ! at the surface |
---|
70 | do ig = 1,ngrid |
---|
71 | do islope = 1,nslope |
---|
72 | do it = 1,timelen |
---|
73 | psv_surf = exp(alpha/tsurf(ig,islope,it) +beta) |
---|
74 | if ((isnan(psv_surf)).or.(isnan(pvapor(ig,it)))) then |
---|
75 | write(*,*) 'pb vapor',ig,islope,it |
---|
76 | stop |
---|
77 | endif |
---|
78 | rhovapor(ig,islope,it) = min(psv_surf,pvapor(ig,it))/tsurf(ig,islope,it) |
---|
79 | enddo |
---|
80 | enddo |
---|
81 | enddo |
---|
82 | deallocate(pvapor) |
---|
83 | |
---|
84 | ! 1.3 Density at the surface yearly averaged |
---|
85 | rhovapor_avg(:,:) = SUM(rhovapor(:,:,:),3)/timelen |
---|
86 | |
---|
87 | deallocate(rhovapor) |
---|
88 | |
---|
89 | ! 2. Compute rho soil vapor |
---|
90 | |
---|
91 | allocate(rho_soil_avg(ngrid,nslope,nsoil_PEM)) |
---|
92 | allocate(rho_soil(ngrid,nslope,nsoil_PEM,timelen)) |
---|
93 | |
---|
94 | do ig = 1,ngrid |
---|
95 | do islope = 1,nslope |
---|
96 | do isoil = 1,nsoil_PEM |
---|
97 | do it = 1,timelen |
---|
98 | rho_soil(ig,islope,isoil,it) = exp(alpha/tsoil(ig,isoil,islope,it) +beta)/tsoil(ig,isoil,islope,it) |
---|
99 | enddo |
---|
100 | enddo |
---|
101 | enddo |
---|
102 | enddo |
---|
103 | |
---|
104 | rho_soil_avg(:,:,:) = SUM( rho_soil(:,:,:,:),4)/timelen |
---|
105 | deallocate(rho_soil) |
---|
106 | |
---|
107 | ! 3. Computing ice table |
---|
108 | |
---|
109 | ice_table (:,:) = -1e4 |
---|
110 | |
---|
111 | allocate(diff_rho(ngrid,nslope,nsoil_PEM)) |
---|
112 | |
---|
113 | do isoil = 1,nsoil_PEM |
---|
114 | diff_rho(:,:,isoil) = rhovapor_avg(:,:) - rho_soil_avg(:,:,isoil) |
---|
115 | ! write(*,*) 'diff =',ig,islope,isoil,diff_rho(ig,islope,isoil),rhovapor_avg(ig,islope) ,rho_soil_avg(ig,islope,isoil) |
---|
116 | enddo |
---|
117 | |
---|
118 | deallocate(rhovapor_avg) |
---|
119 | deallocate(rho_soil_avg) |
---|
120 | |
---|
121 | do ig = 1,ngrid |
---|
122 | do islope = 1,nslope |
---|
123 | if(diff_rho(ig,islope,1) > 0) then |
---|
124 | ice_table(ig,islope) = 0. |
---|
125 | else |
---|
126 | do isoil = 1,nsoil_PEM -1 |
---|
127 | if((diff_rho(ig,islope,isoil).lt.0).and.(diff_rho(ig,islope,isoil+1).gt.0.)) then |
---|
128 | z1 = (diff_rho(ig,islope,isoil) - diff_rho(ig,islope,isoil+1))/(layer_PEM(isoil) - layer_PEM(isoil+1)) |
---|
129 | z2 = -layer_PEM(isoil+1)*z1 + diff_rho(ig,islope,isoil+1) |
---|
130 | ice_table(ig,islope) = -z2/z1 |
---|
131 | exit |
---|
132 | endif |
---|
133 | enddo |
---|
134 | endif |
---|
135 | enddo |
---|
136 | enddo |
---|
137 | |
---|
138 | deallocate(diff_rho) |
---|
139 | |
---|
140 | !======================================================================= |
---|
141 | RETURN |
---|
142 | |
---|
143 | |
---|
144 | END |
---|