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

Last change on this file since 5451 was 5451, checked in by fhourdin, 21 hours ago

Concerns iophys_ini

File size: 8.0 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 subroutine iophys_ecrit
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,nlev)
113      USE dimphy, ONLY: klev
114      USE mod_phys_lmdz_para, ONLY: is_mpi_root
115      USE vertical_layers_mod, ONLY: presnivs
116      USE regular_lonlat_mod, ONLY: lon_reg, lat_reg
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
141integer, intent(in) :: nlev
142real, intent(in) :: timestep
143
144real pi
145INTEGER nlat_eff
146INTEGER jour0,mois0,an0
147REAL t0
148CHARACTER(len=20) :: calendrier
149integer ilev
150real coord_vert(nlev)
151
152!   Arguments:
153!   ----------
154
155
156!$OMP MASTER
157    IF (is_mpi_root) THEN       
158
159! Bidouille pour gerer le fait que lat_reg contient deux latitudes
160! en version uni-dimensionnelle (chose qui pourrait être résolue
161! par ailleurs)
162IF (klon_glo==1) THEN
163   nlat_eff=1
164ELSE
165   nlat_eff=size(lat_reg)
166ENDIF
167pi=2.*asin(1.)
168
169! print*,'day_ini,annee_ref,day_ref',day_ini,annee_ref,day_ref
170! print*,'jD_ref,jH_ref,start_time, calend',jD_ref,jH_ref,start_time, calend
171
172! Attention : les lignes ci dessous supposent un calendrier en 360 jours
173! Pourrait être retravaillé
174
175jour0=day_ref-30*(day_ref/30)
176mois0=day_ref/30+1
177an0=annee_ref
178!FH BIZARE QUAND 1D ...  t0=(day_ini-1)*RDAY
179t0=0.
180calendrier=calend
181if ( calendrier == "earth_360d" ) calendrier="360_day"
182
183print*,'iophys_ini annee_ref day_ref',annee_ref,day_ref,day_ini,calend,t0
184
185if ( nlev == klev ) then
186     coord_vert=presnivs
187print*,'ON EST LA '
188else
189     do ilev=1,nlev
190        coord_vert(ilev)=ilev
191     enddo
192endif
193print*,'nlev=',nlev
194print*,'coord_vert',coord_vert
195call iotd_ini('phys.nc', &
196size(lon_reg),nlat_eff,nlev,lon_reg(:)*180./pi,lat_reg*180./pi,coord_vert,jour0,mois0,an0,t0,timestep,calendrier)
197 !  SUBROUTINE iotd_ini(fichnom,iim,jjm,llm,prlon,prlat,pcoordv,jour0,mois0,an0,t0,timestep,calendrier)
198!   -------
199    ENDIF
200!$OMP END MASTER
201
202      END SUBROUTINE iophys_ini
203
204#ifdef und
205      SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit)
206      IMPLICIT none
207
208!=======================================================================
209      INTEGER nfield,nlon,iim,jjmp1, jjm
210      REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield)
211
212      INTEGER i, n, ig
213
214      jjm = jjmp1 - 1
215      DO n = 1, nfield
216         DO i=1,iim
217            ecrit(i,n) = fi(1,n)
218            ecrit(i+jjm*iim,n) = fi(nlon,n)
219         ENDDO
220         DO ig = 1, nlon - 2
221           ecrit(iim+ig,n) = fi(1+ig,n)
222         ENDDO
223      ENDDO
224      RETURN
225      END SUBROUTINE gr_fi_ecrit
226
227#endif
228!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
229! Interface pour ecrire en netcdf avec les routines d'enseignement
230! iotd de Frederic Hourdin
231!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
232
233      SUBROUTINE iotd_ecrit_seq(nom,lllm,titre,unite,px)
234!call iotd_ecrit_seq('f0',1,'f0 in thermcell_plume_6A',' ',f0(1:ngrid))
235
236        USE iotd_mod_h
237
238      IMPLICIT NONE
239
240
241! Arguments on input:
242      integer lllm
243      character (len=*) :: nom,titre,unite
244      integer imjmax
245      parameter (imjmax=100000)
246      real px(imjmax*lllm)
247      real, allocatable :: zx(:,:,:)
248      integer i,j,l,ijl
249
250      !print*,'iotd_ecrit_seq ,nom,lllm,titre,unite,px',nom,lllm,titre,unite,px
251      allocate(zx(imax,jmax,lllm))
252
253      ijl=0
254      do l=1,lllm
255         ! Pole nord
256         ijl=ijl+1
257         do i=1,imax
258            zx(i,1,l)=px(ijl)
259         enddo
260         ! Grille normale
261         do j=2,jmax-1
262            do i=1,imax
263               ijl=ijl+1
264               zx(i,j,l)=px(ijl)
265            enddo
266         enddo
267         ! Pole sud
268         if ( jmax > 1 ) then
269            ijl=ijl+1
270            do i=1,imax
271               zx(i,jmax,l)=px(ijl)
272            enddo
273         endif
274      enddo
275
276      call iotd_ecrit(nom,lllm,titre,unite,zx)
277      deallocate(zx)
278
279      return
280      end subroutine iotd_ecrit_seq
281
Note: See TracBrowser for help on using the repository browser.