1 | MODULE compute_temperature_mod |
---|
2 | USE earth_const, ONLY : cpp, cppv, kappa, Rd, Rv, preff, Treff, nu, & |
---|
3 | caldyn_thermo, physics_thermo, thermo_fake_moist, & |
---|
4 | thermo_theta, thermo_entropy, thermo_variable_Cp, thermo_moist |
---|
5 | USE grid_param |
---|
6 | IMPLICIT NONE |
---|
7 | PRIVATE |
---|
8 | SAVE |
---|
9 | |
---|
10 | PUBLIC :: temperature, compute_temperature |
---|
11 | |
---|
12 | #include "../unstructured/unstructured.h90" |
---|
13 | |
---|
14 | CONTAINS |
---|
15 | |
---|
16 | SUBROUTINE temperature(f_pmid,f_q,f_temp) |
---|
17 | USE icosa |
---|
18 | TYPE(t_field), POINTER :: f_pmid(:) ! IN |
---|
19 | TYPE(t_field), POINTER :: f_q(:) ! IN |
---|
20 | TYPE(t_field), POINTER :: f_temp(:) ! INOUT |
---|
21 | |
---|
22 | REAL(rstd), POINTER :: pmid(:,:) |
---|
23 | REAL(rstd), POINTER :: q(:,:,:) |
---|
24 | REAL(rstd), POINTER :: temp(:,:) |
---|
25 | INTEGER :: ind |
---|
26 | |
---|
27 | DO ind=1,ndomain |
---|
28 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
29 | CALL swap_dimensions(ind) |
---|
30 | CALL swap_geometry(ind) |
---|
31 | pmid=f_pmid(ind) |
---|
32 | q=f_q(ind) |
---|
33 | temp=f_temp(ind) |
---|
34 | CALL compute_temperature(pmid,q,temp) |
---|
35 | END DO |
---|
36 | END SUBROUTINE temperature |
---|
37 | |
---|
38 | ! Macros loop_compute_temperature_fake_moist() and loop_compute_temperature() |
---|
39 | ! are here to inline the thermodynamical formula inside the innermost loop. |
---|
40 | ! Tests are made outside the outer loop. |
---|
41 | |
---|
42 | #ifdef BEGIN_DYSL |
---|
43 | |
---|
44 | {%- macro loop_compute_temperature() %} |
---|
45 | {%- set comp_temp=caller() %} |
---|
46 | FORALL_CELLS() |
---|
47 | ON_PRIMAL |
---|
48 | p_ik = pmid(CELL) |
---|
49 | theta_ik = temp(CELL) |
---|
50 | qv = q(CELL,1) ! water vapor mixing ratio = mv/md |
---|
51 | {{ comp_temp }} |
---|
52 | temp(CELL) = temp_ik |
---|
53 | END_BLOCK |
---|
54 | END_BLOCK |
---|
55 | {%- endmacro %} |
---|
56 | |
---|
57 | {%- macro loop_compute_temperature_fake_moist() %} |
---|
58 | {%- set comp_temp=caller() %} |
---|
59 | IF(physics_thermo==thermo_fake_moist) THEN |
---|
60 | {% call loop_compute_temperature() %} |
---|
61 | {{ comp_temp }} |
---|
62 | temp_ik = temp_ik/(1+0.608*qv) |
---|
63 | {% endcall %} |
---|
64 | ELSE |
---|
65 | {% call loop_compute_temperature() %} |
---|
66 | {{ comp_temp }} |
---|
67 | {% endcall %} |
---|
68 | END IF |
---|
69 | {%- endmacro %} |
---|
70 | |
---|
71 | KERNEL(compute_temperature) |
---|
72 | SELECT CASE(caldyn_thermo) |
---|
73 | CASE(thermo_theta) |
---|
74 | {% call loop_compute_temperature_fake_moist() %} |
---|
75 | temp_ik = theta_ik*((p_ik/preff)**kappa) |
---|
76 | {% endcall %} |
---|
77 | CASE(thermo_entropy) |
---|
78 | {% call loop_compute_temperature_fake_moist() %} |
---|
79 | temp_ik = Treff*exp((theta_ik + Rd*log(p_ik/preff))/cpp) |
---|
80 | {% endcall %} |
---|
81 | CASE(thermo_variable_Cp) |
---|
82 | {% call loop_compute_temperature() %} |
---|
83 | Cp_ik = nu*( theta_ik + Rd*log(p_ik/preff) ) |
---|
84 | temp_ik = Treff* (Cp_ik/cpp)**(1./nu) |
---|
85 | {% endcall %} |
---|
86 | CASE(thermo_moist) |
---|
87 | {% call loop_compute_temperature() %} |
---|
88 | Rmix = Rd+qv*Rv |
---|
89 | chi = ( theta_ik + Rmix*log(p_ik/preff) ) / (cpp + qv*cppv) ! log(T/Treff) |
---|
90 | temp_ik = Treff*exp(chi) |
---|
91 | {% endcall %} |
---|
92 | END SELECT |
---|
93 | END_BLOCK |
---|
94 | |
---|
95 | #endif END_DYSL |
---|
96 | |
---|
97 | SUBROUTINE compute_temperature_unst(pmid, q, temp) |
---|
98 | USE prec |
---|
99 | REAL(rstd),INTENT(IN) :: pmid(llm, primal_num) |
---|
100 | REAL(rstd),INTENT(IN) :: q(llm, primal_num, nqtot) |
---|
101 | REAL(rstd),INTENT(INOUT) :: temp(llm, primal_num) |
---|
102 | REAL(rstd) :: p_ik, theta_ik, temp_ik, qv, chi, Rmix, Cp_ik |
---|
103 | DECLARE_INDICES |
---|
104 | #include "../kernels_unst/compute_temperature.k90" |
---|
105 | END SUBROUTINE compute_temperature_unst |
---|
106 | |
---|
107 | SUBROUTINE compute_temperature_hex(pmid,q,temp) |
---|
108 | USE icosa |
---|
109 | USE omp_para |
---|
110 | REAL(rstd),INTENT(IN) :: pmid(iim*jjm,llm) |
---|
111 | REAL(rstd),INTENT(IN) :: q(iim*jjm,llm,nqtot) |
---|
112 | REAL(rstd),INTENT(INOUT) :: temp(iim*jjm,llm) |
---|
113 | |
---|
114 | REAL(rstd) :: p_ik, theta_ik, temp_ik, qv, chi, Rmix, Cp_ik |
---|
115 | INTEGER :: ij,l |
---|
116 | #include "../kernels_hex/compute_temperature.k90" |
---|
117 | END SUBROUTINE compute_temperature_hex |
---|
118 | |
---|
119 | SUBROUTINE compute_temperature(pmid,q,temp) |
---|
120 | USE icosa |
---|
121 | USE omp_para |
---|
122 | REAL(rstd),INTENT(IN) :: pmid(iim*jjm,llm) |
---|
123 | REAL(rstd),INTENT(IN) :: q(iim*jjm,llm,nqtot) |
---|
124 | REAL(rstd),INTENT(INOUT) :: temp(iim*jjm,llm) |
---|
125 | |
---|
126 | REAL(rstd) :: p_ik, theta_ik, temp_ik, qv, chi, Rmix |
---|
127 | INTEGER :: ij,l |
---|
128 | |
---|
129 | Rd = kappa*cpp |
---|
130 | DO l=ll_begin,ll_end |
---|
131 | DO ij=ij_begin,ij_end |
---|
132 | p_ik = pmid(ij,l) |
---|
133 | theta_ik = temp(ij,l) |
---|
134 | qv = q(ij,l,1) ! water vapor mixing ratio = mv/md |
---|
135 | SELECT CASE(caldyn_thermo) |
---|
136 | CASE(thermo_theta) |
---|
137 | temp_ik = theta_ik*((p_ik/preff)**kappa) |
---|
138 | CASE(thermo_entropy) |
---|
139 | temp_ik = Treff*exp((theta_ik + Rd*log(p_ik/preff))/cpp) |
---|
140 | CASE(thermo_moist) |
---|
141 | Rmix = Rd+qv*Rv |
---|
142 | chi = ( theta_ik + Rmix*log(p_ik/preff) ) / (cpp + qv*cppv) ! log(T/Treff) |
---|
143 | temp_ik = Treff*exp(chi) |
---|
144 | END SELECT |
---|
145 | IF(physics_thermo==thermo_fake_moist) temp_ik=temp_ik/(1+0.608*qv) |
---|
146 | temp(ij,l)=temp_ik |
---|
147 | END DO |
---|
148 | END DO |
---|
149 | END SUBROUTINE compute_temperature |
---|
150 | |
---|
151 | SUBROUTINE Tv2T(f_Tv, f_q, f_T) |
---|
152 | USE icosa |
---|
153 | TYPE(t_field), POINTER :: f_TV(:) |
---|
154 | TYPE(t_field), POINTER :: f_q(:) |
---|
155 | TYPE(t_field), POINTER :: f_T(:) |
---|
156 | |
---|
157 | REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:) |
---|
158 | INTEGER :: ind |
---|
159 | |
---|
160 | DO ind=1,ndomain |
---|
161 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
162 | CALL swap_dimensions(ind) |
---|
163 | CALL swap_geometry(ind) |
---|
164 | Tv=f_Tv(ind) |
---|
165 | T=f_T(ind) |
---|
166 | SELECT CASE(physics_thermo) |
---|
167 | CASE(thermo_dry) |
---|
168 | T=Tv |
---|
169 | CASE(thermo_fake_moist) |
---|
170 | q=f_q(ind) |
---|
171 | T=Tv/(1+0.608*q(:,:,1)) |
---|
172 | END SELECT |
---|
173 | END DO |
---|
174 | END SUBROUTINE Tv2T |
---|
175 | |
---|
176 | END MODULE compute_temperature_mod |
---|