open_nc Subroutine

public subroutine open_nc(fname, ncid, nrec, n1, n2, n3, ns, nfcts, nlyrs)

Uses

  • proc~~open_nc~~UsesGraph proc~open_nc modstat_nc::open_nc module~modglobal modglobal proc~open_nc->module~modglobal

Arguments

Type IntentOptional Attributes Name
character(len=40), intent(in) :: fname
integer, intent(out) :: ncid
integer, intent(out) :: nrec
integer, intent(in), optional :: n1
integer, intent(in), optional :: n2
integer, intent(in), optional :: n3
integer, intent(in), optional :: ns
integer, intent(in), optional :: nfcts
integer, intent(in), optional :: nlyrs

Calls

proc~~open_nc~~CallsGraph proc~open_nc modstat_nc::open_nc nf90_create nf90_create proc~open_nc->nf90_create nf90_def_dim nf90_def_dim proc~open_nc->nf90_def_dim nf90_def_var nf90_def_var proc~open_nc->nf90_def_var nf90_enddef nf90_enddef proc~open_nc->nf90_enddef nf90_get_var nf90_get_var proc~open_nc->nf90_get_var nf90_inq_dimid nf90_inq_dimid proc~open_nc->nf90_inq_dimid nf90_inq_varid nf90_inq_varid proc~open_nc->nf90_inq_varid nf90_inquire nf90_inquire proc~open_nc->nf90_inquire nf90_inquire_dimension nf90_inquire_dimension proc~open_nc->nf90_inquire_dimension nf90_open nf90_open proc~open_nc->nf90_open nf90_put_att nf90_put_att proc~open_nc->nf90_put_att nf90_sync nf90_sync proc~open_nc->nf90_sync

Called by

proc~~open_nc~~CalledByGraph proc~open_nc modstat_nc::open_nc proc~initeb modEB::initEB proc~initeb->proc~open_nc proc~initfielddump modfielddump::initfielddump proc~initfielddump->proc~open_nc proc~initstatsdump modstatsdump::initstatsdump proc~initstatsdump->proc~open_nc program~dalesurban DALESURBAN program~dalesurban->proc~initeb program~dalesurban->proc~initfielddump program~dalesurban->proc~initstatsdump

Contents

Source Code


