Changeset 374


Ignore:
Timestamp:
Nov 10, 2011, 3:08:51 PM (13 years ago)
Author:
emillour
Message:

Generic GCM: Upgrade: The location of the 'datagcm' directory can now be given in the callphys.def file ( datadir = /absolute/path/to/datagcm ). Changed "datafile.h" into a F90 module "datafile_mod.F90" and spread this change to all routines that used to use "datafile.h".
EM

Location:
trunk/LMDZ.GENERIC
Files:
13 edited
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/README

    r371 r374  
    541541  'create_make_gcm'.
    542542 
     543== 10/11/2011 == EM
     544- Upgrade: The location of the 'datagcm' directory can now be given in the
     545  callphys.def file ( datadir = /absolute/path/to/datagcm ). Changed
     546  "datafile.h" into a F90 module "datafile_mod.F90" and spread this change
     547  to all routines that used to use "datafile.h".
  • trunk/LMDZ.GENERIC/libf/dyn3d/newstart.F

    r367 r374  
    18361836!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    18371837      subroutine load_MONS_data(MONS_Hdn,MONS_d21)
     1838      use datafile_mod, only: datadir
    18381839      implicit none
    18391840      ! routine to load Benedicte Diez MONS dataset, fill in date in southern
     
    18411842#include"dimensions.h"
    18421843#include"paramet.h"
    1843 #include"datafile.h"
    18441844#include"comgeom.h"
    18451845      ! arguments:
     
    18711871
    18721872      ! Open MONS datafile:
    1873       open(42,file=trim(datafile)//"/"//trim(filename),
     1873      open(42,file=trim(datadir)//"/"//trim(filename),
    18741874     &     status="old",iostat=ierr)
    18751875      if (ierr/=0) then
    18761876        write(*,*) "Error in load_MONS_data:"
    18771877        write(*,*) "Failed opening file ",
    1878      &             trim(datafile)//"/"//trim(filename)
    1879         write(*,*)'1) You can change the path to the file in '
    1880         write(*,*)'   file phymars/datafile.h'
    1881         write(*,*)'2) If necessary ',trim(filename),
     1878     &             trim(datadir)//"/"//trim(filename)
     1879        write(*,*)'Check that your path to datagcm:',trim(datadir)
     1880        write(*,*)' is correct. You can change it in callphys.def with:'
     1881        write(*,*)' datadir = /absolute/path/to/datagcm'
     1882        write(*,*)'If necessary ',trim(filename),
    18821883     &                 ' (and other datafiles)'
    18831884        write(*,*)'   can be obtained online at:'
  • trunk/LMDZ.GENERIC/libf/phystd/ave_stelspec.F90

    r253 r374  
    2424      use radinc_h, only: L_NSPECTV
    2525      use radcommon_h, only: BWNV, DWNV, tstellar
     26      use datafile_mod, only: datadir
    2627
    2728      implicit none
    2829
    29 #include "datafile.h"
    3030#include "callkeys.h"
    3131
     
    4747      real lam_temp
    4848      double precision stel_temp
     49     
     50      integer :: ios ! file opening/reading status
    4951
    5052      STELLAR(:)=0.0
     
    5254      ! load high resolution wavenumber data
    5355      file_id='/stellar_spectra/lam.txt'
    54       file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     56      file_path=TRIM(datadir)//TRIM(file_id)
    5557
    56       open(110,file=file_path,form='formatted')
    57       do ifine=1,Nfine
    58          read(110,*) lam(ifine)
    59       enddo
    60       close(110)
     58      open(110,file=file_path,form='formatted',status='old',iostat=ios)
     59      if (ios.ne.0) then        ! file not found
     60        write(*,*) 'Error from ave_stelspec'
     61        write(*,*) 'Data file ',trim(file_id),' not found.'
     62        write(*,*)'Check that your path to datagcm:',trim(datadir)
     63        write(*,*)' is correct. You can change it in callphys.def with:'
     64        write(*,*)' datadir = /absolute/path/to/datagcm'
     65        write(*,*)' Also check that there is a ',trim(file_id),' there.'
     66        call abort
     67      else
     68        do ifine=1,Nfine
     69          read(110,*) lam(ifine)
     70        enddo
     71        close(110)
     72      endif
    6173
    6274      dl=lam(2)-lam(1)
     
    92104            call abort
    93105         endif
    94          file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
    95 
    96          open(111,file=file_path,form='formatted')
    97          do ifine=1,Nfine
    98             read(111,*) stel_f(ifine)
    99          enddo
    100          close(111)
     106         
     107         file_path=TRIM(datadir)//TRIM(file_id)
     108         open(111,file=file_path,form='formatted',status='old',iostat=ios)
     109         if (ios.ne.0) then        ! file not found
     110           write(*,*) 'Error from ave_stelspec'
     111           write(*,*) 'Data file ',trim(file_id),' not found.'
     112           write(*,*)'Check that your path to datagcm:',trim(datadir)
     113           write(*,*)' is correct. You can change it in callphys.def with:'
     114           write(*,*)' datadir = /absolute/path/to/datagcm'
     115           write(*,*)' Also check that there is a ',trim(file_id),' there.'
     116           call abort
     117         else
     118           do ifine=1,Nfine
     119             read(111,*) stel_f(ifine)
     120           enddo
     121           close(111)
     122         endif
    101123      endif
    102124
  • trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90

    r366 r374  
    1111      use radcommon_h
    1212      use watercommon_h
     13      use datafile_mod, only: datadir
    1314      use ioipsl_getincom
    1415
     
    3132
    3233#include "dimphys.h"
    33 #include "datafile.h"
    3434#include "comcstfi.h"
    3535#include "callkeys.h"
     
    233233         endif
    234234
    235          print*, "Correlated-k data base folder:",datafile
     235         print*, "callcorrk: Correlated-k data base folder:",trim(datadir)
    236236         call getin("corrkdir",corrkdir)
    237237         print*, "corrkdir = ",corrkdir
  • trunk/LMDZ.GENERIC/libf/phystd/datafile_mod.F90

    r371 r374  
    11!-----------------------------------------------------------------------
    2 ! INCLUDE datafile.h
     2      module datafile_mod
     3!  Address of the directory containing tables of data needed by the UCM
     4      implicit none
    35
    4 !  Address of the directory containing tables of data needed by the UCM
     6      ! Default for Gnome Idataplex:
     7      character(len=300) :: datadir='/san/home/rdword/gcm/datagcm'
     8      ! Default for LMD machines:
     9!      character(len=300) :: datadir='/u/rwlmd/datagcm'
    510
    6       character (len=300) :: datafile
    7 !      data datafile /'/home/rdword/gcm/datagcm_lite'/  ! laptop
    8       data datafile /'/san/home/rdword/gcm/datagcm'/  ! gnome iDataPlex
    9 !      data datafile /'/u/rwlmd/datagcm'/  ! planetary cluster
     11      end module datafile_mod
    1012!-----------------------------------------------------------------------
  • trunk/LMDZ.GENERIC/libf/phystd/datareadnc.F

    r135 r374  
    4242c=======================================================================
    4343
     44      use datafile_mod, only: datadir
    4445      implicit none
    4546
     
    4950#include "comconst.h"
    5051#include "netcdf.inc"
    51 #include "datafile.h"
    5252
    5353c=======================================================================
     
    112112
    113113      write(*,*) 'Choice of surface data is:'
    114       filestring='ls '//datafile(1:lnblnk(datafile))//'/*.nc'
     114      filestring='ls '//trim(datadir)//'/*.nc'
    115115      call system(filestring)
    116116
     
    120120      read(*,fmt='(a20)') filename
    121121
    122 !      ierr = NF_OPEN (datafile(1:lnblnk(datafile))//'/surface.nc',
    123 !     &  NF_NOWRITE,unit)
    124       ierr = NF_OPEN (datafile(1:lnblnk(datafile))//'/'//filename,
     122      ierr = NF_OPEN (trim(datadir)//'/'//trim(adjustl(filename)),
    125123     &  NF_NOWRITE,unit)
    126124      IF (ierr.NE.NF_NOERR) THEN
    127         write(*,*)'Error : cannot open file '//filename
     125        write(*,*)'Error : cannot open file '//trim(filename)
    128126        write(*,*)'(in phystd/datareadnc.F)'
    129         write(*,*)'It should be in :',datafile(1:lnblnk(datafile)),'/'
    130         write(*,*)'1) You can change this directory address in '
    131         write(*,*)'   file phystd/datafile.h'
    132         write(*,*)'2) If necessary surface.nc (and other datafiles)'
    133         write(*,*)'   can be obtained online on:'
     127        write(*,*)'It should be in :',trim(datadir),'/'
     128        write(*,*)'Check that your path to datagcm:',trim(datadir)
     129        write(*,*)' is correct. You can change it in callphys.def with:'
     130        write(*,*)' datadir = /absolute/path/to/datagcm'
     131        write(*,*)'If necessary surface.nc (and other datafiles)'
     132        write(*,*)' can be obtained online on:'
    134133        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
    135134        CALL ABORT
  • trunk/LMDZ.GENERIC/libf/phystd/inifis.F

    r305 r374  
    55
    66      use radinc_h, only : naerkind
    7 
     7      use datafile_mod, only: datadir
    88
    99!=======================================================================
     
    4242!   declarations:
    4343!   -------------
     44      use datafile_mod, only: datadir
    4445! to use  'getin'
    4546      USE ioipsl_getincom
     
    121122         PRINT*,'--------------------------------------------'
    122123
     124         write(*,*) "Directory where external input files are:"
     125         ! default 'datadir' is set in "datadir_mod"
     126         call getin("datadir",datadir) ! default path
     127         write(*,*) " datadir = ",trim(datadir)
     128
    123129         write(*,*) "Run with or without tracer transport ?"
    124130         tracer=.false. ! default value
  • trunk/LMDZ.GENERIC/libf/phystd/interpolateH2H2.F90

    r305 r374  
    11     subroutine interpolateH2H2(wn,temp,pres,abcoef,firstcall)
    2       implicit none
    32
    43!==================================================================
     
    1514!==================================================================
    1615
    17 #include "datafile.h"
     16      use datafile_mod, only: datadir
     17      implicit none
    1818
    1919      ! input
     
    5050
    5151         ! cold
    52          dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2H2_CIA_LT.dat'
     52         dt_file=TRIM(datadir)//'/continuum_data/H2H2_CIA_LT.dat'
    5353         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
    5454         if (ios.ne.0) then        ! file not found
    55             write(*,*) 'Error from interpolateH2H2.for'
    56             write(*,*) 'Data file could not be found:'
    57             write(*,*) dt_file
    58             call abort
     55           write(*,*) 'Error from interpolateH2H2'
     56           write(*,*) 'Data file ',trim(dt_file),' not found.'
     57           write(*,*)'Check that your path to datagcm:',trim(datadir)
     58           write(*,*)' is correct. You can change it in callphys.def with:'
     59           write(*,*)' datadir = /absolute/path/to/datagcm'
     60           write(*,*)' Also check that there is a continuum_data/H2H2_CIA_LT.dat there.'
     61           call abort
    5962         else
    6063            do k=1,nS
     
    6770
    6871         ! hot
    69          dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2H2_CIA_HT.dat'
     72         dt_file=TRIM(datadir)//'/continuum_data/H2H2_CIA_HT.dat'
    7073         open(34,file=dt_file,form='formatted',status='old',iostat=ios)
    7174         if (ios.ne.0) then        ! file not found
    72             write(*,*) 'Error from interpolateH2H2.for'
    73             write(*,*) 'Data file could not be found:'
    74             write(*,*) dt_file
    75             call abort
     75           write(*,*) 'Error from interpolateH2H2'
     76           write(*,*) 'Data file ',trim(dt_file),' not found.'
     77           write(*,*)'Check that your path to datagcm:',trim(datadir)
     78           write(*,*)' is correct. You can change it in callphys.def with:'
     79           write(*,*)' datadir = /absolute/path/to/datagcm'
     80           write(*,*)' Also check that there is a continuum_data/H2H2_CIA_HT.dat there.'
     81           call abort
    7682         else
    7783            do k=1,nS
     
    100106
    101107
    102          print*,'At wavenumber ',wn,' cm^-1'
     108         print*,'interpolateH2H2: At wavenumber ',wn,' cm^-1'
    103109         print*,'   temperature ',temp,' K'
    104110         print*,'   pressure ',pres,' Pa'
     
    172178     
    173179      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
    174          write(*,*) 'Warning from bilinear.for:'
     180         write(*,*) 'Warning from bilinearH2H2:'
    175181         write(*,*) 'Outside continuum temperature range!'
    176182         if(y.lt.y_arr(1))then
  • trunk/LMDZ.GENERIC/libf/phystd/interpolateH2Ocont.F90

    r305 r374  
    11     subroutine interpolateH2Ocont(wn,temp,presS,presF,abcoef,firstcall)
    2       implicit none
    32
    43!==================================================================
     
    1514!==================================================================
    1615
    17 #include "datafile.h"
     16      use datafile_mod, only: datadir
     17      implicit none
    1818
    1919      ! input
     
    5656
    5757         ! nu array
    58          dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2O_CONT_NU.dat'
     58         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_NU.dat'
    5959         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
    6060         if (ios.ne.0) then        ! file not found
    61             write(*,*) 'Error from interpolateH2O_cont.for'
    62             write(*,*) 'Data file could not be found:'
    63             write(*,*) dt_file
    64             call abort
     61           write(*,*) 'Error from interpolateH2O_cont'
     62           write(*,*) 'Data file ',trim(dt_file),' not found.'
     63           write(*,*)'Check that your path to datagcm:',trim(datadir)
     64           write(*,*)' is correct. You can change it in callphys.def with:'
     65           write(*,*)' datadir = /absolute/path/to/datagcm'
     66           write(*,*)' Also check that there is a continuum_data/H2O_CONT_NU.dat there.'
     67           call abort
    6568         else
    6669            do k=1,nS
     
    7174
    7275         ! self broadening
    73          dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2O_CONT_SELF.dat'
     76         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_SELF.dat'
    7477         open(34,file=dt_file,form='formatted',status='old',iostat=ios)
    7578         if (ios.ne.0) then        ! file not found
    76             write(*,*) 'Error from interpolateH2O_cont.for'
    77             write(*,*) 'Data file could not be found:'
    78             write(*,*) dt_file
    79             call abort
     79           write(*,*) 'Error from interpolateH2O_cont'
     80           write(*,*) 'Data file ',trim(dt_file),' not found.'
     81           write(*,*)'Check that your path to datagcm:',trim(datadir)
     82           write(*,*)' is correct. You can change it in callphys.def with:'
     83           write(*,*)' datadir = /absolute/path/to/datagcm'
     84           write(*,*)' Also check that there is a continuum_data/H2O_CONT_SELF.dat there.'
     85           call abort
    8086         else
    8187            do k=1,nS
     
    8793
    8894         ! foreign (N2+O2+Ar) broadening
    89          dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2O_CONT_FOREIGN.dat'
     95         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_FOREIGN.dat'
    9096         open(35,file=dt_file,form='formatted',status='old',iostat=ios)
    9197         if (ios.ne.0) then        ! file not found
    92             write(*,*) 'Error from interpolateH2O_cont.for'
    93             write(*,*) 'Data file could not be found:'
    94             write(*,*) dt_file
    95             call abort
     98           write(*,*) 'Error from interpolateH2O_cont'
     99           write(*,*) 'Data file ',trim(dt_file),' not found.'
     100           write(*,*)'Check that your path to datagcm:',trim(datadir)
     101           write(*,*)' is correct. You can change it in callphys.def with:'
     102           write(*,*)' datadir = /absolute/path/to/datagcm'
     103           write(*,*)' Also check that there is a continuum_data/H2O_CONT_FOREIGN.dat there.'
     104           call abort
    96105         else
    97106            do k=1,nS
     
    114123         temp_arr(11) = 700.
    115124
    116          print*,'At wavenumber ',wn,' cm^-1'
     125         print*,'interpolateH2Ocont: At wavenumber ',wn,' cm^-1'
    117126         print*,'   temperature ',temp,' K'
    118127         print*,'   H2O pressure ',presS,' Pa'
     
    210219     
    211220      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
    212          write(*,*) 'Warning from bilinear.for:'
     221         write(*,*) 'Warning from bilinearH2Ocont:'
    213222         write(*,*) 'Outside continuum temperature range!'
    214223         if(y.lt.y_arr(1))then
  • trunk/LMDZ.GENERIC/libf/phystd/setspi.F90

    r253 r374  
    2424      use radinc_h,    only: L_NSPECTI,corrkdir,banddir,NTstar,NTstop
    2525      use radcommon_h, only: BWNI,BLAMI,WNOI,DWNI,WAVEI,planckir,sigma
     26      use datafile_mod, only: datadir
    2627
    2728      implicit none
     
    2930#include "callkeys.h"
    3031#include "comcstfi.h"
    31 #include "datafile.h"
    3232
    3333      logical file_ok
     
    7474      !file_id='/corrk_data/' // corrkdir(1:LEN_TRIM(corrkdir)) // '/narrowbands_IR.in'
    7575      file_id='/corrk_data/'//trim(adjustl(banddir))//'/narrowbands_IR.in'
    76       file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     76      file_path=TRIM(datadir)//TRIM(file_id)
    7777
    7878      ! check that the file exists
    7979      inquire(FILE=file_path,EXIST=file_ok)
    8080      if(.not.file_ok) then
    81          write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     81         write(*,*)'The file ',TRIM(file_path)
    8282         write(*,*)'was not found by setspi.F90, exiting.'
     83         write(*,*)'Check that your path to datagcm:',trim(datadir)
     84         write(*,*)' is correct. You can change it in callphys.def with:'
     85         write(*,*)' datadir = /absolute/path/to/datagcm'
     86         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
    8387         call abort
    8488      endif
     
    9296
    9397      if(file_entries-1.ne.L_NSPECTI) then
    94          write(*,*) 'L_NSPECTI = ',L_NSPECTI, 'in the model, but there are '
     98         write(*,*) 'setspi: L_NSPECTI = ',L_NSPECTI, 'in the model, but there are '
    9599         write(*,*) file_entries-1,'entries in ', &
    96               file_path(1:LEN_TRIM(file_path)),', exiting.'
     100              TRIM(file_path),', exiting.'
    97101         call abort
    98102      endif
     
    110114
    111115      print*,''
    112       print*,'IR band limits:'
     116      print*,'setspi: IR band limits:'
    113117      do M=1,L_NSPECTI+1
    114118         print*,m,'-->',BWNI(M),' cm^-1'
     
    134138
    135139      print*,''
    136       print*,'Current Planck integration range:'
     140      print*,'setspi: Current Planck integration range:'
    137141      print*,'T = ',dble(NTstar)/1.0D+1, ' to ',dble(NTstop)/1.0D+1,' K.'
    138142
     
    157161      ! force planck=sigma*eps*T^4 for each temperature in array
    158162      if(forceEC)then
    159          print*,'Force F=sigma*eps*T^4 for all values of T!'
     163         print*,'setspi: Force F=sigma*eps*T^4 for all values of T!'
    160164         do nt=NTstar,NTstop
    161165            plancksum=0.0
     
    182186            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
    183187         end do
    184          print*,'At lower limit:'
     188         print*,'setspi: At lower limit:'
    185189         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
    186190         print*,'actual sig*T^4   = ',sigma*(dble(nt)/1.0D+1)**4,' W m^-2'
     
    192196            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
    193197         end do
    194          print*,'At upper limit:'
     198         print*,'setspi: At upper limit:'
    195199         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
    196200         print*,'actual sig*T^4   = ',sigma*(dble(nt)/1.0D+1)**4,' W m^-2'
  • trunk/LMDZ.GENERIC/libf/phystd/setspv.F90

    r253 r374  
    2626      use radcommon_h, only: BWNV,BLAMV,WNOV,DWNV,WAVEV, &
    2727                             STELLARF,TAURAY,gastype
     28      use datafile_mod, only: datadir
    2829
    2930      implicit none
     
    3132#include "comcstfi.h"
    3233#include "callkeys.h"
    33 #include "datafile.h"
    3434
    3535      logical file_ok
     
    5252      write(temp1,'(i2.2)') L_NSPECTV
    5353      file_id='/corrk_data/'//trim(adjustl(banddir))//'/narrowbands_VI.in'
    54       file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     54      file_path=TRIM(datadir)//TRIM(file_id)
    5555
    5656      ! check that the file exists
    5757      inquire(FILE=file_path,EXIST=file_ok)
    5858      if(.not.file_ok) then
    59          write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     59         write(*,*)'The file ',TRIM(file_path)
    6060         write(*,*)'was not found by setspv.F90, exiting.'
     61         write(*,*)'Check that your path to datagcm:',trim(datadir)
     62         write(*,*)' is correct. You can change it in callphys.def with:'
     63         write(*,*)' datadir = /absolute/path/to/datagcm'
     64         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
    6165         call abort
    6266      endif
     
    7074
    7175      if(file_entries-1.ne.L_NSPECTV) then
    72          write(*,*) 'L_NSPECTV = ',L_NSPECTV, 'in the model, but there are '
     76         write(*,*) 'setspv: L_NSPECTV = ',L_NSPECTV, 'in the model, but there are '
    7377         write(*,*) file_entries-1,'entries in ', &
    74               file_path(1:LEN_TRIM(file_path)),', exiting.'
     78              TRIM(file_path),', exiting.'
    7579         call abort
    7680      endif
     
    8892
    8993
    90       print*,'VI band limits:'
     94      print*,'setspv: VI band limits:'
    9195      do M=1,L_NSPECTV+1
    9296         print*,m,'-->',BWNV(M),' cm^-1'
     
    109113!     Set up stellar spectrum
    110114
    111       write(*,*)'Interpolating stellar spectrum from the hires data...'
     115      write(*,*)'setspv: Interpolating stellar spectrum from the hires data...'
    112116      call ave_stelspec(STELLAR)
    113117
     
    118122         sum         = sum+STELLARF(N)
    119123      end do
    120       write(6,'(" Stellar flux at 1 AU = ",f7.2," W m-2")') sum
     124      write(6,'("setspv: Stellar flux at 1 AU = ",f7.2," W m-2")') sum
    121125      print*,' '
    122126
     
    130134         call calc_rayleigh
    131135      else
    132          print*,'No Rayleigh scattering, check for NaN in output!'
     136         print*,'setspv: No Rayleigh scattering, check for NaN in output!'
    133137         do N=1,L_NSPECTV
    134138            TAURAY(N) = 1E-16
  • trunk/LMDZ.GENERIC/libf/phystd/suaer_corrk.F90

    r253 r374  
    44      use radinc_h,    only: L_NSPECTI,L_NSPECTV,naerkind,nsizemax,iim,jjm
    55      use radcommon_h, only: blamv,blami,lamrefir,lamrefvis
     6      use datafile_mod, only: datadir
    67
    78      ! outputs
     
    4344
    4445#include "callkeys.h"
    45 #include "datafile.h"
    4646
    4747!     Optical properties (read in external ASCII files)
     
    173173!     1.1 Open the ASCII file
    174174
    175             INQUIRE(FILE=datafile(1:LEN_TRIM(datafile))//&
    176          '/'//file_id(iaer,idomain)(1:LEN_TRIM(file_id(iaer,idomain))),&
     175            INQUIRE(FILE=TRIM(datadir)//&
     176         '/'//TRIM(file_id(iaer,idomain)),&
    177177            EXIST=file_ok)
    178178            IF(.NOT.file_ok) THEN
    179                write(*,*)'Problem opening ',&
    180                file_id(iaer,idomain)(1:LEN_TRIM(file_id(iaer,idomain)))
    181                write(*,*)'It should be in: ',&
    182                datafile(1:LEN_TRIM(datafile))
    183                write(*,*)'1) You can change this directory address in '
    184                write(*,*)' file phymars/datafile.h'
     179               write(*,*)'suaer_corrk: Problem opening ',&
     180               TRIM(file_id(iaer,idomain))
     181               write(*,*)'It should be in: ',TRIM(datadir)
     182               write(*,*)'1) You can set this directory address ',&
     183               'in callphys.def with:'
     184               write(*,*)' datadir = /absolute/path/to/datagcm'
    185185               write(*,*)'2) If ',&
    186               file_id(iaer,idomain)(1:LEN_TRIM(file_id(iaer,idomain))),&
     186              TRIM(file_id(iaer,idomain)),&
    187187               ' is a LMD reference datafile, it'
    188188               write(*,*)' can be obtained online at:'
    189189               write(*,*)' http://www.lmd.jussieu.fr/',&
    190190               '~forget/datagcm/datafile'
    191                write(*,*)'3) If the name of the file is wrong, you can'
    192                write(*,*)' change it in file phyuniv/suaer_corrk.F90. Just'
    193                write(*,*)' modify the variable called file_id.'
    194191               CALL ABORT
    195192            ENDIF
    196193            OPEN(UNIT=file_unit,&
    197             FILE=datafile(1:LEN_TRIM(datafile))//&
    198          '/'//file_id(iaer,idomain)(1:LEN_TRIM(file_id(iaer,idomain))),&
     194            FILE=TRIM(datadir)//&
     195         '/'//TRIM(file_id(iaer,idomain)),&
    199196            FORM='formatted')
    200197
  • trunk/LMDZ.GENERIC/libf/phystd/sugas_corrk.F90

    r368 r374  
    2323      use radcommon_h, only : gasv,gasi,FZEROI,FZEROV,gweight
    2424      use radcommon_h, only : wrefvar,gastype
    25 
     25      use datafile_mod, only: datadir
    2626      implicit none
    2727
    28 #include "datafile.h"
    2928#include "callkeys.h"
    3029#include "gases.h"
     
    3837
    3938      character(len=100) :: file_id
    40       character(len=100) :: file_path
     39      character(len=300) :: file_path
    4140
    4241      real*8 gasi8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTI,L_NGAUSS)
     
    5150!=======================================================================
    5251!     Load variable species data, exit if we have wrong database
    53       file_id='/corrk_data/' // corrkdir(1:LEN_TRIM(corrkdir)) // '/Q.dat'
    54       file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     52      file_id='/corrk_data/' // TRIM(corrkdir) // '/Q.dat'
     53      file_path=TRIM(datadir)//TRIM(file_id)
    5554
    5655
     
    5857      inquire(FILE=file_path,EXIST=file_ok)
    5958      if(.not.file_ok) then
    60          write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     59         write(*,*)'The file ',TRIM(file_path)
    6160         write(*,*)'was not found by sugas_corrk.F90, exiting.'
    62          write(*,*)'Check that your path to datagcm in datafile.h is correct and that'
    63          write(*,*)'the corrkdir you chose in callphys.def exists.'
     61         write(*,*)'Check that your path to datagcm:',trim(datadir)
     62         write(*,*)' is correct. You can change it in callphys.def with:'
     63         write(*,*)' datadir = /absolute/path/to/datagcm'
     64         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
    6465         call abort
    6566      endif
    6667
    6768      ! check that database matches varactive toggle
    68       open(111,file=file_path(1:LEN_TRIM(file_path)),form='formatted')
     69      open(111,file=TRIM(file_path),form='formatted')
    6970      read(111,*) ngas
    7071
     
    139140!     Set the weighting in g-space
    140141
    141       file_id='/corrk_data/' // corrkdir(1:LEN_TRIM(corrkdir)) // '/g.dat'
    142       file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     142      file_id='/corrk_data/' // TRIM(corrkdir) // '/g.dat'
     143      file_path=TRIM(datadir)//TRIM(file_id)
    143144
    144145      ! check that the file exists
    145146      inquire(FILE=file_path,EXIST=file_ok)
    146147      if(.not.file_ok) then
    147          write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     148         write(*,*)'The file ',TRIM(file_path)
    148149         write(*,*)'was not found by sugas_corrk.F90, exiting.'
    149          write(*,*)'Check that your path to datagcm in datafile.h is correct and that'
    150          write(*,*)'the corrkdir you chose in callphys.def exists.'
     150         write(*,*)'Check that your path to datagcm:',trim(datadir)
     151         write(*,*)' is correct. You can change it in callphys.def with:'
     152         write(*,*)' datadir = /absolute/path/to/datagcm'
     153         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
    151154         call abort
    152155      endif
    153156     
    154157      ! check the array size is correct, load the coefficients
    155       open(111,file=file_path(1:LEN_TRIM(file_path)),form='formatted')
     158      open(111,file=TRIM(file_path),form='formatted')
    156159      read(111,*) L_NGAUSScheck
    157160      if(.not.(L_NGAUSScheck.eq.L_NGAUSS)) then
     
    177180! pressure
    178181
    179       file_id='/corrk_data/' // corrkdir(1:LEN_TRIM(corrkdir)) // '/p.dat'
    180       file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     182      file_id='/corrk_data/' // TRIM(corrkdir) // '/p.dat'
     183      file_path=TRIM(datadir)//TRIM(file_id)
    181184
    182185      ! check that the file exists
    183186      inquire(FILE=file_path,EXIST=file_ok)
    184187      if(.not.file_ok) then
    185          write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     188         write(*,*)'The file ',TRIM(file_path)
    186189         write(*,*)'was not found by sugas_corrk.F90, exiting.'
     190         write(*,*)'Check that your path to datagcm:',trim(datadir)
     191         write(*,*)' is correct. You can change it in callphys.def with:'
     192         write(*,*)' datadir = /absolute/path/to/datagcm'
     193         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
    187194         call abort
    188195      endif
    189196     
    190197      ! check the array size is correct, load the coefficients
    191       open(111,file=file_path(1:LEN_TRIM(file_path)),form='formatted')
     198      open(111,file=TRIM(file_path),form='formatted')
    192199      read(111,*) L_NPREFcheck
    193200      if(.not.(L_NPREFcheck.eq.L_NPREF)) then
     
    222229! temperature
    223230
    224       file_id='/corrk_data/' // corrkdir(1:LEN_TRIM(corrkdir)) // '/T.dat'
    225       file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     231      file_id='/corrk_data/' // TRIM(corrkdir) // '/T.dat'
     232      file_path=TRIM(datadir)//TRIM(file_id)
    226233
    227234      ! check that the file exists
    228235      inquire(FILE=file_path,EXIST=file_ok)
    229236      if(.not.file_ok) then
    230          write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     237         write(*,*)'The file ',TRIM(file_path)
    231238         write(*,*)'was not found by sugas_corrk.F90, exiting.'
     239         write(*,*)'Check that your path to datagcm:',trim(datadir)
     240         write(*,*)' is correct. You can change it in callphys.def with:'
     241         write(*,*)' datadir = /absolute/path/to/datagcm'
     242         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
    232243         call abort
    233244      endif
    234245
    235246      ! check the array size is correct, load the coefficients
    236       open(111,file=file_path(1:LEN_TRIM(file_path)),form='formatted')
     247      open(111,file=TRIM(file_path),form='formatted')
    237248      read(111,*) L_NTREFcheck
    238249      if(.not.(L_NTREFcheck.eq.L_NTREF)) then
     
    258269
    259270      ! VISIBLE
    260       if (callgasvis.and..not.corrkdir(1:LEN_TRIM(corrkdir)).eq.'null') then
     271      if (callgasvis.and..not.TRIM(corrkdir).eq.'null') then
    261272         file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_VI.dat'
    262          file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     273         file_path=TRIM(datadir)//TRIM(file_id)
    263274         
    264275         ! check that the file exists
    265276         inquire(FILE=file_path,EXIST=file_ok)
    266277         if(.not.file_ok) then
    267             write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     278            write(*,*)'The file ',TRIM(file_path)
    268279            write(*,*)'was not found by sugas_corrk.F90.'
    269280            write(*,*)'Are you sure you have absorption data for these bands?'
     
    271282         endif
    272283         
    273          open(111,file=file_path(1:LEN_TRIM(file_path)),form='formatted')
     284         open(111,file=TRIM(file_path),form='formatted')
    274285         read(111,*) gasv8
    275286         close(111)
     
    281292
    282293      ! INFRA-RED
    283       if (.not.corrkdir(1:LEN_TRIM(corrkdir)).eq.'null') then
     294      if (.not.TRIM(corrkdir).eq.'null') then
    284295         file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_IR.dat'
    285          file_path=datafile(1:LEN_TRIM(datafile))//file_id(1:LEN_TRIM(file_id))
     296         file_path=TRIM(datadir)//TRIM(file_id)
    286297     
    287298         ! check that the file exists
    288299         inquire(FILE=file_path,EXIST=file_ok)
    289300         if(.not.file_ok) then
    290             write(*,*)'The file ',file_path(1:LEN_TRIM(file_path))
     301            write(*,*)'The file ',TRIM(file_path)
    291302            write(*,*)'was not found by sugas_corrk.F90.'
    292303            write(*,*)'Are you sure you have absorption data for these bands?'
     
    294305         endif
    295306         
    296          open(111,file=file_path(1:LEN_TRIM(file_path)),form='formatted')
     307         open(111,file=TRIM(file_path),form='formatted')
    297308         read(111,*) gasi8
    298309         close(111)
  • trunk/LMDZ.GENERIC/tmp

    r253 r374  
    33mklsequential
    44align
    5 common
     5all
     6static
Note: See TracChangeset for help on using the changeset viewer.