source: codes/icosagcm/devel/src/diagnostics/compute_pression.F90 @ 1027

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

devel : towards conformity to F2008 standard

File size: 6.7 KB
Line 
1MODULE compute_pression_mod
2  USE compute_diagnostics_mod
3  USE icosa
4  USE omp_para
5  USE disvert_mod, ONLY : ap, bp, ap_bp_present, ptop, caldyn_eta, eta_lag
6  IMPLICIT NONE
7  PRIVATE
8
9#include "../unstructured/unstructured.h90"
10
11  PUBLIC :: pression, compute_pression_hex, compute_pression_unst, &
12       pression_mid, compute_pression_mid_hex, compute_pression_mid_unst, &
13       hydrostatic_pressure, compute_hydrostatic_pressure_unst, compute_hydrostatic_pressure_hex
14
15CONTAINS
16
17#ifdef BEGIN_DYSL
18
19{%- macro compute_pression(llmax)%} 
20{%- set inner_loop=caller() %}
21{%- set llmax="'%s'"%llmax %}
22IF(ap_bp_present) THEN
23   IF(offset>0) THEN
24     FORALL_CELLS_EXT('1',{{ llmax }})
25       ON_PRIMAL
26         {{ inner_loop }}
27       END_BLOCK
28     END_BLOCK
29  ELSE
30     FORALL_CELLS('1',{{ llmax }})
31       ON_PRIMAL
32         {{ inner_loop }}
33       END_BLOCK
34     END_BLOCK
35  END IF
36END IF
37{%- endmacro %}
38
39KERNEL(compute_pression)
40{% call compute_pression('llm+1') %}
41p(CELL) = AP(CELL) + BP(CELL) * ps(HIDX(CELL)) 
42{% endcall %}
43END_BLOCK
44
45KERNEL(compute_pmid)
46{% call compute_pression('llm') %}
47pmid(CELL) = .5*(AP(CELL)+AP(UP(CELL)) + (BP(CELL)+BP(UP(CELL))) * ps(HIDX(CELL)) ) 
48{% endcall %}
49END_BLOCK
50
51#endif END_DYSL
52
53  SUBROUTINE pression(f_ps,f_p)
54    TYPE(t_field), POINTER :: f_ps(:)
55    TYPE(t_field), POINTER :: f_p(:)
56 
57    REAL(rstd), POINTER :: ps(:)
58    REAL(rstd), POINTER :: p(:,:)
59    INTEGER :: ind
60
61!$OMP BARRIER
62    DO ind=1,ndomain
63      IF (.NOT. assigned_domain(ind)) CYCLE
64      CALL swap_dimensions(ind)
65      CALL swap_geometry(ind)
66      ps=f_ps(ind)
67      p=f_p(ind)
68      CALL compute_pression(ps, p,0)
69    ENDDO
70!$OMP BARRIER
71 
72  END SUBROUTINE pression
73
74  SUBROUTINE pression_mid(f_ps,f_pmid)
75    TYPE(t_field), POINTER :: f_ps(:)
76    TYPE(t_field), POINTER :: f_pmid(:)
77 
78    REAL(rstd), POINTER :: ps(:)
79    REAL(rstd), POINTER :: pmid(:,:)
80    INTEGER :: ind
81
82!$OMP BARRIER
83    DO ind=1,ndomain
84      IF (.NOT. assigned_domain(ind)) CYCLE
85      CALL swap_dimensions(ind)
86      CALL swap_geometry(ind)
87      ps=f_ps(ind)
88      pmid=f_pmid(ind)
89      CALL compute_pression_mid(ps, pmid,0)
90    ENDDO
91!$OMP BARRIER
92 
93  END SUBROUTINE pression_mid
94
95  SUBROUTINE hydrostatic_pressure(f_rhodz, f_theta_rhodz, f_ps, f_p)
96    TYPE(t_field), POINTER :: f_rhodz(:), f_theta_rhodz(:), f_ps(:), f_p(:)
97    REAL(rstd), POINTER :: ps(:), rhodz(:,:), p(:,:), theta_rhodz(:,:,:)
98    INTEGER :: ind
99    DO ind=1,ndomain
100      IF (.NOT. assigned_domain(ind)) CYCLE
101      CALL swap_dimensions(ind)
102      CALL swap_geometry(ind)
103      rhodz=f_rhodz(ind)
104      theta_rhodz=f_theta_rhodz(ind)
105      ps=f_ps(ind)
106      p=f_p(ind)
107      CALL compute_hydrostatic_pressure(rhodz, theta_rhodz, ps, p)
108    ENDDO
109  END SUBROUTINE hydrostatic_pressure
110
111!-------------- Wrappers for F2008 conformity -----------------
112
113  SUBROUTINE compute_pression_hex(ps,p,offset)
114    REAL(rstd),INTENT(IN) :: ps(:)
115    REAL(rstd),INTENT(OUT) :: p(:,:)
116    INTEGER,INTENT(IN) :: offset
117    CALL compute_pression_hex_(ps,p,offset)
118  END SUBROUTINE compute_pression_hex
119
120  SUBROUTINE compute_pression_unst(ps,p,offset)
121    REAL(rstd),INTENT(IN) :: ps(:)
122    REAL(rstd),INTENT(OUT) :: p(:,:)
123    INTEGER,INTENT(IN) :: offset
124    CALL compute_pression_unst_(ps,p,offset)
125  END SUBROUTINE compute_pression_unst
126
127  SUBROUTINE compute_pression_mid_hex(ps,p,offset)
128    REAL(rstd),INTENT(IN) :: ps(:)
129    REAL(rstd),INTENT(OUT) :: p(:,:)
130    INTEGER,INTENT(IN) :: offset
131    CALL compute_pression_mid_hex_(ps,p,offset)
132  END SUBROUTINE compute_pression_mid_hex
133
134  SUBROUTINE compute_pression_mid_unst(ps,p,offset)
135    REAL(rstd),INTENT(IN) :: ps(:)
136    REAL(rstd),INTENT(OUT) :: p(:,:)
137    INTEGER,INTENT(IN) :: offset
138    CALL compute_pression_mid_unst_(ps,p,offset)
139  END SUBROUTINE compute_pression_mid_unst
140
141  SUBROUTINE compute_hydrostatic_pressure_hex(rhodz, theta_rhodz, ps, p)
142    REAL(rstd),INTENT(IN)  :: rhodz(:,:), theta_rhodz(:,:,:)
143    REAL(rstd),INTENT(OUT) :: ps(:), p(:,:)
144    CALL compute_hydrostatic_pressure_hex_(rhodz, theta_rhodz, ps, p)
145  END SUBROUTINE compute_hydrostatic_pressure_hex
146
147  SUBROUTINE compute_hydrostatic_pressure_unst(rhodz, theta_rhodz, ps, p)
148    REAL(rstd),INTENT(IN)  :: rhodz(:,:), theta_rhodz(:,:,:)
149    REAL(rstd),INTENT(OUT) :: ps(:), p(:,:)
150    CALL compute_hydrostatic_pressure_unst_(rhodz, theta_rhodz, ps, p)
151  END SUBROUTINE compute_hydrostatic_pressure_unst
152 
153!------------- hexagonal-mesh compute kernels --------
154
155#define AP(ij,l) ap(l)
156#define BP(ij,l) bp(l)
157
158  SUBROUTINE compute_pression_hex_(ps,p,offset)
159    REAL(rstd),INTENT(IN) :: ps(iim*jjm)
160    REAL(rstd),INTENT(OUT) :: p(iim*jjm,llm+1)
161    INTEGER,INTENT(IN) :: offset
162    INTEGER :: ij,l
163#include "../kernels_hex/compute_pression.k90"
164  END SUBROUTINE compute_pression_hex_
165 
166  SUBROUTINE compute_pression_mid_hex_(ps,pmid,offset)
167    REAL(rstd),INTENT(IN) :: ps(iim*jjm)
168    REAL(rstd),INTENT(OUT) :: pmid(iim*jjm,llm)
169    INTEGER,INTENT(IN) :: offset
170    INTEGER :: ij,l
171#include "../kernels_hex/compute_pmid.k90"
172  END SUBROUTINE compute_pression_mid_hex_
173
174#undef AP
175#undef BP
176
177  SUBROUTINE compute_hydrostatic_pressure_hex_(rhodz, theta_rhodz, ps, pk)
178    REAL(rstd),INTENT(IN)  :: rhodz(iim*jjm,llm) ! mass per unit surface in each model level
179    REAL(rstd),INTENT(IN)  :: theta_rhodz(iim*jjm,llm, nqdyn) ! dynamical tracers (theta/entropy)
180    REAL(rstd),INTENT(OUT) :: ps(iim*jjm)        ! surface pressure, diagnosed if Lagrangian vertical coordinate
181    REAL(rstd),INTENT(OUT) :: pk(iim*jjm,llm)    ! pressure at full levels
182    INTEGER :: ij,l, ij_omp_begin_ext, ij_omp_end_ext
183    !$OMP BARRIER
184    CALL distrib_level(ij_begin_ext,ij_end_ext, ij_omp_begin_ext,ij_omp_end_ext)
185#include "../kernels_hex/compute_hydrostatic_pressure.k90"
186    !$OMP BARRIER
187  END SUBROUTINE compute_hydrostatic_pressure_hex_
188
189!----------- unstructured-mesh compute kernels --------
190
191#define AP(l,ij) ap(l)
192#define BP(l,ij) bp(l)
193 
194  SUBROUTINE compute_pression_unst_(ps, p, offset)
195    FIELD_PS,     INTENT(IN)  :: ps
196    FIELD_GEOPOT, INTENT(OUT) :: p
197    INTEGER,      INTENT(IN)  :: offset
198    DECLARE_INDICES
199#include "../kernels_unst/compute_pression.k90"
200  END SUBROUTINE compute_pression_unst_
201
202  SUBROUTINE compute_pression_mid_unst_(ps, pmid, offset)
203    FIELD_PS,   INTENT(IN)  :: ps
204    FIELD_MASS, INTENT(OUT) :: pmid
205    INTEGER,    INTENT(IN)  :: offset
206    DECLARE_INDICES
207#include "../kernels_unst/compute_pmid.k90"
208  END SUBROUTINE compute_pression_mid_unst_
209
210#undef AP
211#undef BP
212
213  SUBROUTINE compute_hydrostatic_pressure_unst_(rhodz, theta_rhodz, ps, pk)
214    FIELD_MASS,  INTENT(IN)  :: rhodz
215    FIELD_THETA, INTENT(IN)  :: theta_rhodz
216    FIELD_PS,    INTENT(OUT) :: ps
217    FIELD_MASS,  INTENT(OUT) :: pk
218    DECLARE_INDICES
219#include "../kernels_unst/compute_hydrostatic_pressure.k90"
220  END SUBROUTINE compute_hydrostatic_pressure_unst_
221 
222
223END MODULE compute_pression_mod
Note: See TracBrowser for help on using the repository browser.