source: codes/icosagcm/devel/src/diagnostics/compute_temperature.F90

Last change on this file was 1027, checked in by dubos, 4 years ago

devel : towards conformity to F2008 standard

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