#define TRACEBACK write (gol,'("in ",a," (",a,i6,")")') rname, __FILE__, __LINE__ ; call goErr #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if ! #include "tm5.inc" ! ! ! Read Fortuin-Kelder ozone climatology ! module file_fortkeld implicit none ! --- in/out -------------------------- private public :: FortuinKelder_Info public :: FortuinKelder_Read ! --- const ----------------------------- character(len=*), parameter :: mname = 'file_fortkeld' ! latitudes -80,-70,...,70,80 integer, parameter :: nlat = 17 ! fixed pressure levels: integer, parameter :: nlev = 19 real, parameter :: pclim(nlev) = & (/1000.0, 700.0, 500.0, 300.0, 200.0, 150.0, & 100.0, 70.0, 50.0, 30.0, 20.0, 10.0, & 7.0, 5.0, 3.0, 2.0, 1.0, 0.5, 0.3 /) ! number of months: integer, parameter :: nmonth = 12 contains ! ===================================================================== subroutine FortuinKelder_Info( status, nlats, nlevs ) ! --- in/out ------------------------- integer, intent(inout) :: status integer, intent(out), optional :: nlats integer, intent(out), optional :: nlevs ! --- const --------------------------------- character(len=*), parameter :: rname = mname//'/FortuinKelder_Info' ! --- begin --------------------------- ! trap previous errors: if (status>0) return ! return info: if ( present(nlats) ) nlats = nlat if ( present(nlevs) ) nlevs = nlev ! ok status = 0 end subroutine FortuinKelder_Info ! *** subroutine FortuinKelder_Read( fname, imonth, lats_deg, pres_Pa, o3_ppmv, status ) use GO, only : gol, goErr, goPr use GO, only : goGetFU ! --- in/out -------------------------------------- character(len=*), intent(in) :: fname integer, intent(in) :: imonth real, intent(out) :: lats_deg(nlat) ! deg real, intent(out) :: pres_Pa(nlev) ! Pa real, intent(out) :: o3_ppmv(nlat,nlev) ! o3 volume mixing ratio integer, intent(inout) :: status ! --- const ------------------------------------- character(len=*), parameter :: rname = mname//'/FortuinKelder_Read' ! --- local -------------------------------------- integer :: j, l, imon logical :: exist integer :: fu real :: row(nlat) ! --- begin -------------------------------------- ! trap previous errors: if (status>0) return write (gol,'("read ozone climat for month ",i2)') imonth; call goPr ! fill latitudes do j = 1, nlat lats_deg(j) = -80.0 + (j-1)*10.0 ! deg end do ! pressure levels pres_Pa = pclim * 100.0 ! hPa -> Pa ! file exist ? inquire( file=fname, exist=exist ) if ( .not. exist ) then write (gol,'("file not found : ")'); call goErr write (gol,'(" ",a)') fname; call goErr write (gol,'("in ",a)') rname; call goErr; status=1; return end if ! select free file unit: call goGetFU( fu , status ) ! open file: open( fu, file=fname, iostat=status ) if ( status /= 0 ) then write (gol,'("while opening Fortuin climatology file:")'); call goErr write (gol,'(" ",a)') fname; call goErr write (gol,'("in ",a)') rname; call goErr; status=1; return end if ! read original fortuin data - ozone mixing ratios ! loop over months: do imon = 1, imonth ! skip empty line: read (fu,*) ! loop over pressure levels: do l = 1, nlev ! read row for this pressure level: read (fu,'(17f9.4)',iostat=status) row if ( status /= 0 ) then write (gol,'("reading from climatology file:")'); call goErr write (gol,'(" name : ",a)') fname; call goErr write (gol,'(" month : ",i2)') imon; call goErr write (gol,'(" level : ",i2)') l; call goErr write (gol,'("in ",a)') rname; call goErr; status=1; return end if ! store ? if ( imon == imonth ) o3_ppmv(:,l) = row end do end do ! close file: close( fu, iostat=status ) if ( status /= 0 ) then write (gol,'("while closing Fortuin climatology file:")'); call goErr write (gol,'(" ",a)') fname; call goErr write (gol,'("in ",a)') rname; call goErr; status=1; return end if ! ok status = 0 end subroutine FortuinKelder_Read end module file_fortkeld