c======================================================================= *=*=*=*= datareaddrs.html =*=*=*=*
SUBROUTINE datareaddrs

SUBROUTINE datareaddrs


      SUBROUTINE datareaddrs(relief,phisinit,alb,ith,
     .                    zmea,zstd,zsig,zgam,zthe)
c=======================================================================
c
c
c   Author: F. Hourdin      01/1997
c   -------
c
c   Object: To read data from Martian surface to use in a GCM
c   ------				from DRS file "surface.dat"
c
c
c   Arguments:
c   ----------
c
c     Inputs:
c     ------
c
c     Outputs:
c     --------
c
c=======================================================================
c   donnees ALBEDO, INERTIE THERMIQUE, RELIEF:
c
c       Ces donnees sont au format DRS dans le fichier "surface.dat"
c				(/dw/chourdin/MARS/data_mars/surface.dat)
c
c   360 valeurs en longitude (de -179.5 a 179.5)
c   180 valeurs en latitudes (de 89.5 a -89.5)
c
c   Pour les passer au format de la grille, on utilise "interp_horiz.F"
c
c   Il faut donc que ces donnees soient au format grille scalaire
c               (imold+1 jmold+1)
c       avec passage des coordonnees de la "boite" (rlonu, rlatv)
c
c   On prend imd (d pour donnees!)
c           imd = 360 avec copie de la 1ere valeur sur la imd+1
c                   (rlonud de -179 a -181)
c           jmd = 179
c                   (rlatvd de 89 a -89)
c=======================================================================

      implicit none

#include "dimensions.h"
#include "paramet.h"
#include "comgeom.h"
#include "comconst.h"
#include "drsdef.h"

c=======================================================================
c   Declarations:
C=======================================================================

      INTEGER	imd,jmd,imdp1,jmdp1
      parameter	(imd=360,jmd=179,imdp1=361,jmdp1=180)

      INTEGER	iimp1
      parameter	(iimp1=iim+1-1/iim)

      CHARACTER relief*3

      REAL		zdata(imd*jmdp1)
      REAL		zdataS(imdp1*jmdp1)
      REAL		pfield(iimp1*jjp1)

      REAL		alb(iimp1*jjp1)
      REAL		ith(iimp1*jjp1)
      REAL		phisinit(iimp1*jjp1)

      REAL		zmea(imdp1*jmdp1)
      REAL		zstd(imdp1*jmdp1)
      REAL		zsig(imdp1*jmdp1)
      REAL		zgam(imdp1*jmdp1)
      REAL		zthe(imdp1*jmdp1)

      INTEGER   r4, size
      parameter (r4 = IDRS_BYTES_PER_WORD)

      INTEGER   aslun,setname,getdat,cluvdb
      INTEGER   lnblnk, ierr
      EXTERNAL	lnblnk

      INTEGER   unit
      data		unit/20/

      INTEGER	i,j,k

      INTEGER klatdat,klongdat
      PARAMETER (klatdat=180,klongdat=360)

c	on passe une grille en rlonu rlatv et im+1 jm a interp_horiz)

      REAL longitude(imd),latitude(jmdp1) ! Pour lecture des donnees
      REAL rlonud(imdp1),rlatvd(jmd)

      CHARACTER*20 string
      DIMENSION string(4)

      CHARACTER*80	file

#include "lmdstd.h"
#include "fxyprim.h"

      pi=2.*ASIN(1.)

c-----------------------------------------------------------------------
c Chemin pour trouver les donnees de surface (albedo, relief, th.inertia...)
c-----------------------------------------------------------------------

#ifdef CRAY
      file = '/u/rech/vnz/rvnz001/gcm/data_mars/surface'
#else
      file = '/dw/chourdin/gcm/MARS/data_mars/surface'
#endif

c=======================================================================
c	rlonud, rlatvd
c=======================================================================

