source: LMDZ4/branches/unlabeled-1.1.1/libf/dyn3d/startvar.F @ 3722

Last change on this file since 3722 was 524, checked in by lmdzadmin, 20 years ago

Initial revision

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 39.5 KB
Line 
1!
2! $Header$
3!
4C
5C
6      MODULE startvar
7    !
8    !
9    !      There are three ways to access data from the database of atmospheric data which
10    !       can be used to initialize the model. This depends on the type of field which needs
11    !       to be extracted. In any case the call should come after a restget and should be of the type :
12    !                CALL startget(...)
13    !
14    !       We will details the possible arguments to startget here :
15    !
16    !        - A 2D variable on the dynamical grid :
17    !           CALL startget(varname, iml, jml, lon_in, lat_in, champ, val_ex, jml2, lon_in2, lat_in2, interbar )             
18    !
19    !        - A 1D variable on the physical grid :
20    !            CALL startget(varname, iml, jml, lon_in, lat_in, nbindex, champ, val_exp, jml2, lon_in2, lat_in2, interbar )
21    !
22    !
23    !         - A 3D variable on the dynamical grid :
24    !            CALL startget(varname, iml, jml, lon_in, lat_in, lml, pls, workvar, champ, val_exp, jml2, lon_in2, lat_in2, interbar )
25    !
26    !
27    !         There is special constraint on the atmospheric data base except that the
28    !         the data needs to be in netCDF and the variables should have the the following
29    !        names in the file :
30    !
31    !      'RELIEF'  : High resolution orography
32    !       'ST'            : Surface temperature
33    !       'CDSW'     : Soil moisture
34    !       'Z'               : Surface geopotential
35    !       'SP'            : Surface pressure
36    !        'U'              : East ward wind
37    !        'V'              : Northward wind
38    !        'TEMP'             : Temperature
39    !        'R'             : Relative humidity
40    !     
41      USE ioipsl
42    !
43    !
44      IMPLICIT NONE
45    !
46    !
47      PRIVATE
48      PUBLIC startget
49    !
50    !
51      INTERFACE startget
52        MODULE PROCEDURE startget_phys2d, startget_phys1d, startget_dyn
53      END INTERFACE
54    !
55      INTEGER, SAVE :: fid_phys, fid_dyn
56      INTEGER, SAVE  :: iml_phys, iml_rel, iml_dyn
57      INTEGER, SAVE :: jml_phys,  jml_rel, jml_dyn
58      INTEGER, SAVE ::  llm_dyn, ttm_dyn
59      REAL, ALLOCATABLE, SAVE, DIMENSION (:,:)  :: lon_phys, lon_rug,
60     . lon_alb, lon_rel, lon_dyn
61      REAL, ALLOCATABLE, SAVE, DIMENSION (:,:)  :: lat_phys, lat_rug,
62     . lat_alb, lat_rel, lat_dyn
63      REAL, ALLOCATABLE, SAVE, DIMENSION (:)  :: levdyn_ini
64      REAL, ALLOCATABLE, SAVE, DIMENSION (:,:)  :: relief, zstd, zsig,
65     . zgam, zthe, zpic, zval
66      REAL, ALLOCATABLE, SAVE, DIMENSION (:,:)  :: rugo, masque, phis
67    !
68      REAL, ALLOCATABLE, SAVE, DIMENSION (:,:)  :: tsol, qsol, psol_dyn
69      REAL, ALLOCATABLE, SAVE, DIMENSION (:,:,:)  ::   var_ana3d
70    !
71      CONTAINS
72    !
73    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
74    !
75      SUBROUTINE startget_phys2d(varname, iml, jml, lon_in, lat_in,
76     . champ, val_exp, jml2, lon_in2, lat_in2 , interbar, masque_lu )
77    !
78    !    There is a big mess with the size in logitude, should it be iml or iml+1.
79    !    I have chosen to use the iml+1 as an argument to this routine and we declare
80    !   internaly smaler fields when needed. This needs to be cleared once and for all in LMDZ.
81    !  A convention is required.
82    !
83    !
84      CHARACTER*(*), INTENT(in) :: varname
85      INTEGER, INTENT(in) :: iml, jml ,jml2
86      REAL, INTENT(in) :: lon_in(iml), lat_in(jml)
87      REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2)
88      REAL, INTENT(inout) :: champ(iml,jml)
89      REAL, INTENT(in) :: val_exp
90      REAL, INTENT(in), optional :: masque_lu(iml,jml)
91      LOGICAL interbar
92    !
93    !   This routine only works if the variable does not exist or is constant
94    !
95      IF ( MINVAL(champ(:,:)).EQ.MAXVAL(champ(:,:)) .AND.
96     .MINVAL(champ(:,:)).EQ.val_exp ) THEN
97          !
98          SELECTCASE(varname)
99              !
100              CASE ('relief')
101                  !
102                  !  If we do not have the orography we need to get it
103                  !
104                  IF ( .NOT.ALLOCATED(relief)) THEN
105                      !
106                    if (present(masque_lu)) then
107                      CALL start_init_orog( iml, jml, lon_in, lat_in,
108     .                    jml2,lon_in2,lat_in2, interbar, masque_lu )
109                    else
110                      CALL start_init_orog( iml, jml, lon_in, lat_in,
111     .                    jml2,lon_in2,lat_in2, interbar)
112                    endif
113                      !
114                  ENDIF
115                  !
116                  IF (SIZE(relief) .NE. SIZE(champ)) THEN
117                      !
118                      WRITE(*,*) 'STARTVAR module has been',
119     .' initialized to the wrong size'
120                      STOP
121                      !
122                  ENDIF
123                  !
124                  champ(:,:) = relief(:,:)
125                  !
126              CASE ('rugosite')
127                  !
128                  !  If we do not have the orography we need to get it
129                  !
130                  IF ( .NOT.ALLOCATED(rugo)) THEN
131                      !
132                      CALL start_init_orog( iml, jml, lon_in, lat_in,
133     .                    jml2,lon_in2,lat_in2 , interbar )
134                      !
135                  ENDIF
136                  !
137                  IF (SIZE(rugo) .NE. SIZE(champ)) THEN
138                      !
139                      WRITE(*,*)
140     .  'STARTVAR module has been initialized to the wrong size'
141                      STOP
142                      !
143                  ENDIF
144                  !
145                  champ(:,:) = rugo(:,:)
146                  !
147              CASE ('masque')
148                  !
149                  !  If we do not have the orography we need to get it
150                  !
151                  IF ( .NOT.ALLOCATED(masque)) THEN
152                      !
153                      CALL start_init_orog( iml, jml, lon_in, lat_in,
154     .                     jml2,lon_in2,lat_in2 , interbar )
155                      !
156                  ENDIF
157                  !
158                  IF (SIZE(masque) .NE. SIZE(champ)) THEN
159                      !
160                      WRITE(*,*)
161     .   'STARTVAR module has been initialized to the wrong size'
162                      STOP
163                      !
164                  ENDIF
165                  !
166                  champ(:,:) = masque(:,:)
167                  !
168              CASE ('surfgeo')
169                  !
170                  !  If we do not have the orography we need to get it
171                  !
172                  IF ( .NOT.ALLOCATED(phis)) THEN
173                      !
174                      CALL start_init_orog( iml, jml, lon_in, lat_in,
175     .                   jml2,lon_in2, lat_in2 , interbar )
176                      !
177                  ENDIF
178                  !
179                  IF (SIZE(phis) .NE. SIZE(champ)) THEN
180                      !
181                      WRITE(*,*)
182     . 'STARTVAR module has been initialized to the wrong size'
183                      STOP
184                      !
185                  ENDIF
186                  !
187                  champ(:,:) = phis(:,:)
188                  !
189              CASE ('psol')
190                  !
191                  !  If we do not have the orography we need to get it
192                  !
193                  IF ( .NOT.ALLOCATED(psol_dyn)) THEN
194                      !
195                      CALL start_init_dyn( iml, jml, lon_in, lat_in,
196     .                   jml2,lon_in2, lat_in2 , interbar )
197                      !
198                  ENDIF
199                  !
200                  IF (SIZE(psol_dyn) .NE. SIZE(champ)) THEN
201                      !
202                      WRITE(*,*)
203     . 'STARTVAR module has been initialized to the wrong size'
204                      STOP
205                      !
206                  ENDIF
207                  !
208                  champ(:,:) = psol_dyn(:,:)
209                  !
210              CASE DEFAULT
211                  !
212                  WRITE(*,*) 'startget_phys2d'
213                  WRITE(*,*) 'No rule is present to extract variable',
214     .                 varname(:LEN_TRIM(varname)),' from any data set'
215                  STOP
216                  !
217          END SELECT
218          !
219      ELSE
220          !
221          ! There are a few fields we might need if we need to interpolate 3D filed. Thus if they come through here we
222          ! will catch them
223          !
224          SELECTCASE(varname)
225              !
226              CASE ('surfgeo')
227                  !
228                  IF ( .NOT.ALLOCATED(phis)) THEN
229                      ALLOCATE(phis(iml,jml))
230                  ENDIF
231                  !
232                  IF (SIZE(phis) .NE. SIZE(champ)) THEN
233                      !
234                      WRITE(*,*)
235     .  'STARTVAR module has been initialized to the wrong size'
236                      STOP
237                      !
238                  ENDIF
239                  !
240                  phis(:,:) = champ(:,:)
241                  !
242          END SELECT
243          !
244      ENDIF
245    !
246      END SUBROUTINE startget_phys2d
247    !
248    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
249    !
250      SUBROUTINE start_init_orog ( iml,jml,lon_in, lat_in,jml2,lon_in2 ,
251     ,   lat_in2 , interbar, masque_lu )
252    !
253      INTEGER, INTENT(in) :: iml, jml, jml2
254      REAL, INTENT(in) :: lon_in(iml), lat_in(jml)
255      REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2)
256      REAL, intent(in), optional :: masque_lu(iml,jml)
257      LOGICAL interbar
258    !
259    !  LOCAL
260    !
261      LOGICAL interbar2
262      REAL :: lev(1), date, dt,chmin,chmax
263      INTEGER :: itau(1), fid
264      INTEGER ::  llm_tmp, ttm_tmp
265      INTEGER :: i, j
266      INTEGER :: iret
267      CHARACTER*25 title
268      REAL, ALLOCATABLE :: relief_hi(:,:)
269      REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:)
270      REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:)
271      REAL, ALLOCATABLE :: tmp_var(:,:)
272      INTEGER, ALLOCATABLE :: tmp_int(:,:)
273    !
274      CHARACTER*120 :: orogfname
275      LOGICAL :: check=.TRUE.
276    !
277    !
278      orogfname = 'Relief.nc'
279    !
280      IF ( check ) WRITE(*,*) 'Reading the high resolution orography'
281    !
282      CALL flininfo(orogfname,iml_rel, jml_rel, llm_tmp, ttm_tmp, fid)
283    !
284      ALLOCATE (lat_rel(iml_rel,jml_rel), stat=iret)
285      ALLOCATE (lon_rel(iml_rel,jml_rel), stat=iret)
286      ALLOCATE (relief_hi(iml_rel,jml_rel), stat=iret)
287    !
288      CALL flinopen(orogfname, .FALSE., iml_rel, jml_rel,
289     .llm_tmp, lon_rel, lat_rel, lev, ttm_tmp,
290     .      itau, date, dt, fid)
291    !
292      CALL flinget(fid, 'RELIEF', iml_rel, jml_rel, llm_tmp,
293     . ttm_tmp, 1, 1, relief_hi)
294    !
295      CALL flinclo(fid)
296    !
297    !   In case we have a file which is in degrees we do the transformation
298    !
299      ALLOCATE(lon_rad(iml_rel))
300      ALLOCATE(lon_ini(iml_rel))
301
302      IF ( MAXVAL(lon_rel(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
303          lon_ini(:) = lon_rel(:,1) * 2.0 * ASIN(1.0) / 180.0
304      ELSE
305          lon_ini(:) = lon_rel(:,1)
306      ENDIF
307
308      ALLOCATE(lat_rad(jml_rel))
309      ALLOCATE(lat_ini(jml_rel))
310
311      IF ( MAXVAL(lat_rel(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
312          lat_ini(:) = lat_rel(1,:) * 2.0 * ASIN(1.0) / 180.0
313      ELSE
314          lat_ini(:) = lat_rel(1,:)
315      ENDIF
316    !
317    !
318
319      title='RELIEF'
320
321      interbar2 = .FALSE.
322      CALL conf_dat2d(title,iml_rel, jml_rel, lon_ini, lat_ini,
323     . lon_rad, lat_rad, relief_hi , interbar2  )
324
325      IF ( check ) WRITE(*,*) 'Computes all the parameters needed',
326     .' for the gravity wave drag code'
327    !
328    !    Allocate the data we need to put in the interpolated fields
329    !
330    !            RELIEF:  orographie moyenne
331      ALLOCATE(relief(iml,jml))
332    !            zphi :  orographie moyenne
333      ALLOCATE(phis(iml,jml))
334    !             zstd:  deviation standard de l'orographie sous-maille
335      ALLOCATE(zstd(iml,jml))
336    !             zsig:  pente de l'orographie sous-maille
337      ALLOCATE(zsig(iml,jml))
338    !             zgam:  anisotropy de l'orographie sous maille
339      ALLOCATE(zgam(iml,jml))
340    !             zthe:  orientation de l'axe oriente dans la direction
341    !                    de plus grande pente de l'orographie sous maille
342      ALLOCATE(zthe(iml,jml))
343    !             zpic:  hauteur pics de la SSO
344      ALLOCATE(zpic(iml,jml))
345    !             zval:  hauteur vallees de la SSO
346      ALLOCATE(zval(iml,jml))
347    !             masque : Masque terre ocean
348      ALLOCATE(tmp_int(iml,jml))
349      ALLOCATE(masque(iml,jml))
350
351      masque = -99999.
352      if (present(masque_lu)) then
353        masque = masque_lu
354      endif
355    !
356      CALL grid_noro(iml_rel, jml_rel, lon_rad, lat_rad, relief_hi,
357     . iml-1, jml, lon_in, lat_in,
358     . phis, relief, zstd, zsig, zgam, zthe, zpic, zval, masque)
359      phis = phis * 9.81
360    !
361!      masque(:,:) = FLOAT(tmp_int(:,:))
362    !
363    !  Compute surface roughness
364    !
365      IF ( check ) WRITE(*,*)
366     .'Compute surface roughness induced by the orography'
367    !
368      ALLOCATE(rugo(iml,jml))
369      ALLOCATE(tmp_var(iml-1,jml))
370    !
371      CALL rugsoro(iml_rel, jml_rel, lon_rad, lat_rad, relief_hi,
372     . iml-1, jml, lon_in, lat_in, tmp_var)
373    !
374      DO j = 1, jml
375        DO i = 1, iml-1
376          rugo(i,j) = tmp_var(i,j)
377        ENDDO
378        rugo(iml,j) = tmp_var(1,j)
379      ENDDO
380c
381cc   ***   rugo  n'est pas utilise pour l'instant  ******
382    !
383    !   Build land-sea mask
384    !
385    !
386      RETURN
387    !
388      END SUBROUTINE start_init_orog
389    !
390    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
391    !
392      SUBROUTINE startget_phys1d(varname, iml, jml, lon_in,
393     .lat_in, nbindex, champ, val_exp ,jml2, lon_in2, lat_in2,interbar)
394    !
395      CHARACTER*(*), INTENT(in) :: varname
396      INTEGER, INTENT(in) :: iml, jml, nbindex, jml2
397      REAL, INTENT(in) :: lon_in(iml), lat_in(jml)
398      REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2)
399      REAL, INTENT(inout) :: champ(nbindex)
400      REAL, INTENT(in) :: val_exp
401      LOGICAL interbar
402    !
403    !
404    !   This routine only works if the variable does not exist or is constant
405    !
406      IF ( MINVAL(champ(:)).EQ.MAXVAL(champ(:)) .AND.
407     .MINVAL(champ(:)).EQ.val_exp ) THEN
408          SELECTCASE(varname)
409            CASE ('tsol')
410              IF ( .NOT.ALLOCATED(tsol)) THEN
411                CALL start_init_phys( iml, jml, lon_in, lat_in,
412     .              jml2, lon_in2, lat_in2, interbar )
413              ENDIF
414              IF ( SIZE(tsol) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
415                WRITE(*,*)
416     . 'STARTVAR module has been initialized to the wrong size'
417                 STOP
418              ENDIF
419              CALL gr_dyn_fi(1, iml, jml, nbindex, tsol, champ)
420            CASE ('qsol')
421              IF ( .NOT.ALLOCATED(qsol)) THEN
422                CALL start_init_phys( iml, jml, lon_in, lat_in,
423     .              jml2, lon_in2,lat_in2 , interbar )
424              ENDIF
425              IF ( SIZE(qsol) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
426                WRITE(*,*)
427     . 'STARTVAR module has been initialized to the wrong size'
428                STOP
429              ENDIF
430              CALL gr_dyn_fi(1, iml, jml, nbindex, qsol, champ)
431            CASE ('psol')
432              IF ( .NOT.ALLOCATED(psol_dyn)) THEN
433                CALL start_init_dyn( iml, jml, lon_in, lat_in,
434     .              jml2, lon_in2,lat_in2 , interbar )
435              ENDIF
436              IF (SIZE(psol_dyn) .NE. SIZE(lon_in)*SIZE(lat_in)) THEN
437                WRITE(*,*)
438     . 'STARTVAR module has been initialized to the wrong size'
439                STOP
440              ENDIF
441              CALL gr_dyn_fi(1, iml, jml, nbindex, psol_dyn, champ)
442            CASE ('zmea')
443              IF ( .NOT.ALLOCATED(relief)) THEN
444                CALL start_init_orog( iml, jml, lon_in, lat_in,
445     .            jml2, lon_in2,lat_in2 , interbar )
446              ENDIF
447              IF ( SIZE(relief) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
448                WRITE(*,*)
449     . 'STARTVAR module has been initialized to the wrong size'
450                 STOP
451              ENDIF
452              CALL gr_dyn_fi(1, iml, jml, nbindex, relief, champ)
453            CASE ('zstd')
454              IF ( .NOT.ALLOCATED(zstd)) THEN
455                CALL start_init_orog( iml, jml, lon_in, lat_in,
456     .              jml2, lon_in2,lat_in2 , interbar )
457              ENDIF
458              IF ( SIZE(zstd) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
459                WRITE(*,*)
460     . 'STARTVAR module has been initialized to the wrong size'
461                 STOP
462              ENDIF
463              CALL gr_dyn_fi(1, iml, jml, nbindex,zstd, champ)
464            CASE ('zsig')
465              IF ( .NOT.ALLOCATED(zsig)) THEN
466                CALL start_init_orog( iml, jml, lon_in, lat_in,
467     .               jml2, lon_in2,lat_in2 , interbar )
468              ENDIF
469              IF ( SIZE(zsig) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
470                WRITE(*,*)
471     . 'STARTVAR module has been initialized to the wrong size'
472                 STOP
473              ENDIF
474              CALL gr_dyn_fi(1, iml, jml, nbindex,zsig, champ)
475            CASE ('zgam')
476              IF ( .NOT.ALLOCATED(zgam)) THEN
477                CALL start_init_orog( iml, jml, lon_in, lat_in,
478     .            jml2, lon_in2,lat_in2 , interbar )
479              ENDIF
480              IF ( SIZE(zgam) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
481                WRITE(*,*)
482     . 'STARTVAR module has been initialized to the wrong size'
483                 STOP
484              ENDIF
485              CALL gr_dyn_fi(1, iml, jml, nbindex,zgam, champ)
486            CASE ('zthe')
487              IF ( .NOT.ALLOCATED(zthe)) THEN
488                CALL start_init_orog( iml, jml, lon_in, lat_in,
489     .            jml2, lon_in2,lat_in2 , interbar )
490              ENDIF
491              IF ( SIZE(zthe) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
492                WRITE(*,*)
493     . 'STARTVAR module has been initialized to the wrong size'
494                 STOP
495              ENDIF
496              CALL gr_dyn_fi(1, iml, jml, nbindex,zthe, champ)
497            CASE ('zpic')
498              IF ( .NOT.ALLOCATED(zpic)) THEN
499                CALL start_init_orog( iml, jml, lon_in, lat_in,
500     .            jml2, lon_in2,lat_in2 , interbar )
501              ENDIF
502              IF ( SIZE(zpic) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
503                WRITE(*,*)
504     . 'STARTVAR module has been initialized to the wrong size'
505                 STOP
506              ENDIF
507              CALL gr_dyn_fi(1, iml, jml, nbindex,zpic, champ)
508            CASE ('zval')
509              IF ( .NOT.ALLOCATED(zval)) THEN
510                CALL start_init_orog( iml, jml, lon_in, lat_in,
511     .            jml2, lon_in2,lat_in2 , interbar )
512              ENDIF
513              IF ( SIZE(zval) .NE. SIZE(lon_in)*SIZE(lat_in) ) THEN
514                WRITE(*,*)
515     . 'STARTVAR module has been initialized to the wrong size'
516                 STOP
517              ENDIF
518              CALL gr_dyn_fi(1, iml, jml, nbindex,zval, champ)
519            CASE ('rads')
520                  champ(:) = 0.0
521            CASE ('snow')
522                  champ(:) = 0.0
523            CASE ('deltat')
524                  champ(:) = 0.0
525            CASE ('rugmer')
526                  champ(:) = 0.001
527            CASE ('agsno')
528                  champ(:) = 50.0
529            CASE DEFAULT
530              WRITE(*,*) 'startget_phys1d'
531              WRITE(*,*) 'No rule is present to extract variable  ',
532     . varname(:LEN_TRIM(varname)),' from any data set'
533              STOP
534          END SELECT
535      ELSE
536        !
537        ! If we see tsol we catch it as we may need it for a 3D interpolation
538        !
539        SELECTCASE(varname)
540          CASE ('tsol')
541            IF ( .NOT.ALLOCATED(tsol)) THEN
542              ALLOCATE(tsol(SIZE(lon_in),SIZE(lat_in) ))
543            ENDIF
544            CALL gr_fi_dyn(1, iml, jml, nbindex, champ, tsol)
545        END SELECT
546      ENDIF
547      END SUBROUTINE startget_phys1d
548    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
549    !
550      SUBROUTINE start_init_phys( iml, jml, lon_in, lat_in, jml2,
551     .                 lon_in2, lat_in2 , interbar )
552    !
553      INTEGER, INTENT(in) :: iml, jml ,jml2
554      REAL, INTENT(in) :: lon_in(iml), lat_in(jml)
555      REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2)
556      LOGICAL interbar
557    !
558    !  LOCAL
559    !
560      REAL :: lev(1), date, dt
561      INTEGER :: itau(1)
562      INTEGER ::  llm_tmp, ttm_tmp
563      INTEGER :: i, j
564      INTEGER :: iret
565    !
566      CHARACTER*25 title
567      CHARACTER*120 :: physfname
568      LOGICAL :: check=.TRUE.
569    !
570      REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:)
571      REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:)
572      REAL, ALLOCATABLE :: var_ana(:,:), tmp_var(:,:)
573    !
574      physfname = 'ECPHY.nc'
575    !
576      IF ( check ) WRITE(*,*) 'Opening the surface analysis'
577    !
578      CALL flininfo(physfname, iml_phys, jml_phys, llm_tmp,
579     . ttm_tmp, fid_phys)
580    !
581    !
582      ALLOCATE (lat_phys(iml_phys,jml_phys))
583      ALLOCATE (lon_phys(iml_phys,jml_phys))
584    !
585      CALL flinopen(physfname, .FALSE., iml_phys, jml_phys,
586     . llm_tmp, lon_phys, lat_phys, lev, ttm_tmp,
587     . itau, date, dt, fid_phys)
588    !
589    ! Allocate the space we will need to get the data out of this file
590    !
591      ALLOCATE(var_ana(iml_phys, jml_phys))
592    !
593    !   In case we have a file which is in degrees we do the transformation
594    !
595      if (allocated(lon_rad)) DEALLOCATE(lon_rad, stat=iret)
596      if (allocated(lon_ini)) DEALLOCATE(lon_ini, stat=iret)
597
598      ALLOCATE(lon_rad(iml_phys))
599      ALLOCATE(lon_ini(iml_phys))
600
601      IF ( MAXVAL(lon_phys(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
602          lon_ini(:) = lon_phys(:,1) * 2.0 * ASIN(1.0) / 180.0
603      ELSE
604          lon_ini(:) = lon_phys(:,1)
605      ENDIF
606
607      ALLOCATE(lat_rad(jml_phys))
608      ALLOCATE(lat_ini(jml_phys))
609
610      IF ( MAXVAL(lat_phys(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
611          lat_ini(:) = lat_phys(1,:) * 2.0 * ASIN(1.0) / 180.0
612      ELSE
613          lat_ini(:) = lat_phys(1,:)
614      ENDIF
615
616
617    !
618    !   We get the two standard varibales
619    !   Surface temperature
620    !
621      ALLOCATE(tsol(iml,jml))
622      ALLOCATE(tmp_var(iml-1,jml))
623    !
624    !
625
626      CALL flinget(fid_phys, 'ST', iml_phys, jml_phys,
627     .llm_tmp, ttm_tmp, 1, 1, var_ana)
628
629      title='ST'
630      CALL conf_dat2d(title,iml_phys, jml_phys, lon_ini, lat_ini,
631     . lon_rad, lat_rad, var_ana , interbar  )
632
633      IF ( interbar )   THEN
634        WRITE(6,*) '-------------------------------------------------',
635     ,'--------------'
636        WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ',
637     , ' pour  ST $$$ '
638        WRITE(6,*) '-------------------------------------------------',
639     ,'--------------'
640        CALL inter_barxy ( iml_phys,jml_phys -1,lon_rad,lat_rad ,
641     ,   var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var   )
642      ELSE
643        CALL grille_m(iml_phys, jml_phys, lon_rad, lat_rad,
644     .    var_ana, iml-1, jml, lon_in, lat_in, tmp_var     )
645      ENDIF
646
647      CALL gr_int_dyn(tmp_var, tsol, iml-1, jml)
648    !
649    ! Soil moisture
650    !
651      ALLOCATE(qsol(iml,jml))
652      CALL flinget(fid_phys, 'CDSW', iml_phys, jml_phys,
653     . llm_tmp, ttm_tmp, 1, 1, var_ana)
654
655      title='CDSW'
656      CALL conf_dat2d(title,iml_phys, jml_phys, lon_ini, lat_ini,
657     . lon_rad, lat_rad, var_ana, interbar  )
658
659      IF ( interbar )   THEN
660        WRITE(6,*) '-------------------------------------------------',
661     ,'--------------'
662        WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ',
663     , ' pour  CDSW $$$ '
664        WRITE(6,*) '-------------------------------------------------',
665     ,'--------------'
666        CALL inter_barxy ( iml_phys,jml_phys -1,lon_rad,lat_rad ,
667     ,   var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var  )
668      ELSE
669        CALL grille_m(iml_phys, jml_phys, lon_rad, lat_rad,
670     .    var_ana, iml-1, jml, lon_in, lat_in, tmp_var     )
671      ENDIF
672c
673        CALL gr_int_dyn(tmp_var, qsol, iml-1, jml)
674    !
675       CALL flinclo(fid_phys)
676    !
677      END SUBROUTINE start_init_phys
678    !
679    !
680    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
681    !
682    !
683      SUBROUTINE startget_dyn(varname, iml, jml, lon_in, lat_in,
684     . lml, pls, workvar, champ, val_exp,jml2, lon_in2, lat_in2 ,
685     ,  interbar )
686    !
687    !   ARGUMENTS
688    !
689      CHARACTER*(*), INTENT(in) :: varname
690      INTEGER, INTENT(in) :: iml, jml, lml, jml2
691      REAL, INTENT(in) :: lon_in(iml), lat_in(jml)
692      REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2)
693      REAL, INTENT(in) :: pls(iml, jml, lml)
694      REAL, INTENT(in) :: workvar(iml, jml, lml)
695      REAL, INTENT(inout) :: champ(iml, jml, lml)
696      REAL, INTENT(in) :: val_exp
697      LOGICAL interbar
698    !
699    !    LOCAL
700    !
701      INTEGER :: il, ij, ii
702      REAL :: xppn, xpps
703    !
704#include "dimensions.h"
705#include "paramet.h"
706#include "comgeom2.h"
707#include "comconst.h"
708    !
709    !   This routine only works if the variable does not exist or is constant
710    !
711      IF ( MINVAL(champ(:,:,:)).EQ.MAXVAL(champ(:,:,:)) .AND.
712     . MINVAL(champ(:,:,:)).EQ.val_exp ) THEN
713        !
714        SELECTCASE(varname)
715          CASE ('u')
716            IF ( .NOT.ALLOCATED(psol_dyn)) THEN
717              CALL start_init_dyn( iml, jml, lon_in, lat_in, jml2 ,
718     .          lon_in2,lat_in2 , interbar )
719            ENDIF
720            CALL start_inter_3d('U', iml, jml, lml, lon_in,
721     .       lat_in, jml2, lon_in2, lat_in2,  pls, champ,interbar )
722            DO il=1,lml
723              DO ij=1,jml
724                DO ii=1,iml-1
725                  champ(ii,ij,il) = champ(ii,ij,il) * cu(ii,ij)
726                ENDDO
727                champ(iml,ij, il) = champ(1,ij, il)
728              ENDDO
729            ENDDO
730          CASE ('v')
731            IF ( .NOT.ALLOCATED(psol_dyn)) THEN
732              CALL start_init_dyn( iml, jml, lon_in, lat_in , jml2,
733     .           lon_in2, lat_in2 , interbar )
734            ENDIF
735            CALL start_inter_3d('V', iml, jml, lml, lon_in,
736     .       lat_in, jml2, lon_in2, lat_in2,  pls, champ, interbar )
737            DO il=1,lml
738              DO ij=1,jml
739                DO ii=1,iml-1
740                  champ(ii,ij,il) = champ(ii,ij,il) * cv(ii,ij)
741                ENDDO
742                champ(iml,ij, il) = champ(1,ij, il)
743              ENDDO
744            ENDDO
745          CASE ('t')
746            IF ( .NOT.ALLOCATED(psol_dyn)) THEN
747              CALL start_init_dyn( iml, jml, lon_in, lat_in, jml2 ,
748     .           lon_in2, lat_in2 ,interbar )
749            ENDIF
750            CALL start_inter_3d('TEMP', iml, jml, lml, lon_in,
751     .       lat_in, jml2, lon_in2, lat_in2,  pls, champ, interbar )
752 
753          CASE ('tpot')
754            IF ( .NOT.ALLOCATED(psol_dyn)) THEN
755              CALL start_init_dyn( iml, jml, lon_in, lat_in , jml2 ,
756     .            lon_in2, lat_in2 , interbar )
757            ENDIF
758            CALL start_inter_3d('TEMP', iml, jml, lml, lon_in,
759     .       lat_in, jml2, lon_in2, lat_in2,  pls, champ, interbar )
760            IF ( MINVAL(workvar(:,:,:)) .NE. MAXVAL(workvar(:,:,:)) )
761     .                                    THEN
762              DO il=1,lml
763                DO ij=1,jml
764                  DO ii=1,iml-1
765                    champ(ii,ij,il) = champ(ii,ij,il) * cpp
766     .                                 / workvar(ii,ij,il)
767                  ENDDO
768                  champ(iml,ij,il) = champ(1,ij,il)
769                ENDDO
770              ENDDO
771              DO il=1,lml
772                xppn = SUM(aire(:,1)*champ(:,1,il))/apoln
773                xpps = SUM(aire(:,jml)*champ(:,jml,il))/apols
774                champ(:,1,il) = xppn
775                champ(:,jml,il) = xpps
776              ENDDO
777            ELSE
778              WRITE(*,*)'Could not compute potential temperature as the'
779              WRITE(*,*)'Exner function is missing or constant.'
780              STOP
781            ENDIF
782          CASE ('q')
783            IF ( .NOT.ALLOCATED(psol_dyn)) THEN
784              CALL start_init_dyn( iml, jml, lon_in, lat_in, jml2 ,
785     .           lon_in2, lat_in2 , interbar )
786            ENDIF
787            CALL start_inter_3d('R', iml, jml, lml, lon_in, lat_in,
788     .        jml2, lon_in2, lat_in2,  pls, champ, interbar )
789            IF ( MINVAL(workvar(:,:,:)) .NE. MAXVAL(workvar(:,:,:)) )
790     .                                     THEN
791              DO il=1,lml
792                DO ij=1,jml
793                  DO ii=1,iml-1
794                    champ(ii,ij,il) = 0.01 * champ(ii,ij,il) *
795     .                                       workvar(ii,ij,il)
796                  ENDDO
797                  champ(iml,ij,il) = champ(1,ij,il)
798                ENDDO
799              ENDDO
800              WHERE ( champ .LT. 0.) champ = 1.0E-10
801              DO il=1,lml
802                xppn = SUM(aire(:,1)*champ(:,1,il))/apoln
803                xpps = SUM(aire(:,jml)*champ(:,jml,il))/apols
804                champ(:,1,il) = xppn
805                champ(:,jml,il) = xpps
806              ENDDO
807            ELSE
808              WRITE(*,*)'Could not compute specific humidity as the'
809              WRITE(*,*)'saturated humidity is missing or constant.'
810              STOP
811            ENDIF
812          CASE DEFAULT
813            WRITE(*,*) 'startget_dyn'
814            WRITE(*,*) 'No rule is present to extract variable  ',
815     . varname(:LEN_TRIM(varname)),' from any data set'
816            STOP
817          END SELECT
818      ENDIF
819      END SUBROUTINE startget_dyn
820    !
821    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
822    !
823      SUBROUTINE start_init_dyn( iml, jml, lon_in, lat_in,jml2,lon_in2 ,
824     ,             lat_in2 , interbar )
825    !
826      INTEGER, INTENT(in) :: iml, jml, jml2
827      REAL, INTENT(in) :: lon_in(iml), lat_in(jml)
828      REAL, INTENT(in) :: lon_in2(iml), lat_in2(jml2)
829      LOGICAL interbar
830    !
831    !  LOCAL
832    !
833      REAL :: lev(1), date, dt
834      INTEGER :: itau(1)
835      INTEGER :: i, j
836      integer :: iret
837    !
838      CHARACTER*120 :: physfname
839      LOGICAL :: check=.TRUE.
840    !
841      REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:)
842      REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:)
843      REAL, ALLOCATABLE :: var_ana(:,:), tmp_var(:,:), z(:,:)
844      REAL, ALLOCATABLE :: xppn(:), xpps(:)
845      LOGICAL :: allo
846    !
847    !
848#include "dimensions.h"
849#include "paramet.h"
850#include "comgeom2.h"
851
852      CHARACTER*25 title
853
854    !
855      physfname = 'ECDYN.nc'
856    !
857      IF ( check ) WRITE(*,*) 'Opening the surface analysis'
858    !
859      CALL flininfo(physfname, iml_dyn, jml_dyn, llm_dyn,
860     .                            ttm_dyn, fid_dyn)
861      IF ( check ) WRITE(*,*) 'Values read: ', iml_dyn, jml_dyn,
862     .                                         llm_dyn, ttm_dyn
863    !
864      ALLOCATE (lat_dyn(iml_dyn,jml_dyn), stat=iret)
865      ALLOCATE (lon_dyn(iml_dyn,jml_dyn), stat=iret)
866      ALLOCATE (levdyn_ini(llm_dyn), stat=iret)
867    !
868      CALL flinopen(physfname, .FALSE., iml_dyn, jml_dyn, llm_dyn,
869     . lon_dyn, lat_dyn, levdyn_ini, ttm_dyn,
870     . itau, date, dt, fid_dyn)
871    !
872
873      allo = allocated (var_ana)
874      if (allo) then
875        DEALLOCATE(var_ana, stat=iret)
876      endif
877      ALLOCATE(var_ana(iml_dyn, jml_dyn), stat=iret)
878
879      allo = allocated (lon_rad)
880      if (allo) then
881        DEALLOCATE(lon_rad, stat=iret)
882      endif
883
884      ALLOCATE(lon_rad(iml_dyn), stat=iret)
885      ALLOCATE(lon_ini(iml_dyn))
886       
887      IF ( MAXVAL(lon_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
888          lon_ini(:) = lon_dyn(:,1) * 2.0 * ASIN(1.0) / 180.0
889      ELSE
890          lon_ini(:) = lon_dyn(:,1)
891      ENDIF
892
893      ALLOCATE(lat_rad(jml_dyn))
894      ALLOCATE(lat_ini(jml_dyn))
895
896      IF ( MAXVAL(lat_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
897          lat_ini(:) = lat_dyn(1,:) * 2.0 * ASIN(1.0) / 180.0
898      ELSE
899          lat_ini(:) = lat_dyn(1,:)
900      ENDIF
901    !
902
903
904      ALLOCATE(z(iml, jml))
905      ALLOCATE(tmp_var(iml-1,jml))
906    !
907      CALL flinget(fid_dyn, 'Z', iml_dyn, jml_dyn, 0, ttm_dyn,
908     .              1, 1, var_ana)
909c
910      title='Z'
911      CALL conf_dat2d( title,iml_dyn, jml_dyn,lon_ini, lat_ini,
912     . lon_rad, lat_rad, var_ana, interbar  )
913c
914      IF ( interbar )   THEN
915        WRITE(6,*) '-------------------------------------------------',
916     ,'--------------'
917        WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ',
918     , ' pour  Z  $$$ '
919        WRITE(6,*) '-------------------------------------------------',
920     ,'--------------'
921        CALL inter_barxy ( iml_dyn,jml_dyn -1,lon_rad,lat_rad ,
922     ,    var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var)
923      ELSE
924        CALL grille_m(iml_dyn, jml_dyn , lon_rad, lat_rad, var_ana,
925     .               iml-1, jml, lon_in, lat_in, tmp_var)
926      ENDIF
927
928      CALL gr_int_dyn(tmp_var, z, iml-1, jml)
929    !
930      ALLOCATE(psol_dyn(iml, jml))
931    !
932      CALL flinget(fid_dyn, 'SP', iml_dyn, jml_dyn, 0, ttm_dyn,
933     .              1, 1, var_ana)
934
935       title='SP'
936      CALL conf_dat2d( title,iml_dyn, jml_dyn,lon_ini, lat_ini,
937     . lon_rad, lat_rad, var_ana, interbar  )
938
939      IF ( interbar )   THEN
940        WRITE(6,*) '-------------------------------------------------',
941     ,'--------------'
942        WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ',
943     , ' pour  SP  $$$ '
944        WRITE(6,*) '-------------------------------------------------',
945     ,'--------------'
946        CALL inter_barxy ( iml_dyn,jml_dyn -1,lon_rad,lat_rad ,
947     ,    var_ana, iml-1, jml-1, lon_in2, lat_in2, jml, tmp_var)
948      ELSE
949        CALL grille_m(iml_dyn, jml_dyn , lon_rad, lat_rad, var_ana,
950     .             iml-1, jml, lon_in, lat_in, tmp_var  )
951      ENDIF
952
953      CALL gr_int_dyn(tmp_var, psol_dyn, iml-1, jml)
954    !
955      IF ( .NOT.ALLOCATED(tsol)) THEN
956    !   These variables may have been allocated by the need to
957    !   create a start field for them or by the varibale
958    !   coming out of the restart file. In case we dor have it we will initialize it.
959    !
960        CALL start_init_phys( iml, jml, lon_in, lat_in,jml2,lon_in2,
961     .                 lat_in2 , interbar )
962      ELSE
963        IF ( SIZE(tsol) .NE. SIZE(psol_dyn) ) THEN
964        WRITE(*,*) 'start_init_dyn :'
965        WRITE(*,*) 'The temperature field we have does not ',
966     .             'have the right size'
967        STOP
968      ENDIF
969      ENDIF
970      IF ( .NOT.ALLOCATED(phis)) THEN
971            !
972            !    These variables may have been allocated by the need to create a start field for them or by the varibale
973            !     coming out of the restart file. In case we dor have it we will initialize it.
974            !
975        CALL start_init_orog( iml, jml, lon_in, lat_in, jml2, lon_in2 ,
976     .      lat_in2 , interbar )
977            !
978      ELSE
979            !
980          IF (SIZE(phis) .NE. SIZE(psol_dyn)) THEN
981                !
982              WRITE(*,*) 'start_init_dyn :'
983              WRITE(*,*) 'The orography field we have does not ',
984     .                   ' have the right size'
985              STOP
986          ENDIF
987            !
988      ENDIF
989    !
990    !     PSOL is computed in Pascals
991    !
992    !
993      DO j = 1, jml
994        DO i = 1, iml-1
995          psol_dyn(i,j) = psol_dyn(i,j)*(1.0+(z(i,j)-phis(i,j))
996     .                    /287.0/tsol(i,j))
997        ENDDO
998        psol_dyn(iml,j) = psol_dyn(1,j)
999      ENDDO
1000    !
1001    !
1002      ALLOCATE(xppn(iml-1))
1003      ALLOCATE(xpps(iml-1))
1004    !
1005      DO  i   = 1, iml-1
1006        xppn(i) = aire( i,1) * psol_dyn( i,1)
1007        xpps(i) = aire( i,jml) * psol_dyn( i,jml)
1008      ENDDO
1009    !
1010      DO i   = 1, iml
1011        psol_dyn(i,1    )  = SUM(xppn)/apoln
1012        psol_dyn(i,jml)  = SUM(xpps)/apols
1013      ENDDO
1014    !
1015      RETURN
1016    !
1017      END SUBROUTINE start_init_dyn
1018    !
1019    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1020    !
1021      SUBROUTINE start_inter_3d(varname, iml, jml, lml, lon_in,
1022     .      lat_in, jml2, lon_in2, lat_in2, pls_in, var3d, interbar )
1023    !
1024    !    This subroutine gets a variables from a 3D file and does the interpolations needed
1025    !
1026    !
1027    !    ARGUMENTS
1028    !
1029      CHARACTER*(*) :: varname
1030      INTEGER :: iml, jml, lml, jml2
1031      REAL :: lon_in(iml), lat_in(jml), pls_in(iml, jml, lml)
1032      REAL :: lon_in2(iml) , lat_in2(jml2)
1033      REAL :: var3d(iml, jml, lml)
1034      LOGICAL interbar
1035      real chmin,chmax
1036    !
1037    !  LOCAL
1038    !
1039      CHARACTER*25 title
1040      INTEGER :: ii, ij, il, jsort,i,j,l
1041      REAL :: bx, by
1042      REAL, ALLOCATABLE :: lon_rad(:), lat_rad(:)
1043      REAL, ALLOCATABLE :: lon_ini(:), lat_ini(:) , lev_dyn(:)
1044      REAL, ALLOCATABLE :: var_tmp2d(:,:), var_tmp3d(:,:,:)
1045      REAL, ALLOCATABLE :: ax(:), ay(:), yder(:)
1046       REAL, ALLOCATABLE :: varrr(:,:,:)
1047      INTEGER, ALLOCATABLE :: lind(:)
1048    !
1049      LOGICAL :: check = .TRUE.
1050    !
1051      IF ( .NOT. ALLOCATED(var_ana3d)) THEN
1052          ALLOCATE(var_ana3d(iml_dyn, jml_dyn, llm_dyn))
1053      ENDIF
1054          ALLOCATE(varrr(iml_dyn, jml_dyn, llm_dyn))
1055    !
1056    !
1057      IF ( check) WRITE(*,*) 'Going into flinget to extract the 3D ',
1058     .  ' field.', fid_dyn
1059      IF ( check) WRITE(*,*) fid_dyn, varname, iml_dyn, jml_dyn,
1060     .                        llm_dyn,ttm_dyn
1061    !
1062      CALL flinget(fid_dyn, varname, iml_dyn, jml_dyn, llm_dyn,
1063     . ttm_dyn, 1, 1, var_ana3d)
1064    !
1065      IF ( check) WRITE(*,*) 'Allocating space for the interpolation',
1066     . iml, jml, llm_dyn
1067    !
1068      ALLOCATE(lon_rad(iml_dyn))
1069      ALLOCATE(lon_ini(iml_dyn))
1070
1071      IF ( MAXVAL(lon_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
1072          lon_ini(:) = lon_dyn(:,1) * 2.0 * ASIN(1.0) / 180.0
1073      ELSE
1074          lon_ini(:) = lon_dyn(:,1)
1075      ENDIF
1076
1077      ALLOCATE(lat_rad(jml_dyn))
1078      ALLOCATE(lat_ini(jml_dyn))
1079
1080      ALLOCATE(lev_dyn(llm_dyn))
1081
1082      IF ( MAXVAL(lat_dyn(:,:)) .GT. 2.0 * ASIN(1.0) ) THEN
1083          lat_ini(:) = lat_dyn(1,:) * 2.0 * ASIN(1.0) / 180.0
1084      ELSE
1085          lat_ini(:) = lat_dyn(1,:)
1086      ENDIF
1087    !
1088
1089      CALL conf_dat3d ( varname,iml_dyn, jml_dyn, llm_dyn, lon_ini,
1090     . lat_ini, levdyn_ini, lon_rad, lat_rad, lev_dyn, var_ana3d  ,
1091     ,  interbar                                                   )
1092
1093      ALLOCATE(var_tmp2d(iml-1, jml))
1094      ALLOCATE(var_tmp3d(iml, jml, llm_dyn))
1095      ALLOCATE(ax(llm_dyn))
1096      ALLOCATE(ay(llm_dyn))
1097      ALLOCATE(yder(llm_dyn))
1098      ALLOCATE(lind(llm_dyn))
1099    !
1100 
1101      DO il=1,llm_dyn
1102        !
1103      IF( interbar )  THEN
1104       IF( il.EQ.1 )  THEN
1105        WRITE(6,*) '-------------------------------------------------',
1106     ,'--------------'
1107        WRITE(6,*) '$$$ Utilisation de l interpolation barycentrique ',
1108     , ' pour ', varname
1109        WRITE(6,*) '-------------------------------------------------',
1110     ,'--------------'
1111       ENDIF
1112       CALL inter_barxy ( iml_dyn, jml_dyn -1,lon_rad, lat_rad,
1113     , var_ana3d(:,:,il),iml-1, jml2, lon_in2, lat_in2,jml,var_tmp2d )
1114      ELSE
1115       CALL grille_m(iml_dyn, jml_dyn, lon_rad, lat_rad,
1116     .  var_ana3d(:,:,il), iml-1, jml, lon_in, lat_in, var_tmp2d )
1117      ENDIF
1118        !
1119        CALL gr_int_dyn(var_tmp2d, var_tmp3d(:,:,il), iml-1, jml)
1120        !
1121       ENDDO
1122       !
1123          DO il=1,llm_dyn
1124            lind(il) = llm_dyn-il+1
1125          ENDDO
1126    !
1127c
1128c  ... Pour l'interpolation verticale ,on interpole du haut de l'atmosphere
1129c                    vers  le  sol  ...
1130c
1131      DO ij=1,jml
1132        DO ii=1,iml-1
1133          !
1134          ax(:) = lev_dyn(lind(:))
1135          ay(:) = var_tmp3d(ii, ij, lind(:))
1136          !
1137         
1138          CALL SPLINE(ax, ay, llm_dyn, 1.e30, 1.e30, yder)
1139          !
1140          DO il=1,lml
1141            bx = pls_in(ii, ij, il)
1142            CALL SPLINT(ax, ay, yder, llm_dyn, bx, by)
1143            var3d(ii, ij, il) = by
1144          ENDDO
1145          !
1146        ENDDO
1147        var3d(iml, ij, :) = var3d(1, ij, :)
1148      ENDDO
1149
1150      do il=1,lml
1151        call minmax(iml*jml,var3d(1,1,il),chmin,chmax)
1152      SELECTCASE(varname)
1153       CASE('U')
1154          WRITE(*,*) ' U  min max l ',il,chmin,chmax
1155       CASE('V')
1156          WRITE(*,*) ' V  min max l ',il,chmin,chmax
1157       CASE('TEMP')
1158          WRITE(*,*) ' TEMP  min max l ',il,chmin,chmax
1159       CASE('R')
1160          WRITE(*,*) ' R  min max l ',il,chmin,chmax
1161      END SELECT
1162           enddo
1163
1164      DEALLOCATE(lon_rad)
1165      DEALLOCATE(lat_rad)
1166      DEALLOCATE(var_tmp2d)
1167      DEALLOCATE(var_tmp3d)
1168      DEALLOCATE(ax)
1169      DEALLOCATE(ay)
1170      DEALLOCATE(yder)
1171      DEALLOCATE(lind)
1172
1173    !
1174      RETURN
1175    !
1176      END SUBROUTINE start_inter_3d
1177    !
1178      END MODULE startvar
Note: See TracBrowser for help on using the repository browser.