Ignore:
Timestamp:
Jun 12, 2013, 9:53:44 AM (12 years ago)
Author:
emillour
Message:

Generic GCM:

  • Moved "newstart" (and related "lect_start_archive.F") to phystd directory
  • Adapted makegcm_* scripts to enable compiling main prog from physics
  • Added in newstart the possibility to not read in any surface.nc file (when loading a start_archive) with keyword "none" (instead of surface file name)
  • Some general cleanup:
    • in bibio: removed unused lmdstd.h readstd.F writestd.F mywrite.F

readcoord.F scatter.F gather.F ini36.F from36.F to36.F
lnblnk.F (F90 len_trim() should be used instead)

  • in dyn3d: removed unused inigrads.F wrgrads.F gradsdef.h

xvik.F (specific to current Mars)

EM

Location:
trunk/LMDZ.GENERIC/libf
Files:
15 deleted
1 edited
2 moved

Legend:

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

    r837 r988  
    11c=======================================================================
    2       SUBROUTINE datareadnc(relief,phisinit,alb,ith,
     2      SUBROUTINE datareadnc(relief,filename,phisinit,alb,ith,
    33     .                    zmea,zstd,zsig,zgam,zthe)
    44c=======================================================================
     
    6363      parameter    (iimp1=iim+1-1/iim)
    6464
    65       CHARACTER relief*3
    66 
     65      character(len=3),intent(inout) :: relief*3
     66      character(len=*),intent(in) :: filename ! surface.nc file
     67      real,intent(out) :: phisinit(iimp1*jjp1) ! surface geopotential
     68      real,intent(out) :: alb(iimp1*jjp1) ! albedo
     69      real,intent(out) :: ith(iimp1*jjp1) ! thermal inertia
     70      real,intent(out) :: zmea(imdp1*jmdp1)
     71      real,intent(out) :: zstd(imdp1*jmdp1)
     72      real,intent(out) :: zsig(imdp1*jmdp1)
     73      real,intent(out) :: zgam(imdp1*jmdp1)
     74      real,intent(out) :: zthe(imdp1*jmdp1)
     75     
    6776      REAL        zdata(imd*jmdp1)
    6877      REAL        zdataS(imdp1*jmdp1)
    6978      REAL        pfield(iimp1*jjp1)
    7079
    71       REAL        alb(iimp1*jjp1)
    72       REAL        ith(iimp1*jjp1)
    73       REAL        phisinit(iimp1*jjp1)
    74 
    75       REAL        zmea(imdp1*jmdp1)
    76       REAL        zstd(imdp1*jmdp1)
    77       REAL        zsig(imdp1*jmdp1)
    78       REAL        zgam(imdp1*jmdp1)
    79       REAL        zthe(imdp1*jmdp1)
    80 
    8180      INTEGER   ierr
    8281
     
    9695      DIMENSION string(4)
    9796
    98       CHARACTER*20 filename
    99       CHARACTER*50 filestring
    100 
    101 #include "lmdstd.h"
    10297#include "fxyprim.h"
    10398
     
    111106c    Lecture NetCDF des donnees latitude et longitude
    112107c-----------------------------------------------------------------------
    113 !    First get the corret value of dtadir, if not already done:
    114       ! default 'datadir' is set in "datadir_mod"
    115       call getin("datadir",datadir)
    116       write(*,*) 'Choice of surface data is:'
    117       filestring='ls '//trim(datadir)//'/*.nc'
    118       call system(filestring)
    119 
    120       write(*,*) ''
    121       write(*,*) 'Please choose the relevant file:'
    122       write(*,*) '(e.g. type "surface_mars.nc")'
    123       read(*,fmt='(a20)') filename
    124 
    125108      ierr = NF_OPEN (trim(datadir)//'/'//trim(adjustl(filename)),
    126109     &  NF_NOWRITE,unit)
     
    135118        write(*,*)' can be obtained online on:'
    136119        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
    137         CALL ABORT
     120        STOP
    138121      ENDIF
    139122
  • trunk/LMDZ.GENERIC/libf/phystd/lect_start_archive.F

    r985 r988  
    3636#include "ener.h"
    3737#include "temps.h"
    38 #include "lmdstd.h"
    3938#include "netcdf.inc"
    4039#include"advtrac.h"
  • trunk/LMDZ.GENERIC/libf/phystd/newstart.F

    r985 r988  
    1919      USE surfdat_h
    2020      USE comgeomfi_h
    21 
     21      use datafile_mod, only: datadir
     22! to use  'getin'
     23      USE ioipsl_getincom, only: getin
    2224      implicit none
    2325
     
    3436#include "ener.h"
    3537#include "temps.h"
    36 #include "lmdstd.h"
    3738#include "comdissnew.h"
    3839#include "clesph0.h"
     
    130131      real choix_1,pp
    131132      character*80      fichnom
     133      character*250  filestring
    132134      integer Lmodif,iq,thermo
    133135      character modif*20
     
    146148      INTEGER :: nq,numvanle
    147149      character(len=20) :: txt ! to store some text
     150      character(len=50) :: surfacefile ! "surface.nc" (or equivalent file)
     151      character(len=150) :: longline
    148152      integer :: count
    149153      real :: profile(llm+1) ! to store an atmospheric profile + surface value
     
    437441c=======================================================================
    438442
    439       if (choix_1.ne.1) then  ! pour ne pas avoir besoin du fichier
    440                               ! surface.dat dans le cas des start
     443      if (choix_1.eq.0) then  ! for start_archive files,
     444                              ! where an external "surface.nc" file is needed
    441445
    442446c do while((relief(1:3).ne.'mol').AND.(relief(1:3).ne.'pla'))
     
    448452c     enddo
    449453
    450       CALL datareadnc(relief,phis,alb,surfith,zmeaS,zstdS,zsigS,zgamS,
    451      .          ztheS)
    452 
    453       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
    454       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,surfith,surfithfi)
    455       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
    456       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zmeaS,zmea)
    457       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zstdS,zstd)
    458       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zsigS,zsig)
    459       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zgamS,zgam)
    460       CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ztheS,zthe)
    461 
    462       endif ! of if (choix_1.ne.1)
     454!    First get the correct value of datadir, if not already done:
     455        ! default 'datadir' is set in "datafile_mod"
     456        call getin("datadir",datadir)
     457        write(*,*) 'Available surface data files are:'
     458        filestring='ls '//trim(datadir)//' | grep .nc'
     459        call system(filestring)
     460
     461        write(*,*) ''
     462        write(*,*) 'Please choose the relevant file',
     463     &  ' (e.g. "surface_mars.nc")'
     464        write(*,*) ' or "none" to not use any (and set planetary'
     465        write(*,*) '  albedo and surface thermal inertia)'
     466        read(*,fmt='(a50)') surfacefile
     467
     468        if (surfacefile.ne."none") then
     469
     470          CALL datareadnc(relief,surfacefile,phis,alb,surfith,
     471     &          zmeaS,zstdS,zsigS,zgamS,ztheS)
     472        else
     473        ! specific case when not using a "surface.nc" file
     474          phis(:,:)=0
     475          zmeaS(:,:)=0
     476          zstdS(:,:)=0
     477          zsigS(:,:)=0
     478          zgamS(:,:)=0
     479          ztheS(:,:)=0
     480         
     481          write(*,*) "Enter value of albedo of the bare ground:"
     482          read(*,*) alb(1,1)
     483          alb(:,:)=alb(1,1)
     484         
     485          write(*,*) "Enter value of thermal inertia of soil:"
     486          read(*,*) surfith(1,1)
     487          surfith(:,:)=surfith(1,1)
     488       
     489        endif ! of if (surfacefile.ne."none")
     490
     491        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
     492        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,surfith,surfithfi)
     493        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
     494        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zmeaS,zmea)
     495        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zstdS,zstd)
     496        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zsigS,zsig)
     497        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zgamS,zgam)
     498        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ztheS,zthe)
     499
     500      endif ! of if (choix_1.eq.0)
    463501
    464502
Note: See TracChangeset for help on using the changeset viewer.