source: LMDZ6/trunk/libf/phylmd/iophys.F90 @ 5322

Last change on this file since 5322 was 5291, checked in by abarral, 3 weeks ago

Move thermcell_old.h iotd.h to module

File size: 7.4 KB
Line 
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
86! Arguments on input:
87    INTEGER lllm
88    CHARACTER (len=*) :: nom,titre,unite
89    REAL px(klon_omp,lllm)
90    INTEGER, INTENT(IN)                  :: knon     ! nomber of points on compressed grid
91    INTEGER, DIMENSION(klon), INTENT(IN) :: knindex  ! grid point number for compressed grid
92    REAL, DIMENSION(klon,lllm) :: varout
93
94    INTEGER :: i,l
95
96    IF (klon/=klon_omp) THEN
97      print*,'klon, klon_omp',klon,klon_omp
98      CALL abort_physic('iophys_ecrit','probleme de dimension parallele',1)
99    ENDIF
100
101    varout(1:klon,1:lllm)=0.
102    DO l = 1, lllm
103    DO i = 1, knon
104       varout(knindex(i),l) = px(i,l)
105    END DO
106    END DO
107    CALL iophys_ecrit(nom,lllm,titre,unite,varout)
108
109  END SUBROUTINE iophys_ecrit_index
110
111!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
112      SUBROUTINE iophys_ini(timestep)
113      USE mod_phys_lmdz_para, ONLY: is_mpi_root
114      USE vertical_layers_mod, ONLY: presnivs
115      USE regular_lonlat_mod, ONLY: lon_reg, lat_reg
116      USE dimphy, ONLY: klev
117      USE mod_grid_phy_lmdz, ONLY: klon_glo
118      USE time_phylmdz_mod, ONLY : annee_ref, day_ref, day_ini
119      USE phys_cal_mod, ONLY : calend
120
121      USE yomcst_mod_h
122IMPLICIT NONE
123
124
125!=======================================================================
126!
127!   Auteur:  L. Fairhead  ,  P. Le Van, Y. Wanherdrick, F. Forget
128!   -------
129!
130!   Objet:
131!   ------
132!
133!   'Initialize' the diagfi.nc file: write down dimensions as well
134!   as time-independent fields (e.g: geopotential, mesh area, ...)
135!
136!=======================================================================
137!-----------------------------------------------------------------------
138!   Declarations:
139!   -------------
140
141real pi
142INTEGER nlat_eff
143INTEGER jour0,mois0,an0
144REAL timestep,t0
145CHARACTER(len=20) :: calendrier
146
147!   Arguments:
148!   ----------
149
150
151!$OMP MASTER
152    IF (is_mpi_root) THEN       
153
154! Bidouille pour gerer le fait que lat_reg contient deux latitudes
155! en version uni-dimensionnelle (chose qui pourrait être résolue
156! par ailleurs)
157IF (klon_glo==1) THEN
158   nlat_eff=1
159ELSE
160   nlat_eff=size(lat_reg)
161ENDIF
162pi=2.*asin(1.)
163
164! print*,'day_ini,annee_ref,day_ref',day_ini,annee_ref,day_ref
165! print*,'jD_ref,jH_ref,start_time, calend',jD_ref,jH_ref,start_time, calend
166
167! Attention : les lignes ci dessous supposent un calendrier en 360 jours
168! Pourrait être retravaillé
169
170jour0=day_ref-30*(day_ref/30)
171mois0=day_ref/30+1
172an0=annee_ref
173!FH BIZARE QUAND 1D ...  t0=(day_ini-1)*RDAY
174t0=0.
175calendrier=calend
176if ( calendrier == "earth_360d" ) calendrier="360_day"
177
178print*,'iophys_ini annee_ref day_ref',annee_ref,day_ref,day_ini,calend,t0
179
180
181call iotd_ini('phys.nc', &
182size(lon_reg),nlat_eff,klev,lon_reg(:)*180./pi,lat_reg*180./pi,presnivs,jour0,mois0,an0,t0,timestep,calendrier)
183    ENDIF
184!$OMP END MASTER
185
186      END
187
188#ifdef und
189      SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit)
190      IMPLICIT none
191
192!=======================================================================
193      INTEGER nfield,nlon,iim,jjmp1, jjm
194      REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield)
195
196      INTEGER i, n, ig
197
198      jjm = jjmp1 - 1
199      DO n = 1, nfield
200         DO i=1,iim
201            ecrit(i,n) = fi(1,n)
202            ecrit(i+jjm*iim,n) = fi(nlon,n)
203         ENDDO
204         DO ig = 1, nlon - 2
205           ecrit(iim+ig,n) = fi(1+ig,n)
206         ENDDO
207      ENDDO
208      RETURN
209      END
210
211#endif
212!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
213! Interface pour ecrire en netcdf avec les routines d'enseignement
214! iotd de Frederic Hourdin
215!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
216
217      SUBROUTINE iotd_ecrit_seq(nom,lllm,titre,unite,px)
218        USE iotd_mod_h
219
220      IMPLICIT NONE
221
222
223! Arguments on input:
224      integer lllm
225      character (len=*) :: nom,titre,unite
226      integer imjmax
227      parameter (imjmax=100000)
228      real px(imjmax*lllm)
229      real, allocatable :: zx(:,:,:)
230      integer i,j,l,ijl
231
232      allocate(zx(imax,jmax,lllm))
233
234      ijl=0
235      do l=1,lllm
236         ! Pole nord
237         ijl=ijl+1
238         do i=1,imax
239            zx(i,1,l)=px(ijl)
240         enddo
241         ! Grille normale
242         do j=2,jmax-1
243            do i=1,imax
244               ijl=ijl+1
245               zx(i,j,l)=px(ijl)
246            enddo
247         enddo
248         ! Pole sud
249         if ( jmax > 1 ) then
250            ijl=ijl+1
251            do i=1,imax
252               zx(i,jmax,l)=px(ijl)
253            enddo
254         endif
255      enddo
256
257      call iotd_ecrit(nom,lllm,titre,unite,zx)
258      deallocate(zx)
259
260      return
261      end
262
Note: See TracBrowser for help on using the repository browser.