source: LMDZ6/branches/Amaury_dev/libf/phylmdiso/phyaqua_mod.F90 @ 5225

Last change on this file since 5225 was 5158, checked in by abarral, 3 months ago

Add missing klon on strataer_emiss_mod.F90
Correct various missing explicit declarations
Replace tabs by spaces (tabs are not part of the fortran charset)
Continue cleaning modules
Removed unused arguments and variables

File size: 30.7 KB
Line 
1! $Id: phyaqua_mod.F90 3579 2019-10-09 13:11:07Z fairhead $
2
3MODULE phyaqua_mod
4  ! Routines complementaires pour la physique planetaire.
5  USE lmdz_abort_physic, ONLY: abort_physic
6  IMPLICIT NONE
7
8CONTAINS
9
10  SUBROUTINE iniaqua(nlon, year_len, iflag_phys)
11
12    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
13    ! Creation d'un etat initial et de conditions aux limites
14    ! (resp startphy.nc et limit.nc) pour des configurations idealisees
15    ! du modele LMDZ dans sa version terrestre.
16    ! iflag_phys est un parametre qui controle
17    ! iflag_phys = N
18    ! de 100 a 199 : aqua planetes avec SST forcees
19    ! N-100 determine le type de SSTs
20    ! de 200 a 299 : terra planetes avec Ts calcule
21
22    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
23
24    USE dimphy, ONLY: klon
25    USE lmdz_geometry, ONLY: latitude
26    USE surface_data, ONLY: type_ocean, ok_veget
27    USE pbl_surface_mod, ONLY: pbl_surface_init
28    USE fonte_neige_mod, ONLY: fonte_neige_init
29#ifdef ISO
30  USE fonte_neige_mod,  ONLY: fonte_neige_init_iso
31  USE pbl_surface_mod,  ONLY: pbl_surface_init_iso
32#endif
33    USE phys_state_var_mod
34    USE time_phylmdz_mod, ONLY: day_ref, ndays, pdtphys, &
35            day_ini, day_end
36    USE indice_sol_mod
37    USE lmdz_physical_constants, ONLY: pi
38    !    USE ioipsl
39    USE lmdz_phys_para, ONLY: is_master
40    USE lmdz_phys_transfert_para, ONLY: bcast
41    USE lmdz_grid_phy
42    USE lmdz_ioipsl_getin_p, ONLY: getin_p
43    USE phys_cal_mod, ONLY: calend, year_len_phy => year_len
44    USE lmdz_clesphys
45#ifdef ISO
46  USE infotrac_phy, ONLY: niso
47#endif
48    USE lmdz_yomcst
49
50    IMPLICIT NONE
51
52    include "dimsoil.h"
53
54    INTEGER, INTENT (IN) :: nlon, year_len, iflag_phys
55    ! IM ajout latfi, lonfi
56    !    REAL, INTENT (IN) :: lonfi(nlon), latfi(nlon)
57
58    INTEGER type_profil, type_aqua
59
60    ! Ajouts initialisation des surfaces
61    REAL :: run_off_lic_0(nlon)
62    REAL :: qsolsrf(nlon, nbsrf), snsrf(nlon, nbsrf)
63    REAL :: tsoil(nlon, nsoilmx, nbsrf)
64    REAL :: tslab(nlon), seaice(nlon)
65    REAL fder(nlon)
66
67#ifdef ISO
68    REAL :: xtrun_off_lic_0(niso,nlon)
69    REAL :: xtsolsrf(niso,nlon, nbsrf), xtsnsrf(niso,nlon, nbsrf)
70    REAL :: Rland_ice(niso,nlon)
71#endif
72
73
74
75    ! Arguments :
76    ! -----------
77
78    ! integer radpas
79    INTEGER it, unit, i, k, itap
80
81    REAL rugos, albedo
82    REAL tsurf
83    REAL time, timestep, day, day0
84    REAL qsol_f
85    REAL rugsrel(nlon)
86    LOGICAL alb_ocean
87
88    CHARACTER *80 ans, file_forctl, file_fordat, file_start
89    CHARACTER *100 file, var
90    CHARACTER *2 cnbl
91
92    REAL phy_nat(nlon, year_len)
93    REAL phy_alb(nlon, year_len)
94    REAL phy_sst(nlon, year_len)
95    REAL phy_bil(nlon, year_len)
96    REAL phy_rug(nlon, year_len)
97    REAL phy_ice(nlon, year_len)
98    REAL phy_fter(nlon, year_len)
99    REAL phy_foce(nlon, year_len)
100    REAL phy_fsic(nlon, year_len)
101    REAL phy_flic(nlon, year_len)
102
103    INTEGER, SAVE :: read_climoz = 0 ! read ozone climatology
104    !$OMP THREADPRIVATE(read_climoz)
105
106    ! -------------------------------------------------------------------------
107    ! declaration pour l'appel a phyredem
108    ! -------------------------------------------------------------------------
109
110    ! real pctsrf(nlon,nbsrf),ftsol(nlon,nbsrf)
111    REAL falbe(nlon, nbsrf), falblw(nlon, nbsrf)
112    ! real pbl_tke(nlon,llm,nbsrf)
113    ! real rain_fall(nlon),snow_fall(nlon)
114    ! real solsw(nlon), sollw(nlon),radsol(nlon)
115    ! real t_ancien(nlon,llm),q_ancien(nlon,llm),rnebcon(nlon,llm)
116    ! real ratqs(nlon,llm)
117    ! real clwcon(nlon,llm)
118
119    INTEGER longcles
120    PARAMETER (longcles = 20)
121    REAL clesphy0(longcles)
122
123
124    ! -----------------------------------------------------------------------
125    ! dynamial tendencies :
126    ! ---------------------
127
128    INTEGER l, ierr, aslun
129
130    REAL paire
131
132    ! Local
133    CHARACTER (LEN = 20) :: modname = 'phyaqua'
134    CHARACTER (LEN = 80) :: abort_message
135
136
137    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
138    ! INITIALISATIONS
139    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
140
141    ! -----------------------------------------------------------------------
142    ! Initialisations  des constantes
143    ! -------------------------------
144
145    !IF (calend .EQ. "earth_360d") Then
146    year_len_phy = year_len
147    !END IF
148
149    IF (year_len/=360) THEN
150      write (*, *) year_len
151      write (*, *) 'iniaqua: 360 day calendar is required !'
152      stop
153    endif
154
155    type_aqua = iflag_phys / 100
156    type_profil = iflag_phys - type_aqua * 100
157    PRINT *, 'iniaqua:type_aqua, type_profil', type_aqua, type_profil
158
159    IF (klon/=nlon) THEN
160      WRITE (*, *) 'iniaqua: klon=', klon, ' nlon=', nlon
161      abort_message = 'probleme de dimensions dans iniaqua'
162      CALL abort_physic(modname, abort_message, 1)
163    END IF
164    CALL phys_state_var_init(read_climoz)
165
166    read_climoz = 0
167    day0 = 217.
168    day = day0
169    it = 0
170    time = 0.
171
172    ! -----------------------------------------------------------------------
173    ! initialisations de la physique
174    ! -----------------------------------------------------------------------
175
176    day_ini = day_ref
177    day_end = day_ini + ndays
178
179    nbapp_rad = 24
180    CALL getin_p('nbapp_rad', nbapp_rad)
181
182    ! ---------------------------------------------------------------------
183    ! Creation des conditions aux limites:
184    ! ------------------------------------
185    ! Initialisations des constantes
186    ! Ajouter les manquants dans planete.def... (albedo etc)
187    co2_ppm = 348.
188    CALL getin_p('co2_ppm', co2_ppm)
189
190    solaire = 1365.
191    CALL getin_p('solaire', solaire)
192
193    ! CALL getin('albedo',albedo) ! albedo is set below, depending on
194    ! type_aqua
195    alb_ocean = .TRUE.
196    CALL getin_p('alb_ocean', alb_ocean)
197
198    WRITE (*, *) 'iniaqua: co2_ppm=', co2_ppm
199    WRITE (*, *) 'iniaqua: solaire=', solaire
200    WRITE (*, *) 'iniaqua: alb_ocean=', alb_ocean
201
202    radsol = 0.
203    qsol_f = 10.
204
205    ! Conditions aux limites:
206    ! -----------------------
207
208    qsol(:) = qsol_f
209    rugsrel = 0.0 ! (rugsrel = rugoro)
210    rugoro = 0.0
211    u_ancien = 0.0
212    v_ancien = 0.0
213    agesno = 50.0
214    ! Relief plat
215    zmea = 0.
216    zstd = 0.
217    zsig = 0.
218    zgam = 0.
219    zthe = 0.
220    zpic = 0.
221    zval = 0.
222
223    ! Une seule surface
224    pctsrf = 0.
225    IF (type_aqua==1) THEN
226      rugos = 1.E-4
227      albedo = 0.19
228      pctsrf(:, is_oce) = 1.
229    ELSE IF (type_aqua==2) THEN
230      rugos = 0.03
231      albedo = 0.1
232      pctsrf(:, is_ter) = 1.
233    END IF
234
235    CALL getin_p('rugos', rugos)
236
237    WRITE (*, *) 'iniaqua: rugos=', rugos
238    zmasq(:) = pctsrf(:, is_ter)
239
240    ! pctsrf_pot(:,is_oce) = 1. - zmasq(:)
241    ! pctsrf_pot(:,is_sic) = 1. - zmasq(:)
242
243    ! Si alb_ocean on calcule un albedo oceanique moyen
244    ! if (alb_ocean) THEN
245    ! Voir pourquoi on avait ca.
246    ! CALL ini_alb_oce(phy_alb)
247    ! else
248    phy_alb(:, :) = albedo ! albedo land only (old value condsurf_jyg=0.3)
249    ! END IF !alb_ocean
250
251    DO i = 1, year_len
252      ! IM Terraplanete   phy_sst(:,i) = 260.+50.*cos(rlatd(:))**2
253      ! IM ajout calcul profil sst selon le cas considere (cf. FBr)
254
255      phy_nat(:, i) = 1.0 ! 0=ocean libre, 1=land, 2=glacier, 3=banquise
256      phy_bil(:, i) = 1.0 ! ne sert que pour les slab_ocean
257      phy_rug(:, i) = rugos ! longueur rugosite utilisee sur land only
258      phy_ice(:, i) = 0.0 ! fraction de glace (?)
259      phy_fter(:, i) = pctsrf(:, is_ter) ! fraction de glace (?)
260      phy_foce(:, i) = pctsrf(:, is_oce) ! fraction de glace (?)
261      phy_fsic(:, i) = pctsrf(:, is_sic) ! fraction de glace (?)
262      phy_flic(:, i) = pctsrf(:, is_lic) ! fraction de glace (?)
263    END DO
264    ! IM calcul profil sst
265    CALL profil_sst(nlon, latitude, type_profil, phy_sst)
266
267    IF (grid_type==unstructured) THEN
268      CALL writelim_unstruct(klon, phy_nat, phy_alb, phy_sst, phy_bil, phy_rug, phy_ice, &
269              phy_fter, phy_foce, phy_flic, phy_fsic)
270    ELSE
271
272      CALL writelim(klon, phy_nat, phy_alb, phy_sst, phy_bil, phy_rug, phy_ice, &
273              phy_fter, phy_foce, phy_flic, phy_fsic)
274    ENDIF
275
276    ! ---------------------------------------------------------------------
277    ! Ecriture de l'etat initial:
278    ! ---------------------------
279
280
281    ! Ecriture etat initial physique
282
283    timestep = pdtphys
284    radpas = nint(rday / timestep / float(nbapp_rad))
285
286    DO i = 1, longcles
287      clesphy0(i) = 0.
288    END DO
289    clesphy0(1) = float(iflag_con)
290    clesphy0(2) = float(nbapp_rad)
291    ! IF( cycle_diurne  ) clesphy0(3) =  1.
292    clesphy0(3) = 1. ! cycle_diurne
293    clesphy0(4) = 1. ! soil_model
294    clesphy0(5) = 1. ! new_oliq
295    clesphy0(6) = 0. ! ok_orodr
296    clesphy0(7) = 0. ! ok_orolf
297    clesphy0(8) = 0. ! ok_limitvrai
298
299
300    ! =======================================================================
301    ! Profils initiaux
302    ! =======================================================================
303
304    ! On initialise les temperatures de surfaces comme les sst
305    DO i = 1, nlon
306      ftsol(i, :) = phy_sst(i, 1)
307      tsoil(i, :, :) = phy_sst(i, 1)
308      tslab(i) = phy_sst(i, 1)
309    END DO
310
311    falbe(:, :) = albedo
312    falblw(:, :) = albedo
313    rain_fall(:) = 0.
314    snow_fall(:) = 0.
315    solsw(:) = 0.
316    sollw(:) = 0.
317    radsol(:) = 0.
318
319    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
320    ! intialisation bidon mais pas grave
321    t_ancien(:, :) = 0.
322    q_ancien(:, :) = 0.
323    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
324    rnebcon = 0.
325    ratqs = 0.
326    clwcon = 0.
327    pbl_tke = 1.E-8
328
329    ! variables supplementaires pour appel a plb_surface_init
330    fder(:) = 0.
331    seaice(:) = 0.
332    run_off_lic_0 = 0.
333    fevap = 0.
334#ifdef ISO
335    xtrun_off_lic_0 = 0.
336#endif
337
338
339    ! Initialisations necessaires avant phyredem
340    type_ocean = 'force'
341    CALL fonte_neige_init(run_off_lic_0)
342#ifdef ISO
343   CALL fonte_neige_init_iso(xtrun_off_lic_0)
344#endif
345    qsolsrf(:, :) = qsol(1) ! humidite du sol des sous surface
346    snsrf(:, :) = 0. ! couverture de neige des sous surface
347    z0m(:, :) = rugos ! couverture de neige des sous surface
348    z0h = z0m
349
350    CALL pbl_surface_init(fder, snsrf, qsolsrf, tsoil)
351#ifdef ISO
352  CALL pbl_surface_init_iso(xtsnsrf,Rland_ice)
353#endif
354
355    PRINT *, 'iniaqua: before phyredem'
356
357    pbl_tke(:, :, :) = 1.e-8
358    falb1 = albedo
359    falb2 = albedo
360    zmax0 = 0.
361    f0 = 0.
362    sig1 = 0.
363    w01 = 0.
364    wake_deltat = 0.
365    wake_deltaq = 0.
366    wake_s = 0.
367    wake_dens = 0.
368    wake_cstar = 0.
369    wake_pe = 0.
370    wake_fip = 0.
371    fm_therm = 0.
372    entr_therm = 0.
373    detr_therm = 0.
374    ale_bl = 0.
375    ale_bl_trig = 0.
376    alp_bl = 0.
377    treedrg(:, :, :) = 0.
378
379    u10m = 0.
380    v10m = 0.
381
382    ql_ancien = 0.
383    qs_ancien = 0.
384    u_ancien = 0.
385    v_ancien = 0.
386    prw_ancien = 0.
387    prlw_ancien = 0.
388    prsw_ancien = 0.
389
390    ale_wake = 0.
391    ale_bl_stat = 0.
392
393
394    !ym error : the sub surface dimension is the third not second : forgotten for iniaqua
395    !    falb_dir(:,is_ter,:)=0.08; falb_dir(:,is_lic,:)=0.6
396    !    falb_dir(:,is_oce,:)=0.5;  falb_dir(:,is_sic,:)=0.6
397    falb_dir(:, :, is_ter) = 0.08; falb_dir(:, :, is_lic) = 0.6
398    falb_dir(:, :, is_oce) = 0.5;  falb_dir(:, :, is_sic) = 0.6
399
400    !ym falb_dif has been forgotten, initialize with defaukt value found in phyetat0 or 0 ?
401    !ym probably the uninitialized value was 0 for standard (regular grid) case
402    falb_dif(:, :, :) = 0
403
404    CALL phyredem('startphy.nc')
405
406    PRINT *, 'iniaqua: after phyredem'
407    CALL phys_state_var_end
408
409  END SUBROUTINE iniaqua
410
411
412  ! ====================================================================
413  ! ====================================================================
414  SUBROUTINE zenang_an(cycle_diurne, gmtime, rlat, rlon, rmu0, fract)
415    USE dimphy
416    USE lmdz_yomcst
417
418    IMPLICIT NONE
419    ! ====================================================================
420    ! =============================================================
421    ! CALL zenang(cycle_diurne,gmtime,rlat,rlon,rmu0,fract)
422    ! Auteur : A. Campoy et F. Hourdin
423    ! Objet  : calculer les valeurs moyennes du cos de l'angle zenithal
424    ! et l'ensoleillement moyen entre gmtime1 et gmtime2
425    ! connaissant la declinaison, la latitude et la longitude.
426
427    ! Dans cette version particuliere, on calcule le rayonnement
428    ! moyen sur l'année à chaque latitude.
429    ! angle zenithal calculé pour obtenir un
430    ! Fit polynomial de  l'ensoleillement moyen au sommet de l'atmosphere
431    ! en moyenne annuelle.
432    ! Spécifique de la terre. Utilisé pour les aqua planetes.
433
434    ! Rque   : Different de la routine angle en ce sens que zenang
435    ! fournit des moyennes de pmu0 et non des valeurs
436    ! instantanees, du coup frac prend toutes les valeurs
437    ! entre 0 et 1.
438    ! Date   : premiere version le 13 decembre 1994
439    ! revu pour  GCM  le 30 septembre 1996
440    ! ===============================================================
441    ! longi----INPUT : la longitude vraie de la terre dans son plan
442    ! solaire a partir de l'equinoxe de printemps (degre)
443    ! gmtime---INPUT : temps universel en fraction de jour
444    ! pdtrad---INPUT : pas de temps du rayonnement (secondes)
445    ! lat------INPUT : latitude en degres
446    ! long-----INPUT : longitude en degres
447    ! pmu0-----OUTPUT: angle zenithal moyen entre gmtime et gmtime+pdtrad
448    ! frac-----OUTPUT: ensoleillement moyen entre gmtime et gmtime+pdtrad
449
450    LOGICAL cycle_diurne
451    REAL gmtime
452    REAL rlat(klon), rlon(klon), rmu0(klon), fract(klon)
453    ! ================================================================
454    INTEGER i
455    REAL gmtime1, gmtime2
456    REAL pi_local
457
458    REAL rmu0m(klon), rmu0a(klon)
459
460    pi_local = 4.0 * atan(1.0)
461
462    ! ================================================================
463    ! Calcul de l'angle zenithal moyen sur la journee
464    ! ================================================================
465
466    DO i = 1, klon
467      fract(i) = 1.
468      ! Calcule du flux moyen
469      IF (abs(rlat(i))<=28.75) THEN
470        rmu0m(i) = (210.1924 + 206.6059 * cos(0.0174533 * rlat(i))**2) / 1365.
471      ELSE IF (abs(rlat(i))<=43.75) THEN
472        rmu0m(i) = (187.4562 + 236.1853 * cos(0.0174533 * rlat(i))**2) / 1365.
473      ELSE IF (abs(rlat(i))<=71.25) THEN
474        rmu0m(i) = (162.4439 + 284.1192 * cos(0.0174533 * rlat(i))**2) / 1365.
475      ELSE
476        rmu0m(i) = (172.8125 + 183.7673 * cos(0.0174533 * rlat(i))**2) / 1365.
477      END IF
478    END DO
479
480    ! ================================================================
481    ! Avec ou sans cycle diurne
482    ! ================================================================
483
484    IF (cycle_diurne) THEN
485
486      ! On redecompose flux  au sommet suivant un cycle diurne idealise
487      ! identique a toutes les latitudes.
488
489      DO i = 1, klon
490        rmu0a(i) = 2. * rmu0m(i) * sqrt(2.) * pi_local / (4. - pi_local)
491        rmu0(i) = rmu0a(i) * abs(sin(pi_local * gmtime + pi_local * rlon(i) / 360.)) - &
492                rmu0a(i) / sqrt(2.)
493      END DO
494
495      DO i = 1, klon
496        IF (rmu0(i)<=0.) THEN
497          rmu0(i) = 0.
498          fract(i) = 0.
499        ELSE
500          fract(i) = 1.
501        END IF
502      END DO
503
504      ! Affichage de l'angel zenitale
505      ! PRINT*,'************************************'
506      ! PRINT*,'************************************'
507      ! PRINT*,'************************************'
508      ! PRINT*,'latitude=',rlat(i),'longitude=',rlon(i)
509      ! PRINT*,'rmu0m=',rmu0m(i)
510      ! PRINT*,'rmu0a=',rmu0a(i)
511      ! PRINT*,'rmu0=',rmu0(i)
512
513    ELSE
514
515      DO i = 1, klon
516        fract(i) = 0.5
517        rmu0(i) = rmu0m(i) * 2.
518      END DO
519
520    END IF
521
522  END SUBROUTINE zenang_an
523
524  ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
525
526  SUBROUTINE writelim_unstruct(klon, phy_nat, phy_alb, phy_sst, phy_bil, phy_rug, &
527          phy_ice, phy_fter, phy_foce, phy_flic, phy_fsic)
528
529    USE lmdz_phys_para, ONLY: is_omp_master, klon_mpi
530    USE lmdz_phys_transfert_para, ONLY: gather_omp
531    USE lmdz_xios
532    IMPLICIT NONE
533
534    INTEGER, INTENT (IN) :: klon
535    REAL, INTENT (IN) :: phy_nat(klon, 360)
536    REAL, INTENT (IN) :: phy_alb(klon, 360)
537    REAL, INTENT (IN) :: phy_sst(klon, 360)
538    REAL, INTENT (IN) :: phy_bil(klon, 360)
539    REAL, INTENT (IN) :: phy_rug(klon, 360)
540    REAL, INTENT (IN) :: phy_ice(klon, 360)
541    REAL, INTENT (IN) :: phy_fter(klon, 360)
542    REAL, INTENT (IN) :: phy_foce(klon, 360)
543    REAL, INTENT (IN) :: phy_flic(klon, 360)
544    REAL, INTENT (IN) :: phy_fsic(klon, 360)
545
546    REAL :: phy_mpi(klon_mpi, 360) ! temporary variable, to store phy_***(:)
547    ! on the whole physics grid
548
549    PRINT *, 'writelim: Ecriture du fichier limit'
550
551    CALL gather_omp(phy_foce, phy_mpi)
552    IF (is_omp_master) CALL xios_send_field('foce_limout', phy_mpi)
553
554    CALL gather_omp(phy_fsic, phy_mpi)
555    IF (is_omp_master) CALL xios_send_field('fsic_limout', phy_mpi)
556
557    CALL gather_omp(phy_fter, phy_mpi)
558    IF (is_omp_master) CALL xios_send_field('fter_limout', phy_mpi)
559
560    CALL gather_omp(phy_flic, phy_mpi)
561    IF (is_omp_master) CALL xios_send_field('flic_limout', phy_mpi)
562
563    CALL gather_omp(phy_sst, phy_mpi)
564    IF (is_omp_master) CALL xios_send_field('sst_limout', phy_mpi)
565
566    CALL gather_omp(phy_bil, phy_mpi)
567    IF (is_omp_master) CALL xios_send_field('bils_limout', phy_mpi)
568
569    CALL gather_omp(phy_alb, phy_mpi)
570    IF (is_omp_master) CALL xios_send_field('alb_limout', phy_mpi)
571
572    CALL gather_omp(phy_rug, phy_mpi)
573    IF (is_omp_master) CALL xios_send_field('rug_limout', phy_mpi)
574
575  END SUBROUTINE writelim_unstruct
576
577
578  SUBROUTINE writelim(klon, phy_nat, phy_alb, phy_sst, phy_bil, phy_rug, &
579          phy_ice, phy_fter, phy_foce, phy_flic, phy_fsic)
580
581    USE lmdz_phys_para, ONLY: is_master
582    USE lmdz_grid_phy, ONLY: klon_glo
583    USE lmdz_phys_transfert_para, ONLY: gather
584    USE phys_cal_mod, ONLY: year_len
585    USE netcdf, ONLY: nf90_clobber, nf90_close, nf90_noerr, nf90_strerror, nf90_put_att, nf90_def_var, &
586            nf90_def_dim, nf90_create, nf90_put_var, nf90_unlimited, nf90_global, nf90_64bit_offset, &
587            nf90_enddef
588    USE lmdz_cppkeys_wrapper, ONLY: nf90_format
589    IMPLICIT NONE
590
591    INTEGER, INTENT (IN) :: klon
592    REAL, INTENT (IN) :: phy_nat(klon, year_len)
593    REAL, INTENT (IN) :: phy_alb(klon, year_len)
594    REAL, INTENT (IN) :: phy_sst(klon, year_len)
595    REAL, INTENT (IN) :: phy_bil(klon, year_len)
596    REAL, INTENT (IN) :: phy_rug(klon, year_len)
597    REAL, INTENT (IN) :: phy_ice(klon, year_len)
598    REAL, INTENT (IN) :: phy_fter(klon, year_len)
599    REAL, INTENT (IN) :: phy_foce(klon, year_len)
600    REAL, INTENT (IN) :: phy_flic(klon, year_len)
601    REAL, INTENT (IN) :: phy_fsic(klon, year_len)
602
603    REAL :: phy_glo(klon_glo, year_len) ! temporary variable, to store phy_***(:)
604    ! on the whole physics grid
605    INTEGER :: k
606    INTEGER ierr
607    INTEGER dimfirst(3)
608    INTEGER dimlast(3)
609
610    INTEGER nid, ndim, ntim
611    INTEGER dims(2), debut(2), epais(2)
612    INTEGER id_tim
613    INTEGER id_nat, id_sst, id_bils, id_rug, id_alb
614    INTEGER id_fter, id_foce, id_fsic, id_flic
615
616    IF (is_master) THEN
617
618      PRINT *, 'writelim: Ecriture du fichier limit'
619
620      ierr = nf90_create('limit.nc', IOR(nf90_clobber, nf90_64bit_offset), nid)
621
622      ierr = nf90_put_att(nid, nf90_global, 'title', &
623              'Fichier conditions aux limites')
624      !        ierr = nf90_def_dim (nid, "points_physiques", klon, ndim)
625      ierr = nf90_def_dim(nid, 'points_physiques', klon_glo, ndim)
626      ierr = nf90_def_dim(nid, 'time', nf90_unlimited, ntim)
627
628      dims(1) = ndim
629      dims(2) = ntim
630
631      ierr = nf90_def_var(nid, 'TEMPS', nf90_format, ntim, id_tim)
632      ierr = nf90_put_att(nid, id_tim, 'title', 'Jour dans l annee')
633
634      ierr = nf90_def_var(nid, 'NAT', nf90_format, dims, id_nat)
635      ierr = nf90_put_att(nid, id_nat, 'title', &
636              'Nature du sol (0,1,2,3)')
637
638      ierr = nf90_def_var(nid, 'SST', nf90_format, dims, id_sst)
639      ierr = nf90_put_att(nid, id_sst, 'title', &
640              'Temperature superficielle de la mer')
641
642      ierr = nf90_def_var(nid, 'BILS', nf90_format, dims, id_bils)
643      ierr = nf90_put_att(nid, id_bils, 'title', &
644              'Reference flux de chaleur au sol')
645
646      ierr = nf90_def_var(nid, 'ALB', nf90_format, dims, id_alb)
647      ierr = nf90_put_att(nid, id_alb, 'title', 'Albedo a la surface')
648
649      ierr = nf90_def_var(nid, 'RUG', nf90_format, dims, id_rug)
650      ierr = nf90_put_att(nid, id_rug, 'title', 'Rugosite')
651
652      ierr = nf90_def_var(nid, 'FTER', nf90_format, dims, id_fter)
653      ierr = nf90_put_att(nid, id_fter, 'title', 'Frac. Land')
654      ierr = nf90_def_var(nid, 'FOCE', nf90_format, dims, id_foce)
655      ierr = nf90_put_att(nid, id_foce, 'title', 'Frac. Ocean')
656      ierr = nf90_def_var(nid, 'FSIC', nf90_format, dims, id_fsic)
657      ierr = nf90_put_att(nid, id_fsic, 'title', 'Frac. Sea Ice')
658      ierr = nf90_def_var(nid, 'FLIC', nf90_format, dims, id_flic)
659      ierr = nf90_put_att(nid, id_flic, 'title', 'Frac. Land Ice')
660
661      ierr = nf90_enddef(nid)
662      IF (ierr/=nf90_noerr) THEN
663        WRITE (*, *) 'writelim error: failed to end define mode'
664        WRITE (*, *) nf90_strerror(ierr)
665      END IF
666
667
668      ! write the 'times'
669      DO k = 1, year_len
670        ierr = nf90_put_var(nid, id_tim, k, [k])
671        IF (ierr/=nf90_noerr) THEN
672          WRITE (*, *) 'writelim error with temps(k),k=', k
673          WRITE (*, *) nf90_strerror(ierr)
674        END IF
675      END DO
676
677    END IF ! of if (is_master)
678
679    ! write the fields, after having collected them on master
680
681    CALL gather(phy_nat, phy_glo)
682    IF (is_master) THEN
683      ierr = nf90_put_var(nid, id_nat, phy_glo)
684      IF (ierr/=nf90_noerr) THEN
685        WRITE (*, *) 'writelim error with phy_nat'
686        WRITE (*, *) nf90_strerror(ierr)
687      END IF
688    END IF
689
690    CALL gather(phy_sst, phy_glo)
691    IF (is_master) THEN
692      ierr = nf90_put_var(nid, id_sst, phy_glo)
693      IF (ierr/=nf90_noerr) THEN
694        WRITE (*, *) 'writelim error with phy_sst'
695        WRITE (*, *) nf90_strerror(ierr)
696      END IF
697    END IF
698
699    CALL gather(phy_bil, phy_glo)
700    IF (is_master) THEN
701      ierr = nf90_put_var(nid, id_bils, phy_glo)
702      IF (ierr/=nf90_noerr) THEN
703        WRITE (*, *) 'writelim error with phy_bil'
704        WRITE (*, *) nf90_strerror(ierr)
705      END IF
706    END IF
707
708    CALL gather(phy_alb, phy_glo)
709    IF (is_master) THEN
710      ierr = nf90_put_var(nid, id_alb, phy_glo)
711      IF (ierr/=nf90_noerr) THEN
712        WRITE (*, *) 'writelim error with phy_alb'
713        WRITE (*, *) nf90_strerror(ierr)
714      END IF
715    END IF
716
717    CALL gather(phy_rug, phy_glo)
718    IF (is_master) THEN
719      ierr = nf90_put_var(nid, id_rug, phy_glo)
720      IF (ierr/=nf90_noerr) THEN
721        WRITE (*, *) 'writelim error with phy_rug'
722        WRITE (*, *) nf90_strerror(ierr)
723      END IF
724    END IF
725
726    CALL gather(phy_fter, phy_glo)
727    IF (is_master) THEN
728      ierr = nf90_put_var(nid, id_fter, phy_glo)
729      IF (ierr/=nf90_noerr) THEN
730        WRITE (*, *) 'writelim error with phy_fter'
731        WRITE (*, *) nf90_strerror(ierr)
732      END IF
733    END IF
734
735    CALL gather(phy_foce, phy_glo)
736    IF (is_master) THEN
737      ierr = nf90_put_var(nid, id_foce, phy_glo)
738      IF (ierr/=nf90_noerr) THEN
739        WRITE (*, *) 'writelim error with phy_foce'
740        WRITE (*, *) nf90_strerror(ierr)
741      END IF
742    END IF
743
744    CALL gather(phy_fsic, phy_glo)
745    IF (is_master) THEN
746      ierr = nf90_put_var(nid, id_fsic, phy_glo)
747      IF (ierr/=nf90_noerr) THEN
748        WRITE (*, *) 'writelim error with phy_fsic'
749        WRITE (*, *) nf90_strerror(ierr)
750      END IF
751    END IF
752
753    CALL gather(phy_flic, phy_glo)
754    IF (is_master) THEN
755      ierr = nf90_put_var(nid, id_flic, phy_glo)
756      IF (ierr/=nf90_noerr) THEN
757        WRITE (*, *) 'writelim error with phy_flic'
758        WRITE (*, *) nf90_strerror(ierr)
759      END IF
760    END IF
761
762    ! close file:
763    IF (is_master) THEN
764      ierr = nf90_close(nid)
765    END IF
766
767  END SUBROUTINE writelim
768
769  ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
770
771  SUBROUTINE profil_sst(nlon, rlatd, type_profil, phy_sst)
772    USE dimphy
773    USE phys_cal_mod, ONLY: year_len
774    IMPLICIT NONE
775
776    INTEGER nlon, type_profil, i, k, j
777    REAL :: rlatd(nlon), phy_sst(nlon, year_len)
778    INTEGER imn, imx, amn, amx, kmn, kmx
779    INTEGER p, pplus, nlat_max
780    PARAMETER (nlat_max = 72)
781    REAL x_anom_sst(nlat_max)
782    CHARACTER (LEN = 20) :: modname = 'profil_sst'
783    CHARACTER (LEN = 80) :: abort_message
784
785    IF (klon/=nlon) THEN
786      abort_message = 'probleme de dimensions dans profil_sst'
787      CALL abort_physic(modname, abort_message, 1)
788    ENDIF
789    WRITE (*, *) ' profil_sst: type_profil=', type_profil
790    DO i = 1, year_len
791      ! phy_sst(:,i) = 260.+50.*cos(rlatd(:))**2
792
793      ! Rajout fbrlmd
794
795      IF (type_profil==1) THEN
796        ! Méthode 1 "Control" faible plateau à l'Equateur
797        DO j = 1, klon
798          phy_sst(j, i) = 273. + 27. * (1 - sin(1.5 * rlatd(j))**2)
799          ! PI/3=1.047197551
800          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
801            phy_sst(j, i) = 273.
802          END IF
803        END DO
804      END IF
805      IF (type_profil==2) THEN
806        ! Méthode 2 "Flat" fort plateau à l'Equateur
807        DO j = 1, klon
808          phy_sst(j, i) = 273. + 27. * (1 - sin(1.5 * rlatd(j))**4)
809          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
810            phy_sst(j, i) = 273.
811          END IF
812        END DO
813      END IF
814
815      IF (type_profil==3) THEN
816        ! Méthode 3 "Qobs" plateau réel à l'Equateur
817        DO j = 1, klon
818          phy_sst(j, i) = 273. + 0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 * &
819                  rlatd(j))**4)
820          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
821            phy_sst(j, i) = 273.
822          END IF
823        END DO
824      END IF
825
826      IF (type_profil==4) THEN
827        ! Méthode 4 : Méthode 3 + SST+2 "Qobs" plateau réel à l'Equateur
828        DO j = 1, klon
829          phy_sst(j, i) = 273. + 0.5 * 29. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 * &
830                  rlatd(j))**4)
831          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
832            phy_sst(j, i) = 273.
833          END IF
834        END DO
835      END IF
836
837      IF (type_profil==5) THEN
838        ! Méthode 5 : Méthode 3 + +2K "Qobs" plateau réel à l'Equateur
839        DO j = 1, klon
840          phy_sst(j, i) = 273. + 2. + 0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 &
841                  * rlatd(j))**4)
842          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
843            phy_sst(j, i) = 273. + 2.
844          END IF
845
846        END DO
847      END IF
848
849      IF (type_profil==6) THEN
850        ! Méthode 6 "cst" valeur constante de SST
851        DO j = 1, klon
852          phy_sst(j, i) = 288.
853        END DO
854      END IF
855
856      IF (type_profil==7) THEN
857        ! Méthode 7 "cst" valeur constante de SST +2
858        DO j = 1, klon
859          phy_sst(j, i) = 288. + 2.
860        END DO
861      END IF
862
863      p = 0
864      IF (type_profil==8) THEN
865        ! Méthode 8 profil anomalies SST du modèle couplé AR4
866        DO j = 1, klon
867          IF (rlatd(j)==rlatd(j - 1)) THEN
868            phy_sst(j, i) = 273. + x_anom_sst(pplus) + &
869                    0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 * rlatd(j))**4)
870            IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
871              phy_sst(j, i) = 273. + x_anom_sst(pplus)
872            END IF
873          ELSE
874            p = p + 1
875            pplus = 73 - p
876            phy_sst(j, i) = 273. + x_anom_sst(pplus) + &
877                    0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 * rlatd(j))**4)
878            IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
879              phy_sst(j, i) = 273. + x_anom_sst(pplus)
880            END IF
881            WRITE (*, *) rlatd(j), x_anom_sst(pplus), phy_sst(j, i)
882          END IF
883        END DO
884      END IF
885
886      IF (type_profil==9) THEN
887        ! Méthode 5 : Méthode 3 + -2K "Qobs" plateau réel à l'Equateur
888        DO j = 1, klon
889          phy_sst(j, i) = 273. - 2. + 0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 &
890                  * rlatd(j))**4)
891          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
892            phy_sst(j, i) = 273. - 2.
893          END IF
894        END DO
895      END IF
896
897      IF (type_profil==10) THEN
898        ! Méthode 10 : Méthode 3 + +4K "Qobs" plateau réel à l'Equateur
899        DO j = 1, klon
900          phy_sst(j, i) = 273. + 4. + 0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 &
901                  * rlatd(j))**4)
902          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
903            phy_sst(j, i) = 273. + 4.
904          END IF
905        END DO
906      END IF
907
908      IF (type_profil==11) THEN
909        ! Méthode 11 : Méthode 3 + 4CO2 "Qobs" plateau réel à l'Equateur
910        DO j = 1, klon
911          phy_sst(j, i) = 273. + 0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 * &
912                  rlatd(j))**4)
913          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
914            phy_sst(j, i) = 273.
915          END IF
916        END DO
917      END IF
918
919      IF (type_profil==12) THEN
920        ! Méthode 12 : Méthode 10 + 4CO2 "Qobs" plateau réel à l'Equateur
921        DO j = 1, klon
922          phy_sst(j, i) = 273. + 4. + 0.5 * 27. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 &
923                  * rlatd(j))**4)
924          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
925            phy_sst(j, i) = 273. + 4.
926          END IF
927        END DO
928      END IF
929
930      IF (type_profil==13) THEN
931        ! Méthode 13 "Qmax" plateau réel à l'Equateur augmenté !
932        DO j = 1, klon
933          phy_sst(j, i) = 273. + 0.5 * 29. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 * &
934                  rlatd(j))**4)
935          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
936            phy_sst(j, i) = 273.
937          END IF
938        END DO
939      END IF
940
941      IF (type_profil==14) THEN
942        ! Méthode 13 "Qmax2K" plateau réel à l'Equateur augmenté +2K !
943        DO j = 1, klon
944          phy_sst(j, i) = 273. + 2. + 0.5 * 29. * (2 - sin(1.5 * rlatd(j))**2 - sin(1.5 &
945                  * rlatd(j))**4)
946          IF ((rlatd(j)>1.0471975) .OR. (rlatd(j)<-1.0471975)) THEN
947            phy_sst(j, i) = 273.
948          END IF
949        END DO
950      END IF
951
952      IF (type_profil==20) THEN
953        PRINT*, 'Profile SST 20'
954        !     Méthode 13 "Qmax2K" plateau réel é|  l'Equateur augmenté +2K
955
956        DO j = 1, klon
957          phy_sst(j, i) = 248. + 55. * (1 - sin(rlatd(j))**2)
958        enddo
959      endif
960
961      IF (type_profil==21) THEN
962        PRINT*, 'Profile SST 21'
963        !     Méthode 13 "Qmax2K" plateau réel é|  l'Equateur augmenté +2K
964        DO j = 1, klon
965          phy_sst(j, i) = 252. + 55. * (1 - sin(rlatd(j))**2)
966        enddo
967      endif
968
969    END DO
970
971    ! IM beg : verif profil SST: phy_sst
972    amn = min(phy_sst(1, 1), 1000.)
973    amx = max(phy_sst(1, 1), -1000.)
974    imn = 1
975    kmn = 1
976    imx = 1
977    kmx = 1
978    DO k = 1, year_len
979      DO i = 2, nlon
980        IF (phy_sst(i, k)<amn) THEN
981          amn = phy_sst(i, k)
982          imn = i
983          kmn = k
984        END IF
985        IF (phy_sst(i, k)>amx) THEN
986          amx = phy_sst(i, k)
987          imx = i
988          kmx = k
989        END IF
990      END DO
991    END DO
992
993    PRINT *, 'profil_sst: imn, kmn, phy_sst(imn,kmn) ', imn, kmn, amn
994    PRINT *, 'profil_sst: imx, kmx, phy_sst(imx,kmx) ', imx, kmx, amx
995    ! IM end : verif profil SST: phy_sst
996
997  END SUBROUTINE profil_sst
998
999END MODULE phyaqua_mod
Note: See TracBrowser for help on using the repository browser.