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

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

devel : DYSL for compute_temperature

File size: 4.9 KB
Line 
1MODULE compute_temperature_mod
2  USE earth_const, ONLY : cpp, cppv, kappa, Rd, Rv, preff, Treff, &
3       caldyn_thermo, physics_thermo, &
4       thermo_theta, thermo_entropy, thermo_moist, thermo_fake_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
14CONTAINS
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() %}
59IF(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 %}
64ELSE
65   {% call loop_compute_temperature() %}
66     {{ comp_temp }}
67   {% endcall %}
68END IF
69{%- endmacro %}
70
71KERNEL(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_moist)
82     {% call loop_compute_temperature() %}
83       Rmix = Rd+qv*Rv
84       chi = ( theta_ik + Rmix*log(p_ik/preff) ) / (cpp + qv*cppv) ! log(T/Treff)
85       temp_ik = Treff*exp(chi)
86     {% endcall %}
87  END SELECT
88END_BLOCK
89
90#endif END_DYSL
91
92  SUBROUTINE compute_temperature_unst(pmid, q, temp)
93    USE prec
94    REAL(rstd),INTENT(IN)    :: pmid(llm, primal_num)
95    REAL(rstd),INTENT(IN)    :: q(llm, primal_num, nqtot)
96    REAL(rstd),INTENT(INOUT) :: temp(llm, primal_num)
97    REAL(rstd) :: p_ik, theta_ik, temp_ik, qv, chi, Rmix 
98    DECLARE_INDICES
99#include "../kernels_unst/compute_temperature.k90"
100  END SUBROUTINE compute_temperature_unst
101
102  SUBROUTINE compute_temperature_hex(pmid,q,temp)
103    USE icosa
104    USE omp_para
105    REAL(rstd),INTENT(IN)    :: pmid(iim*jjm,llm)
106    REAL(rstd),INTENT(IN)    :: q(iim*jjm,llm,nqtot)
107    REAL(rstd),INTENT(INOUT) :: temp(iim*jjm,llm)
108
109    REAL(rstd) :: p_ik, theta_ik, temp_ik, qv, chi, Rmix
110    INTEGER :: ij,l
111#include "../kernels_hex/compute_temperature.k90"
112  END SUBROUTINE compute_temperature_hex
113
114  SUBROUTINE compute_temperature(pmid,q,temp)
115    USE icosa
116    USE omp_para
117    REAL(rstd),INTENT(IN)    :: pmid(iim*jjm,llm)
118    REAL(rstd),INTENT(IN)    :: q(iim*jjm,llm,nqtot)
119    REAL(rstd),INTENT(INOUT) :: temp(iim*jjm,llm)
120
121    REAL(rstd) :: p_ik, theta_ik, temp_ik, qv, chi, Rmix
122    INTEGER :: ij,l
123
124    Rd = kappa*cpp
125    DO l=ll_begin,ll_end
126       DO ij=ij_begin,ij_end
127          p_ik = pmid(ij,l)
128          theta_ik = temp(ij,l)
129          qv = q(ij,l,1) ! water vapor mixing ratio = mv/md
130          SELECT CASE(caldyn_thermo)
131          CASE(thermo_theta)
132             temp_ik = theta_ik*((p_ik/preff)**kappa)
133          CASE(thermo_entropy)
134             temp_ik = Treff*exp((theta_ik + Rd*log(p_ik/preff))/cpp)
135          CASE(thermo_moist)
136             Rmix = Rd+qv*Rv
137             chi = ( theta_ik + Rmix*log(p_ik/preff) ) / (cpp + qv*cppv) ! log(T/Treff)
138             temp_ik = Treff*exp(chi)
139          END SELECT
140          IF(physics_thermo==thermo_fake_moist) temp_ik=temp_ik/(1+0.608*qv)
141          temp(ij,l)=temp_ik
142       END DO
143    END DO
144  END SUBROUTINE compute_temperature
145
146  SUBROUTINE Tv2T(f_Tv, f_q, f_T)
147    USE icosa
148    TYPE(t_field), POINTER :: f_TV(:)
149    TYPE(t_field), POINTER :: f_q(:)
150    TYPE(t_field), POINTER :: f_T(:)
151   
152    REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:)
153    INTEGER :: ind
154   
155    DO ind=1,ndomain
156       IF (.NOT. assigned_domain(ind)) CYCLE
157       CALL swap_dimensions(ind)
158       CALL swap_geometry(ind)
159       Tv=f_Tv(ind)
160       T=f_T(ind)
161       SELECT CASE(physics_thermo)
162       CASE(thermo_dry)
163          T=Tv
164       CASE(thermo_fake_moist)
165          q=f_q(ind)
166          T=Tv/(1+0.608*q(:,:,1))
167       END SELECT
168    END DO
169  END SUBROUTINE Tv2T
170
171END MODULE compute_temperature_mod
Note: See TracBrowser for help on using the repository browser.