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

Last change on this file since 977 was 952, checked in by dubos, 5 years ago

devel : temperature diagnostics for unstructured mesh

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