Ignore:
Timestamp:
Nov 25, 2020, 7:24:05 PM (4 years ago)
Author:
yjaziri
Message:

Generic GCM:
Implementation of an option for a new reading process of "traceur.def"
Use "#ModernTrac-v1" flag as first line of "traceur.def" to use this option
Further details in "LMDZ.GENERIC/deftank/traceur.def.modern"
YJ + JVO

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/dyn3d_common/infotrac.F90

    r2307 r2436  
    115115    INTEGER :: ifils,ipere,generation ! CRisi
    116116    LOGICAL :: continu,nouveau_traceurdef
     117    LOGICAL :: moderntracdef=.false. ! JVO, YJ : For planets modern traceur.def
    117118    INTEGER :: IOstatus ! gestion de la retrocompatibilite de traceur.def
    118119    CHARACTER(len=30) :: tchaine   
    119120   
    120     character(len=80) :: line ! to store a line of text
     121    character(len=500) :: line ! to store a line of text
    121122 
    122123    character(len=*),parameter :: modname="infotrac_init"
     
    240241    ELSE  ! not Earth
    241242       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
    242        IF(ierr.EQ.0) THEN
     243       IF (ierr.EQ.0) THEN
    243244          WRITE(lunout,*) 'Open traceur.def : ok'
    244           READ(90,*) nqtrue
     245!! - Modif. by JVO and YJ for modern planetary traceur.def ---------------
     246          READ(90,'(A)') line
     247          IF (trim(line).ne.'#ModernTrac-v1') THEN ! Test modern traceur.def
     248             READ(line,*) nqtrue ! Try standard traceur.def
     249          ELSE
     250             moderntracdef = .true.
     251             DO
     252               READ(90,'(A)',iostat=ierr1) line
     253               IF (ierr1==0) THEN
     254                 IF (index(line,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header
     255                   READ(line,*) nbtr
     256                   WRITE(lunout,*) trim(modname),': Number of tracers in the physics: nbtr=',nbtr
     257                   EXIT
     258                   ! JVO 20 : Here nbtr corresponds to ALL what we consider as modern 'tracers' (species)
     259                   ! from physics, including non-dynamic species (radiative, compo..), and nqtrue will be set
     260                   ! below, according to the flags in traceur.def, by read and rewind.
     261                 ENDIF
     262               ELSE ! If pb, or if reached EOF without having found nbtr
     263                 CALL abort_gcm('infotrac_init','Unable to read numbers of tracers in traceur.def',1)
     264               ENDIF
     265             ENDDO
     266          ENDIF ! if modern or standard traceur.def
     267!! -----------------------------------------------------------------------
    245268       ELSE
    246269          WRITE(lunout,*) 'Problem in opening traceur.def'
     
    257280! Allocate variables depending on nqtrue and nbtr
    258281!
    259     ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue),tnom_transp(nqtrue))
     282    IF (.NOT.moderntracdef) ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue),tnom_transp(nqtrue))
    260283!
    261284!jyg<
     
    526549
    527550    ELSE  ! not Earth
    528        ! Other planets (for now); we have the same number of tracers
    529        ! in the dynamics than in the physics
    530        nbtr=nqtrue
    531        ! NB: Reading a traceur.def with isotopes remains to be done...
    532        IF(ierr.EQ.0) THEN
    533           ! Continue to read tracer.def
    534           DO iq=1,nqtrue
    535              !READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
    536             ! try to be smart when reading traceur.def
    537             read(90,'(80a)') line ! store the line from traceur.def
    538             ! if format is hadv, vadv, tnom_0, tnom_transp
    539             read(line,*,iostat=ierr1) hadv(iq),vadv(iq),tnom_0(iq),tnom_transp(iq)
    540             if (ierr1.ne.0) then
    541             ! assume format is hadv,vadv,tnom_0
    542              read(line,*,iostat=ierr2) hadv(iq),vadv(iq),tnom_0(iq)
    543              if (ierr2.ne.0) then
    544               ! maybe format is tnom0,hadv,vadv
    545               read(line,*,iostat=ierr3) tnom_0(iq),hadv(iq),vadv(iq)
    546               if (ierr3.ne.0) then
    547                 !assuming values of hadv et vadv
    548                 hadv(iq)=10
    549                 vadv(iq)=10
    550                 read(line,*, iostat=ierr4) tnom_0(iq), tnom_transp(iq)
    551                 if (ierr4.ne.0) then
    552                 ! assume only tnom0 is provided (hadv and vadv default to 10)
    553                   read(line,*) tnom_0(iq)
    554                   tnom_transp(iq)='air'
    555                 endif
    556               else
    557                 !format is tnom0,hadv,vadv
    558                 tnom_transp(iq)='air' ! no isotopes
    559               endif ! of if (ierr3.ne.0)
    560              else
    561                ! format is hadv,vadv,tnom_0
    562                tnom_transp(iq)='air' ! no isotopes
    563              endif ! of if (ierr2.ne.0)
    564             endif ! of if(ierr1.ne.0)
    565           END DO ! of DO iq=1,nqtrue
    566           CLOSE(90) 
    567        ELSE ! Without tracer.def
    568           hadv(1) = 10
    569           vadv(1) = 10
    570           tnom_0(1) = 'dummy'
    571           tnom_transp(1)='air'
    572        END IF
     551
     552       IF (.not.moderntracdef) THEN ! Standard traceur.def or no traceur.def file
     553
     554         ! Other planets (for now); we have the same number of tracers
     555         ! in the dynamics than in the physics
     556         nbtr=nqtrue
     557         ! NB: Reading a traceur.def with isotopes remains to be done...
     558         IF(ierr.EQ.0) THEN
     559            ! Continue to read tracer.def
     560            DO iq=1,nqtrue
     561               !READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
     562              ! try to be smart when reading traceur.def
     563              read(90,'(80a)') line ! store the line from traceur.def
     564              ! if format is hadv, vadv, tnom_0, tnom_transp
     565              read(line,*,iostat=ierr1) hadv(iq),vadv(iq),tnom_0(iq),tnom_transp(iq)
     566              if (ierr1.ne.0) then
     567              ! assume format is hadv,vadv,tnom_0
     568               read(line,*,iostat=ierr2) hadv(iq),vadv(iq),tnom_0(iq)
     569               if (ierr2.ne.0) then
     570                ! maybe format is tnom0,hadv,vadv
     571                read(line,*,iostat=ierr3) tnom_0(iq),hadv(iq),vadv(iq)
     572                if (ierr3.ne.0) then
     573                  !assuming values of hadv et vadv
     574                  hadv(iq)=10
     575                  vadv(iq)=10
     576                  read(line,*, iostat=ierr4) tnom_0(iq), tnom_transp(iq)
     577                  if (ierr4.ne.0) then
     578                  ! assume only tnom0 is provided (hadv and vadv default to 10)
     579                    read(line,*) tnom_0(iq)
     580                    tnom_transp(iq)='air'
     581                  endif
     582                else
     583                  !format is tnom0,hadv,vadv
     584                  tnom_transp(iq)='air' ! no isotopes
     585                endif ! of if (ierr3.ne.0)
     586               else
     587                 ! format is hadv,vadv,tnom_0
     588                 tnom_transp(iq)='air' ! no isotopes
     589               endif ! of if (ierr2.ne.0)
     590              endif ! of if(ierr1.ne.0)
     591            END DO ! of DO iq=1,nqtrue
     592            CLOSE(90) 
     593         ELSE ! Without tracer.def
     594            hadv(1) = 10
     595            vadv(1) = 10
     596            tnom_0(1) = 'dummy'
     597            tnom_transp(1)='air'
     598         END IF
    573599       
     600       ELSE ! Modern traceur.def (moderntracdef=true) - Authors : JVO, YJ - 2020
     601
     602         nqtrue = 0
     603         DO iq=1,nbtr ! This loops on ALL the 'tracers' including the non-dyn ones.
     604           READ(90,'(A)') line ! store the line from traceur.def
     605           IF (index(line,'is_dyn=0').ne.0) CYCLE ! Only if explicitly not dyn., otherwise default behaviour
     606           nqtrue = nqtrue + 1
     607         ENDDO
     608         WRITE(lunout,*) trim(modname),': Number of tracers in the dynamics: nqtrue=',nqtrue
     609         ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue),tnom_transp(nqtrue))
     610
     611         ! Now that we know the dyn. tracers and allocated variables, let's read again the file.
     612         REWIND(90)
     613         DO
     614           READ(90,'(A)') line ! no need to check on iostat, already done at first read
     615           IF (index(line,'#').ne.1) THEN ! Skip arbitary number of commented lines in the header
     616             READ(line,*) ! This is the line of nbtr (see above)
     617             EXIT
     618           ENDIF
     619         ENDDO
     620         iiq=0
     621         DO iq=1,nbtr ! This loops on ALL the 'tracers' including the non-dyn ones.
     622           
     623           READ(90,'(A)') line ! store the line from traceur.def
     624           IF (index(line,'is_dyn=0').ne.0) CYCLE ! Skip the non-dyn ones
     625           iiq = iiq +  1 ! This ensure that iiq is on 1:nqtrue and not 1:nbtr
     626           ! Name must be first parameter in this version, but all others can be in any order
     627           read(line,*) tnom_0(iiq)
     628           write(*,*)"infotrac_init: iq=",iiq,"noms(iq)=",trim(tnom_0(iiq))
     629
     630           if (index(line,'vadv=') /= 0) then
     631             read(line(index(line,'vadv='   )+len('vadv='):),*) vadv(iiq)
     632             write(*,*) ' Parameter value (traceur.def) : vadv=', vadv(iiq)
     633           else
     634             vadv(iiq)=10
     635             write(*,*) ' Parameter value (default)     : vadv=',vadv(iiq)
     636           end if
     637
     638           if (index(line,'hadv=') /= 0) then
     639             read(line(index(line,'hadv='   )+len('hadv='):),*) hadv(iiq)
     640             write(*,*) ' Parameter value (traceur.def) : hadv=', hadv(iiq)
     641           else
     642             hadv(iiq)=10
     643             write(*,*) ' Parameter value (default)     : hadv=',hadv(iiq)
     644           end if
     645
     646           if (index(line,'tnom_transp=' ) /= 0) then
     647             read(line(index(line,'tnom_transp='   )+len('tnom_transp='):),*) tnom_transp(iiq)
     648             write(*,*) ' Parameter value (traceur.def) : tnom_transp=', tnom_transp(iiq)
     649           else
     650             tnom_transp(iiq)='air'
     651             write(*,*) ' Parameter value (default)     : tnom_transp=',tnom_transp(iiq)
     652           end if
     653
     654         ENDDO
     655         CLOSE(90)
     656
     657       ENDIF ! if (moderntracdef)
     658
    574659       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
    575660       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
Note: See TracChangeset for help on using the changeset viewer.