1 | MODULE diagflux_mod |
---|
2 | USE icosa |
---|
3 | IMPLICIT NONE |
---|
4 | SAVE |
---|
5 | |
---|
6 | TYPE(t_field),POINTER :: & |
---|
7 | f_masst(:), f_qmasst(:), & ! time-integrated mass, tracer mass, |
---|
8 | f_massfluxt(:), f_qfluxt(:), & ! mass flux and tracer flux |
---|
9 | f_qfluxt_lon(:), f_qfluxt_lat(:) ! scalar flux reconstructed cell centers |
---|
10 | LOGICAL :: diagflux_on |
---|
11 | !$OMP THREADPRIVATE(diagflux_on) |
---|
12 | |
---|
13 | CONTAINS |
---|
14 | |
---|
15 | SUBROUTINE init_diagflux |
---|
16 | USE getin_mod |
---|
17 | diagflux_on = .FALSE. |
---|
18 | CALL getin("diagflux", diagflux_on) |
---|
19 | IF(diagflux_on) THEN |
---|
20 | CALL allocate_field(f_masst, field_t,type_real,llm, name="masst") |
---|
21 | CALL allocate_field(f_qmasst, field_t,type_real,llm,nqtot, name="qmasst") |
---|
22 | CALL allocate_field(f_massfluxt, field_u,type_real,llm, name="massfluxt") |
---|
23 | CALL allocate_field(f_qfluxt, field_u,type_real,llm,nqtot, name="qfluxt") |
---|
24 | CALL allocate_field(f_qfluxt_lon, field_t,type_real,llm,nqtot, name="qfluxt_lon") |
---|
25 | CALL allocate_field(f_qfluxt_lat, field_t,type_real,llm,nqtot, name="qfluxt_lat") |
---|
26 | CALL zero_qfluxt |
---|
27 | ELSE |
---|
28 | CALL allocate_field(f_masst, field_t,type_real,0, name="masst") |
---|
29 | CALL allocate_field(f_qmasst, field_t,type_real,llm,0, name="qmasst") |
---|
30 | CALL allocate_field(f_massfluxt, field_u,type_real,0, name="massfluxt") |
---|
31 | CALL allocate_field(f_qfluxt, field_u,type_real,llm,0, name="qfluxt") |
---|
32 | CALL allocate_field(f_qfluxt_lon, field_t,type_real,llm,0, name="qfluxt_lon") |
---|
33 | CALL allocate_field(f_qfluxt_lat, field_t,type_real,llm,0, name="qfluxt_lat") |
---|
34 | END IF |
---|
35 | END SUBROUTINE init_diagflux |
---|
36 | |
---|
37 | SUBROUTINE zero_qfluxt |
---|
38 | USE mpipara |
---|
39 | USE omp_para |
---|
40 | INTEGER :: ind |
---|
41 | REAL(rstd), POINTER :: buf2(:,:),buf3(:,:,:) |
---|
42 | DO ind=1,ndomain |
---|
43 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
44 | CALL swap_dimensions(ind) |
---|
45 | buf2=f_masst(ind) |
---|
46 | buf2(:,ll_begin:ll_end)=0. |
---|
47 | buf2=f_massfluxt(ind) |
---|
48 | buf2(:,ll_begin:ll_end)=0. |
---|
49 | buf3=f_qmasst(ind) |
---|
50 | buf3(:,ll_begin:ll_end,:)=0. |
---|
51 | buf3=f_qfluxt(ind) |
---|
52 | buf3(:,ll_begin:ll_end,:)=0. |
---|
53 | END DO |
---|
54 | END SUBROUTINE zero_qfluxt |
---|
55 | |
---|
56 | SUBROUTINE flux_centered_lonlat(scale, f_massflux, f_flux, f_massflux_lon, f_massflux_lat, f_flux_lon, f_flux_lat) |
---|
57 | REAL(rstd), INTENT(IN) :: scale |
---|
58 | TYPE(t_field),POINTER :: f_flux(:), f_flux_lon(:), f_flux_lat(:), & |
---|
59 | f_massflux(:), f_massflux_lon(:), f_massflux_lat(:) |
---|
60 | REAL(rstd), POINTER :: flux(:,:,:), flux_lon(:,:,:), flux_lat(:,:,:), & |
---|
61 | massflux(:,:), massflux_lon(:,:), massflux_lat(:,:) |
---|
62 | INTEGER :: ind, itrac |
---|
63 | DO ind=1,ndomain |
---|
64 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
65 | CALL swap_dimensions(ind) |
---|
66 | CALL swap_geometry(ind) |
---|
67 | flux=f_flux(ind) |
---|
68 | flux_lon=f_flux_lon(ind) |
---|
69 | flux_lat=f_flux_lat(ind) |
---|
70 | DO itrac=1,nqtot |
---|
71 | CALL compute_flux_centered_lonlat(scale, flux(:,:,itrac), flux_lon(:,:,itrac), flux_lat(:,:,itrac)) |
---|
72 | END DO |
---|
73 | massflux=f_massflux(ind) |
---|
74 | massflux_lon=f_massflux_lon(ind) |
---|
75 | massflux_lat=f_massflux_lat(ind) |
---|
76 | CALL compute_flux_centered_lonlat(scale, massflux, massflux_lon, massflux_lat) |
---|
77 | END DO |
---|
78 | END SUBROUTINE flux_centered_lonlat |
---|
79 | |
---|
80 | SUBROUTINE compute_flux_centered_lonlat(scale, flux, flux_lon, flux_lat) |
---|
81 | USE wind_mod |
---|
82 | REAL(rstd), INTENT(IN) :: scale |
---|
83 | REAL(rstd), INTENT(IN) :: flux(3*iim*jjm,llm) |
---|
84 | REAL(rstd), INTENT(OUT) :: flux_lon(iim*jjm,llm), flux_lat(iim*jjm,llm) |
---|
85 | REAL(rstd) :: flux_3d(iim*jjm,llm,3) |
---|
86 | CALL compute_flux_centered(scale, flux, flux_3d) |
---|
87 | CALL compute_wind_centered_lonlat_compound(flux_3d, flux_lon, flux_lat) |
---|
88 | END SUBROUTINE compute_flux_centered_lonlat |
---|
89 | |
---|
90 | END MODULE diagflux_mod |
---|