Changeset 988 for trunk/LMDZ.GENERIC/libf
- Timestamp:
- Jun 12, 2013, 9:53:44 AM (12 years ago)
- 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 1 1 c======================================================================= 2 SUBROUTINE datareadnc(relief, phisinit,alb,ith,2 SUBROUTINE datareadnc(relief,filename,phisinit,alb,ith, 3 3 . zmea,zstd,zsig,zgam,zthe) 4 4 c======================================================================= … … 63 63 parameter (iimp1=iim+1-1/iim) 64 64 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 67 76 REAL zdata(imd*jmdp1) 68 77 REAL zdataS(imdp1*jmdp1) 69 78 REAL pfield(iimp1*jjp1) 70 79 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 81 80 INTEGER ierr 82 81 … … 96 95 DIMENSION string(4) 97 96 98 CHARACTER*20 filename99 CHARACTER*50 filestring100 101 #include "lmdstd.h"102 97 #include "fxyprim.h" 103 98 … … 111 106 c Lecture NetCDF des donnees latitude et longitude 112 107 c----------------------------------------------------------------------- 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)') filename124 125 108 ierr = NF_OPEN (trim(datadir)//'/'//trim(adjustl(filename)), 126 109 & NF_NOWRITE,unit) … … 135 118 write(*,*)' can be obtained online on:' 136 119 write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile' 137 CALL ABORT120 STOP 138 121 ENDIF 139 122 -
trunk/LMDZ.GENERIC/libf/phystd/lect_start_archive.F
r985 r988 36 36 #include "ener.h" 37 37 #include "temps.h" 38 #include "lmdstd.h"39 38 #include "netcdf.inc" 40 39 #include"advtrac.h" -
trunk/LMDZ.GENERIC/libf/phystd/newstart.F
r985 r988 19 19 USE surfdat_h 20 20 USE comgeomfi_h 21 21 use datafile_mod, only: datadir 22 ! to use 'getin' 23 USE ioipsl_getincom, only: getin 22 24 implicit none 23 25 … … 34 36 #include "ener.h" 35 37 #include "temps.h" 36 #include "lmdstd.h"37 38 #include "comdissnew.h" 38 39 #include "clesph0.h" … … 130 131 real choix_1,pp 131 132 character*80 fichnom 133 character*250 filestring 132 134 integer Lmodif,iq,thermo 133 135 character modif*20 … … 146 148 INTEGER :: nq,numvanle 147 149 character(len=20) :: txt ! to store some text 150 character(len=50) :: surfacefile ! "surface.nc" (or equivalent file) 151 character(len=150) :: longline 148 152 integer :: count 149 153 real :: profile(llm+1) ! to store an atmospheric profile + surface value … … 437 441 c======================================================================= 438 442 439 if (choix_1. ne.1) then ! pour ne pas avoir besoin du fichier440 ! surface.dat dans le cas des start443 if (choix_1.eq.0) then ! for start_archive files, 444 ! where an external "surface.nc" file is needed 441 445 442 446 c do while((relief(1:3).ne.'mol').AND.(relief(1:3).ne.'pla')) … … 448 452 c enddo 449 453 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) 463 501 464 502
Note: See TracChangeset
for help on using the changeset viewer.