Changeset 3387


Ignore:
Timestamp:
Sep 4, 2018, 7:42:50 PM (6 years ago)
Author:
oboucher
Message:

adding infocfields_init routine and
its call from physiq. This routines reads
which fields need to be transferred between
model components in ESM configuration.
No impact whatsoever in LMDZ mode.

Location:
LMDZ6/trunk/libf/phylmd
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • LMDZ6/trunk/libf/phylmd/carbon_cycle_mod.F90

    r3385 r3387  
    11MODULE carbon_cycle_mod
    2 ! Controle module for the carbon CO2 tracers :
     2!=======================================================================
     3!   Authors: Patricia Cadule and Laurent Fairhead
     4!            base sur un travail anterieur mene par Patricia Cadule et Josefine Ghattas
     5!
     6!  Purpose and description:
     7!  -----------------------
     8! Control module for the carbon CO2 tracers :
    39!   - Identification
    410!   - Get concentrations comming from coupled model or read from file to tracers
     
    612!   - Calculate new carbon flux for sending to coupled models (PISCES and ORCHIDEE)
    713!
    8 ! Author : Josefine GHATTAS, Patricia CADULE
     14! Module permettant de mettre a jour les champs (puits et sources) pour le
     15! transport de CO2 en online (IPSL-CM et LMDZOR) et offline (lecture de carte)
     16!
     17! Le cas online/offline est defini par le flag carbon_cycle_cpl (y/n)
     18! Le transport du traceur CO2 est defini par le flag carbon_cycle_tr (y/n)
     19! la provenance des champs (termes de puits) est denini par le flag level_coupling_esm
     20!
     21! level_coupling_esm : level of coupling of the biogeochemical fields between
     22! LMDZ, ORCHIDEE and NEMO
     23! Definitions of level_coupling_esm in physiq.def
     24! level_coupling_esm = 0  ! No field exchange between LMDZ and ORCHIDEE models
     25!                         ! No field exchange between LMDZ and NEMO
     26! level_coupling_esm = 1  ! Field exchange between LMDZ and ORCHIDEE models
     27!                         ! No field exchange between LMDZ and NEMO models
     28! level_coupling_esm = 2  ! No field exchange between LMDZ and ORCHIDEE models
     29!                         ! Field exchange between LMDZ and NEMO models
     30! level_coupling_esm = 3  ! Field exchange between LMDZ and ORCHIDEE models
     31!                         ! Field exchange between LMDZ and NEMO models
     32!=======================================================================
    933
    1034  IMPLICIT NONE
    1135  SAVE
    1236  PRIVATE
    13   PUBLIC :: carbon_cycle_init, carbon_cycle
     37  PUBLIC :: carbon_cycle_init, carbon_cycle, infocfields_init
    1438
    1539! Variables read from parmeter file physiq.def
     
    1842  LOGICAL, PUBLIC :: carbon_cycle_cpl       ! Coupling of CO2 fluxes between LMDZ/ORCHIDEE and LMDZ/OCEAN(PISCES)
    1943!$OMP THREADPRIVATE(carbon_cycle_cpl)
    20   INTEGER, SAVE, PUBLIC :: level_coupling_esm ! Level of coupling for the ESM - 0, 1, 2, 3
     44  INTEGER, PUBLIC :: level_coupling_esm ! Level of coupling for the ESM - 0, 1, 2, 3
    2145!$OMP THREADPRIVATE(level_coupling_esm)
    2246
     
    6488!$OMP THREADPRIVATE(co2_send)
    6589
     90! nbfields : total number of fields
     91  INTEGER, PUBLIC :: nbcf
     92!$OMP THREADPRIVATE(nbcf)
     93
     94! nbcf_in : number of fields IN
     95  INTEGER, PUBLIC  :: nbcf_in
     96!$OMP THREADPRIVATE(nbcf_in)
     97
     98! nbcf_in_orc : number of fields IN
     99  INTEGER, PUBLIC  :: nbcf_in_orc
     100!$OMP THREADPRIVATE(nbcf_in_orc)
     101
     102! nbcf_in_inca : number of fields IN (from INCA)
     103  INTEGER, PUBLIC  :: nbcf_in_inca
     104!$OMP THREADPRIVATE(nbcf_in_inca)
     105
     106! nbcf_in_nemo : number of fields IN (from nemo)
     107  INTEGER, PUBLIC  :: nbcf_in_nemo
     108!$OMP THREADPRIVATE(nbcf_in_nemo)
     109
     110! nbcf_in_ant : number of fields IN (from anthropogenic sources)
     111  INTEGER, PUBLIC  :: nbcf_in_ant
     112!$OMP THREADPRIVATE(nbcf_in_ant)
     113
     114! nbcf_out : number of fields OUT
     115  INTEGER, PUBLIC :: nbcf_out
     116!$OMP THREADPRIVATE(nbcf_out)
     117
     118! Name of variables
     119  CHARACTER(len=25), ALLOCATABLE, DIMENSION(:), PUBLIC :: cfname     ! coupling field short name for restart (?) and diagnostics
     120!$OMP THREADPRIVATE(cfname)
     121
     122  CHARACTER(len=25), ALLOCATABLE, DIMENSION(:), PUBLIC :: cfname_in  ! coupling field short name for restart (?) and diagnostics
     123!$OMP THREADPRIVATE(cfname_in)
     124
     125  CHARACTER(len=25), ALLOCATABLE, DIMENSION(:), PUBLIC :: cfname_out ! coupling field short name for restart (?) and diagnostics
     126!$OMP THREADPRIVATE(cfname_out)
     127
     128  CHARACTER(len=15), ALLOCATABLE, DIMENSION(:), PUBLIC :: cfunits_in  !  coupling field units for diagnostics
     129!$OMP THREADPRIVATE(cfunits_in)
     130
     131  CHARACTER(len=15), ALLOCATABLE, DIMENSION(:), PUBLIC :: cfunits_out !  coupling field units for diagnostics
     132!$OMP THREADPRIVATE(cfunits_out)
     133
     134  CHARACTER(len=120), ALLOCATABLE, DIMENSION(:), PUBLIC :: cftext_in  ! coupling field long name for diagnostics
     135!$OMP THREADPRIVATE(cftext_in)
     136
     137  CHARACTER(len=120), ALLOCATABLE, DIMENSION(:), PUBLIC :: cftext_out ! coupling field long name for diagnostics
     138!$OMP THREADPRIVATE(cftext_out)
     139
     140  CHARACTER(len=5), ALLOCATABLE, DIMENSION(:), PUBLIC :: cfmod1 ! model 1 (rreference) : LMDz
     141!$OMP THREADPRIVATE(cfmod1)
     142
     143  CHARACTER(len=5), ALLOCATABLE, DIMENSION(:), PUBLIC :: cfmod2 ! model 2
     144!$OMP THREADPRIVATE(cfmod2)
     145
     146  REAL, ALLOCATABLE, DIMENSION(:,:), PUBLIC :: zcfields_in
     147!$OMP THREADPRIVATE(zcfields_in)
     148
     149  REAL, ALLOCATABLE, DIMENSION(:,:), PUBLIC :: zcfields_out
     150!$OMP THREADPRIVATE(zcfields_out)
    66151
    67152  TYPE, PUBLIC ::   co2_trac_type
     
    111196    INTEGER               :: teststop
    112197
    113 
    114 
    115198! 1) Read controle parameters from .def input file
    116199! ------------------------------------------------
     
    125208       WRITE(lunout,*) 'carbon_cycle_fos_fuel = ', fos_fuel_s
    126209    END IF
    127 
    128210
    129211    ! Read parmeter for calculation compatible emission
     
    297379  END SUBROUTINE carbon_cycle_init
    298380
    299 
    300381  SUBROUTINE carbon_cycle(nstep, pdtphys, pctsrf, tr_seri, source)
    301382! Subroutine for injection of co2 in the tracers
     
    483564  END SUBROUTINE carbon_cycle
    484565 
     566  SUBROUTINE infocfields_init
     567
     568    USE control_mod, ONLY: planet_type
     569    USE phys_cal_mod, ONLY : mth_cur
     570    USE mod_synchro_omp
     571    USE mod_phys_lmdz_para, ONLY: is_mpi_root, is_omp_root
     572    USE mod_phys_lmdz_transfert_para
     573    USE mod_phys_lmdz_omp_transfert
     574    USE dimphy, ONLY: klon
     575
     576    IMPLICIT NONE
     577
     578!=======================================================================
     579!
     580!   Authors: Patricia Cadule and Laurent Fairhead 
     581!   -------
     582!
     583!  Purpose and description:
     584!  -----------------------
     585!
     586! Infofields
     587! this routine enables to define the field exchanges in both directions between
     588! the atmospheric circulation model (LMDZ) and ORCHIDEE. In the future this
     589! routing might apply to other models (e.g., NEMO, INCA, ...).
     590! Therefore, currently with this routine, it is possible to define the coupling
     591! fields only between LMDZ and ORCHIDEE.
     592! The coupling_fields.def file enables to define the name of the exchanged
     593! fields at the coupling interface.
     594! field_in_names : the set of names of the exchanged fields in input to ORCHIDEE
     595! (LMDZ to ORCHIDEE)
     596! field_out_names : the set of names of the exchanged fields in output of
     597! ORCHIDEE (ORCHIDEE to LMDZ)
     598! n : the number of exchanged fields at th coupling interface
     599! nb_fields_in : number of inputs fields to ORCHIDEE (LMDZ to ORCHIDEE)
     600! nb_fields_out : number of ouput fields of ORCHIDEE (ORCHIDEE to LMDZ)
     601!
     602! The syntax for coupling_fields.def is as follows:
     603! IMPORTANT: each column entry must be separated from the previous one by 3
     604! spaces and only that
     605! field name         coupling          model 1         model 2         long_name
     606!                    direction
     607!   10char  -3spaces-  3char  -3spaces- 4char -3spaces- 4char -3spaces-  30char
     608!
     609! n
     610! FIELD1 IN LMDZ ORC
     611! ....
     612! FIELD(j) IN LMDZ ORC
     613! FIELD(j+1) OUT LMDZ ORC
     614! ...
     615! FIELDn OUT LMDZ ORC   
     616!
     617!=======================================================================
     618!   ... 22/12/2017 ....
     619!-----------------------------------------------------------------------
     620! Declarations
     621
     622  INCLUDE "clesphys.h"
     623  INCLUDE "dimensions.h"
     624  INCLUDE "iniprint.h"
     625
     626! Local variables
     627
     628  INTEGER :: iq,  ierr, stat, error
     629
     630  CHARACTER(LEN=20), ALLOCATABLE, DIMENSION(:), SAVE  :: cfname_root
     631  CHARACTER(LEN=120), ALLOCATABLE, DIMENSION(:), SAVE :: cftext_root
     632  CHARACTER(LEN=15), ALLOCATABLE, DIMENSION(:), SAVE  :: cfunits_root
     633
     634  CHARACTER(len=3), ALLOCATABLE, DIMENSION(:) :: cfintent_root
     635  CHARACTER(len=5), ALLOCATABLE, DIMENSION(:) :: cfmod1_root
     636  CHARACTER(len=5), ALLOCATABLE, DIMENSION(:) :: cfmod2_root
     637
     638  LOGICAL, ALLOCATABLE, DIMENSION(:), SAVE :: mask_in_root
     639  LOGICAL, ALLOCATABLE, DIMENSION(:), SAVE :: mask_out_root
     640
     641  CHARACTER(len=*),parameter :: modname="infocfields"
     642
     643!-----------------------------------------------------------------------
     644
     645nbcf=0
     646nbcf_in=0
     647nbcf_out=0
     648
     649IF (planet_type=='earth') THEN
     650
     651    IF (is_mpi_root .AND. is_omp_root) THEN
     652
     653       IF (level_coupling_esm.GT.0) THEN
     654
     655          OPEN(200,file='coupling_fields.def',form='formatted',status='old', iostat=ierr)
     656
     657          IF (ierr.EQ.0) THEN
     658
     659             WRITE(lunout,*) trim(modname),': Open coupling_fields.def : ok'
     660             READ(200,*) nbcf
     661             WRITE(lunout,*) 'infocfields_mod.F90 --- nbcf=',nbcf
     662             IF (.NOT. ALLOCATED(cfname_root))   ALLOCATE(cfname_root(nbcf))
     663             IF (.NOT. ALLOCATED(cfintent_root)) ALLOCATE(cfintent_root(nbcf))
     664             IF (.NOT. ALLOCATED(cfmod1_root))   ALLOCATE(cfmod1_root(nbcf))
     665             IF (.NOT. ALLOCATED(cfmod2_root))   ALLOCATE(cfmod2_root(nbcf))
     666             IF (.NOT. ALLOCATED(cftext_root))   ALLOCATE(cftext_root(nbcf))
     667             IF (.NOT. ALLOCATED(cfunits_root))  ALLOCATE(cfunits_root(nbcf))
     668             IF (.NOT. ALLOCATED(mask_in_root))  ALLOCATE(mask_in_root(nbcf))
     669             IF (.NOT. ALLOCATED(mask_out_root)) ALLOCATE(mask_out_root(nbcf))
     670
     671             nbcf_in=0
     672             nbcf_out=0
     673             DO iq=1,nbcf
     674                WRITE(lunout,*) 'infofields : field=',iq
     675                READ(200,'(A15,3X,A3,3X,A5,3X,A5,3X,A120,3X,A15)',IOSTAT=ierr) &
     676                   cfname_root(iq),cfintent_root(iq),cfmod1_root(iq),cfmod2_root(iq),cftext_root(iq),cfunits_root(iq)
     677                cfname_root(iq)=TRIM(cfname_root(iq))
     678                cfintent_root(iq)=TRIM(cfintent_root(iq))
     679                cfmod1_root(iq)=TRIM(cfmod1_root(iq))
     680                cfmod2_root(iq)=TRIM(cfmod2_root(iq))
     681                cftext_root(iq)=TRIM(cftext_root(iq))
     682                cfunits_root(iq)=TRIM(cfunits_root(iq))
     683                WRITE(lunout,*) 'coupling field: ',cfname_root(iq), & 
     684                               ', number: ',iq,', INTENT: ',cfintent_root(iq)
     685                WRITE(lunout,*) 'coupling field: ',cfname_root(iq), &
     686                               ', number: ',iq,', model 1 (ref): ',cfmod1_root(iq),', model 2: ',cfmod2_root(iq)
     687                WRITE(lunout,*) 'coupling field: ',cfname_root(iq), &
     688                               ', number: ',iq,', long name: ',cftext_root(iq),', units ',cfunits_root(iq)
     689                IF ((TRIM(cfintent_root(iq)).NE.'OUT') .AND. (nbcf_in.LE.nbcf)) THEN
     690                  nbcf_in=nbcf_in+1
     691                  mask_in_root(iq)=.TRUE.
     692                  mask_out_root(iq)=.FALSE.
     693                ELSE IF ((TRIM(cfintent_root(iq)).EQ.'OUT') .AND. (nbcf_out.LE.nbcf)) THEN
     694                  nbcf_out=nbcf_out+1
     695                  mask_in_root(iq)=.FALSE.
     696                  mask_out_root(iq)=.TRUE.
     697                ELSE
     698                  WRITE(lunout,*) 'abort_gcm --- nbcf: ',nbcf
     699                  WRITE(lunout,*) 'abort_gcm --- nbcf_in: ',nbcf_in
     700                  WRITE(lunout,*) 'abort_gcm --- nbcf_out: ',nbcf_out
     701                  CALL abort_gcm('infocfields_init','Problem in the definition of the coupling fields',1)
     702               ENDIF
     703             ENDDO !DO iq=1,nbcf
     704          ELSE
     705             WRITE(lunout,*) trim(modname),': infocfields_mod.F90 --- Problem in opening coupling_fields.def'
     706             WRITE(lunout,*) trim(modname),': infocfields_mod.F90 --- WARNING using defaut values'
     707          ENDIF ! ierr
     708          CLOSE(200)
     709
     710       ENDIF ! level_coupling_esm
     711
     712       IF (nbcf_out .EQ. 0) nbcf_out=-1
     713       IF (nbcf_in .EQ. 0) nbcf_in=-1
     714
     715    ENDIF !   (is_mpi_root .AND. is_omp_root)
     716!$OMP BARRIER
     717
     718    CALL bcast(nbcf)
     719    CALL bcast(nbcf_in)
     720    CALL bcast(nbcf_out)
     721
     722    WRITE(lunout,*) 'infocfields_mod.F90 --- nbcf=',nbcf
     723    WRITE(lunout,*) 'infocfields_mod.F90 --- nbcf_in=',nbcf_in
     724    WRITE(lunout,*) 'infocfields_mod.F90 --- nbcf_out=',nbcf_out
     725
     726    IF ((nbcf.GT.0) .AND. (.NOT. ALLOCATED(cfname))) ALLOCATE(cfname(nbcf))
     727    IF ((nbcf_in.GT.0) .AND. (.NOT. ALLOCATED(cfname_in))) ALLOCATE(cfname_in(nbcf_in))
     728    IF ((nbcf_in.GT.0) .AND. (.NOT. ALLOCATED(cftext_in))) ALLOCATE(cftext_in(nbcf_in))
     729    IF ((nbcf_out.GT.0) .AND. (.NOT. ALLOCATED(cfname_out))) ALLOCATE(cfname_out(nbcf_out))
     730    IF ((nbcf_out.GT.0) .AND. (.NOT. ALLOCATED(cftext_out))) ALLOCATE(cftext_out(nbcf_out))
     731    IF ((nbcf.GT.0) .AND. (.NOT. ALLOCATED(cfmod1))) ALLOCATE(cfmod1(nbcf))
     732    IF ((nbcf.GT.0) .AND. (.NOT. ALLOCATED(cfmod2))) ALLOCATE(cfmod2(nbcf))
     733    IF ((nbcf_in.GT.0) .AND. (.NOT. ALLOCATED(cfunits_in))) ALLOCATE(cfunits_in(nbcf_in))
     734    IF ((nbcf_out.GT.0) .AND. (.NOT. ALLOCATED(cfunits_out))) ALLOCATE(cfunits_out(nbcf_out))
     735       
     736    IF (is_mpi_root .AND. is_omp_root) THEN
     737
     738        IF (nbcf.GT.0) cfname=cfname_root
     739        IF (nbcf_in.GT.0) cfname_in=PACK(cfname_root,mask_in_root)
     740        IF (nbcf_out.GT.0) cfname_out=PACK(cfname_root,mask_out_root)
     741        IF (nbcf_in.GT.0) cftext_in=PACK(cftext_root,mask_in_root)
     742        IF (nbcf_out.GT.0) cftext_out=PACK(cftext_root,mask_out_root)
     743        IF (nbcf.GT.0) cfmod1=cfmod1_root
     744        IF (nbcf.GT.0) cfmod2=cfmod2_root
     745        IF (nbcf_in.GT.0) cfunits_in=PACK(cfunits_root,mask_in_root)
     746        IF (nbcf_out.GT.0) cfunits_out=PACK(cfunits_root,mask_out_root)
     747
     748        nbcf_in_orc=0
     749        nbcf_in_nemo=0
     750        nbcf_in_inca=0
     751        nbcf_in_ant=0
     752
     753        DO iq=1,nbcf
     754            IF (cfmod1(iq) == "ORC") nbcf_in_orc=nbcf_in_orc+1 
     755            IF (cfmod1(iq) == "NEMO") nbcf_in_nemo=nbcf_in_nemo+1 
     756            IF (cfmod1(iq) == "INCA") nbcf_in_inca=nbcf_in_inca+1
     757            IF (cfmod1(iq) == "ALL") nbcf_in_orc=nbcf_in_orc+1    ! ALL = ORC/NEMO/INCA
     758            IF (cfmod1(iq) == "ALL") nbcf_in_nemo=nbcf_in_nemo+1  ! ALL = ORC/NEMO/INCA
     759            IF (cfmod1(iq) == "ALL") nbcf_in_inca=nbcf_in_inca+1  ! ALL = ORC/NEMO/INCA
     760            IF (cfmod1(iq) == "ANT") nbcf_in_ant=nbcf_in_ant+1 
     761        ENDDO
     762
     763    ENDIF !   (is_mpi_root .AND. is_omp_root)
     764!$OMP BARRIER
     765
     766    CALL bcast(nbcf_in_orc)
     767    CALL bcast(nbcf_in_nemo)
     768    CALL bcast(nbcf_in_inca)
     769    CALL bcast(nbcf_in_ant)
     770
     771    WRITE(lunout,*) 'nbcf_in_orc ',nbcf_in_orc
     772    WRITE(lunout,*) 'nbcf_in_nemo ',nbcf_in_nemo
     773    WRITE(lunout,*) 'nbcf_in_inca ',nbcf_in_inca
     774    WRITE(lunout,*) 'nbcf_in_ant ',nbcf_in_ant
     775
     776    IF (nbcf_in.GT.0) THEN
     777        DO iq=1,nbcf_in
     778          CALL bcast(cfname_in(iq))
     779          CALL bcast(cftext_in(iq))
     780          CALL bcast(cfunits_in(iq))
     781        ENDDO
     782    ENDIF
     783
     784    IF (nbcf_out.GT.0) THEN
     785        DO iq=1,nbcf_out
     786          CALL bcast(cfname_out(iq))
     787          CALL bcast(cftext_out(iq))
     788          CALL bcast(cfunits_out(iq))
     789        ENDDO
     790    ENDIF
     791
     792    IF (nbcf.GT.0) THEN
     793        DO iq=1,nbcf
     794          CALL bcast(cfmod1(iq))
     795          CALL bcast(cfmod2(iq))
     796        ENDDO
     797    ENDIF
     798
     799    IF (nbcf_in.GT.0) WRITE(lunout,*)'infocfields_mod --- cfname_in: ',cfname_in
     800    IF (nbcf_out.GT.0) WRITE(lunout,*)'infocfields_mod --- cfname_out: ',cfname_out
     801
     802    IF (nbcf_in.GT.0) WRITE(lunout,*)'infocfields_mod --- cftext_in: ',cftext_in
     803    IF (nbcf_out.GT.0) WRITE(lunout,*)'infocfields_mod --- cftext_out: ',cftext_out
     804
     805    IF (nbcf.GT.0) WRITE(lunout,*)'infocfields_mod --- cfmod1: ',cfmod1
     806    IF (nbcf.GT.0) WRITE(lunout,*)'infocfields_mod --- cfmod2: ',cfmod2
     807
     808    IF (nbcf_in.GT.0) WRITE(lunout,*)'infocfunits_mod --- cfunits_in: ',cfunits_in
     809    IF (nbcf_out.GT.0) WRITE(lunout,*)'infocfunits_mod --- cfunits_out: ',cfunits_out
     810
     811    IF (nbcf_in.GT.0) WRITE(*,*)'infocfields_init --- number of fields in to LMDZ: ',nbcf_in
     812    IF (nbcf_out.GT.0) WRITE(*,*)'infocfields_init --- number of fields out of LMDZ: ',nbcf_out
     813
     814 ELSE
     815 ! Default values for other planets
     816    nbcf=0
     817    nbcf_in=0
     818    nbcf_out=0
     819 ENDIF ! planet_type
     820
     821 IF ((nbcf_in.GT.0) .AND. (.NOT. ALLOCATED(zcfields_in))) ALLOCATE(zcfields_in(klon,nbcf_in),stat=error)
     822        IF (error /= 0)  CALL abort_gcm(modname,'Pb in allocation zcfields_in',1)
     823 IF ((nbcf_out.GT.0) .AND. (.NOT. ALLOCATED(zcfields_out))) ALLOCATE(zcfields_out(klon,nbcf_out),stat=error)
     824        IF (error /= 0)  CALL abort_gcm(modname,'Pb in allocation zcfields_out',1)
     825
     826END SUBROUTINE infocfields_init
     827
    485828END MODULE carbon_cycle_mod
  • LMDZ6/trunk/libf/phylmd/physiq_mod.F90

    r3370 r3387  
    1616       d_u, d_v, d_t, d_qx, d_ps)
    1717
    18     use assert_m, only: assert
     18    USE assert_m, only: assert
    1919    USE ioipsl, only: histbeg, histvert, histdef, histend, histsync, &
    2020         histwrite, ju2ymds, ymds2ju, getin
     
    168168       omega,  &
    169169       epmax_diag,  &
    170 
    171170       !    Deep convective variables used in phytrac
    172171       pmflxr, pmflxs,  &
     
    184183       epmlmMm, eplaMm, &
    185184       sij, &
    186 
     185       !
    187186       cldemi,  &
    188187       cldfra, cldtau, fiwc,  &
     
    199198       prfl, psfl, fraca, Vprecip,  &
    200199       zw2,  &
    201        
     200       !
    202201       fluxu, fluxv,  &
    203202       fluxt,  &
    204 
     203       !
    205204       uwriteSTD, vwriteSTD, &                !pour calcul_STDlev.h
    206205       wwriteSTD, phiwriteSTD, &              !pour calcul_STDlev.h
    207206       qwriteSTD, twriteSTD, rhwriteSTD, &    !pour calcul_STDlev.h
    208        
     207       !
    209208       beta_prec,  &
    210209       rneb,  &
     
    213212    USE phys_state_var_mod ! Variables sauvegardees de la physique
    214213#ifdef CPP_Dust
    215   USE phys_output_write_spl_mod
     214    USE phys_output_write_spl_mod
    216215#else
    217216    USE phys_output_var_mod ! Variables pour les ecritures des sorties
     
    222221    USE phys_output_mod
    223222    USE phys_output_ctrlout_mod
    224     use open_climoz_m, only: open_climoz ! ozone climatology from a file
    225     use regr_pr_time_av_m, only: regr_pr_time_av
    226     use netcdf95, only: nf95_close
     223    USE open_climoz_m, only: open_climoz ! ozone climatology from a file
     224    USE regr_pr_time_av_m, only: regr_pr_time_av
     225    USE netcdf95, only: nf95_close
    227226    !IM for NMC files
    228     !     use netcdf, only: nf90_fill_real
    229     use netcdf, only: nf90_fill_real
    230     use mod_phys_lmdz_mpi_data, only: is_mpi_root
     227    USE netcdf, only: nf90_fill_real
     228    USE mod_phys_lmdz_mpi_data, only: is_mpi_root
    231229    USE aero_mod
    232     use ozonecm_m, only: ozonecm ! ozone of J.-F. Royer
    233     use conf_phys_m, only: conf_phys
    234     use radlwsw_m, only: radlwsw
    235     use phyaqua_mod, only: zenang_an
     230    USE ozonecm_m, only: ozonecm ! ozone of J.-F. Royer
     231    USE conf_phys_m, only: conf_phys
     232    USE radlwsw_m, only: radlwsw
     233    USE phyaqua_mod, only: zenang_an
    236234    USE time_phylmdz_mod, only: day_step_phy, annee_ref, day_ref, itau_phy, &
    237235         start_time, pdtphys, day_ini
     
    246244    USE indice_sol_mod
    247245    USE phytrac_mod, ONLY : phytrac
     246    USE carbon_cycle_mod, ONLY : infocfields_init
    248247
    249248#ifdef CPP_RRTM
     
    263262    !IM stations CFMIP
    264263    USE CFMIP_point_locations
    265     use FLOTT_GWD_rando_m, only: FLOTT_GWD_rando
    266     use ACAMA_GWD_rando_m, only: ACAMA_GWD_rando
     264    USE FLOTT_GWD_rando_m, only: FLOTT_GWD_rando
     265    USE ACAMA_GWD_rando_m, only: ACAMA_GWD_rando
    267266    USE VERTICAL_LAYERS_MOD, ONLY: aps,bps
    268267
    269 
    270     IMPLICIT none
     268    IMPLICIT NONE
    271269    !>======================================================================
    272270    !!
     
    13051303       nvm_lmdz = 13
    13061304       CALL getin_p('NVM',nvm_lmdz)
     1305
     1306       !--PC: defining fields to be exchanged between LMDz, ORCHIDEE and NEMO
     1307       WRITE(lunout,*) 'Call to infocfields from physiq'
     1308       CALL infocfields_init
     1309
    13071310    ENDIF
    13081311
    13091312    IF (prt_level.ge.1) print *,'CONVERGENCE PHYSIQUE THERM 1 '
    1310 
    13111313
    13121314    !======================================================================
Note: See TracChangeset for help on using the changeset viewer.