Ignore:
Timestamp:
Apr 3, 2014, 9:09:47 AM (11 years ago)
Author:
emillour
Message:

Generic model:
Major cleanup, in order to ease the use of LMDZ.GENERIC with (parallel) dynamics
in LMDZ.COMMON: (NB: this will break LMDZ.UNIVERSAL, which should be thrashed
in the near future)

  • Updated makegcm_* scripts (and makdim) and added the "-full" (to enforce full recomputation of the model) option
  • In dyn3d: converted control.h to module control_mod.F90 and converted iniadvtrac.F to module infotrac.F90
  • Added module mod_const_mpi.F90 in dyn3d (not used in serial mode)
  • Rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallelism)
  • added created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the max and min of a field over the whole planet. This should be further imroved with computation of means (possibly area weighed), etc.

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/start2archive.F

    r993 r1216  
    1919c=======================================================================
    2020
     21      use infotrac, only: iniadvtrac, nqtot, tname
    2122      USE comsoil_h
    22       USE comgeomfi_h
    23 
     23      USE comgeomfi_h, ONLY: lati, long, area
     24!      use control_mod
     25      use comgeomphy, only: initcomgeomphy
    2426      implicit none
    2527
     
    3234#include "logic.h"
    3335#include "temps.h"
    34 #include "control.h"
     36!#include "control.h"
    3537#include "ener.h"
    3638
    3739#include "dimphys.h"
    3840#include "planete.h"
    39 #include"advtrac.h"
     41!#include"advtrac.h"
    4042#include "netcdf.inc"
    4143
     
    4850      REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm) ! vents covariants
    4951      REAL teta(ip1jmp1,llm)                    ! temperature potentielle
    50       REAL q(ip1jmp1,llm,nqmx)               ! champs advectes
     52      REAL,ALLOCATABLE :: q(:,:,:)   ! champs advectes
    5153      REAL pks(ip1jmp1)                      ! exner (f pour filtre)
    5254      REAL pk(ip1jmp1,llm)
     
    6365      REAL tsoil(ngridmx,nsoilmx) ! Soil temperature
    6466      REAL co2ice(ngridmx)      ! CO2 ice layer
    65       REAL q2(ngridmx,nlayermx+1),qsurf(ngridmx,nqmx)
     67      REAL q2(ngridmx,nlayermx+1)
     68      REAL,ALLOCATABLE :: qsurf(:,:)
    6669      REAL emis(ngridmx)
    6770      INTEGER start,length
     
    8386      REAL ithS(ip1jmp1,nsoilmx) ! Soil Thermal Inertia
    8487      REAL co2iceS(ip1jmp1)
    85       REAL q2S(ip1jmp1,llm+1),qsurfS(ip1jmp1,nqmx)
     88      REAL q2S(ip1jmp1,llm+1)
     89      REAL,ALLOCATABLE :: qsurfS(:,:)
    8690      REAL emisS(ip1jmp1)
    8791
     
    122126c-----------------------------------------------------------------------
    123127
     128      CALL defrun_new(99, .TRUE. )
    124129      grireg   = .TRUE.
    125130      cursor = 1 ! added by AS in dimphys.
     131
     132! initialize "serial/parallel" related stuff
     133      CALL init_phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
     134      call initcomgeomphy
    126135
    127136      ! ALLOCATE ARRAYS IN comgeomfi_h (usually done in inifis)
     
    134143c Lecture des donnees
    135144c=======================================================================
    136 ! Load tracer names:
    137       call iniadvtrac(nq,numvanle)
    138 ! Initialize global tracer indexes (stored in tracer.h)
    139 ! ... this has to be done before phyetat0
    140       call initracer(ngridmx,nqmx,tnom)
     145! Load tracer number and names:
     146      call iniadvtrac(nqtot,numvanle)
     147
     148! allocate arrays:
     149      allocate(q(ip1jmp1,llm,nqtot))
     150      allocate(qsurf(ngridmx,nqtot))
     151      allocate(qsurfS(ip1jmp1,nqtot))
    141152
    142153      fichnom = 'start.nc'
    143       CALL dynetat0(fichnom,nqmx,vcov,ucov,teta,q,masse,
     154      CALL dynetat0(fichnom,nqtot,vcov,ucov,teta,q,masse,
    144155     .       ps,phis,timedyn)
    145156
     
    173184      Lmodif=0
    174185
    175       CALL phyetat0 (ngridmx,fichnom,0,Lmodif,nsoilmx,nqmx,day_ini_fi,
     186      CALL phyetat0 (ngridmx,fichnom,0,Lmodif,nsoilmx,nqtot,day_ini_fi,
    176187     .      timefi,
    177188     .      tsurf,tsoil,emis,q2,qsurf,
     
    286297      call gr_fi_dyn(1,ngridmx,iip1,jjp1,emis,emisS)
    287298      call gr_fi_dyn(llm+1,ngridmx,iip1,jjp1,q2,q2S)
    288       call gr_fi_dyn(nqmx,ngridmx,iip1,jjp1,qsurf,qsurfS)
     299      call gr_fi_dyn(nqtot,ngridmx,iip1,jjp1,qsurf,qsurfS)
    289300      call gr_fi_dyn(llm,ngridmx,iip1,jjp1,cloudfrac,cloudfracS)
    290301      call gr_fi_dyn(1,ngridmx,iip1,jjp1,hice,hiceS)
     
    390401
    391402c-----------------------------------------------------------------------
    392 c Ecriture du champs  q  ( q[1,nqmx] )
    393 c-----------------------------------------------------------------------
    394       do iq=1,nqmx
    395         call write_archive(nid,ntime,tnom(iq),'tracer','kg/kg',
     403c Ecriture du champs  q  ( q[1,nqtot] )
     404c-----------------------------------------------------------------------
     405      do iq=1,nqtot
     406        call write_archive(nid,ntime,tname(iq),'tracer','kg/kg',
    396407     &         3,q(1,1,iq))
    397408      end do
    398409c-----------------------------------------------------------------------
    399 c Ecriture du champs  qsurf  ( qsurf[1,nqmx] )
    400 c-----------------------------------------------------------------------
    401       do iq=1,nqmx
    402         txt=trim(tnom(iq))//"_surf"
     410c Ecriture du champs  qsurf  ( qsurf[1,nqtot] )
     411c-----------------------------------------------------------------------
     412      do iq=1,nqtot
     413        txt=trim(tname(iq))//"_surf"
    403414        call write_archive(nid,ntime,txt,'Tracer on surface',
    404415     &  'kg.m-2',2,qsurfS(1,iq))
     
    448459      ierr=NF_CLOSE(nid)
    449460
     461      write(*,*) "start2archive: All is well that ends well."
     462
    450463      end
Note: See TracChangeset for help on using the changeset viewer.