Source Code

  subroutine open_nc (fname, ncid,nrec,n1, n2, n3, ns, nfcts, nlyrs)
    use modglobal, only : author,version,timee
    implicit none
    integer, intent (out) :: ncid,nrec
    integer, optional, intent (in) :: n1, n2, n3, ns, nfcts, nlyrs
    character (len=40), intent (in) :: fname

    character (len=12):: date='',time=''
    integer :: iret,varid,ncall,RecordDimID
    real, allocatable :: xtimes(:)
    logical :: exans

    inquire(file=trim(fname),exist=exans)
    !write(*,*) 'opennc'
    ncall = 0
    if (.not.exans) then

      call date_and_time(date,time)
      iret = nf90_create(fname,NF90_SHARE,ncid)
      iret = nf90_put_att(ncid,NF90_GLOBAL,'title',fname)
      iret = nf90_put_att(ncid,NF90_GLOBAL,'history','Created on '//trim(date)//' at '//trim(time))
      iret = nf90_put_att(ncid, NF90_GLOBAL, 'Source',trim(version))
      iret = nf90_put_att(ncid, NF90_GLOBAL, 'Author',trim(author))
      iret = nf90_def_dim(ncID, 'time', NF90_UNLIMITED, timeID)
      if (present(n1)) then
        iret = nf90_def_dim(ncID, 'xt', n1, xtID)
        iret = nf90_def_dim(ncID, 'xm', n1, xmID)
        iret = nf90_def_var(ncID,'xt',NF90_FLOAT,(/xtID/) ,VarID)
        iret=nf90_put_att(ncID,VarID,'longname','West-East displacement of cell centers')
        iret=nf90_put_att(ncID,VarID,'units','m')
        iret = nf90_def_var(ncID,'xm',NF90_FLOAT,(/xmID/),VarID)
        iret=nf90_put_att(ncID,VarID,'longname','West-East displacement of cell edges')
        iret=nf90_put_att(ncID,VarID,'units','m')
      end if
      if (present(n2)) then
        iret = nf90_def_dim(ncID, 'yt', n2, ytID)
        iret = nf90_def_dim(ncID, 'ym', n2, ymID)
        iret = nf90_def_var(ncID,'yt',NF90_FLOAT,ytID ,VarID)
        iret=nf90_put_att(ncID,VarID,'longname','South-North displacement of cell centers')
        iret=nf90_put_att(ncID,VarID,'units','m')
        iret = nf90_def_var(ncID,'ym',NF90_FLOAT,ymID,VarID)
        iret=nf90_put_att(ncID,VarID,'longname','South-North displacement of cell edges')
        iret=nf90_put_att(ncID,VarID,'units','m')
      end if
      if (present(n3)) then
        iret = nf90_def_dim(ncID, 'zt', n3, ztID)
        iret = nf90_def_dim(ncID, 'zm', n3, zmID)
        iret = nf90_def_var(ncID,'zt',NF90_FLOAT,(/ztID/) ,VarID)
        iret=nf90_put_att(ncID,VarID,'longname','Vertical displacement of cell centers')
        iret=nf90_put_att(ncID,VarID,'units','m')
        iret = nf90_def_var(ncID,'zm',NF90_FLOAT,(/zmID/),VarID)
        iret=nf90_put_att(ncID,VarID,'longname','Vertical displacement of cell edges')
        iret=nf90_put_att(ncID,VarID,'units','m')
      end if
      if (present(ns)) then
        iret = nf90_def_dim(ncID, 'zts', ns, ztsID)
        iret = nf90_def_var(ncID,'zts',NF90_FLOAT,(/ztsID/) ,VarID)
        iret=nf90_put_att(ncID,VarID,'longname','Soil level depth of cell centers')
        iret=nf90_put_att(ncID,VarID,'units','m')
      end if
      if (present(nfcts)) then
        iret = nf90_def_dim(ncID, 'fct', nfcts, fctID)
        iret = nf90_def_var(ncID, 'fct',NF90_INT,(/fctID/) ,VarID)
        iret=nf90_put_att(ncID,VarID,'longname','Facet number')
      end if
      if (present(nlyrs)) then
        iret = nf90_def_dim(ncID, 'lyr', nlyrs, lyrID)
        iret = nf90_def_var(ncID, 'lyr',NF90_INT,(/lyrID/) ,VarID)
        iret=nf90_put_att(ncID,VarID,'longname','Number of wall layers')
      end if

    else
       nrec = 0
       ncall= 0
       iret = nf90_open (trim(fname), NF90_WRITE, ncid)
       iret = nf90_inquire(ncid, unlimitedDimId = RecordDimID)
       iret = nf90_inquire_dimension(ncid, RecordDimID, len=nrec)
       if (nrec>0) then
        iret = nf90_inq_varid(ncid,'time',timeID)
        allocate (xtimes(nrec))
        iret = nf90_get_var(ncid, timeId, xtimes(1:nrec))

        do while(xtimes(ncall+1) < timee - spacing(1.))
            ncall=ncall+1
            if (ncall >= nrec) exit 
        end do
        deallocate(xtimes)
       end if
       if (present(n1)) then
         iret = nf90_inq_dimid(ncid,'xt',xtId)
         iret = nf90_inq_dimid(ncid,'xm',xmId)
       end if
       if (present(n2)) then
         iret = nf90_inq_dimid(ncid,'yt',ytId)
         iret = nf90_inq_dimid(ncid,'ym',ymId)
       end if
       if (present(n3)) then
         iret = nf90_inq_dimid(ncid,'zt',ztId)
         iret = nf90_inq_dimid(ncid,'zm',zmId)
       end if
       if (present(ns)) then
         iret = nf90_inq_dimid(ncid,'zts',ztsId)
       end if
       if (present(nfcts)) then
         iret = nf90_inq_dimid(ncid,'fct',fctId)
       end if
    end if
    nrec = ncall
    iret = nf90_sync(ncid)

    iret= nf90_enddef(ncID)

  end subroutine open_nc