NCEPLIBS-sfcio  1.4.0
All Data Structures Files
sfcio_module.f
Go to the documentation of this file.
1 
397 module sfcio_module
398  implicit none
399  private
400 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
401 ! Public Variables
402  integer,parameter,public:: sfcio_lhead1=32
403  integer,parameter,public:: sfcio_intkind=4,sfcio_realkind=4,sfcio_dblekind=8
404  real(sfcio_realkind),parameter,public:: sfcio_realfill=-9999.
405  real(sfcio_dblekind),parameter,public:: sfcio_dblefill=sfcio_realfill
406 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
407 ! Public Types
408  type,public:: sfcio_head
409  character(sfcio_lhead1):: clabsfc=' '
410  real(sfcio_realkind):: fhour=0.
411  integer(sfcio_intkind):: idate(4)=(/0,0,0,0/),latb=0,lonb=0,lsoil=0,ivs=0
412  integer(sfcio_intkind):: irealf=1
413  integer(sfcio_intkind),allocatable:: lpl(:)
414  real(sfcio_realkind),allocatable:: zsoil(:)
415  end type
416  type,public:: sfcio_data
417  real(sfcio_realkind),pointer:: tsea(:,:)=>null()
418  real(sfcio_realkind),pointer:: smc(:,:,:)=>null()
419  real(sfcio_realkind),pointer:: sheleg(:,:)=>null()
420  real(sfcio_realkind),pointer:: stc(:,:,:)=>null()
421  real(sfcio_realkind),pointer:: tg3(:,:)=>null()
422  real(sfcio_realkind),pointer:: zorl(:,:)=>null()
423  real(sfcio_realkind),pointer:: cv(:,:)=>null()
424  real(sfcio_realkind),pointer:: cvb(:,:)=>null()
425  real(sfcio_realkind),pointer:: cvt(:,:)=>null()
426  real(sfcio_realkind),pointer:: alvsf(:,:)=>null()
427  real(sfcio_realkind),pointer:: alvwf(:,:)=>null()
428  real(sfcio_realkind),pointer:: alnsf(:,:)=>null()
429  real(sfcio_realkind),pointer:: alnwf(:,:)=>null()
430  real(sfcio_realkind),pointer:: slmsk(:,:)=>null()
431  real(sfcio_realkind),pointer:: vfrac(:,:)=>null()
432  real(sfcio_realkind),pointer:: canopy(:,:)=>null()
433  real(sfcio_realkind),pointer:: f10m(:,:)=>null()
434  real(sfcio_realkind),pointer:: t2m(:,:)=>null()
435  real(sfcio_realkind),pointer:: q2m(:,:)=>null()
436  real(sfcio_realkind),pointer:: vtype(:,:)=>null()
437  real(sfcio_realkind),pointer:: stype(:,:)=>null()
438  real(sfcio_realkind),pointer:: facsf(:,:)=>null()
439  real(sfcio_realkind),pointer:: facwf(:,:)=>null()
440  real(sfcio_realkind),pointer:: uustar(:,:)=>null()
441  real(sfcio_realkind),pointer:: ffmm(:,:)=>null()
442  real(sfcio_realkind),pointer:: ffhh(:,:)=>null()
443  real(sfcio_realkind),pointer:: hice(:,:)=>null()
444  real(sfcio_realkind),pointer:: fice(:,:)=>null()
445  real(sfcio_realkind),pointer:: tisfc(:,:)=>null()
446  real(sfcio_realkind),pointer:: tprcp(:,:)=>null()
447  real(sfcio_realkind),pointer:: srflag(:,:)=>null()
448  real(sfcio_realkind),pointer:: snwdph(:,:)=>null()
449  real(sfcio_realkind),pointer:: slc(:,:,:)=>null()
450  real(sfcio_realkind),pointer:: shdmin(:,:)=>null()
451  real(sfcio_realkind),pointer:: shdmax(:,:)=>null()
452  real(sfcio_realkind),pointer:: slope(:,:)=>null()
453  real(sfcio_realkind),pointer:: snoalb(:,:)=>null()
454  real(sfcio_realkind),pointer:: orog(:,:)=>null()
455  end type
456  type,public:: sfcio_dbta
457  real(sfcio_dblekind),pointer:: tsea(:,:)=>null()
458  real(sfcio_dblekind),pointer:: smc(:,:,:)=>null()
459  real(sfcio_dblekind),pointer:: sheleg(:,:)=>null()
460  real(sfcio_dblekind),pointer:: stc(:,:,:)=>null()
461  real(sfcio_dblekind),pointer:: tg3(:,:)=>null()
462  real(sfcio_dblekind),pointer:: zorl(:,:)=>null()
463  real(sfcio_dblekind),pointer:: cv(:,:)=>null()
464  real(sfcio_dblekind),pointer:: cvb(:,:)=>null()
465  real(sfcio_dblekind),pointer:: cvt(:,:)=>null()
466  real(sfcio_dblekind),pointer:: alvsf(:,:)=>null()
467  real(sfcio_dblekind),pointer:: alvwf(:,:)=>null()
468  real(sfcio_dblekind),pointer:: alnsf(:,:)=>null()
469  real(sfcio_dblekind),pointer:: alnwf(:,:)=>null()
470  real(sfcio_dblekind),pointer:: slmsk(:,:)=>null()
471  real(sfcio_dblekind),pointer:: vfrac(:,:)=>null()
472  real(sfcio_dblekind),pointer:: canopy(:,:)=>null()
473  real(sfcio_dblekind),pointer:: f10m(:,:)=>null()
474  real(sfcio_dblekind),pointer:: t2m(:,:)=>null()
475  real(sfcio_dblekind),pointer:: q2m(:,:)=>null()
476  real(sfcio_dblekind),pointer:: vtype(:,:)=>null()
477  real(sfcio_dblekind),pointer:: stype(:,:)=>null()
478  real(sfcio_dblekind),pointer:: facsf(:,:)=>null()
479  real(sfcio_dblekind),pointer:: facwf(:,:)=>null()
480  real(sfcio_dblekind),pointer:: uustar(:,:)=>null()
481  real(sfcio_dblekind),pointer:: ffmm(:,:)=>null()
482  real(sfcio_dblekind),pointer:: ffhh(:,:)=>null()
483  real(sfcio_dblekind),pointer:: hice(:,:)=>null()
484  real(sfcio_dblekind),pointer:: fice(:,:)=>null()
485  real(sfcio_dblekind),pointer:: tisfc(:,:)=>null()
486  real(sfcio_dblekind),pointer:: tprcp(:,:)=>null()
487  real(sfcio_dblekind),pointer:: srflag(:,:)=>null()
488  real(sfcio_dblekind),pointer:: snwdph(:,:)=>null()
489  real(sfcio_dblekind),pointer:: slc(:,:,:)=>null()
490  real(sfcio_dblekind),pointer:: shdmin(:,:)=>null()
491  real(sfcio_dblekind),pointer:: shdmax(:,:)=>null()
492  real(sfcio_dblekind),pointer:: slope(:,:)=>null()
493  real(sfcio_dblekind),pointer:: snoalb(:,:)=>null()
494  real(sfcio_dblekind),pointer:: orog(:,:)=>null()
495  end type
496 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
497 ! Public Subprograms
498  public sfcio_sropen,sfcio_swopen,sfcio_sclose,sfcio_srhead,sfcio_swhead
499  public sfcio_alhead,sfcio_aldata,sfcio_axdata,sfcio_srdata,sfcio_swdata
500  public sfcio_aldbta,sfcio_axdbta,sfcio_srdbta,sfcio_swdbta
501  public sfcio_srohdc,sfcio_swohdc
502  interface sfcio_srohdc
503  module procedure sfcio_srohdca,sfcio_srohdcb
504  end interface
505  interface sfcio_swohdc
506  module procedure sfcio_swohdca,sfcio_swohdcb
507  end interface
508 contains
509 !-------------------------------------------------------------------------------
510  subroutine sfcio_sropen(lu,cfname,iret)
511  implicit none
512  integer(sfcio_intkind),intent(in):: lu
513  character*(*),intent(in):: cfname
514  integer(sfcio_intkind),intent(out):: iret
515  integer ios
516 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
517  open(lu,file=cfname,form='unformatted',&
518  status='old',action='read',iostat=ios)
519  iret=ios
520  if(iret.ne.0) iret=-1
521 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
522  end subroutine
523 !-------------------------------------------------------------------------------
524  subroutine sfcio_swopen(lu,cfname,iret)
525  implicit none
526  integer(sfcio_intkind),intent(in):: lu
527  character*(*),intent(in):: cfname
528  integer(sfcio_intkind),intent(out):: iret
529  integer ios
530 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
531  open(lu,file=cfname,form='unformatted',&
532  status='unknown',action='readwrite',iostat=ios)
533  iret=ios
534  if(iret.ne.0) iret=-1
535 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
536  end subroutine
537 !-------------------------------------------------------------------------------
538  subroutine sfcio_sclose(lu,iret)
539  implicit none
540  integer(sfcio_intkind),intent(in):: lu
541  integer(sfcio_intkind),intent(out):: iret
542  integer ios
543 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
544  close(lu,iostat=ios)
545  iret=ios
546  if(iret.ne.0) iret=-1
547 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
548  end subroutine
549 !-------------------------------------------------------------------------------
550  subroutine sfcio_srhead(lu,head,iret)
551  implicit none
552  integer(sfcio_intkind),intent(in):: lu
553  type(sfcio_head),intent(out):: head
554  integer(sfcio_intkind),intent(out):: iret
555  integer:: ios
556  character(4):: cgfs,csfc
557  integer(sfcio_intkind):: nhead,ndata,nresv(3)
558 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
559  iret=-2
560  rewind lu
561  read(lu,iostat=ios) head%clabsfc
562  if(ios.ne.0) return
563  if(head%clabsfc(1:8).eq.'GFS SFC ' .or. head%clabsfc(1:8).eq.' SFG CFS') then ! modern surface file
564 !port if(head%clabsfc(1:8).eq.'GFS SFC ') then ! modern surface file
565 ! Commented by Moorthi
566 ! if(head%clabsfc(1:8).eq.'GFS SFC ') then
567 ! print *,' MODERN BIG ENDIAN SURFACE FILE'
568 ! else
569 ! print *,' MODERN LITTLE ENDIAN SURFACE FILE'
570 ! endif
571 
572  rewind lu
573  read(lu,iostat=ios) cgfs,csfc,head%ivs,nhead,ndata,nresv
574  if(ios.ne.0) return
575  if(head%ivs.eq.200509) then
576  read(lu,iostat=ios)
577  if(ios.ne.0) return
578  read(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,&
579  head%lsoil,head%irealf
580  if(ios.ne.0) return
581  call sfcio_alhead(head,ios)
582  if(ios.ne.0) return
583  read(lu,iostat=ios) head%lpl
584  if(ios.ne.0) return
585  read(lu,iostat=ios) head%zsoil
586  if(ios.ne.0) return
587  elseif(head%ivs.eq.200501) then
588  read(lu,iostat=ios)
589  if(ios.ne.0) return
590  read(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,head%lsoil
591  if(ios.ne.0) return
592  call sfcio_alhead(head,ios)
593  if(ios.ne.0) return
594  read(lu,iostat=ios) head%lpl
595  if(ios.ne.0) return
596  read(lu,iostat=ios) head%zsoil
597  if(ios.ne.0) return
598  else
599  return
600  endif
601  else
602  read(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,head%ivs
603  if(ios.ne.0) return
604  if(head%ivs.eq.199802) then
605  head%lsoil=2
606  call sfcio_alhead(head,ios)
607  if(ios.ne.0) return
608  head%lpl=head%lonb
609  head%zsoil=(/-0.1,-2.0/)
610  elseif(head%ivs.eq.200004) then
611  head%lsoil=2
612  call sfcio_alhead(head,ios)
613  if(ios.ne.0) return
614  rewind lu
615  read(lu) head%clabsfc
616  read(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,head%ivs,&
617  head%lpl
618  if(ios.ne.0) return
619  head%zsoil=(/-0.1,-2.0/)
620  elseif(head%ivs.eq.200412) then
621  head%lsoil=4
622  call sfcio_alhead(head,ios)
623  if(ios.ne.0) return
624  rewind lu
625  read(lu) head%clabsfc
626  read(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,head%ivs,&
627  head%lpl
628  if(ios.ne.0) return
629  head%zsoil=(/-0.1,-0.4,-1.0,-2.0/)
630  else
631  return
632  endif
633  iret=0
634  endif
635  iret=0
636 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
637  end subroutine
638 !-------------------------------------------------------------------------------
639  subroutine sfcio_swhead(lu,head,iret)
640  implicit none
641  integer(sfcio_intkind),intent(in):: lu
642  type(sfcio_head),intent(in):: head
643  integer(sfcio_intkind),intent(out):: iret
644  integer:: ios
645  integer i
646 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
647  iret=-2
648  if(head%ivs.eq.200509) then
649  rewind lu
650  write(lu,iostat=ios) 'GFS SFC ',head%ivs,5,29+3*head%lsoil,0,0,0
651  if(ios.ne.0) return
652  write(lu,iostat=ios) 4*(/8,5+29+3*head%lsoil,25,head%latb/2,head%lsoil/),&
653  4*head%irealf*(/(head%lonb*head%latb,&
654  i=1,29+3*head%lsoil)/)
655  if(ios.ne.0) return
656  write(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,&
657  head%lsoil,head%irealf,(0,i=1,16)
658  if(ios.ne.0) return
659  write(lu,iostat=ios) head%lpl
660  if(ios.ne.0) return
661  write(lu,iostat=ios) head%zsoil
662  if(ios.ne.0) return
663  iret=0
664  elseif(head%ivs.eq.200501) then
665  rewind lu
666  write(lu,iostat=ios) 'GFS SFC ',head%ivs,5,29+3*head%lsoil,0,0,0
667  if(ios.ne.0) return
668  write(lu,iostat=ios) 4*(/8,5+29+3*head%lsoil,8,head%latb/2,head%lsoil/),&
669  4*(/(head%lonb*head%latb,i=1,29+3*head%lsoil)/)
670  if(ios.ne.0) return
671  write(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,head%lsoil
672  if(ios.ne.0) return
673  write(lu,iostat=ios) head%lpl
674  if(ios.ne.0) return
675  write(lu,iostat=ios) head%zsoil
676  if(ios.ne.0) return
677  iret=0
678  elseif(head%ivs.eq.200004.and.head%lsoil.eq.2) then
679  rewind lu
680  write(lu,iostat=ios) head%clabsfc
681  if(ios.ne.0) return
682  write(lu,iostat=ios) head%fhour,head%idate,head%lonb,head%latb,head%ivs,&
683  head%lpl
684  if(ios.ne.0) return
685  iret=0
686  endif
687 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
688  end subroutine
689 !-------------------------------------------------------------------------------
690  subroutine sfcio_alhead(head,iret,latb,lsoil)
691  implicit none
692  type(sfcio_head),intent(inout):: head
693  integer(sfcio_intkind),intent(out):: iret
694  integer(sfcio_intkind),optional,intent(in):: latb,lsoil
695  integer dim1l,dim1z
696 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
697  if(present(latb)) then
698  dim1l=(latb+1)/2
699  else
700  dim1l=(head%latb+1)/2
701  endif
702  if(present(lsoil)) then
703  dim1z=lsoil
704  else
705  dim1z=head%lsoil
706  endif
707  if(allocated(head%lpl)) deallocate(head%lpl)
708  if(allocated(head%zsoil)) deallocate(head%zsoil)
709  allocate(head%lpl(dim1l),head%zsoil(dim1z),stat=iret)
710  if(iret.eq.0) then
711  head%lpl=0
712  head%zsoil=sfcio_realfill
713  endif
714  if(iret.ne.0) iret=-3
715 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
716  end subroutine
717 !-------------------------------------------------------------------------------
718  subroutine sfcio_aldata(head,data,iret)
719  implicit none
720  type(sfcio_head),intent(in):: head
721  type(sfcio_data),intent(inout):: data
722  integer(sfcio_intkind),intent(out):: iret
723  integer dim1,dim2,dim3
724 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
725  call sfcio_axdata(data,iret)
726  dim1=head%lonb
727  dim2=head%latb
728  dim3=head%lsoil
729  allocate(&
730  data%tsea(dim1,dim2),&
731  data%smc(dim1,dim2,dim3),&
732  data%sheleg(dim1,dim2),&
733  data%stc(dim1,dim2,dim3),&
734  data%tg3(dim1,dim2),&
735  data%zorl(dim1,dim2),&
736  data%cv(dim1,dim2),&
737  data%cvb(dim1,dim2),&
738  data%cvt(dim1,dim2),&
739  data%alvsf(dim1,dim2),&
740  data%alvwf(dim1,dim2),&
741  data%alnsf(dim1,dim2),&
742  data%alnwf(dim1,dim2),&
743  data%slmsk(dim1,dim2),&
744  data%vfrac(dim1,dim2),&
745  data%canopy(dim1,dim2),&
746  data%f10m(dim1,dim2),&
747  data%t2m(dim1,dim2),&
748  data%q2m(dim1,dim2),&
749  data%vtype(dim1,dim2),&
750  data%stype(dim1,dim2),&
751  data%facsf(dim1,dim2),&
752  data%facwf(dim1,dim2),&
753  data%uustar(dim1,dim2),&
754  data%ffmm(dim1,dim2),&
755  data%ffhh(dim1,dim2),&
756  data%hice(dim1,dim2),&
757  data%fice(dim1,dim2),&
758  data%tisfc(dim1,dim2),&
759  data%tprcp(dim1,dim2),&
760  data%srflag(dim1,dim2),&
761  data%snwdph(dim1,dim2),&
762  data%slc(dim1,dim2,dim3),&
763  data%shdmin(dim1,dim2),&
764  data%shdmax(dim1,dim2),&
765  data%slope(dim1,dim2),&
766  data%snoalb(dim1,dim2),&
767  data%orog(dim1,dim2),&
768  stat=iret)
769  if(iret.ne.0) iret=-3
770 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
771  end subroutine
772 !-------------------------------------------------------------------------------
773  subroutine sfcio_axdata(data,iret)
774  implicit none
775  type(sfcio_data),intent(inout):: data
776  integer(sfcio_intkind),intent(out):: iret
777 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
778  deallocate(&
779  data%tsea,&
780  data%smc,&
781  data%sheleg,&
782  data%stc,&
783  data%tg3,&
784  data%zorl,&
785  data%cv,&
786  data%cvb,&
787  data%cvt,&
788  data%alvsf,&
789  data%alvwf,&
790  data%alnsf,&
791  data%alnwf,&
792  data%slmsk,&
793  data%vfrac,&
794  data%canopy,&
795  data%f10m,&
796  data%t2m,&
797  data%q2m,&
798  data%vtype,&
799  data%stype,&
800  data%facsf,&
801  data%facwf,&
802  data%uustar,&
803  data%ffmm,&
804  data%ffhh,&
805  data%hice,&
806  data%fice,&
807  data%tisfc,&
808  data%tprcp,&
809  data%srflag,&
810  data%snwdph,&
811  data%slc,&
812  data%shdmin,&
813  data%shdmax,&
814  data%slope,&
815  data%snoalb,&
816  data%orog,&
817  stat=iret)
818  nullify(&
819  data%tsea,&
820  data%smc,&
821  data%sheleg,&
822  data%stc,&
823  data%tg3,&
824  data%zorl,&
825  data%cv,&
826  data%cvb,&
827  data%cvt,&
828  data%alvsf,&
829  data%alvwf,&
830  data%alnsf,&
831  data%alnwf,&
832  data%slmsk,&
833  data%vfrac,&
834  data%canopy,&
835  data%f10m,&
836  data%t2m,&
837  data%q2m,&
838  data%vtype,&
839  data%stype,&
840  data%facsf,&
841  data%facwf,&
842  data%uustar,&
843  data%ffmm,&
844  data%ffhh,&
845  data%hice,&
846  data%fice,&
847  data%tisfc,&
848  data%tprcp,&
849  data%srflag,&
850  data%snwdph,&
851  data%slc,&
852  data%shdmin,&
853  data%shdmax,&
854  data%slope,&
855  data%snoalb,&
856  data%orog)
857  if(iret.ne.0) iret=-3
858 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
859  end subroutine
860 !-------------------------------------------------------------------------------
861  subroutine sfcio_srdata(lu,head,data,iret)
862  implicit none
863  integer(sfcio_intkind),intent(in):: lu
864  type(sfcio_head),intent(in):: head
865  type(sfcio_data),intent(inout):: data
866  integer(sfcio_intkind),intent(out):: iret
867  integer:: dim1,dim2,dim3,mdim1,mdim2,mdim3
868  integer:: ios
869  integer i
870  type(sfcio_dbta) dbta
871 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
872  dim1=head%lonb
873  dim2=head%latb
874  dim3=head%lsoil
875  mdim1=min(&
876  size(data%tsea,1),&
877  size(data%smc,1),&
878  size(data%sheleg,1),&
879  size(data%stc,1),&
880  size(data%tg3,1),&
881  size(data%zorl,1),&
882  size(data%alvsf,1),&
883  size(data%alvwf,1),&
884  size(data%alnsf,1),&
885  size(data%alnwf,1),&
886  size(data%slmsk,1),&
887  size(data%vfrac,1),&
888  size(data%canopy,1),&
889  size(data%f10m,1),&
890  size(data%t2m,1),&
891  size(data%q2m,1),&
892  size(data%vtype,1),&
893  size(data%stype,1),&
894  size(data%facsf,1),&
895  size(data%facwf,1),&
896  size(data%uustar,1),&
897  size(data%ffmm,1),&
898  size(data%ffhh,1),&
899  size(data%hice,1),&
900  size(data%fice,1),&
901  size(data%tisfc,1),&
902  size(data%tprcp,1),&
903  size(data%srflag,1),&
904  size(data%snwdph,1),&
905  size(data%slc,1),&
906  size(data%shdmin,1),&
907  size(data%shdmax,1),&
908  size(data%slope,1),&
909  size(data%snoalb,1),&
910  size(data%orog,1))
911  mdim2=min(&
912  size(data%tsea,2),&
913  size(data%smc,2),&
914  size(data%sheleg,2),&
915  size(data%stc,2),&
916  size(data%tg3,2),&
917  size(data%zorl,2),&
918  size(data%alvsf,2),&
919  size(data%alvwf,2),&
920  size(data%alnsf,2),&
921  size(data%alnwf,2),&
922  size(data%slmsk,2),&
923  size(data%vfrac,2),&
924  size(data%canopy,2),&
925  size(data%f10m,2),&
926  size(data%t2m,2),&
927  size(data%q2m,2),&
928  size(data%vtype,2),&
929  size(data%stype,2),&
930  size(data%facsf,2),&
931  size(data%facwf,2),&
932  size(data%uustar,2),&
933  size(data%ffmm,2),&
934  size(data%ffhh,2),&
935  size(data%hice,2),&
936  size(data%fice,2),&
937  size(data%tisfc,2),&
938  size(data%tprcp,2),&
939  size(data%srflag,2),&
940  size(data%snwdph,2),&
941  size(data%slc,2),&
942  size(data%shdmin,2),&
943  size(data%shdmax,2),&
944  size(data%slope,2),&
945  size(data%snoalb,2),&
946  size(data%orog,2))
947  mdim3=min(&
948  size(data%smc,3),&
949  size(data%stc,3),&
950  size(data%slc,3))
951  iret=-5
952  if(mdim1.lt.dim1.or.&
953  mdim2.lt.dim2.or.&
954  mdim3.lt.dim3) return
955  iret=-4
956 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
957  data%t2m(:dim1,:dim2)=sfcio_realfill
958  data%q2m(:dim1,:dim2)=sfcio_realfill
959  data%tisfc(:dim1,:dim2)=sfcio_realfill
960  if(head%ivs.eq.200509) then
961  if(head%irealf.ne.2) then
962  read(lu,iostat=ios) data%slmsk(:dim1,:dim2)
963  if(ios.ne.0) return
964  read(lu,iostat=ios) data%orog(:dim1,:dim2)
965  if(ios.ne.0) return
966  read(lu,iostat=ios) data%tsea(:dim1,:dim2)
967  if(ios.ne.0) return
968  read(lu,iostat=ios) data%sheleg(:dim1,:dim2)
969  if(ios.ne.0) return
970  read(lu,iostat=ios) data%tg3(:dim1,:dim2)
971  if(ios.ne.0) return
972  read(lu,iostat=ios) data%zorl(:dim1,:dim2)
973  if(ios.ne.0) return
974  read(lu,iostat=ios) data%alvsf(:dim1,:dim2)
975  if(ios.ne.0) return
976  read(lu,iostat=ios) data%alvwf(:dim1,:dim2)
977  if(ios.ne.0) return
978  read(lu,iostat=ios) data%alnsf(:dim1,:dim2)
979  if(ios.ne.0) return
980  read(lu,iostat=ios) data%alnwf(:dim1,:dim2)
981  if(ios.ne.0) return
982  read(lu,iostat=ios) data%vfrac(:dim1,:dim2)
983  if(ios.ne.0) return
984  read(lu,iostat=ios) data%canopy(:dim1,:dim2)
985  if(ios.ne.0) return
986  read(lu,iostat=ios) data%f10m(:dim1,:dim2)
987  if(ios.ne.0) return
988  read(lu,iostat=ios) data%t2m(:dim1,:dim2)
989  if(ios.ne.0) return
990  read(lu,iostat=ios) data%q2m(:dim1,:dim2)
991  if(ios.ne.0) return
992  read(lu,iostat=ios) data%vtype(:dim1,:dim2)
993  if(ios.ne.0) return
994  read(lu,iostat=ios) data%stype(:dim1,:dim2)
995  if(ios.ne.0) return
996  read(lu,iostat=ios) data%facsf(:dim1,:dim2)
997  if(ios.ne.0) return
998  read(lu,iostat=ios) data%facwf(:dim1,:dim2)
999  if(ios.ne.0) return
1000  read(lu,iostat=ios) data%uustar(:dim1,:dim2)
1001  if(ios.ne.0) return
1002  read(lu,iostat=ios) data%ffmm(:dim1,:dim2)
1003  if(ios.ne.0) return
1004  read(lu,iostat=ios) data%ffhh(:dim1,:dim2)
1005  if(ios.ne.0) return
1006  read(lu,iostat=ios) data%hice(:dim1,:dim2)
1007  if(ios.ne.0) return
1008  read(lu,iostat=ios) data%fice(:dim1,:dim2)
1009  if(ios.ne.0) return
1010  read(lu,iostat=ios) data%tisfc(:dim1,:dim2)
1011  if(ios.ne.0) return
1012  read(lu,iostat=ios) data%tprcp(:dim1,:dim2)
1013  if(ios.ne.0) return
1014  read(lu,iostat=ios) data%srflag(:dim1,:dim2)
1015  if(ios.ne.0) return
1016  read(lu,iostat=ios) data%snwdph(:dim1,:dim2)
1017  if(ios.ne.0) return
1018  read(lu,iostat=ios) data%shdmin(:dim1,:dim2)
1019  if(ios.ne.0) return
1020  read(lu,iostat=ios) data%shdmax(:dim1,:dim2)
1021  if(ios.ne.0) return
1022  read(lu,iostat=ios) data%slope(:dim1,:dim2)
1023  if(ios.ne.0) return
1024  read(lu,iostat=ios) data%snoalb(:dim1,:dim2)
1025  if(ios.ne.0) return
1026  do i=1,head%lsoil
1027  read(lu,iostat=ios) data%stc(:dim1,:dim2,i)
1028  if(ios.ne.0) return
1029  enddo
1030  do i=1,head%lsoil
1031  read(lu,iostat=ios) data%smc(:dim1,:dim2,i)
1032  if(ios.ne.0) return
1033  enddo
1034  do i=1,head%lsoil
1035  read(lu,iostat=ios) data%slc(:dim1,:dim2,i)
1036  if(ios.ne.0) return
1037  enddo
1038  data%cv(:dim1,:dim2)=sfcio_realfill
1039  data%cvb(:dim1,:dim2)=sfcio_realfill
1040  data%cvt(:dim1,:dim2)=sfcio_realfill
1041  else
1042  call sfcio_aldbta(head,dbta,iret)
1043  if(iret.ne.0) return
1044  call sfcio_srdbta(lu,head,dbta,iret)
1045  if(iret.ne.0) return
1046  data%tsea(:dim1,:dim2)=dbta%tsea(:dim1,:dim2)
1047  data%smc(:dim1,:dim2,:dim3)=dbta%smc(:dim1,:dim2,:dim3)
1048  data%sheleg(:dim1,:dim2)=dbta%sheleg(:dim1,:dim2)
1049  data%stc(:dim1,:dim2,:dim3)=dbta%stc(:dim1,:dim2,:dim3)
1050  data%tg3(:dim1,:dim2)=dbta%tg3(:dim1,:dim2)
1051  data%zorl(:dim1,:dim2)=dbta%zorl(:dim1,:dim2)
1052  data%cv(:dim1,:dim2)=dbta%cv(:dim1,:dim2)
1053  data%cvb(:dim1,:dim2)=dbta%cvb(:dim1,:dim2)
1054  data%cvt(:dim1,:dim2)=dbta%cvt(:dim1,:dim2)
1055  data%alvsf(:dim1,:dim2)=dbta%alvsf(:dim1,:dim2)
1056  data%alvwf(:dim1,:dim2)=dbta%alvwf(:dim1,:dim2)
1057  data%alnsf(:dim1,:dim2)=dbta%alnsf(:dim1,:dim2)
1058  data%alnwf(:dim1,:dim2)=dbta%alnwf(:dim1,:dim2)
1059  data%slmsk(:dim1,:dim2)=dbta%slmsk(:dim1,:dim2)
1060  data%vfrac(:dim1,:dim2)=dbta%vfrac(:dim1,:dim2)
1061  data%canopy(:dim1,:dim2)=dbta%canopy(:dim1,:dim2)
1062  data%f10m(:dim1,:dim2)=dbta%f10m(:dim1,:dim2)
1063  data%t2m(:dim1,:dim2)=dbta%t2m(:dim1,:dim2)
1064  data%q2m(:dim1,:dim2)=dbta%q2m(:dim1,:dim2)
1065  data%vtype(:dim1,:dim2)=dbta%vtype(:dim1,:dim2)
1066  data%stype(:dim1,:dim2)=dbta%stype(:dim1,:dim2)
1067  data%facsf(:dim1,:dim2)=dbta%facsf(:dim1,:dim2)
1068  data%facwf(:dim1,:dim2)=dbta%facwf(:dim1,:dim2)
1069  data%uustar(:dim1,:dim2)=dbta%uustar(:dim1,:dim2)
1070  data%ffmm(:dim1,:dim2)=dbta%ffmm(:dim1,:dim2)
1071  data%ffhh(:dim1,:dim2)=dbta%ffhh(:dim1,:dim2)
1072  data%hice(:dim1,:dim2)=dbta%hice(:dim1,:dim2)
1073  data%fice(:dim1,:dim2)=dbta%fice(:dim1,:dim2)
1074  data%tisfc(:dim1,:dim2)=dbta%tisfc(:dim1,:dim2)
1075  data%tprcp(:dim1,:dim2)=dbta%tprcp(:dim1,:dim2)
1076  data%srflag(:dim1,:dim2)=dbta%srflag(:dim1,:dim2)
1077  data%snwdph(:dim1,:dim2)=dbta%snwdph(:dim1,:dim2)
1078  data%slc(:dim1,:dim2,:dim3)=dbta%slc(:dim1,:dim2,:dim3)
1079  data%shdmin(:dim1,:dim2)=dbta%shdmin(:dim1,:dim2)
1080  data%shdmax(:dim1,:dim2)=dbta%shdmax(:dim1,:dim2)
1081  data%slope(:dim1,:dim2)=dbta%slope(:dim1,:dim2)
1082  data%snoalb(:dim1,:dim2)=dbta%snoalb(:dim1,:dim2)
1083  data%orog(:dim1,:dim2)=dbta%orog(:dim1,:dim2)
1084  call sfcio_axdbta(dbta,iret)
1085  endif
1086 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1087  elseif(head%ivs.eq.200501.and.head%irealf.ne.2) then
1088  read(lu,iostat=ios) data%slmsk(:dim1,:dim2)
1089  if(ios.ne.0) return
1090  read(lu,iostat=ios) data%orog(:dim1,:dim2)
1091  if(ios.ne.0) return
1092  read(lu,iostat=ios) data%tsea(:dim1,:dim2)
1093  if(ios.ne.0) return
1094  read(lu,iostat=ios) data%sheleg(:dim1,:dim2)
1095  if(ios.ne.0) return
1096  read(lu,iostat=ios) data%tg3(:dim1,:dim2)
1097  if(ios.ne.0) return
1098  read(lu,iostat=ios) data%zorl(:dim1,:dim2)
1099  if(ios.ne.0) return
1100  read(lu,iostat=ios) data%alvsf(:dim1,:dim2)
1101  if(ios.ne.0) return
1102  read(lu,iostat=ios) data%alvwf(:dim1,:dim2)
1103  if(ios.ne.0) return
1104  read(lu,iostat=ios) data%alnsf(:dim1,:dim2)
1105  if(ios.ne.0) return
1106  read(lu,iostat=ios) data%alnwf(:dim1,:dim2)
1107  if(ios.ne.0) return
1108  read(lu,iostat=ios) data%vfrac(:dim1,:dim2)
1109  if(ios.ne.0) return
1110  read(lu,iostat=ios) data%canopy(:dim1,:dim2)
1111  if(ios.ne.0) return
1112  read(lu,iostat=ios) data%f10m(:dim1,:dim2)
1113  if(ios.ne.0) return
1114  read(lu,iostat=ios) data%vtype(:dim1,:dim2)
1115  if(ios.ne.0) return
1116  read(lu,iostat=ios) data%stype(:dim1,:dim2)
1117  if(ios.ne.0) return
1118  read(lu,iostat=ios) data%facsf(:dim1,:dim2)
1119  if(ios.ne.0) return
1120  read(lu,iostat=ios) data%facwf(:dim1,:dim2)
1121  if(ios.ne.0) return
1122  read(lu,iostat=ios) data%uustar(:dim1,:dim2)
1123  if(ios.ne.0) return
1124  read(lu,iostat=ios) data%ffmm(:dim1,:dim2)
1125  if(ios.ne.0) return
1126  read(lu,iostat=ios) data%ffhh(:dim1,:dim2)
1127  if(ios.ne.0) return
1128  read(lu,iostat=ios) data%hice(:dim1,:dim2)
1129  if(ios.ne.0) return
1130  read(lu,iostat=ios) data%fice(:dim1,:dim2)
1131  if(ios.ne.0) return
1132  read(lu,iostat=ios) data%tprcp(:dim1,:dim2)
1133  if(ios.ne.0) return
1134  read(lu,iostat=ios) data%srflag(:dim1,:dim2)
1135  if(ios.ne.0) return
1136  read(lu,iostat=ios) data%snwdph(:dim1,:dim2)
1137  if(ios.ne.0) return
1138  read(lu,iostat=ios) data%shdmin(:dim1,:dim2)
1139  if(ios.ne.0) return
1140  read(lu,iostat=ios) data%shdmax(:dim1,:dim2)
1141  if(ios.ne.0) return
1142  read(lu,iostat=ios) data%slope(:dim1,:dim2)
1143  if(ios.ne.0) return
1144  read(lu,iostat=ios) data%snoalb(:dim1,:dim2)
1145  if(ios.ne.0) return
1146  do i=1,head%lsoil
1147  read(lu,iostat=ios) data%stc(:dim1,:dim2,i)
1148  if(ios.ne.0) return
1149  enddo
1150  do i=1,head%lsoil
1151  read(lu,iostat=ios) data%smc(:dim1,:dim2,i)
1152  if(ios.ne.0) return
1153  enddo
1154  do i=1,head%lsoil
1155  read(lu,iostat=ios) data%slc(:dim1,:dim2,i)
1156  if(ios.ne.0) return
1157  enddo
1158  data%cv(:dim1,:dim2)=sfcio_realfill
1159  data%cvb(:dim1,:dim2)=sfcio_realfill
1160  data%cvt(:dim1,:dim2)=sfcio_realfill
1161 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1162  else
1163  read(lu,iostat=ios) data%tsea(:dim1,:dim2)
1164  if(ios.ne.0) return
1165  read(lu,iostat=ios) data%smc(:dim1,:dim2,:dim3)
1166  if(ios.ne.0) return
1167  read(lu,iostat=ios) data%sheleg(:dim1,:dim2)
1168  if(ios.ne.0) return
1169  read(lu,iostat=ios) data%stc(:dim1,:dim2,:dim3)
1170  if(ios.ne.0) return
1171  read(lu,iostat=ios) data%tg3(:dim1,:dim2)
1172  if(ios.ne.0) return
1173  read(lu,iostat=ios) data%zorl(:dim1,:dim2)
1174  if(ios.ne.0) return
1175  read(lu,iostat=ios) data%cv(:dim1,:dim2)
1176  if(ios.ne.0) return
1177  read(lu,iostat=ios) data%cvb(:dim1,:dim2)
1178  if(ios.ne.0) return
1179  read(lu,iostat=ios) data%cvt(:dim1,:dim2)
1180  if(ios.ne.0) return
1181  read(lu,iostat=ios) data%alvsf(:dim1,:dim2),&
1182  data%alvwf(:dim1,:dim2),&
1183  data%alnsf(:dim1,:dim2),&
1184  data%alnwf(:dim1,:dim2)
1185  if(ios.ne.0) return
1186  read(lu,iostat=ios) data%slmsk(:dim1,:dim2)
1187  if(ios.ne.0) return
1188  read(lu,iostat=ios) data%vfrac(:dim1,:dim2)
1189  if(ios.ne.0) return
1190  read(lu,iostat=ios) data%canopy(:dim1,:dim2)
1191  if(ios.ne.0) return
1192  read(lu,iostat=ios) data%f10m(:dim1,:dim2)
1193  if(ios.ne.0) return
1194  read(lu,iostat=ios) data%vtype(:dim1,:dim2)
1195  if(ios.ne.0) return
1196  read(lu,iostat=ios) data%stype(:dim1,:dim2)
1197  if(ios.ne.0) return
1198  read(lu,iostat=ios) data%facsf(:dim1,:dim2),&
1199  data%facwf(:dim1,:dim2)
1200  if(ios.ne.0) return
1201  if(head%ivs.ge.200004) then
1202  read(lu,iostat=ios) data%uustar(:dim1,:dim2)
1203  if(ios.ne.0) return
1204  read(lu,iostat=ios) data%ffmm(:dim1,:dim2)
1205  if(ios.ne.0) return
1206  read(lu,iostat=ios) data%ffhh(:dim1,:dim2)
1207  if(ios.ne.0) return
1208  else
1209  data%uustar(:dim1,:dim2)=sfcio_realfill
1210  data%ffmm(:dim1,:dim2)=sfcio_realfill
1211  data%ffhh(:dim1,:dim2)=sfcio_realfill
1212  endif
1213  if(head%ivs.eq.200412) then
1214  read(lu,iostat=ios) data%hice(:dim1,:dim2)
1215  if(ios.ne.0) return
1216  read(lu,iostat=ios) data%fice(:dim1,:dim2)
1217  if(ios.ne.0) return
1218  read(lu,iostat=ios) data%tprcp(:dim1,:dim2)
1219  if(ios.ne.0) return
1220  read(lu,iostat=ios) data%srflag(:dim1,:dim2)
1221  if(ios.ne.0) return
1222  read(lu,iostat=ios) data%snwdph(:dim1,:dim2)
1223  if(ios.ne.0) return
1224  read(lu,iostat=ios) data%slc(:dim1,:dim2,:dim3)
1225  if(ios.ne.0) return
1226  read(lu,iostat=ios) data%shdmin(:dim1,:dim2)
1227  if(ios.ne.0) return
1228  read(lu,iostat=ios) data%shdmax(:dim1,:dim2)
1229  if(ios.ne.0) return
1230  read(lu,iostat=ios) data%slope(:dim1,:dim2)
1231  if(ios.ne.0) return
1232  read(lu,iostat=ios) data%snoalb(:dim1,:dim2)
1233  if(ios.ne.0) return
1234  data%orog(:dim1,:dim2)=sfcio_realfill
1235  else
1236  data%hice(:dim1,:dim2)=sfcio_realfill
1237  data%fice(:dim1,:dim2)=sfcio_realfill
1238  data%tprcp(:dim1,:dim2)=sfcio_realfill
1239  data%srflag(:dim1,:dim2)=sfcio_realfill
1240  data%snwdph(:dim1,:dim2)=sfcio_realfill
1241  data%slc(:dim1,:dim2,:dim3)=sfcio_realfill
1242  data%shdmin(:dim1,:dim2)=sfcio_realfill
1243  data%shdmax(:dim1,:dim2)=sfcio_realfill
1244  data%slope(:dim1,:dim2)=sfcio_realfill
1245  data%snoalb(:dim1,:dim2)=sfcio_realfill
1246  data%orog(:dim1,:dim2)=sfcio_realfill
1247  endif
1248  endif
1249 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1250  iret=0
1251  end subroutine
1252 !-------------------------------------------------------------------------------
1253  subroutine sfcio_swdata(lu,head,data,iret)
1254  implicit none
1255  integer(sfcio_intkind),intent(in):: lu
1256  type(sfcio_head),intent(in):: head
1257  type(sfcio_data),intent(in):: data
1258  integer(sfcio_intkind),intent(out):: iret
1259  integer:: dim1,dim2,dim3,mdim1,mdim2,mdim3
1260  integer:: ios
1261  integer i
1262  type(sfcio_dbta) dbta
1263 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1264  dim1=head%lonb
1265  dim2=head%latb
1266  dim3=head%lsoil
1267  mdim1=min(&
1268  size(data%tsea,1),&
1269  size(data%smc,1),&
1270  size(data%sheleg,1),&
1271  size(data%stc,1),&
1272  size(data%tg3,1),&
1273  size(data%zorl,1),&
1274  size(data%alvsf,1),&
1275  size(data%alvwf,1),&
1276  size(data%alnsf,1),&
1277  size(data%alnwf,1),&
1278  size(data%slmsk,1),&
1279  size(data%vfrac,1),&
1280  size(data%canopy,1),&
1281  size(data%f10m,1),&
1282  size(data%t2m,1),&
1283  size(data%q2m,1),&
1284  size(data%vtype,1),&
1285  size(data%stype,1),&
1286  size(data%facsf,1),&
1287  size(data%facwf,1),&
1288  size(data%uustar,1),&
1289  size(data%ffmm,1),&
1290  size(data%ffhh,1),&
1291  size(data%hice,1),&
1292  size(data%fice,1),&
1293  size(data%tisfc,1),&
1294  size(data%tprcp,1),&
1295  size(data%srflag,1),&
1296  size(data%snwdph,1),&
1297  size(data%slc,1),&
1298  size(data%shdmin,1),&
1299  size(data%shdmax,1),&
1300  size(data%slope,1),&
1301  size(data%snoalb,1),&
1302  size(data%orog,1))
1303  mdim2=min(&
1304  size(data%tsea,2),&
1305  size(data%smc,2),&
1306  size(data%sheleg,2),&
1307  size(data%stc,2),&
1308  size(data%tg3,2),&
1309  size(data%zorl,2),&
1310  size(data%alvsf,2),&
1311  size(data%alvwf,2),&
1312  size(data%alnsf,2),&
1313  size(data%alnwf,2),&
1314  size(data%slmsk,2),&
1315  size(data%vfrac,2),&
1316  size(data%canopy,2),&
1317  size(data%f10m,2),&
1318  size(data%t2m,2),&
1319  size(data%q2m,2),&
1320  size(data%vtype,2),&
1321  size(data%stype,2),&
1322  size(data%facsf,2),&
1323  size(data%facwf,2),&
1324  size(data%uustar,2),&
1325  size(data%ffmm,2),&
1326  size(data%ffhh,2),&
1327  size(data%hice,2),&
1328  size(data%fice,2),&
1329  size(data%tisfc,2),&
1330  size(data%tprcp,2),&
1331  size(data%srflag,2),&
1332  size(data%snwdph,2),&
1333  size(data%slc,2),&
1334  size(data%shdmin,2),&
1335  size(data%shdmax,2),&
1336  size(data%slope,2),&
1337  size(data%snoalb,2),&
1338  size(data%orog,2))
1339  mdim3=min(&
1340  size(data%smc,3),&
1341  size(data%stc,3),&
1342  size(data%slc,3))
1343  iret=-5
1344  if(mdim1.lt.dim1.or.&
1345  mdim2.lt.dim2.or.&
1346  mdim3.lt.dim3) return
1347  iret=-4
1348 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1349  if(head%ivs.eq.200509) then
1350  if(head%irealf.ne.2) then
1351  write(lu,iostat=ios) data%slmsk(:dim1,:dim2)
1352  if(ios.ne.0) return
1353  write(lu,iostat=ios) data%orog(:dim1,:dim2)
1354  if(ios.ne.0) return
1355  write(lu,iostat=ios) data%tsea(:dim1,:dim2)
1356  if(ios.ne.0) return
1357  write(lu,iostat=ios) data%sheleg(:dim1,:dim2)
1358  if(ios.ne.0) return
1359  write(lu,iostat=ios) data%tg3(:dim1,:dim2)
1360  if(ios.ne.0) return
1361  write(lu,iostat=ios) data%zorl(:dim1,:dim2)
1362  if(ios.ne.0) return
1363  write(lu,iostat=ios) data%alvsf(:dim1,:dim2)
1364  if(ios.ne.0) return
1365  write(lu,iostat=ios) data%alvwf(:dim1,:dim2)
1366  if(ios.ne.0) return
1367  write(lu,iostat=ios) data%alnsf(:dim1,:dim2)
1368  if(ios.ne.0) return
1369  write(lu,iostat=ios) data%alnwf(:dim1,:dim2)
1370  if(ios.ne.0) return
1371  write(lu,iostat=ios) data%vfrac(:dim1,:dim2)
1372  if(ios.ne.0) return
1373  write(lu,iostat=ios) data%canopy(:dim1,:dim2)
1374  if(ios.ne.0) return
1375  write(lu,iostat=ios) data%f10m(:dim1,:dim2)
1376  if(ios.ne.0) return
1377  write(lu,iostat=ios) data%t2m(:dim1,:dim2)
1378  if(ios.ne.0) return
1379  write(lu,iostat=ios) data%q2m(:dim1,:dim2)
1380  if(ios.ne.0) return
1381  write(lu,iostat=ios) data%vtype(:dim1,:dim2)
1382  if(ios.ne.0) return
1383  write(lu,iostat=ios) data%stype(:dim1,:dim2)
1384  if(ios.ne.0) return
1385  write(lu,iostat=ios) data%facsf(:dim1,:dim2)
1386  if(ios.ne.0) return
1387  write(lu,iostat=ios) data%facwf(:dim1,:dim2)
1388  if(ios.ne.0) return
1389  write(lu,iostat=ios) data%uustar(:dim1,:dim2)
1390  if(ios.ne.0) return
1391  write(lu,iostat=ios) data%ffmm(:dim1,:dim2)
1392  if(ios.ne.0) return
1393  write(lu,iostat=ios) data%ffhh(:dim1,:dim2)
1394  if(ios.ne.0) return
1395  write(lu,iostat=ios) data%hice(:dim1,:dim2)
1396  if(ios.ne.0) return
1397  write(lu,iostat=ios) data%fice(:dim1,:dim2)
1398  if(ios.ne.0) return
1399  write(lu,iostat=ios) data%tisfc(:dim1,:dim2)
1400  if(ios.ne.0) return
1401  write(lu,iostat=ios) data%tprcp(:dim1,:dim2)
1402  if(ios.ne.0) return
1403  write(lu,iostat=ios) data%srflag(:dim1,:dim2)
1404  if(ios.ne.0) return
1405  write(lu,iostat=ios) data%snwdph(:dim1,:dim2)
1406  if(ios.ne.0) return
1407  write(lu,iostat=ios) data%shdmin(:dim1,:dim2)
1408  if(ios.ne.0) return
1409  write(lu,iostat=ios) data%shdmax(:dim1,:dim2)
1410  if(ios.ne.0) return
1411  write(lu,iostat=ios) data%slope(:dim1,:dim2)
1412  if(ios.ne.0) return
1413  write(lu,iostat=ios) data%snoalb(:dim1,:dim2)
1414  if(ios.ne.0) return
1415  do i=1,head%lsoil
1416  write(lu,iostat=ios) data%stc(:dim1,:dim2,i)
1417  if(ios.ne.0) return
1418  enddo
1419  do i=1,head%lsoil
1420  write(lu,iostat=ios) data%smc(:dim1,:dim2,i)
1421  if(ios.ne.0) return
1422  enddo
1423  do i=1,head%lsoil
1424  write(lu,iostat=ios) data%slc(:dim1,:dim2,i)
1425  if(ios.ne.0) return
1426  enddo
1427  else
1428  call sfcio_aldbta(head,dbta,iret)
1429  if(iret.ne.0) return
1430  dbta%tsea(:dim1,:dim2)=data%tsea(:dim1,:dim2)
1431  dbta%smc(:dim1,:dim2,:dim3)=data%smc(:dim1,:dim2,:dim3)
1432  dbta%sheleg(:dim1,:dim2)=data%sheleg(:dim1,:dim2)
1433  dbta%stc(:dim1,:dim2,:dim3)=data%stc(:dim1,:dim2,:dim3)
1434  dbta%tg3(:dim1,:dim2)=data%tg3(:dim1,:dim2)
1435  dbta%zorl(:dim1,:dim2)=data%zorl(:dim1,:dim2)
1436  dbta%cv(:dim1,:dim2)=data%cv(:dim1,:dim2)
1437  dbta%cvb(:dim1,:dim2)=data%cvb(:dim1,:dim2)
1438  dbta%cvt(:dim1,:dim2)=data%cvt(:dim1,:dim2)
1439  dbta%alvsf(:dim1,:dim2)=data%alvsf(:dim1,:dim2)
1440  dbta%alvwf(:dim1,:dim2)=data%alvwf(:dim1,:dim2)
1441  dbta%alnsf(:dim1,:dim2)=data%alnsf(:dim1,:dim2)
1442  dbta%alnwf(:dim1,:dim2)=data%alnwf(:dim1,:dim2)
1443  dbta%slmsk(:dim1,:dim2)=data%slmsk(:dim1,:dim2)
1444  dbta%vfrac(:dim1,:dim2)=data%vfrac(:dim1,:dim2)
1445  dbta%canopy(:dim1,:dim2)=data%canopy(:dim1,:dim2)
1446  dbta%f10m(:dim1,:dim2)=data%f10m(:dim1,:dim2)
1447  dbta%t2m(:dim1,:dim2)=data%t2m(:dim1,:dim2)
1448  dbta%q2m(:dim1,:dim2)=data%q2m(:dim1,:dim2)
1449  dbta%vtype(:dim1,:dim2)=data%vtype(:dim1,:dim2)
1450  dbta%stype(:dim1,:dim2)=data%stype(:dim1,:dim2)
1451  dbta%facsf(:dim1,:dim2)=data%facsf(:dim1,:dim2)
1452  dbta%facwf(:dim1,:dim2)=data%facwf(:dim1,:dim2)
1453  dbta%uustar(:dim1,:dim2)=data%uustar(:dim1,:dim2)
1454  dbta%ffmm(:dim1,:dim2)=data%ffmm(:dim1,:dim2)
1455  dbta%ffhh(:dim1,:dim2)=data%ffhh(:dim1,:dim2)
1456  dbta%hice(:dim1,:dim2)=data%hice(:dim1,:dim2)
1457  dbta%fice(:dim1,:dim2)=data%fice(:dim1,:dim2)
1458  dbta%tisfc(:dim1,:dim2)=data%tisfc(:dim1,:dim2)
1459  dbta%tprcp(:dim1,:dim2)=data%tprcp(:dim1,:dim2)
1460  dbta%srflag(:dim1,:dim2)=data%srflag(:dim1,:dim2)
1461  dbta%snwdph(:dim1,:dim2)=data%snwdph(:dim1,:dim2)
1462  dbta%slc(:dim1,:dim2,:dim3)=data%slc(:dim1,:dim2,:dim3)
1463  dbta%shdmin(:dim1,:dim2)=data%shdmin(:dim1,:dim2)
1464  dbta%shdmax(:dim1,:dim2)=data%shdmax(:dim1,:dim2)
1465  dbta%slope(:dim1,:dim2)=data%slope(:dim1,:dim2)
1466  dbta%snoalb(:dim1,:dim2)=data%snoalb(:dim1,:dim2)
1467  dbta%orog(:dim1,:dim2)=data%orog(:dim1,:dim2)
1468  call sfcio_swdbta(lu,head,dbta,iret)
1469  if(iret.ne.0) return
1470  call sfcio_axdbta(dbta,iret)
1471  endif
1472 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1473  elseif(head%ivs.eq.200501.and.head%irealf.ne.2) then
1474  write(lu,iostat=ios) data%slmsk(:dim1,:dim2)
1475  if(ios.ne.0) return
1476  write(lu,iostat=ios) data%orog(:dim1,:dim2)
1477  if(ios.ne.0) return
1478  write(lu,iostat=ios) data%tsea(:dim1,:dim2)
1479  if(ios.ne.0) return
1480  write(lu,iostat=ios) data%sheleg(:dim1,:dim2)
1481  if(ios.ne.0) return
1482  write(lu,iostat=ios) data%tg3(:dim1,:dim2)
1483  if(ios.ne.0) return
1484  write(lu,iostat=ios) data%zorl(:dim1,:dim2)
1485  if(ios.ne.0) return
1486  write(lu,iostat=ios) data%alvsf(:dim1,:dim2)
1487  if(ios.ne.0) return
1488  write(lu,iostat=ios) data%alvwf(:dim1,:dim2)
1489  if(ios.ne.0) return
1490  write(lu,iostat=ios) data%alnsf(:dim1,:dim2)
1491  if(ios.ne.0) return
1492  write(lu,iostat=ios) data%alnwf(:dim1,:dim2)
1493  if(ios.ne.0) return
1494  write(lu,iostat=ios) data%vfrac(:dim1,:dim2)
1495  if(ios.ne.0) return
1496  write(lu,iostat=ios) data%canopy(:dim1,:dim2)
1497  if(ios.ne.0) return
1498  write(lu,iostat=ios) data%f10m(:dim1,:dim2)
1499  if(ios.ne.0) return
1500  write(lu,iostat=ios) data%vtype(:dim1,:dim2)
1501  if(ios.ne.0) return
1502  write(lu,iostat=ios) data%stype(:dim1,:dim2)
1503  if(ios.ne.0) return
1504  write(lu,iostat=ios) data%facsf(:dim1,:dim2)
1505  if(ios.ne.0) return
1506  write(lu,iostat=ios) data%facwf(:dim1,:dim2)
1507  if(ios.ne.0) return
1508  write(lu,iostat=ios) data%uustar(:dim1,:dim2)
1509  if(ios.ne.0) return
1510  write(lu,iostat=ios) data%ffmm(:dim1,:dim2)
1511  if(ios.ne.0) return
1512  write(lu,iostat=ios) data%ffhh(:dim1,:dim2)
1513  if(ios.ne.0) return
1514  write(lu,iostat=ios) data%hice(:dim1,:dim2)
1515  if(ios.ne.0) return
1516  write(lu,iostat=ios) data%fice(:dim1,:dim2)
1517  if(ios.ne.0) return
1518  write(lu,iostat=ios) data%tprcp(:dim1,:dim2)
1519  if(ios.ne.0) return
1520  write(lu,iostat=ios) data%srflag(:dim1,:dim2)
1521  if(ios.ne.0) return
1522  write(lu,iostat=ios) data%snwdph(:dim1,:dim2)
1523  if(ios.ne.0) return
1524  write(lu,iostat=ios) data%shdmin(:dim1,:dim2)
1525  if(ios.ne.0) return
1526  write(lu,iostat=ios) data%shdmax(:dim1,:dim2)
1527  if(ios.ne.0) return
1528  write(lu,iostat=ios) data%slope(:dim1,:dim2)
1529  if(ios.ne.0) return
1530  write(lu,iostat=ios) data%snoalb(:dim1,:dim2)
1531  if(ios.ne.0) return
1532  do i=1,head%lsoil
1533  write(lu,iostat=ios) data%stc(:dim1,:dim2,i)
1534  if(ios.ne.0) return
1535  enddo
1536  do i=1,head%lsoil
1537  write(lu,iostat=ios) data%smc(:dim1,:dim2,i)
1538  if(ios.ne.0) return
1539  enddo
1540  do i=1,head%lsoil
1541  write(lu,iostat=ios) data%slc(:dim1,:dim2,i)
1542  if(ios.ne.0) return
1543  enddo
1544 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1545  elseif(head%ivs.eq.200004.and.head%lsoil.eq.2) then
1546  write(lu,iostat=ios) data%tsea(:dim1,:dim2)
1547  if(ios.ne.0) return
1548  write(lu,iostat=ios) data%smc(:dim1,:dim2,:dim3)
1549  if(ios.ne.0) return
1550  write(lu,iostat=ios) data%sheleg(:dim1,:dim2)
1551  if(ios.ne.0) return
1552  write(lu,iostat=ios) data%stc(:dim1,:dim2,:dim3)
1553  if(ios.ne.0) return
1554  write(lu,iostat=ios) data%tg3(:dim1,:dim2)
1555  if(ios.ne.0) return
1556  write(lu,iostat=ios) data%zorl(:dim1,:dim2)
1557  if(ios.ne.0) return
1558  write(lu,iostat=ios) data%cv(:dim1,:dim2)
1559  if(ios.ne.0) return
1560  write(lu,iostat=ios) data%cvb(:dim1,:dim2)
1561  if(ios.ne.0) return
1562  write(lu,iostat=ios) data%cvt(:dim1,:dim2)
1563  if(ios.ne.0) return
1564  write(lu,iostat=ios) data%alvsf(:dim1,:dim2),&
1565  data%alvwf(:dim1,:dim2),&
1566  data%alnsf(:dim1,:dim2),&
1567  data%alnwf(:dim1,:dim2)
1568  if(ios.ne.0) return
1569  write(lu,iostat=ios) data%slmsk(:dim1,:dim2)
1570  if(ios.ne.0) return
1571  write(lu,iostat=ios) data%vfrac(:dim1,:dim2)
1572  if(ios.ne.0) return
1573  write(lu,iostat=ios) data%canopy(:dim1,:dim2)
1574  if(ios.ne.0) return
1575  write(lu,iostat=ios) data%f10m(:dim1,:dim2)
1576  if(ios.ne.0) return
1577  write(lu,iostat=ios) data%vtype(:dim1,:dim2)
1578  if(ios.ne.0) return
1579  write(lu,iostat=ios) data%stype(:dim1,:dim2)
1580  if(ios.ne.0) return
1581  write(lu,iostat=ios) data%facsf(:dim1,:dim2),&
1582  data%facwf(:dim1,:dim2)
1583  if(ios.ne.0) return
1584  write(lu,iostat=ios) data%uustar(:dim1,:dim2)
1585  if(ios.ne.0) return
1586  write(lu,iostat=ios) data%ffmm(:dim1,:dim2)
1587  if(ios.ne.0) return
1588  write(lu,iostat=ios) data%ffhh(:dim1,:dim2)
1589  if(ios.ne.0) return
1590  endif
1591 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1592  iret=0
1593  end subroutine
1594 !-------------------------------------------------------------------------------
1595  subroutine sfcio_srohdca(lu,cfname,head,data,iret)
1596  implicit none
1597  integer(sfcio_intkind),intent(in):: lu
1598  character*(*),intent(in):: cfname
1599  type(sfcio_head),intent(inout):: head
1600  type(sfcio_data),intent(inout):: data
1601  integer(sfcio_intkind),intent(out):: iret
1602 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1603  call sfcio_sropen(lu,cfname,iret)
1604  if(iret.ne.0) return
1605 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1606  call sfcio_srhead(lu,head,iret)
1607  if(iret.ne.0) return
1608 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1609  call sfcio_aldata(head,data,iret)
1610  if(iret.ne.0) return
1611 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1612  call sfcio_srdata(lu,head,data,iret)
1613  if(iret.ne.0) return
1614 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1615  call sfcio_sclose(lu,iret)
1616  if(iret.ne.0) return
1617 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1618  end subroutine
1619 !-------------------------------------------------------------------------------
1620  subroutine sfcio_swohdca(lu,cfname,head,data,iret)
1621  implicit none
1622  integer(sfcio_intkind),intent(in):: lu
1623  character*(*),intent(in):: cfname
1624  type(sfcio_head),intent(in):: head
1625  type(sfcio_data),intent(in):: data
1626  integer(sfcio_intkind),intent(out):: iret
1627 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1628  call sfcio_swopen(lu,cfname,iret)
1629  if(iret.ne.0) return
1630 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1631  call sfcio_swhead(lu,head,iret)
1632  if(iret.ne.0) return
1633 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1634  call sfcio_swdata(lu,head,data,iret)
1635  if(iret.ne.0) return
1636 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1637  call sfcio_sclose(lu,iret)
1638  if(iret.ne.0) return
1639 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1640  end subroutine
1641 !-------------------------------------------------------------------------------
1642  subroutine sfcio_aldbta(head,dbta,iret)
1643  implicit none
1644  type(sfcio_head),intent(in):: head
1645  type(sfcio_dbta),intent(inout):: dbta
1646  integer(sfcio_intkind),intent(out):: iret
1647  integer dim1,dim2,dim3
1648 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1649  call sfcio_axdbta(dbta,iret)
1650  dim1=head%lonb
1651  dim2=head%latb
1652  dim3=head%lsoil
1653  allocate(&
1654  dbta%tsea(dim1,dim2),&
1655  dbta%smc(dim1,dim2,dim3),&
1656  dbta%sheleg(dim1,dim2),&
1657  dbta%stc(dim1,dim2,dim3),&
1658  dbta%tg3(dim1,dim2),&
1659  dbta%zorl(dim1,dim2),&
1660  dbta%cv(dim1,dim2),&
1661  dbta%cvb(dim1,dim2),&
1662  dbta%cvt(dim1,dim2),&
1663  dbta%alvsf(dim1,dim2),&
1664  dbta%alvwf(dim1,dim2),&
1665  dbta%alnsf(dim1,dim2),&
1666  dbta%alnwf(dim1,dim2),&
1667  dbta%slmsk(dim1,dim2),&
1668  dbta%vfrac(dim1,dim2),&
1669  dbta%canopy(dim1,dim2),&
1670  dbta%f10m(dim1,dim2),&
1671  dbta%t2m(dim1,dim2),&
1672  dbta%q2m(dim1,dim2),&
1673  dbta%vtype(dim1,dim2),&
1674  dbta%stype(dim1,dim2),&
1675  dbta%facsf(dim1,dim2),&
1676  dbta%facwf(dim1,dim2),&
1677  dbta%uustar(dim1,dim2),&
1678  dbta%ffmm(dim1,dim2),&
1679  dbta%ffhh(dim1,dim2),&
1680  dbta%hice(dim1,dim2),&
1681  dbta%fice(dim1,dim2),&
1682  dbta%tisfc(dim1,dim2),&
1683  dbta%tprcp(dim1,dim2),&
1684  dbta%srflag(dim1,dim2),&
1685  dbta%snwdph(dim1,dim2),&
1686  dbta%slc(dim1,dim2,dim3),&
1687  dbta%shdmin(dim1,dim2),&
1688  dbta%shdmax(dim1,dim2),&
1689  dbta%slope(dim1,dim2),&
1690  dbta%snoalb(dim1,dim2),&
1691  dbta%orog(dim1,dim2),&
1692  stat=iret)
1693  if(iret.ne.0) iret=-3
1694 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1695  end subroutine
1696 !-------------------------------------------------------------------------------
1697  subroutine sfcio_axdbta(dbta,iret)
1698  implicit none
1699  type(sfcio_dbta),intent(inout):: dbta
1700  integer(sfcio_intkind),intent(out):: iret
1701 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1702  deallocate(&
1703  dbta%tsea,&
1704  dbta%smc,&
1705  dbta%sheleg,&
1706  dbta%stc,&
1707  dbta%tg3,&
1708  dbta%zorl,&
1709  dbta%cv,&
1710  dbta%cvb,&
1711  dbta%cvt,&
1712  dbta%alvsf,&
1713  dbta%alvwf,&
1714  dbta%alnsf,&
1715  dbta%alnwf,&
1716  dbta%slmsk,&
1717  dbta%vfrac,&
1718  dbta%canopy,&
1719  dbta%f10m,&
1720  dbta%t2m,&
1721  dbta%q2m,&
1722  dbta%vtype,&
1723  dbta%stype,&
1724  dbta%facsf,&
1725  dbta%facwf,&
1726  dbta%uustar,&
1727  dbta%ffmm,&
1728  dbta%ffhh,&
1729  dbta%hice,&
1730  dbta%fice,&
1731  dbta%tisfc,&
1732  dbta%tprcp,&
1733  dbta%srflag,&
1734  dbta%snwdph,&
1735  dbta%slc,&
1736  dbta%shdmin,&
1737  dbta%shdmax,&
1738  dbta%slope,&
1739  dbta%snoalb,&
1740  dbta%orog,&
1741  stat=iret)
1742  nullify(&
1743  dbta%tsea,&
1744  dbta%smc,&
1745  dbta%sheleg,&
1746  dbta%stc,&
1747  dbta%tg3,&
1748  dbta%zorl,&
1749  dbta%cv,&
1750  dbta%cvb,&
1751  dbta%cvt,&
1752  dbta%alvsf,&
1753  dbta%alvwf,&
1754  dbta%alnsf,&
1755  dbta%alnwf,&
1756  dbta%slmsk,&
1757  dbta%vfrac,&
1758  dbta%canopy,&
1759  dbta%f10m,&
1760  dbta%t2m,&
1761  dbta%q2m,&
1762  dbta%vtype,&
1763  dbta%stype,&
1764  dbta%facsf,&
1765  dbta%facwf,&
1766  dbta%uustar,&
1767  dbta%ffmm,&
1768  dbta%ffhh,&
1769  dbta%hice,&
1770  dbta%fice,&
1771  dbta%tisfc,&
1772  dbta%tprcp,&
1773  dbta%srflag,&
1774  dbta%snwdph,&
1775  dbta%slc,&
1776  dbta%shdmin,&
1777  dbta%shdmax,&
1778  dbta%slope,&
1779  dbta%snoalb,&
1780  dbta%orog)
1781  if(iret.ne.0) iret=-3
1782 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1783  end subroutine
1784 !-------------------------------------------------------------------------------
1785  subroutine sfcio_srdbta(lu,head,dbta,iret)
1786  implicit none
1787  integer(sfcio_intkind),intent(in):: lu
1788  type(sfcio_head),intent(in):: head
1789  type(sfcio_dbta),intent(inout):: dbta
1790  integer(sfcio_intkind),intent(out):: iret
1791  integer:: dim1,dim2,dim3,mdim1,mdim2,mdim3
1792  integer:: ios
1793  integer i
1794  type(sfcio_data):: data
1795 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1796  dim1=head%lonb
1797  dim2=head%latb
1798  dim3=head%lsoil
1799  mdim1=min(&
1800  size(dbta%tsea,1),&
1801  size(dbta%smc,1),&
1802  size(dbta%sheleg,1),&
1803  size(dbta%stc,1),&
1804  size(dbta%tg3,1),&
1805  size(dbta%zorl,1),&
1806  size(dbta%alvsf,1),&
1807  size(dbta%alvwf,1),&
1808  size(dbta%alnsf,1),&
1809  size(dbta%alnwf,1),&
1810  size(dbta%slmsk,1),&
1811  size(dbta%vfrac,1),&
1812  size(dbta%canopy,1),&
1813  size(dbta%f10m,1),&
1814  size(dbta%t2m,1),&
1815  size(dbta%q2m,1),&
1816  size(dbta%vtype,1),&
1817  size(dbta%stype,1),&
1818  size(dbta%facsf,1),&
1819  size(dbta%facwf,1),&
1820  size(dbta%uustar,1),&
1821  size(dbta%ffmm,1),&
1822  size(dbta%ffhh,1),&
1823  size(dbta%hice,1),&
1824  size(dbta%fice,1),&
1825  size(dbta%tisfc,1),&
1826  size(dbta%tprcp,1),&
1827  size(dbta%srflag,1),&
1828  size(dbta%snwdph,1),&
1829  size(dbta%slc,1),&
1830  size(dbta%shdmin,1),&
1831  size(dbta%shdmax,1),&
1832  size(dbta%slope,1),&
1833  size(dbta%snoalb,1),&
1834  size(dbta%orog,1))
1835  mdim2=min(&
1836  size(dbta%tsea,2),&
1837  size(dbta%smc,2),&
1838  size(dbta%sheleg,2),&
1839  size(dbta%stc,2),&
1840  size(dbta%tg3,2),&
1841  size(dbta%zorl,2),&
1842  size(dbta%alvsf,2),&
1843  size(dbta%alvwf,2),&
1844  size(dbta%alnsf,2),&
1845  size(dbta%alnwf,2),&
1846  size(dbta%slmsk,2),&
1847  size(dbta%vfrac,2),&
1848  size(dbta%canopy,2),&
1849  size(dbta%f10m,2),&
1850  size(dbta%t2m,2),&
1851  size(dbta%q2m,2),&
1852  size(dbta%vtype,2),&
1853  size(dbta%stype,2),&
1854  size(dbta%facsf,2),&
1855  size(dbta%facwf,2),&
1856  size(dbta%uustar,2),&
1857  size(dbta%ffmm,2),&
1858  size(dbta%ffhh,2),&
1859  size(dbta%hice,2),&
1860  size(dbta%fice,2),&
1861  size(dbta%tisfc,2),&
1862  size(dbta%tprcp,2),&
1863  size(dbta%srflag,2),&
1864  size(dbta%snwdph,2),&
1865  size(dbta%slc,2),&
1866  size(dbta%shdmin,2),&
1867  size(dbta%shdmax,2),&
1868  size(dbta%slope,2),&
1869  size(dbta%snoalb,2),&
1870  size(dbta%orog,2))
1871  mdim3=min(&
1872  size(dbta%smc,3),&
1873  size(dbta%stc,3),&
1874  size(dbta%slc,3))
1875  iret=-5
1876  if(mdim1.lt.dim1.or.&
1877  mdim2.lt.dim2.or.&
1878  mdim3.lt.dim3) return
1879  iret=-4
1880 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1881  if(head%irealf.ne.2) then
1882  call sfcio_aldata(head,data,iret)
1883  if(iret.ne.0) return
1884  call sfcio_srdata(lu,head,data,iret)
1885  if(iret.ne.0) return
1886  dbta%tsea(:dim1,:dim2)=data%tsea(:dim1,:dim2)
1887  dbta%smc(:dim1,:dim2,:dim3)=data%smc(:dim1,:dim2,:dim3)
1888  dbta%sheleg(:dim1,:dim2)=data%sheleg(:dim1,:dim2)
1889  dbta%stc(:dim1,:dim2,:dim3)=data%stc(:dim1,:dim2,:dim3)
1890  dbta%tg3(:dim1,:dim2)=data%tg3(:dim1,:dim2)
1891  dbta%zorl(:dim1,:dim2)=data%zorl(:dim1,:dim2)
1892  dbta%cv(:dim1,:dim2)=data%cv(:dim1,:dim2)
1893  dbta%cvb(:dim1,:dim2)=data%cvb(:dim1,:dim2)
1894  dbta%cvt(:dim1,:dim2)=data%cvt(:dim1,:dim2)
1895  dbta%alvsf(:dim1,:dim2)=data%alvsf(:dim1,:dim2)
1896  dbta%alvwf(:dim1,:dim2)=data%alvwf(:dim1,:dim2)
1897  dbta%alnsf(:dim1,:dim2)=data%alnsf(:dim1,:dim2)
1898  dbta%alnwf(:dim1,:dim2)=data%alnwf(:dim1,:dim2)
1899  dbta%slmsk(:dim1,:dim2)=data%slmsk(:dim1,:dim2)
1900  dbta%vfrac(:dim1,:dim2)=data%vfrac(:dim1,:dim2)
1901  dbta%canopy(:dim1,:dim2)=data%canopy(:dim1,:dim2)
1902  dbta%f10m(:dim1,:dim2)=data%f10m(:dim1,:dim2)
1903  dbta%t2m(:dim1,:dim2)=data%t2m(:dim1,:dim2)
1904  dbta%q2m(:dim1,:dim2)=data%q2m(:dim1,:dim2)
1905  dbta%vtype(:dim1,:dim2)=data%vtype(:dim1,:dim2)
1906  dbta%stype(:dim1,:dim2)=data%stype(:dim1,:dim2)
1907  dbta%facsf(:dim1,:dim2)=data%facsf(:dim1,:dim2)
1908  dbta%facwf(:dim1,:dim2)=data%facwf(:dim1,:dim2)
1909  dbta%uustar(:dim1,:dim2)=data%uustar(:dim1,:dim2)
1910  dbta%ffmm(:dim1,:dim2)=data%ffmm(:dim1,:dim2)
1911  dbta%ffhh(:dim1,:dim2)=data%ffhh(:dim1,:dim2)
1912  dbta%hice(:dim1,:dim2)=data%hice(:dim1,:dim2)
1913  dbta%fice(:dim1,:dim2)=data%fice(:dim1,:dim2)
1914  dbta%tisfc(:dim1,:dim2)=data%tisfc(:dim1,:dim2)
1915  dbta%tprcp(:dim1,:dim2)=data%tprcp(:dim1,:dim2)
1916  dbta%srflag(:dim1,:dim2)=data%srflag(:dim1,:dim2)
1917  dbta%snwdph(:dim1,:dim2)=data%snwdph(:dim1,:dim2)
1918  dbta%slc(:dim1,:dim2,:dim3)=data%slc(:dim1,:dim2,:dim3)
1919  dbta%shdmin(:dim1,:dim2)=data%shdmin(:dim1,:dim2)
1920  dbta%shdmax(:dim1,:dim2)=data%shdmax(:dim1,:dim2)
1921  dbta%slope(:dim1,:dim2)=data%slope(:dim1,:dim2)
1922  dbta%snoalb(:dim1,:dim2)=data%snoalb(:dim1,:dim2)
1923  dbta%orog(:dim1,:dim2)=data%orog(:dim1,:dim2)
1924  call sfcio_axdata(data,iret)
1925 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1926  elseif(head%ivs == 200509) then
1927  read(lu,iostat=ios) dbta%slmsk(:dim1,:dim2)
1928  if(ios.ne.0) return
1929  read(lu,iostat=ios) dbta%orog(:dim1,:dim2)
1930  if(ios.ne.0) return
1931  read(lu,iostat=ios) dbta%tsea(:dim1,:dim2)
1932  if(ios.ne.0) return
1933  read(lu,iostat=ios) dbta%sheleg(:dim1,:dim2)
1934  if(ios.ne.0) return
1935  read(lu,iostat=ios) dbta%tg3(:dim1,:dim2)
1936  if(ios.ne.0) return
1937  read(lu,iostat=ios) dbta%zorl(:dim1,:dim2)
1938  if(ios.ne.0) return
1939  read(lu,iostat=ios) dbta%alvsf(:dim1,:dim2)
1940  if(ios.ne.0) return
1941  read(lu,iostat=ios) dbta%alvwf(:dim1,:dim2)
1942  if(ios.ne.0) return
1943  read(lu,iostat=ios) dbta%alnsf(:dim1,:dim2)
1944  if(ios.ne.0) return
1945  read(lu,iostat=ios) dbta%alnwf(:dim1,:dim2)
1946  if(ios.ne.0) return
1947  read(lu,iostat=ios) dbta%vfrac(:dim1,:dim2)
1948  if(ios.ne.0) return
1949  read(lu,iostat=ios) dbta%canopy(:dim1,:dim2)
1950  if(ios.ne.0) return
1951  read(lu,iostat=ios) dbta%f10m(:dim1,:dim2)
1952  if(ios.ne.0) return
1953  read(lu,iostat=ios) dbta%t2m(:dim1,:dim2)
1954  if(ios.ne.0) return
1955  read(lu,iostat=ios) dbta%q2m(:dim1,:dim2)
1956  if(ios.ne.0) return
1957  read(lu,iostat=ios) dbta%vtype(:dim1,:dim2)
1958  if(ios.ne.0) return
1959  read(lu,iostat=ios) dbta%stype(:dim1,:dim2)
1960  if(ios.ne.0) return
1961  read(lu,iostat=ios) dbta%facsf(:dim1,:dim2)
1962  if(ios.ne.0) return
1963  read(lu,iostat=ios) dbta%facwf(:dim1,:dim2)
1964  if(ios.ne.0) return
1965  read(lu,iostat=ios) dbta%uustar(:dim1,:dim2)
1966  if(ios.ne.0) return
1967  read(lu,iostat=ios) dbta%ffmm(:dim1,:dim2)
1968  if(ios.ne.0) return
1969  read(lu,iostat=ios) dbta%ffhh(:dim1,:dim2)
1970  if(ios.ne.0) return
1971  read(lu,iostat=ios) dbta%hice(:dim1,:dim2)
1972  if(ios.ne.0) return
1973  read(lu,iostat=ios) dbta%fice(:dim1,:dim2)
1974  if(ios.ne.0) return
1975  read(lu,iostat=ios) dbta%tisfc(:dim1,:dim2)
1976  if(ios.ne.0) return
1977  read(lu,iostat=ios) dbta%tprcp(:dim1,:dim2)
1978  if(ios.ne.0) return
1979  read(lu,iostat=ios) dbta%srflag(:dim1,:dim2)
1980  if(ios.ne.0) return
1981  read(lu,iostat=ios) dbta%snwdph(:dim1,:dim2)
1982  if(ios.ne.0) return
1983  read(lu,iostat=ios) dbta%shdmin(:dim1,:dim2)
1984  if(ios.ne.0) return
1985  read(lu,iostat=ios) dbta%shdmax(:dim1,:dim2)
1986  if(ios.ne.0) return
1987  read(lu,iostat=ios) dbta%slope(:dim1,:dim2)
1988  if(ios.ne.0) return
1989  read(lu,iostat=ios) dbta%snoalb(:dim1,:dim2)
1990  if(ios.ne.0) return
1991  do i=1,head%lsoil
1992  read(lu,iostat=ios) dbta%stc(:dim1,:dim2,i)
1993  if(ios.ne.0) return
1994  enddo
1995  do i=1,head%lsoil
1996  read(lu,iostat=ios) dbta%smc(:dim1,:dim2,i)
1997  if(ios.ne.0) return
1998  enddo
1999  do i=1,head%lsoil
2000  read(lu,iostat=ios) dbta%slc(:dim1,:dim2,i)
2001  if(ios.ne.0) return
2002  enddo
2003  dbta%cv(:dim1,:dim2)=sfcio_realfill
2004  dbta%cvb(:dim1,:dim2)=sfcio_realfill
2005  dbta%cvt(:dim1,:dim2)=sfcio_realfill
2006  endif
2007 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2008  iret=0
2009  end subroutine
2010 !-------------------------------------------------------------------------------
2011  subroutine sfcio_swdbta(lu,head,dbta,iret)
2012  implicit none
2013  integer(sfcio_intkind),intent(in):: lu
2014  type(sfcio_head),intent(in):: head
2015  type(sfcio_dbta),intent(in):: dbta
2016  integer(sfcio_intkind),intent(out):: iret
2017  integer:: dim1,dim2,dim3,mdim1,mdim2,mdim3
2018  integer:: ios
2019  integer i
2020  type(sfcio_data):: data
2021 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2022  dim1=head%lonb
2023  dim2=head%latb
2024  dim3=head%lsoil
2025  mdim1=min(&
2026  size(dbta%tsea,1),&
2027  size(dbta%smc,1),&
2028  size(dbta%sheleg,1),&
2029  size(dbta%stc,1),&
2030  size(dbta%tg3,1),&
2031  size(dbta%zorl,1),&
2032  size(dbta%alvsf,1),&
2033  size(dbta%alvwf,1),&
2034  size(dbta%alnsf,1),&
2035  size(dbta%alnwf,1),&
2036  size(dbta%slmsk,1),&
2037  size(dbta%vfrac,1),&
2038  size(dbta%canopy,1),&
2039  size(dbta%f10m,1),&
2040  size(dbta%t2m,1),&
2041  size(dbta%q2m,1),&
2042  size(dbta%vtype,1),&
2043  size(dbta%stype,1),&
2044  size(dbta%facsf,1),&
2045  size(dbta%facwf,1),&
2046  size(dbta%uustar,1),&
2047  size(dbta%ffmm,1),&
2048  size(dbta%ffhh,1),&
2049  size(dbta%hice,1),&
2050  size(dbta%fice,1),&
2051  size(dbta%tisfc,1),&
2052  size(dbta%tprcp,1),&
2053  size(dbta%srflag,1),&
2054  size(dbta%snwdph,1),&
2055  size(dbta%slc,1),&
2056  size(dbta%shdmin,1),&
2057  size(dbta%shdmax,1),&
2058  size(dbta%slope,1),&
2059  size(dbta%snoalb,1),&
2060  size(dbta%orog,1))
2061  mdim2=min(&
2062  size(dbta%tsea,2),&
2063  size(dbta%smc,2),&
2064  size(dbta%sheleg,2),&
2065  size(dbta%stc,2),&
2066  size(dbta%tg3,2),&
2067  size(dbta%zorl,2),&
2068  size(dbta%alvsf,2),&
2069  size(dbta%alvwf,2),&
2070  size(dbta%alnsf,2),&
2071  size(dbta%alnwf,2),&
2072  size(dbta%slmsk,2),&
2073  size(dbta%vfrac,2),&
2074  size(dbta%canopy,2),&
2075  size(dbta%f10m,2),&
2076  size(dbta%t2m,2),&
2077  size(dbta%q2m,2),&
2078  size(dbta%vtype,2),&
2079  size(dbta%stype,2),&
2080  size(dbta%facsf,2),&
2081  size(dbta%facwf,2),&
2082  size(dbta%uustar,2),&
2083  size(dbta%ffmm,2),&
2084  size(dbta%ffhh,2),&
2085  size(dbta%hice,2),&
2086  size(dbta%fice,2),&
2087  size(dbta%tisfc,2),&
2088  size(dbta%tprcp,2),&
2089  size(dbta%srflag,2),&
2090  size(dbta%snwdph,2),&
2091  size(dbta%slc,2),&
2092  size(dbta%shdmin,2),&
2093  size(dbta%shdmax,2),&
2094  size(dbta%slope,2),&
2095  size(dbta%snoalb,2),&
2096  size(dbta%orog,2))
2097  mdim3=min(&
2098  size(dbta%smc,3),&
2099  size(dbta%stc,3),&
2100  size(dbta%slc,3))
2101  iret=-5
2102  if(mdim1.lt.dim1.or.&
2103  mdim2.lt.dim2.or.&
2104  mdim3.lt.dim3) return
2105  iret=-4
2106 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2107  if(head%irealf.ne.2) then
2108  call sfcio_aldata(head,data,iret)
2109  if(iret.ne.0) return
2110  data%tsea(:dim1,:dim2)=dbta%tsea(:dim1,:dim2)
2111  data%smc(:dim1,:dim2,:dim3)=dbta%smc(:dim1,:dim2,:dim3)
2112  data%sheleg(:dim1,:dim2)=dbta%sheleg(:dim1,:dim2)
2113  data%stc(:dim1,:dim2,:dim3)=dbta%stc(:dim1,:dim2,:dim3)
2114  data%tg3(:dim1,:dim2)=dbta%tg3(:dim1,:dim2)
2115  data%zorl(:dim1,:dim2)=dbta%zorl(:dim1,:dim2)
2116  data%cv(:dim1,:dim2)=dbta%cv(:dim1,:dim2)
2117  data%cvb(:dim1,:dim2)=dbta%cvb(:dim1,:dim2)
2118  data%cvt(:dim1,:dim2)=dbta%cvt(:dim1,:dim2)
2119  data%alvsf(:dim1,:dim2)=dbta%alvsf(:dim1,:dim2)
2120  data%alvwf(:dim1,:dim2)=dbta%alvwf(:dim1,:dim2)
2121  data%alnsf(:dim1,:dim2)=dbta%alnsf(:dim1,:dim2)
2122  data%alnwf(:dim1,:dim2)=dbta%alnwf(:dim1,:dim2)
2123  data%slmsk(:dim1,:dim2)=dbta%slmsk(:dim1,:dim2)
2124  data%vfrac(:dim1,:dim2)=dbta%vfrac(:dim1,:dim2)
2125  data%canopy(:dim1,:dim2)=dbta%canopy(:dim1,:dim2)
2126  data%f10m(:dim1,:dim2)=dbta%f10m(:dim1,:dim2)
2127  data%t2m(:dim1,:dim2)=dbta%t2m(:dim1,:dim2)
2128  data%q2m(:dim1,:dim2)=dbta%q2m(:dim1,:dim2)
2129  data%vtype(:dim1,:dim2)=dbta%vtype(:dim1,:dim2)
2130  data%stype(:dim1,:dim2)=dbta%stype(:dim1,:dim2)
2131  data%facsf(:dim1,:dim2)=dbta%facsf(:dim1,:dim2)
2132  data%facwf(:dim1,:dim2)=dbta%facwf(:dim1,:dim2)
2133  data%uustar(:dim1,:dim2)=dbta%uustar(:dim1,:dim2)
2134  data%ffmm(:dim1,:dim2)=dbta%ffmm(:dim1,:dim2)
2135  data%ffhh(:dim1,:dim2)=dbta%ffhh(:dim1,:dim2)
2136  data%hice(:dim1,:dim2)=dbta%hice(:dim1,:dim2)
2137  data%fice(:dim1,:dim2)=dbta%fice(:dim1,:dim2)
2138  data%tisfc(:dim1,:dim2)=dbta%tisfc(:dim1,:dim2)
2139  data%tprcp(:dim1,:dim2)=dbta%tprcp(:dim1,:dim2)
2140  data%srflag(:dim1,:dim2)=dbta%srflag(:dim1,:dim2)
2141  data%snwdph(:dim1,:dim2)=dbta%snwdph(:dim1,:dim2)
2142  data%slc(:dim1,:dim2,:dim3)=dbta%slc(:dim1,:dim2,:dim3)
2143  data%shdmin(:dim1,:dim2)=dbta%shdmin(:dim1,:dim2)
2144  data%shdmax(:dim1,:dim2)=dbta%shdmax(:dim1,:dim2)
2145  data%slope(:dim1,:dim2)=dbta%slope(:dim1,:dim2)
2146  data%snoalb(:dim1,:dim2)=dbta%snoalb(:dim1,:dim2)
2147  data%orog(:dim1,:dim2)=dbta%orog(:dim1,:dim2)
2148  call sfcio_swdata(lu,head,data,iret)
2149  if(iret.ne.0) return
2150  call sfcio_axdata(data,iret)
2151 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2152  elseif(head%ivs == 200509) then
2153  write(lu,iostat=ios) dbta%slmsk(:dim1,:dim2)
2154  if(ios.ne.0) return
2155  write(lu,iostat=ios) dbta%orog(:dim1,:dim2)
2156  if(ios.ne.0) return
2157  write(lu,iostat=ios) dbta%tsea(:dim1,:dim2)
2158  if(ios.ne.0) return
2159  write(lu,iostat=ios) dbta%sheleg(:dim1,:dim2)
2160  if(ios.ne.0) return
2161  write(lu,iostat=ios) dbta%tg3(:dim1,:dim2)
2162  if(ios.ne.0) return
2163  write(lu,iostat=ios) dbta%zorl(:dim1,:dim2)
2164  if(ios.ne.0) return
2165  write(lu,iostat=ios) dbta%alvsf(:dim1,:dim2)
2166  if(ios.ne.0) return
2167  write(lu,iostat=ios) dbta%alvwf(:dim1,:dim2)
2168  if(ios.ne.0) return
2169  write(lu,iostat=ios) dbta%alnsf(:dim1,:dim2)
2170  if(ios.ne.0) return
2171  write(lu,iostat=ios) dbta%alnwf(:dim1,:dim2)
2172  if(ios.ne.0) return
2173  write(lu,iostat=ios) dbta%vfrac(:dim1,:dim2)
2174  if(ios.ne.0) return
2175  write(lu,iostat=ios) dbta%canopy(:dim1,:dim2)
2176  if(ios.ne.0) return
2177  write(lu,iostat=ios) dbta%f10m(:dim1,:dim2)
2178  if(ios.ne.0) return
2179  write(lu,iostat=ios) dbta%t2m(:dim1,:dim2)
2180  if(ios.ne.0) return
2181  write(lu,iostat=ios) dbta%q2m(:dim1,:dim2)
2182  if(ios.ne.0) return
2183  write(lu,iostat=ios) dbta%vtype(:dim1,:dim2)
2184  if(ios.ne.0) return
2185  write(lu,iostat=ios) dbta%stype(:dim1,:dim2)
2186  if(ios.ne.0) return
2187  write(lu,iostat=ios) dbta%facsf(:dim1,:dim2)
2188  if(ios.ne.0) return
2189  write(lu,iostat=ios) dbta%facwf(:dim1,:dim2)
2190  if(ios.ne.0) return
2191  write(lu,iostat=ios) dbta%uustar(:dim1,:dim2)
2192  if(ios.ne.0) return
2193  write(lu,iostat=ios) dbta%ffmm(:dim1,:dim2)
2194  if(ios.ne.0) return
2195  write(lu,iostat=ios) dbta%ffhh(:dim1,:dim2)
2196  if(ios.ne.0) return
2197  write(lu,iostat=ios) dbta%hice(:dim1,:dim2)
2198  if(ios.ne.0) return
2199  write(lu,iostat=ios) dbta%fice(:dim1,:dim2)
2200  if(ios.ne.0) return
2201  write(lu,iostat=ios) dbta%tisfc(:dim1,:dim2)
2202  if(ios.ne.0) return
2203  write(lu,iostat=ios) dbta%tprcp(:dim1,:dim2)
2204  if(ios.ne.0) return
2205  write(lu,iostat=ios) dbta%srflag(:dim1,:dim2)
2206  if(ios.ne.0) return
2207  write(lu,iostat=ios) dbta%snwdph(:dim1,:dim2)
2208  if(ios.ne.0) return
2209  write(lu,iostat=ios) dbta%shdmin(:dim1,:dim2)
2210  if(ios.ne.0) return
2211  write(lu,iostat=ios) dbta%shdmax(:dim1,:dim2)
2212  if(ios.ne.0) return
2213  write(lu,iostat=ios) dbta%slope(:dim1,:dim2)
2214  if(ios.ne.0) return
2215  write(lu,iostat=ios) dbta%snoalb(:dim1,:dim2)
2216  if(ios.ne.0) return
2217  do i=1,head%lsoil
2218  write(lu,iostat=ios) dbta%stc(:dim1,:dim2,i)
2219  if(ios.ne.0) return
2220  enddo
2221  do i=1,head%lsoil
2222  write(lu,iostat=ios) dbta%smc(:dim1,:dim2,i)
2223  if(ios.ne.0) return
2224  enddo
2225  do i=1,head%lsoil
2226  write(lu,iostat=ios) dbta%slc(:dim1,:dim2,i)
2227  if(ios.ne.0) return
2228  enddo
2229  endif
2230 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2231  iret=0
2232  end subroutine
2233 !-------------------------------------------------------------------------------
2234  subroutine sfcio_srohdcb(lu,cfname,head,dbta,iret)
2235  implicit none
2236  integer(sfcio_intkind),intent(in):: lu
2237  character*(*),intent(in):: cfname
2238  type(sfcio_head),intent(inout):: head
2239  type(sfcio_dbta),intent(inout):: dbta
2240  integer(sfcio_intkind),intent(out):: iret
2241 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2242  call sfcio_sropen(lu,cfname,iret)
2243  if(iret.ne.0) return
2244 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2245  call sfcio_srhead(lu,head,iret)
2246  if(iret.ne.0) return
2247 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2248  call sfcio_aldbta(head,dbta,iret)
2249  if(iret.ne.0) return
2250 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2251  call sfcio_srdbta(lu,head,dbta,iret)
2252  if(iret.ne.0) return
2253 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2254  call sfcio_sclose(lu,iret)
2255  if(iret.ne.0) return
2256 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2257  end subroutine
2258 !-------------------------------------------------------------------------------
2259  subroutine sfcio_swohdcb(lu,cfname,head,dbta,iret)
2260  implicit none
2261  integer(sfcio_intkind),intent(in):: lu
2262  character*(*),intent(in):: cfname
2263  type(sfcio_head),intent(in):: head
2264  type(sfcio_dbta),intent(in):: dbta
2265  integer(sfcio_intkind),intent(out):: iret
2266 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2267  call sfcio_swopen(lu,cfname,iret)
2268  if(iret.ne.0) return
2269 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2270  call sfcio_swhead(lu,head,iret)
2271  if(iret.ne.0) return
2272 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2273  call sfcio_swdbta(lu,head,dbta,iret)
2274  if(iret.ne.0) return
2275 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2276  call sfcio_sclose(lu,iret)
2277  if(iret.ne.0) return
2278 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2279  end subroutine
2280 !-------------------------------------------------------------------------------
2281 end module
sfcio_module::sfcio_head
Definition: sfcio_module.f:408
sfcio_module::sfcio_srohdc
Definition: sfcio_module.f:502
sfcio_module::sfcio_dbta
Definition: sfcio_module.f:456
sfcio_module::sfcio_data
Definition: sfcio_module.f:416