source: LMDZ6/trunk/libf/phylmd/cospv2/lmdz_cosp_interface.F90 @ 3981

Last change on this file since 3981 was 3723, checked in by idelkadi, 4 years ago

Debugging COSP v2 for simulators Calipso, Parasol, Cloudsat

File size: 24.1 KB
Line 
1!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2!
3! Nouveau code d'interface entre LMDZ et COSPv2 (version 2)
4! L'ancienne interface s'appelait "phys_cosp2" et avait ete concue pour COSPv1.4.
5! Dans cette nouvelle version de COSP, le code a ete restructure pour optimiser les calculs
6! des differents simulateurs et pour proposer de nouvelles fonctionnalites (par exemple,
7! intervenir sur les profils sous-maille, ou subcolumns, donnes en entre a COSP afin que
8! leur definition soit coherente avec les parametrisations du modele hote).
9! Cette version de COSP propose aussi de nombreux nouveaux diagnostics, notamment pour
10! le simulateur lidar (diagnostics CALIPSO-OPAQ, lidar sol 532nm et lidar ATLID 355nm).
11!
12! Interface reecrite par R.Guzman (01/2019), a partir de l'interface initiale concue
13! et ecrite par A.Idelkadi
14!
15!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
16!  subroutine phys_cosp2( itap,dtime,freq_cosp, &
17  subroutine lmdz_cosp_interface(itap, dtime, freq_cosp, ok_mensuelCOSP, ok_journeCOSP,   &
18                         ok_hfCOSP, ecrit_mth, ecrit_day, ecrit_hf, ok_all_xml,   &
19                         missing_val, Nptslmdz, Nlevlmdz, lon, lat, presnivs,     &
20                         overlaplmdz, sunlit, ref_liq, ref_ice, fracTerLic,       &
21                         u_wind, v_wind, phis, phi, ph, p, skt, t, sh, rh,        &
22                         tca, cca, mr_lsliq, mr_lsice, fl_lsrainI, fl_lssnowI,    &
23                         fl_ccrainI, fl_ccsnowI, mr_ozone, dtau_s, dem_s)
24
25!--------------  Inputs  ---------------
26! itap,                                 !Increment de la physiq
27! dtime,                                !Pas de temps physiq
28! overlaplmdz,                          !Type Overlap venant de LMDZ
29! Npoints,                              !Nb de points de la grille physiq
30! Nlevels,                              !Nb de niveaux verticaux
31! Ncolumns,                             !Number of subcolumns
32! lon,lat,                              !Longitudes et latitudes de la grille LMDZ
33! ref_liq, ref_ice,                     !Rayons effectifs des particules liq et ice (en micron)
34! fracTerLic,                           !Fraction terre a convertir en masque
35! u_wind, v_wind,                       !Vents a 10m ???
36! phi,                                  !Geopotentiel
37! phis,                                 !Geopotentiel sol
38! ph,                                   !pression pour chaque inter-couche
39! p,                                    !Pression aux milieux des couches
40! skt, t,                               !Temp au sol et temp 3D
41! sh,                                   !Humidite specifique
42! rh,                                   !Humidite relative
43! tca,                                  !Fraction nuageuse
44! cca                                   !Fraction nuageuse convective
45! mr_lsliq,                             !Liq Cloud water content
46! mr_lsice,                             !Ice Cloud water content
47! mr_ccliq,                             !Convective Cloud Liquid water content 
48! mr_ccice,                             !Cloud ice water content
49! fl_lsrain,                            !Large scale precipitation lic
50! fl_lssnow,                            !Large scale precipitation ice
51! fl_ccrain,                            !Convective precipitation lic
52! fl_ccsnow,                            !Convective precipitation ice
53! mr_ozone,                             !Concentration ozone (Kg/Kg)
54! dem_s                                 !Cloud optical emissivity
55! dtau_s                                !Cloud optical thickness
56! emsfc_lw = 1.                         !Surface emissivity dans radlwsw.F90
57
58
59!--------------  Outputs  --------------
60! La liste complete des diagnostics de sortie (observables simulees) que l'on peut
61! avoir avec COSPv2 se trouve au debut du fichier : cosp_read_otputkeys.F90
62
63
64!!! Modules specifiques a l'interface LMDZ-COSP
65  use mod_phys_lmdz_para
66  use mod_grid_phy_lmdz
67  use ioipsl
68  use iophy
69  use lmdz_cosp_output_mod
70  use lmdz_cosp_output_write_mod
71  use lmdz_cosp_read_outputkeys
72  use lmdz_cosp_subsample_and_optics_mod, only : subsample_and_optics
73  use lmdz_cosp_construct_destroy_mod
74
75!!! Modules faisant partie du code source de COSPv2
76  use cosp_kinds,  only: wp                         
77  use MOD_COSP_CONFIG,       only: N_HYDRO,RTTOV_MAX_CHANNELS, &
78                                    niv_sorties, vgrid_z_in
79  use mod_quickbeam_optics,  only: size_distribution,hydro_class_init, &
80                                   quickbeam_optics_init
81  use quickbeam,            only: radar_cfg
82  use mod_cosp,             only: cosp_init,cosp_optical_inputs, &
83                                  cosp_column_inputs,cosp_outputs, &
84                                  cosp_simulator
85
86
87!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
88!
89! Declaration des variables
90!
91!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
92
93  IMPLICIT NONE
94
95  ! Local variables
96  character(len=64),PARAMETER   :: cosp_input_nl  = 'cospv2_input_nl.txt'
97  character(len=64),PARAMETER   :: cosp_output_nl = 'cospv2_output_nl.txt'
98
99  integer, save                 :: isccp_topheight, isccp_topheight_direction, overlap
100  integer, save                 :: Ncolumns             ! Number of subcolumns in SCOPS
101  integer, save                 :: Npoints              ! Number of gridpoints
102!$OMP THREADPRIVATE(Npoints)
103  integer, save                 :: Nlevels              ! Number of model vertical levels
104  integer                       :: Nptslmdz, Nlevlmdz   ! Nb de points issus de physiq.F
105  integer, save                 :: Npoints_it           ! Max number of gridpoints to be
106                                                        ! processed in one iteration
107  type(cosp_config), save       :: cfg                  ! Variable qui contient les cles
108                                                        ! logiques des simulateurs et des
109                                                        ! diagnostics, definie dans:
110                                                        ! lmdz_cosp_construct_destroy_mod
111!$OMP THREADPRIVATE(cfg)
112
113  integer                       :: t0, t1, count_rate, count_max
114  real(wp), save                :: cloudsat_radar_freq, cloudsat_k2, rttov_ZenAng, co2, &
115                                   ch4, n2o, co, emsfc_lw
116!$OMP THREADPRIVATE(emsfc_lw)
117
118  integer, dimension(RTTOV_MAX_CHANNELS), save   :: rttov_Channels
119  real(wp), dimension(RTTOV_MAX_CHANNELS), save  :: rttov_Surfem
120  integer, save                                  :: surface_radar, use_mie_tables, &
121                                                    cloudsat_use_gas_abs, cloudsat_do_ray, &
122                                                    melt_lay
123  integer, save                                  :: lidar_ice_type
124  integer, save                                  :: rttov_platform, rttov_satellite, &
125                                                    rttov_Instrument, rttov_Nchannels
126  logical, save                                  :: use_vgrid_in, csat_vgrid_in, &
127                                                    use_precipitation_fluxes
128
129! Declaration necessaires pour les sorties IOIPSL
130  real          :: ecrit_day, ecrit_hf, ecrit_mth, missing_val
131  logical       :: ok_mensuelCOSP, ok_journeCOSP, ok_hfCOSP, ok_all_xml
132  logical, save :: debut_cosp=.true.
133!$OMP THREADPRIVATE(debut_cosp)
134
135  logical, save :: first_write=.true.
136!$OMP THREADPRIVATE(first_write)
137
138  integer, save :: cosp_init_flag = 0
139!$OMP THREADPRIVATE(cosp_init_flag)
140
141
142!-----------------------------  Input variables from LMDZ-GCM  -------------------------------
143  integer                                    :: overlaplmdz   ! overlap type: 1=max,
144                                                              ! 2=rand, 3=max/rand
145  real, dimension(Nptslmdz,Nlevlmdz)         :: phi, p, ph, T, sh, rh, tca, cca, mr_lsliq,   &
146                                                mr_lsice, mr_ccliq, mr_ccice, fl_lsrain,     &
147                                                fl_lssnow, fl_ccrain, fl_ccsnow, fl_lsgrpl,  &
148                                                zlev, zlev_half, mr_ozone, radliq, radice,   &
149                                                dtau_s, dem_s, dtau_c, dem_c, ref_liq, ref_ice
150  real, dimension(Nptslmdz,Nlevlmdz)         :: fl_lsrainI, fl_lssnowI, fl_ccrainI, fl_ccsnowI
151  real, dimension(Nptslmdz)                  :: lon, lat, skt, fracTerLic, u_wind, v_wind, &
152                                                phis, sunlit
153  real, dimension(Nptslmdz)                  :: land ! variables intermediaire pour masque TerLic
154  real, dimension(Nlevlmdz)                  :: presnivs
155  integer                                    :: itap, k, ip
156  real                                       :: dtime, freq_cosp
157  real, dimension(2)                         :: time_bnds
158
159  double precision                           :: d_dtime
160  double precision, dimension(2)             :: d_time_bnds
161
162
163  ! ######################################################################################
164  ! Declarations specific to COSP2
165  ! ######################################################################################
166
167  ! Local variables
168  logical :: &
169       Lsingle     = .true.,  & ! True if using MMF_v3_single_moment CLOUDSAT
170                                ! microphysical scheme (default)
171       Ldouble     = .false.    ! True if using MMF_v3.5_two_moment CLOUDSAT
172                                ! microphysical scheme
173  type(size_distribution), save              :: sd            ! Hydrometeor description
174!$OMP THREADPRIVATE(sd)
175  type(radar_cfg), save                      :: rcfg_cloudsat ! Radar configuration
176!$OMP THREADPRIVATE(rcfg_cloudsat)
177  real, dimension(Nptslmdz,Nlevlmdz,N_HYDRO) :: Reff          ! Liquid and Ice particles
178                                                              ! effective radius
179  type(cosp_outputs)                         :: cospOUT       ! COSP simulator outputs
180  type(cosp_optical_inputs)                  :: cospIN        ! COSP optical (or derived?)
181                                                              ! fields needed by simulators
182  type(cosp_column_inputs)                   :: cospstateIN   ! COSP model fields needed
183                                                              ! by simulators
184  character(len=256), dimension(100)         :: cosp_status
185  character(len=64), save                    :: cloudsat_micro_scheme
186
187  ! Indices to address arrays of LS and CONV hydrometeors
188  integer,parameter :: &
189       I_LSCLIQ = 1, & ! Large-scale (stratiform) liquid
190       I_LSCICE = 2, & ! Large-scale (stratiform) ice
191       I_LSRAIN = 3, & ! Large-scale (stratiform) rain
192       I_LSSNOW = 4, & ! Large-scale (stratiform) snow
193       I_CVCLIQ = 5, & ! Convective liquid
194       I_CVCICE = 6, & ! Convective ice
195       I_CVRAIN = 7, & ! Convective rain
196       I_CVSNOW = 8, & ! Convective snow
197       I_LSGRPL = 9    ! Large-scale (stratiform) groupel
198
199! Parametres qui sont lus a partir du fichier "cosp_input_nl.txt"
200   namelist/COSP_INPUT/overlap, isccp_topheight, isccp_topheight_direction,  &
201              npoints_it, ncolumns, use_vgrid_in, csat_vgrid_in,  &
202              cloudsat_radar_freq, surface_radar, use_mie_tables,  &
203              cloudsat_use_gas_abs, cloudsat_do_ray, melt_lay, cloudsat_k2,  &
204              cloudsat_micro_scheme, lidar_ice_type, use_precipitation_fluxes,   &
205              rttov_platform, rttov_satellite, rttov_Instrument, rttov_Nchannels, &
206              rttov_Channels, rttov_Surfem, rttov_ZenAng, co2, ch4, n2o, co
207
208!------------------------   Fin declaration des variables   ------------------------
209
210
211!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
212!
213! 1) Lecture du fichier "cosp_input_nl.txt", parametres d'entree pour COSP
214!
215!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
216
217print*,'Entree lmdz_cosp_interface' !phys_cosp2'
218if (debut_cosp) then
219  NPoints=Nptslmdz
220  Nlevels=Nlevlmdz
221! Surface emissivity
222        emsfc_lw = 1.
223 
224! Lecture du namelist input
225!  CALL read_cosp_input
226   IF (is_master) THEN
227      OPEN(10,file=cosp_input_nl,status='old')
228      READ(10,nml=cosp_input)
229      CLOSE(10)
230   ENDIF
231
232
233!$OMP BARRIER
234    CALL bcast(overlap)
235    CALL bcast(isccp_topheight)
236    CALL bcast(isccp_topheight_direction)
237    CALL bcast(npoints_it)
238    CALL bcast(ncolumns)
239    CALL bcast(use_vgrid_in)
240    CALL bcast(csat_vgrid_in)
241    CALL bcast(cloudsat_radar_freq)
242    CALL bcast(surface_radar)
243    CALL bcast(cloudsat_use_gas_abs)
244    CALL bcast(cloudsat_do_ray)
245    CALL bcast(cloudsat_k2)
246    CALL bcast(lidar_ice_type)
247    CALL bcast(use_precipitation_fluxes)
248    CALL bcast(rttov_platform)
249    CALL bcast(rttov_satellite)
250    CALL bcast(rttov_Instrument)
251    CALL bcast(rttov_Nchannels)
252    CALL bcast(rttov_Channels)
253    CALL bcast(rttov_Surfem)
254    CALL bcast(rttov_ZenAng)
255    CALL bcast(co2)
256    CALL bcast(ch4)
257    CALL bcast(n2o)
258    CALL bcast(co)
259    CALL bcast(cloudsat_micro_scheme)
260
261    print*,'ok read  cosp_input_nl'
262
263! Clefs Outputs initialisation
264#ifdef CPP_XIOS
265  call cosp_outputkeys_init(cfg)
266#else
267   call read_cosp_output_nl(itap,cosp_output_nl,cfg)
268#endif
269
270  print*,' Cles des differents simulateurs cosp a itap :',itap
271  print*,'cfg%Lcloudsat, cfg%Lcalipso, cfg%LgrLidar532, cfg%Latlid, cfg%Lparasol, &
272        cfg%Lisccp, cfg%Lmisr, cfg%Lmodis, cfg%Lrttov', &
273        cfg%Lcloudsat, cfg%Lcalipso, cfg%LgrLidar532, cfg%Latlid, cfg%Lparasol, &
274        cfg%Lisccp, cfg%Lmisr, cfg%Lmodis, cfg%Lrttov
275
276    if (overlaplmdz.ne.overlap) then
277       print*,'Attention overlaplmdz different de overlap lu dans namelist '
278    endif
279
280#ifdef CPP_XIOS
281   print*,'On passe par ifdef CPP_XIOS'
282#else
283if (cosp_init_flag .eq. 0) then
284
285    ! Initialize the distributional parameters for hydrometeors in radar simulator.
286    ! In COSPv1.4, this was declared in cosp_defs.f.
287    if (cloudsat_micro_scheme == 'MMF_v3.5_two_moment')  then
288       ldouble = .true.
289       lsingle = .false.
290    endif
291    call hydro_class_init(lsingle,ldouble,sd)
292    call quickbeam_optics_init()
293
294  print*,' just before call COSP_INIT, cosp_init_flag =', cosp_init_flag
295  call COSP_INIT(cfg%Lisccp, cfg%Lmodis, cfg%Lmisr, cfg%Lcloudsat, cfg%Lcalipso, &
296        cfg%LgrLidar532, cfg%Latlid, cfg%Lparasol, cfg%Lrttov,          &
297        cloudsat_radar_freq, cloudsat_k2, cloudsat_use_gas_abs,         &
298        cloudsat_do_ray, isccp_topheight, isccp_topheight_direction,    &
299        surface_radar, rcfg_cloudsat, use_vgrid_in, csat_vgrid_in,      &
300        niv_sorties, Nlevels, cloudsat_micro_scheme)
301  cosp_init_flag = 1
302 print*,' just after call COSP_INIT, cosp_init_flag =', cosp_init_flag
303endif
304#endif
305
306  print*,'Fin lecture Namelists, debut_cosp =',debut_cosp
307
308  endif ! debut_cosp
309
310
311!!! Ici on modifie les cles logiques pour les outputs selon les champs actives dans les .xml
312  if ((itap.ge.1).and.(first_write))then
313#ifdef CPP_XIOS
314    call read_xiosfieldactive(cfg)
315#endif
316    first_write=.false.
317
318if (cosp_init_flag .eq. 0) then
319
320    ! Initialize the distributional parameters for hydrometeors in radar simulator.
321    ! In COSPv1.4, this was declared in cosp_defs.f.
322    if (cloudsat_micro_scheme == 'MMF_v3.5_two_moment')  then
323       ldouble = .true.
324       lsingle = .false.
325    endif
326    call hydro_class_init(lsingle,ldouble,sd)
327    call quickbeam_optics_init()
328
329  print*,' just before call COSP_INIT, cosp_init_flag =', cosp_init_flag
330  call COSP_INIT(cfg%Lisccp, cfg%Lmodis, cfg%Lmisr, cfg%Lcloudsat, cfg%Lcalipso, &
331        cfg%LgrLidar532, cfg%Latlid, cfg%Lparasol, cfg%Lrttov,          &
332        cloudsat_radar_freq, cloudsat_k2, cloudsat_use_gas_abs,         &
333        cloudsat_do_ray, isccp_topheight, isccp_topheight_direction,    &
334        surface_radar, rcfg_cloudsat, use_vgrid_in, csat_vgrid_in,      &
335        niv_sorties, Nlevels, cloudsat_micro_scheme)
336  cosp_init_flag = 1
337 print*,' just after call COSP_INIT, cosp_init_flag =', cosp_init_flag
338endif ! cosp_init_flag
339
340
341print*,' Cles des differents simulateurs cosp a itap :',itap
342print*,'cfg%Lcloudsat, cfg%Lcalipso, cfg%LgrLidar532, cfg%Latlid, cfg%Lparasol, &
343        cfg%Lisccp, cfg%Lmisr, cfg%Lmodis, cfg%Lrttov', &
344        cfg%Lcloudsat, cfg%Lcalipso, cfg%LgrLidar532, cfg%Latlid, cfg%Lparasol, &
345        cfg%Lisccp, cfg%Lmisr, cfg%Lmodis, cfg%Lrttov
346
347  endif !(itap.gt.1).and.(first_write)
348
349  time_bnds(1) = dtime-dtime/2.
350  time_bnds(2) = dtime+dtime/2.
351
352  d_time_bnds=time_bnds
353  d_dtime=dtime
354
355!-------------------------   Fin initialisation de COSP   --------------------------
356
357
358!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
359!
360! 3) Calculs des champs d'entree COSP a partir des variables LMDZ
361!
362!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
363
364    ! 0) Create ptop/ztop for gbx%pf and gbx%zlev are for the the interface,
365    !    also reverse CAM height/pressure values for input into CSOP
366    !    CAM state%pint from top to surface, COSP wants surface to top.
367
368! 0) Altitudes du modele calculees a partir de la variable geopotentiel phi et phis   
369        zlev = phi/9.81
370
371        zlev_half(:,1) = phis(:)/9.81
372        do k = 2, Nlevels
373          do ip = 1, Npoints
374           zlev_half(ip,k) = phi(ip,k)/9.81 + &
375               (phi(ip,k)-phi(ip,k-1))/9.81 * (ph(ip,k)-p(ip,k)) / (p(ip,k)-p(ip,k-1))
376          enddo
377        enddo
378
379! 1) Quantite de nuages (couverture?), convectif (=0) et total
380        cca = 0._wp ! convective_cloud_amount (1)
381        tca = tca   ! total_cloud_amount (1)
382
383! 2) Humidite relative est donnee tel quel (variable rh)
384
385! 3) Masque terre/mer a partir de la variable fracTerLic
386        do ip = 1, Npoints
387          if (fracTerLic(ip).ge.0.5) then
388             land(ip) = 1.
389          else
390             land(ip) = 0.
391          endif
392        enddo
393
394
395! A voir l equivalent LMDZ
396  mr_ccliq = 0.0
397  mr_ccice = 0.0
398!!!        gbx%mr_hydro(:,:,I_LSCLIQ) = mr_lsliq !mixing_ratio_large_scale_cloud_liquid (kg/kg)
399!!!        gbx%mr_hydro(:,:,I_LSCICE) = mr_lsice !mixing_ratio_large_scale_cloud_ic
400!!!        gbx%mr_hydro(:,:,I_CVCLIQ) = mr_ccliq !mixing_ratio_convective_cloud_liquid
401!!!        gbx%mr_hydro(:,:,I_CVCICE) = mr_ccice !mixing_ratio_convective_cloud_ice
402! A revoir
403        fl_lsrain = fl_lsrainI + fl_ccrainI
404        fl_lssnow = fl_lssnowI + fl_ccsnowI
405!!!        gbx%rain_ls = fl_lsrain !flux_large_scale_cloud_rain (kg m^-2 s^-1)
406!!!        gbx%snow_ls = fl_lssnow !flux_large_scale_cloud_snow
407!  A voir l equivalent LMDZ
408        fl_lsgrpl = 0.
409        fl_ccsnow = 0.
410        fl_ccrain = 0.
411!!!        gbx%grpl_ls = fl_lsgrpl  !flux_large_scale_cloud_graupel
412!!!        gbx%rain_cv = fl_ccrain  !flux_convective_cloud_rain
413!!!        gbx%snow_cv = fl_ccsnow  !flux_convective_cloud_snow
414
415        ! ISCCP simulator
416        dtau_c  = 0.
417        dem_c   = 0.
418
419! note: reff_cosp dimensions should be same as cosp (reff_cosp has 9 hydrometeor dimension)
420     Reff(1:Npoints,1:Nlevels,1:N_HYDRO) = 0.
421     Reff(:,:,I_LSCLIQ) = ref_liq*1e-6
422     Reff(:,:,I_LSCICE) = ref_ice*1e-6
423     Reff(:,:,I_CVCLIQ) = ref_liq*1e-6
424     Reff(:,:,I_CVCICE) = ref_ice*1e-6
425
426
427if (cosp_init_flag .eq. 1) then      ! cosp_init_flag = 1
428
429!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
430!
431! 4a) On construit la variable cospOUT qui contient tous les diagnostics de sortie.
432! Elle sera remplie lors de l'appel du simulateur COSP
433!
434!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
435
436    call construct_cosp_outputs(cfg,Npoints,Ncolumns,Nlevels,niv_sorties,0,cospOUT)
437
438
439!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
440!
441! 4b) On construit la variable cospstateIN que l'on va remplir avec les champs LMDZ
442! Les champ verticaux doivent etre donnes a l'envers, c-a-d : (Nlevels:1) = (TOA:SFC)
443!
444!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
445
446    call construct_cospstateIN(Npoints,Nlevels,0,cospstateIN)
447
448    cospstateIN%lat                                     = lat(1:Npoints)
449    cospstateIN%lon                                     = lon(1:Npoints)
450    cospstateIN%at                                      = t(1:Npoints,Nlevels:1:-1)
451    cospstateIN%qv                                      = sh(1:Npoints,Nlevels:1:-1)
452    cospstateIN%o3                                      = mr_ozone(1:Npoints,Nlevels:1:-1) 
453    cospstateIN%sunlit                                  = sunlit(1:Npoints)
454    cospstateIN%skt                                     = skt(1:Npoints)
455    cospstateIN%land                                    = land(1:Npoints)
456    cospstateIN%surfelev                                = zlev_half(1:Npoints,1)
457    cospstateIN%pfull                                   = p(1:Npoints,Nlevels:1:-1)
458    cospstateIN%phalf(1:Npoints,1)                      = 0._wp
459    cospstateIN%phalf(1:Npoints,2:Nlevels+1)            = ph(1:Npoints,Nlevels:1:-1) 
460    cospstateIN%hgt_matrix                              = zlev(1:Npoints,Nlevels:1:-1)
461    cospstateIN%hgt_matrix_half(1:Npoints,Nlevels+1)    = 0._wp
462    cospstateIN%hgt_matrix_half(1:Npoints,1:Nlevels)    = zlev_half(1:Npoints,Nlevels:1:-1)
463
464
465!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
466!
467! 4c) On construit la variable cospIN qui contient les proprietes optiques subcolumn
468! pour COSP. Elle sera essentiellement remplie dans la subroutine subsample_and_optics
469! ou sont appeles SCOPS, PREC_SCOPS et les subroutines qui calculent les signaux
470! simules pour chaque simulateur actif.
471!
472!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
473
474    call construct_cospIN(cfg,Npoints,Ncolumns,Nlevels,cospIN)
475    cospIN%emsfc_lw  = emsfc_lw
476    if (cfg%Lcloudsat) cospIN%rcfg_cloudsat = rcfg_cloudsat
477
478
479!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
480!
481! 5) Appel de subsample_and_optics : Les champs verticaux doivent etre donnes a
482! l'envers comme pour le remplissage de cospstateIN, c-a-d : (Nlevels:1) = (TOA:SFC)
483!
484!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
485
486    call subsample_and_optics(cfg, Npoints, Nlevels, Ncolumns, N_HYDRO,overlap, &
487                              use_precipitation_fluxes, lidar_ice_type, sd, &
488                              tca(1:Npoints,Nlevels:1:-1), cca(1:Npoints,Nlevels:1:-1), &
489                              fl_lsrain(1:Npoints,Nlevels:1:-1), &
490                              fl_lssnow(1:Npoints,Nlevels:1:-1), &
491                              fl_lsgrpl(1:Npoints,Nlevels:1:-1), &
492                              fl_ccrain(1:Npoints,Nlevels:1:-1), &
493                              fl_ccsnow(1:Npoints,Nlevels:1:-1), &
494                              mr_lsliq(1:Npoints,Nlevels:1:-1), &
495                              mr_lsice(1:Npoints,Nlevels:1:-1), &
496                              mr_ccliq(1:Npoints,Nlevels:1:-1), &
497                              mr_ccice(1:Npoints,Nlevels:1:-1), &
498                              Reff(1:Npoints,Nlevels:1:-1,:), &
499                              dtau_c(1:Npoints,Nlevels:1:-1), &
500                              dtau_s(1:Npoints,Nlevels:1:-1), &
501                              dem_c(1:Npoints,Nlevels:1:-1), &
502                              dem_s(1:Npoints,Nlevels:1:-1), cospstateIN, cospIN)
503
504
505!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
506!
507! 6) On appelle le simulateur COSPv2
508!
509!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
510
511      print*,'call simulateur'
512
513cosp_status = COSP_SIMULATOR(cospIN, cospstateIN, cospOUT, 1, Npoints, debug=.false.)
514
515endif ! cosp_init_flag = 1
516
517
518!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
519!
520! 7a) Ecriture des sorties 1: on cree d'abord les fichiers NCDF pour ecrire les sorties
521! en appelant lmdz_cosp_output_open (lors du premier appel de cette interface pour les
522! 2 options d'ecriture), ou sont definis les axes et les caracteristiques
523! des fichiers de sortie avec les diagnostics.
524!
525!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
526
527    if (debut_cosp) then
528      !$OMP MASTER
529
530        print *, ' Open outpts files and define axis'
531        call lmdz_cosp_output_open(Nlevlmdz, Ncolumns, presnivs, dtime, freq_cosp,       &
532                                   ok_mensuelCOSP, ok_journeCOSP, ok_hfCOSP, ok_all_xml, &
533                                   ecrit_mth, ecrit_day, ecrit_hf, use_vgrid_in,         &
534                                   niv_sorties, vgrid_z_in, zlev(1,:))
535
536      !$OMP END MASTER
537      !$OMP BARRIER
538       debut_cosp=.false.
539      endif ! debut_cosp
540
541if (cosp_init_flag .eq. 1) then
542
543!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
544!
545! 7b) Ecriture des sorties 2: le remplissage des fichiers de sortie se fait a chaque
546! appel de cette interface avec une difference entre les 2 options d'ecriture:
547!
548! AVEC xios, on commence a remplir les fichiers de sortie a partir du DEUXIEME
549! appel de cette interface (lorsque cosp_init_flag = 1).
550!
551! SANS xios, on commence a remplir les fichiers de sortie a partir du PREMIER
552! appel de cette interface (lorsque cosp_init_flag = 1).
553!
554!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
555
556       print *, 'Calling write output'
557        call lmdz_cosp_output_write(Nlevlmdz, Npoints, Ncolumns, itap, dtime, freq_COSP, &
558                                    missing_val, cfg, niv_sorties, cospOUT)
559
560!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
561!
562! 8) On libere la memoire allouee lors de cet appel a l'interface
563!
564!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
565
566    call destroy_cospIN(cospIN)
567    call destroy_cospstateIN(cospstateIN)
568    call destroy_cosp_outputs(cospOUT)
569
570endif ! cosp_init_flag = 1
571 
572end subroutine lmdz_cosp_interface
Note: See TracBrowser for help on using the repository browser.