c-----------------------------------------------------------------------
c	Lecture DRS des donnees latitude et longitude
c-----------------------------------------------------------------------
      write(*,*) 'ouverture du fichier surface'

      ierr = ASLUN(unit,file(1:lnblnk(file))//'.dic',
     .             unit+1,file(1:lnblnk(file))//'.dat',IDRS_READ)

      ierr = CLUVDB()
      ierr = SETNAME(' ','longitude',' ',' ',' ')
      size = (imd) * r4
      ierr = GETDAT(unit,longitude,size)

      ierr = CLUVDB()
      ierr = SETNAME(' ','latitude',' ',' ',' ')
      size = (jmdp1) * r4
      ierr = GETDAT(unit,latitude,size)

c-----------------------------------------------------------------------
c	Passage au format boites scalaires
c-----------------------------------------------------------------------

c-----------------------------------------------------------------------
c   	longitude(imd)		-->  	rlonud(imdp1)
c-----------------------------------------------------------------------

c Passage en coordonnees boites scalaires et en radian
      do i=1,imd
      	rlonud(i)=(longitude(i)+.5)*pi/180.
      enddo

c Repetition de la valeur im+1
      rlonud(imdp1)=rlonud(1) + 2*pi

c-----------------------------------------------------------------------
c		latitude(jmdp1)	 	-->		rlonvd(jmd)
c-----------------------------------------------------------------------

c Passage en coordonnees boites scalaires et en radian
      do j=1,jmd
      	rlatvd(j)=(latitude(j)-.5)*pi/180.
      enddo

c=======================================================================
c   lecture DRS de albedo, thermal, relief, zdtm (pour francois Lott)
c=======================================================================

      string(1) = 'albedo'
      string(2) = 'thermal'
      if (relief.ne.'pla') then
      	string(3) = 'z'//relief
      else
      	string(3) = 'zdtm'	! pour qu'il lise qqchose sur DRS
      						! remise a 0 derriere
      endif
      string(4) = 'zdtm'	! lecture pour francois Lott


      DO k=1,4
      	write(*,*) 'string',k,string(k)
      	
c-----------------------------------------------------------------------
c	initialisation
c-----------------------------------------------------------------------
      call initial0(iimp1*jjp1,pfield)
      call initial0(imd*jmdp1,zdata)
      call initial0(imdp1*jmdp1,zdataS)

c-----------------------------------------------------------------------
c	Lecture DRS
c-----------------------------------------------------------------------
      ierr = CLUVDB()
      ierr = SETNAME(' ',string(k),' ',' ',' ')
      size = imd * jmdp1 * r4
      ierr = GETDAT(unit,zdata, size)

c-----------------------------------------------------------------------
c		Cas particulier Francois Lott ( k=4 )
c-----------------------------------------------------------------------
      if (k.eq.4) then

      	call multscal(imd*jmdp1,zdata,1000.,zdata)
      	call multscal(imd,longitude,pi/180.,longitude)
      	call multscal(jmdp1,latitude,pi/180.,latitude)

      	call grid_noro(360, 180, longitude, latitude, zdata,
     .         iim, jjp1, rlonv, rlatu, zmea,zstd,zsig,zgam,zthe)

      CALL dump2d(iip1,jjp1,zmea,'zmea')
      CALL dump2d(iip1,jjp1,zstd,'zstd')
      CALL dump2d(iip1,jjp1,zsig,'zsig')
      CALL dump2d(iip1,jjp1,zgam,'zgam')
      CALL dump2d(iip1,jjp1,zthe,'zthe')

      endif

c-----------------------------------------------------------------------
c   Passage de zdata en grille (imdp1 jmdp1)
c-----------------------------------------------------------------------
      do j=1,jmdp1
      	do i=1,imd
      		zdataS(i+imdp1*(j-1)) = zdata(i+klongdat*(j-1))
      	enddo
      	zdataS(imdp1+imdp1*(j-1)) = zdata(1+klongdat*(j-1))
      enddo

c-----------------------------------------------------------------------
c	Interpolation
c-----------------------------------------------------------------------
      call interp_horiz(zdataS,pfield,imd,jmd,
     .	iim, jjm,1,rlonud,rlatvd,rlonu,rlatv)

c-----------------------------------------------------------------------
c	Periodicite	
c-----------------------------------------------------------------------

      do j=1,jjp1
         pfield(iimp1*j) =  pfield(1+iimp1*(j-1))
      enddo

c-----------------------------------------------------------------------
c	Sauvegarde des champs	
c-----------------------------------------------------------------------

      if (k.eq.1) then					! albedo
     	do i=1,iimp1*jjp1
      		alb(i) = pfield(i)
      	enddo
      elseif (k.eq.2) then				! thermal
     	do i=1,iimp1*jjp1
      		ith(i) = pfield(i)
      	enddo
      elseif (k.eq.3) then				! relief
      	if (relief.eq.'pla') then
      		call initial0(iimp1*jjp1,phisinit)
      	else
     		do i=1,iimp1*jjp1
      			phisinit(i) = pfield(i)
      		enddo
      	endif
      endif

      ENDDO

c-----------------------------------------------------------------------
c	Traitement Phisinit
c-----------------------------------------------------------------------

      DO i=1,iimp1*jjp1
            phisinit(i)=1000.*phisinit(i)
      ENDDO
      CALL dump2d(iimp1,jjp1,phisinit,'Altitude en m')
      CALL multscal(iimp1*jjp1,phisinit,g,phisinit)

c-----------------------------------------------------------------------
c	FIN
c-----------------------------------------------------------------------

      END