source: LMDZ6/branches/Amaury_dev/libf/phylmd/iophys.F90 @ 5135

Last change on this file since 5135 was 5135, checked in by abarral, 8 weeks ago

Put iotd* into lmdz_iotd.f90

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