1 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
2 | ! Interface pour ecrire en netcdf avec les routines d'enseignement |
---|
3 | ! iotd de Frederic Hourdin |
---|
4 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
5 | |
---|
6 | subroutine iophys_ecrit(nom,lllm,titre,unite,px) |
---|
7 | |
---|
8 | USE mod_phys_lmdz_para, ONLY: klon_omp, is_mpi_root |
---|
9 | USE mod_phys_lmdz_transfert_para, ONLY: gather |
---|
10 | USE mod_grid_phy_lmdz, ONLY: klon_glo, nbp_lon, nbp_lat, grid1dto2d_glo |
---|
11 | IMPLICIT NONE |
---|
12 | |
---|
13 | |
---|
14 | |
---|
15 | ! Ecriture de variables diagnostiques au choix dans la physique |
---|
16 | ! dans un fichier NetCDF nomme 'diagfi'. Ces variables peuvent etre |
---|
17 | ! 3d (ex : temperature), 2d (ex : temperature de surface), ou |
---|
18 | ! 0d (pour un scalaire qui ne depend que du temps : ex : la longitude |
---|
19 | ! solaire) |
---|
20 | ! (ou encore 1d, dans le cas de testphys1d, pour sortir une colonne) |
---|
21 | ! La periode d'ecriture est donnee par |
---|
22 | ! "ecritphy " regle dans le fichier de controle de run : run.def |
---|
23 | |
---|
24 | ! writediagfi peut etre appele de n'importe quelle subroutine |
---|
25 | ! de la physique, plusieurs fois. L'initialisation et la creation du |
---|
26 | ! fichier se fait au tout premier appel. |
---|
27 | |
---|
28 | ! WARNING : les variables dynamique (u,v,t,q,ps) |
---|
29 | ! sauvees par writediagfi avec une |
---|
30 | ! date donnee sont legerement differentes que dans le fichier histoire car |
---|
31 | ! on ne leur a pas encore ajoute de la dissipation et de la physique !!! |
---|
32 | ! IL est RECOMMANDE d'ajouter les tendance physique a ces variables |
---|
33 | ! avant l'ecriture dans diagfi (cf. physiq.F) |
---|
34 | |
---|
35 | ! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4 |
---|
36 | |
---|
37 | ! parametres (input) : |
---|
38 | ! ---------- |
---|
39 | ! unit : unite logique du fichier de sortie (toujours la meme) |
---|
40 | ! nom : nom de la variable a sortir (chaine de caracteres) |
---|
41 | ! titre: titre de la variable (chaine de caracteres) |
---|
42 | ! unite : unite de la variable (chaine de caracteres) |
---|
43 | ! px : variable a sortir (real 0, 1, 2, ou 3d) |
---|
44 | |
---|
45 | !================================================================= |
---|
46 | |
---|
47 | |
---|
48 | ! Arguments on input: |
---|
49 | integer lllm |
---|
50 | character (len=*) :: nom,titre,unite |
---|
51 | integer imjmax |
---|
52 | parameter (imjmax=100000) |
---|
53 | real px(klon_omp,lllm) |
---|
54 | real xglo(klon_glo,lllm) |
---|
55 | real zx(nbp_lon,nbp_lat,lllm) |
---|
56 | |
---|
57 | |
---|
58 | |
---|
59 | CALL Gather(px,xglo) |
---|
60 | !$OMP MASTER |
---|
61 | IF (is_mpi_root) THEN |
---|
62 | CALL Grid1Dto2D_glo(xglo,zx) |
---|
63 | call iotd_ecrit(nom,lllm,titre,unite,zx) |
---|
64 | ENDIF |
---|
65 | !$OMP END MASTER |
---|
66 | |
---|
67 | return |
---|
68 | end |
---|
69 | |
---|
70 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
71 | ! Version avec reindexation pour appeler depuis les routines internes |
---|
72 | ! à la sous surface |
---|
73 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
74 | |
---|
75 | |
---|
76 | subroutine iophys_ecrit_index(nom,lllm,titre,unite,knon,knindex,px) |
---|
77 | |
---|
78 | USE mod_phys_lmdz_para, ONLY: klon_omp |
---|
79 | USE dimphy, ONLY: klon |
---|
80 | USE mod_grid_phy_lmdz, ONLY: klon_glo |
---|
81 | IMPLICIT NONE |
---|
82 | |
---|
83 | ! This subroutine returns the sea surface temperature already read from limit.nc |
---|
84 | |
---|
85 | ! Arguments on input: |
---|
86 | INTEGER lllm |
---|
87 | CHARACTER (len=*) :: nom,titre,unite |
---|
88 | REAL px(klon_omp,lllm) |
---|
89 | INTEGER, INTENT(IN) :: knon ! nomber of points on compressed grid |
---|
90 | INTEGER, DIMENSION(klon), INTENT(IN) :: knindex ! grid point number for compressed grid |
---|
91 | REAL, DIMENSION(klon,lllm) :: varout |
---|
92 | |
---|
93 | INTEGER :: i,l |
---|
94 | |
---|
95 | IF (klon/=klon_omp) THEN |
---|
96 | print*,'klon, klon_omp',klon,klon_omp |
---|
97 | CALL abort_physic('iophys_ecrit','probleme de dimension parallele',1) |
---|
98 | ENDIF |
---|
99 | |
---|
100 | varout(1:klon,1:lllm)=0. |
---|
101 | DO l = 1, lllm |
---|
102 | DO i = 1, knon |
---|
103 | varout(knindex(i),l) = px(i,l) |
---|
104 | END DO |
---|
105 | END DO |
---|
106 | CALL iophys_ecrit(nom,lllm,titre,unite,varout) |
---|
107 | |
---|
108 | END SUBROUTINE iophys_ecrit_index |
---|
109 | |
---|
110 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
111 | SUBROUTINE iophys_ini(timestep) |
---|
112 | USE mod_phys_lmdz_para, ONLY: is_mpi_root |
---|
113 | USE vertical_layers_mod, ONLY: presnivs |
---|
114 | USE regular_lonlat_mod, ONLY: lon_reg, lat_reg |
---|
115 | USE dimphy, ONLY: klev |
---|
116 | USE mod_grid_phy_lmdz, ONLY: klon_glo |
---|
117 | USE time_phylmdz_mod, ONLY: annee_ref, day_ref, day_ini |
---|
118 | USE phys_cal_mod, ONLY: calend |
---|
119 | |
---|
120 | IMPLICIT NONE |
---|
121 | |
---|
122 | include "YOMCST.h" |
---|
123 | !======================================================================= |
---|
124 | |
---|
125 | ! Auteur: L. Fairhead , P. Le Van, Y. Wanherdrick, F. Forget |
---|
126 | ! ------- |
---|
127 | |
---|
128 | ! Objet: |
---|
129 | ! ------ |
---|
130 | |
---|
131 | ! 'Initialize' the diagfi.nc file: write down dimensions as well |
---|
132 | ! as time-independent fields (e.g: geopotential, mesh area, ...) |
---|
133 | |
---|
134 | !======================================================================= |
---|
135 | !----------------------------------------------------------------------- |
---|
136 | ! Declarations: |
---|
137 | ! ------------- |
---|
138 | |
---|
139 | real pi |
---|
140 | INTEGER nlat_eff |
---|
141 | INTEGER jour0,mois0,an0 |
---|
142 | REAL timestep,t0 |
---|
143 | CHARACTER(len=20) :: calendrier |
---|
144 | |
---|
145 | ! Arguments: |
---|
146 | ! ---------- |
---|
147 | |
---|
148 | |
---|
149 | !$OMP MASTER |
---|
150 | IF (is_mpi_root) THEN |
---|
151 | |
---|
152 | ! Bidouille pour gerer le fait que lat_reg contient deux latitudes |
---|
153 | ! en version uni-dimensionnelle (chose qui pourrait être résolue |
---|
154 | ! par ailleurs) |
---|
155 | IF (klon_glo==1) THEN |
---|
156 | nlat_eff=1 |
---|
157 | ELSE |
---|
158 | nlat_eff=size(lat_reg) |
---|
159 | ENDIF |
---|
160 | pi=2.*asin(1.) |
---|
161 | |
---|
162 | ! print*,'day_ini,annee_ref,day_ref',day_ini,annee_ref,day_ref |
---|
163 | ! print*,'jD_ref,jH_ref,start_time, calend',jD_ref,jH_ref,start_time, calend |
---|
164 | |
---|
165 | ! Attention : les lignes ci dessous supposent un calendrier en 360 jours |
---|
166 | ! Pourrait être retravaillé |
---|
167 | |
---|
168 | jour0=day_ref-30*(day_ref/30) |
---|
169 | mois0=day_ref/30+1 |
---|
170 | an0=annee_ref |
---|
171 | !FH BIZARE QUAND 1D ... t0=(day_ini-1)*RDAY |
---|
172 | t0=0. |
---|
173 | calendrier=calend |
---|
174 | if ( calendrier == "earth_360d" ) calendrier="360_day" |
---|
175 | |
---|
176 | print*,'iophys_ini annee_ref day_ref',annee_ref,day_ref,day_ini,calend,t0 |
---|
177 | |
---|
178 | |
---|
179 | call iotd_ini('phys.nc', & |
---|
180 | size(lon_reg),nlat_eff,klev,lon_reg(:)*180./pi,lat_reg*180./pi,presnivs,jour0,mois0,an0,t0,timestep,calendrier) |
---|
181 | ENDIF |
---|
182 | !$OMP END MASTER |
---|
183 | |
---|
184 | END |
---|
185 | |
---|
186 | #ifdef und |
---|
187 | SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit) |
---|
188 | IMPLICIT none |
---|
189 | |
---|
190 | !======================================================================= |
---|
191 | INTEGER nfield,nlon,iim,jjmp1, jjm |
---|
192 | REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield) |
---|
193 | |
---|
194 | INTEGER i, n, ig |
---|
195 | |
---|
196 | jjm = jjmp1 - 1 |
---|
197 | DO n = 1, nfield |
---|
198 | DO i=1,iim |
---|
199 | ecrit(i,n) = fi(1,n) |
---|
200 | ecrit(i+jjm*iim,n) = fi(nlon,n) |
---|
201 | ENDDO |
---|
202 | DO ig = 1, nlon - 2 |
---|
203 | ecrit(iim+ig,n) = fi(1+ig,n) |
---|
204 | ENDDO |
---|
205 | ENDDO |
---|
206 | RETURN |
---|
207 | END |
---|
208 | |
---|
209 | #endif |
---|
210 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
211 | ! Interface pour ecrire en netcdf avec les routines d'enseignement |
---|
212 | ! iotd de Frederic Hourdin |
---|
213 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
214 | |
---|
215 | SUBROUTINE iotd_ecrit_seq(nom,lllm,titre,unite,px) |
---|
216 | |
---|
217 | IMPLICIT NONE |
---|
218 | |
---|
219 | ! px arrive |
---|
220 | |
---|
221 | INCLUDE "iotd.h" |
---|
222 | |
---|
223 | |
---|
224 | ! Arguments on input: |
---|
225 | integer lllm |
---|
226 | character (len=*) :: nom,titre,unite |
---|
227 | integer imjmax |
---|
228 | parameter (imjmax=100000) |
---|
229 | real px(imjmax*lllm) |
---|
230 | real, allocatable :: zx(:,:,:) |
---|
231 | integer i,j,l,ijl |
---|
232 | |
---|
233 | allocate(zx(imax,jmax,lllm)) |
---|
234 | |
---|
235 | ijl=0 |
---|
236 | do l=1,lllm |
---|
237 | ! Pole nord |
---|
238 | ijl=ijl+1 |
---|
239 | do i=1,imax |
---|
240 | zx(i,1,l)=px(ijl) |
---|
241 | enddo |
---|
242 | ! Grille normale |
---|
243 | do j=2,jmax-1 |
---|
244 | do i=1,imax |
---|
245 | ijl=ijl+1 |
---|
246 | zx(i,j,l)=px(ijl) |
---|
247 | enddo |
---|
248 | enddo |
---|
249 | ! Pole sud |
---|
250 | if ( jmax > 1 ) then |
---|
251 | ijl=ijl+1 |
---|
252 | do i=1,imax |
---|
253 | zx(i,jmax,l)=px(ijl) |
---|
254 | enddo |
---|
255 | endif |
---|
256 | enddo |
---|
257 | |
---|
258 | call iotd_ecrit(nom,lllm,titre,unite,zx) |
---|
259 | deallocate(zx) |
---|
260 | |
---|
261 | return |
---|
262 | end |
---|
263 | |
---|