Ignore:
Timestamp:
Apr 14, 2010, 4:03:19 PM (14 years ago)
Author:
Ehouarn Millour
Message:

Some cleanup and fixing the possibility to output fields in the dynamics, on the dynamical grids.

CLEANUPS:

  • arch-PW6_VARGAS.fcm : add potentially benefic compiling options
  • removed obsolete "control.h" in dyn3d/dyn3dpar (module control_mod.F90 is used instead)

OUTPUTS in the dynamics (3 sets of files, one for each grid: scalar, u, v):

  • removed "com_io_dyn.h" common; use module "com_io_dyn_mod.F90" instead
  • updated "initdynav.F","inithist.F","writehist.F" and "writedynav.F" in bibio: which field will be written is hard coded there.
  • flags "ok_dyn_ins" and "ok_dyn_ave" (loaded via conf_gcm.F) trigger output of fields in the dynamics: if ok_dyn_ins is true, then files "dyn_hist.nc", "dyn_histu.nc" and "dyn_histv.nc" are written (the frequency of the outputs is given by 'iecri' in run.def; values are written every 'iecri' dynamical step). if ok_dyn_ave is true then files "dyn_hist_ave.nc", "dyn_histu_ave.nc" and "dyn_histv_ave.nc" are written (the rate at which averages and made/written, in days, is given by 'periodav' in run.def).

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ4/branches/LMDZ4V5.0-dev/libf/bibio/writedynav.F

    r1279 r1357  
    22! $Id$
    33!
    4       subroutine writedynav( histid, time, vcov,
    5      ,                          ucov,teta,ppk,phi,q,masse,ps,phis)
     4      subroutine writedynav(time, vcov,
     5     ,                ucov,teta,ppk,phi,q,masse,ps,phis)
    66
    77#ifdef CPP_IOIPSL
     
    99#endif
    1010      USE infotrac, ONLY : nqtot, ttext
     11      use com_io_dyn_mod, only : histaveid,histvaveid,histuaveid
    1112      implicit none
    1213
     
    1718C
    1819C   Entree:
    19 C      histid: ID du fichier histoire
    2020C      time: temps de l'ecriture
    2121C      vcov: vents v covariants
     
    2929C     
    3030C
    31 C   Sortie:
    32 C      fileid: ID du fichier netcdf cree
    3331C
    3432C   L. Fairhead, LMD, 03/99
     
    5351C
    5452
    55       INTEGER histid
    5653      REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm)
    57       REAL teta(ip1jmp1*llm),phi(ip1jmp1,llm),ppk(ip1jmp1*llm)                  
     54      REAL teta(ip1jmp1*llm),phi(ip1jmp1,llm),ppk(ip1jmp1*llm)     
    5855      REAL ps(ip1jmp1),masse(ip1jmp1,llm)                   
    5956      REAL phis(ip1jmp1)                 
     
    6663C   Variables locales
    6764C
    68       integer ndex2d(iip1*jjp1),ndex3d(iip1*jjp1*llm),iq, ii, ll
    69       real us(ip1jmp1*llm), vs(ip1jmp1*llm)
     65      integer ndex2d(ip1jmp1),ndexu(ip1jmp1*llm),ndexv(ip1jm*llm)
     66      INTEGER iq, ii, ll
    7067      real tm(ip1jmp1*llm)
    7168      REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
     
    7572C  Initialisations
    7673C
    77       ndex3d = 0
     74      ndexu = 0
     75      ndexv = 0
    7876      ndex2d = 0
    7977      ok_sync = .TRUE.
    80       us = 999.999
    81       vs = 999.999
    8278      tm = 999.999
    8379      vnat = 999.999
     
    9187C  Appels a histwrite pour l'ecriture des variables a sauvegarder
    9288C
    93 C  Vents U scalaire
     89C  Vents U
    9490C
    95       call gr_u_scal(llm, unat, us)
    96       call histwrite(histid, 'u', itau_w, us,
    97      .               iip1*jjp1*llm, ndex3d)
     91      call histwrite(histuaveid, 'u', itau_w, unat,
     92     .               iip1*jjp1*llm, ndexu)
    9893C
    99 C  Vents V scalaire
     94C  Vents V
    10095C
    101       call gr_v_scal(llm, vnat, vs)
    102       call histwrite(histid, 'v', itau_w, vs,
    103      .               iip1*jjp1*llm, ndex3d)
     96      call histwrite(histvaveid, 'v', itau_w, vnat,
     97     .               iip1*jjm*llm, ndexv)
    10498C
    10599C  Temperature potentielle moyennee
    106100C
    107       call histwrite(histid, 'theta', itau_w, teta,
    108      .                iip1*jjp1*llm, ndex3d)
     101      call histwrite(histaveid, 'theta', itau_w, teta,
     102     .                iip1*jjp1*llm, ndexu)
    109103C
    110104C  Temperature moyennee
     
    113107        tm(ii) = teta(ii) * ppk(ii)/cpp
    114108      enddo
    115       call histwrite(histid, 'temp', itau_w, tm,
    116      .                iip1*jjp1*llm, ndex3d)
     109      call histwrite(histaveid, 'temp', itau_w, tm,
     110     .                iip1*jjp1*llm, ndexu)
    117111C
    118112C  Geopotentiel
    119113C
    120       call histwrite(histid, 'phi', itau_w, phi,
    121      .                iip1*jjp1*llm, ndex3d)
     114      call histwrite(histaveid, 'phi', itau_w, phi,
     115     .                iip1*jjp1*llm, ndexu)
    122116C
    123117C  Traceurs
    124118C
    125         DO iq=1,nqtot
    126           call histwrite(histid, ttext(iq), itau_w, q(:,:,iq),
    127      .                   iip1*jjp1*llm, ndex3d)
    128         enddo
     119!        DO iq=1,nqtot
     120!          call histwrite(histaveid, ttext(iq), itau_w, q(:,:,iq),
     121!     .                   iip1*jjp1*llm, ndexu)
     122!        enddo
    129123C
    130124C  Masse
    131125C
    132        call histwrite(histid, 'masse', itau_w, masse, iip1*jjp1, ndex2d)
     126       call histwrite(histaveid, 'masse', itau_w, masse,
     127     $                   iip1*jjp1*llm, ndexu)
    133128C
    134129C  Pression au sol
    135130C
    136        call histwrite(histid, 'ps', itau_w, ps, iip1*jjp1, ndex2d)
     131       call histwrite(histaveid, 'ps', itau_w, ps, iip1*jjp1, ndex2d)
    137132C
    138133C  Geopotentiel au sol
    139134C
    140        call histwrite(histid, 'phis', itau_w, phis, iip1*jjp1, ndex2d)
     135!       call histwrite(histaveid,'phis',itau_w, phis,iip1*jjp1, ndex2d)
    141136C
    142137C  Fin
    143138C
    144       if (ok_sync) call histsync(histid)
     139      if (ok_sync) then
     140          call histsync(histaveid)
     141          call histsync(histvaveid)
     142          call histsync(histuaveid)
     143      ENDIF
    145144
    146145#else
Note: See TracChangeset for help on using the changeset viewer.