define_nc Subroutine

public subroutine define_nc(ncID, nVar, sx)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ncID
integer, intent(in) :: nVar
character(len=*), intent(in) :: sx(nVar,4)

Calls

proc~~define_nc~~CallsGraph proc~define_nc define_nc nf90_def_var nf90_def_var proc~define_nc->nf90_def_var nf90_enddef nf90_enddef proc~define_nc->nf90_enddef nf90_inq_dimid nf90_inq_dimid proc~define_nc->nf90_inq_dimid nf90_inq_varid nf90_inq_varid proc~define_nc->nf90_inq_varid nf90_put_att nf90_put_att proc~define_nc->nf90_put_att nf90_redef nf90_redef proc~define_nc->nf90_redef proc~nchandle_error nchandle_error proc~define_nc->proc~nchandle_error nf90_strerror nf90_strerror proc~nchandle_error->nf90_strerror

Called by

proc~~define_nc~~CalledByGraph proc~define_nc define_nc proc~initeb initEB proc~initeb->proc~define_nc proc~initfielddump initfielddump proc~initfielddump->proc~define_nc proc~initibm initibm proc~initibm->proc~define_nc proc~initstatsdump initstatsdump proc~initstatsdump->proc~define_nc program~dalesurban DALESURBAN program~dalesurban->proc~initeb program~dalesurban->proc~initfielddump program~dalesurban->proc~initibm program~dalesurban->proc~initstatsdump

Source Code

  subroutine define_nc(ncID, nVar, sx)
    implicit none
    integer, intent (in) :: nVar, ncID
    character (*), intent (in) :: sx(nVar,4)

    integer, save ::  dim_mttt(4) = 0, dim_tmtt(4) = 0, dim_ttmt(4) = 0, dim_tttt(4) = 0, &
                      dim_tt(2)= 0, dim_mt(2)= 0,dim_t0tt(3)=0,dim_m0tt(3)=0,dim_t0mt(3)=0,&
                      dim_m0mt(3)=0, dim_tt0t(3)=0, &
                      dim_mt0t(3)=0,dim_tm0t(3)=0,dim_0ttt(3)=0,dim_0mtt(3)=0,dim_0tmt(3)=0,&
                      dim_tts(2)=0,dim_t0tts(3)=0,dim_0ttts(3)=0,dim_tttts(4)=0,dim_ttt0(3)=0,& !tg3315 added last one
                      dim_mtmt(4),dim_tmmt(4),dim_mmtt(4),& !bss116
                      dim_ft(2), dim_flt(3)!SO

    integer :: iret, n, VarID
    !write(*,*) 'definenc'
    iret = nf90_inq_dimid(ncid,'time',timeId)
    iret = nf90_inq_dimid(ncid,'xt',xtId)
    iret = nf90_inq_dimid(ncid,'xm',xmId)
    iret = nf90_inq_dimid(ncid,'yt',ytId)
    iret = nf90_inq_dimid(ncid,'ym',ymId)
    iret = nf90_inq_dimid(ncid,'zt',ztId)
    iret = nf90_inq_dimid(ncid,'zm',zmId)
    iret = nf90_inq_dimid(ncid,'zts',ztsId)
    iret = nf90_inq_dimid(ncid,'fct',fctId) ! so4718 for energy balance output
    iret = nf90_inq_dimid(ncid,'lyr',lyrId) ! so4718 for energy balance output

    iret = nf90_redef(ncid)
    dim_tt = (/ztId,timeId/)
    dim_mt = (/zmId,timeId/)

    dim_t0tt= (/xtID,ztID,timeId/)! thermo point
    dim_t0mt= (/xtID,zmID,timeId/)! zpoint
    dim_m0tt= (/xmID,ztID,timeId/)! upoint
    dim_m0mt= (/xmID,ztID,timeId/)! uw stats point
    dim_tt0t= (/xtID,ytID,timeId/)! thermo point
    dim_tm0t= (/xtID,ymID,timeId/)! vpoint
    dim_mt0t= (/xmID,ytID,timeId/)! upoint
    dim_0ttt= (/ytID,ztID,timeId/)! thermo point
    dim_0tmt= (/ytID,zmID,timeId/)! wpoint
    dim_0mtt= (/ymID,ztID,timeId/)! vpoint

    dim_tttt= (/xtID,ytID,ztID,timeId/)! thermo point
    dim_ttmt= (/xtID,ytID,zmID,timeId/)! zpoint
    dim_mttt= (/xmID,ytID,ztID,timeId/)! upoint
    dim_tmtt= (/xtID,ymID,ztId,timeId/)! ypoint
    dim_mtmt= (/xmID,ytID,zmId,timeId/)! uw stats point bss116
    dim_tmmt= (/xtID,ymID,zmId,timeId/)! vw stats point bss116
    dim_mmtt= (/xmID,ymID,ztId,timeId/)! uv stats point bss116

    dim_ttt0= (/xtID,ytID,ztID/)! stats point tg3315

    dim_tts = (/ztsId,timeId/)
    dim_t0tts= (/xtID,ztsID,timeId/)! thermo soil point
    dim_0ttts= (/ytID,ztsID,timeId/)! thermo point
    dim_tttts= (/xtID,ytID,ztsID,timeId/)! thermo point

    dim_ft = (/fctID,timeId/)
    dim_flt = (/fctID,lyrID,timeId/)

    do n=1,nVar
!      write(*,*) 'n', n
!      write(*,*) "dummyline1"
!      write(*,*) 'sx1', sx(1,:)
!      write(*,*) 'sx2', sx(2,:)
!      write(*,*) "dummyline2"
!      write(*,*) 'trim(sx(n,1))', trim(sx(n,1))
!      write(*,*) 'trim(sx(n,2))', trim(sx(n,2))
!      write(*,*) 'trim(sx(n,3))', trim(sx(n,3))
!      write(*,*) 'trim(sx(n,4))', trim(sx(n,4))
!      write (*,*) 'ncID', ncID
      iret = nf90_inq_varid(ncid, trim(sx(n,1)), VarID)
      if (iret == 0) cycle
      select case(trim(sx(n,4)))
        case ('time')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,(/timeID/) ,VarID)
        case ('tt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tt ,VarID)
        case ('mt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_mt,VarID)
  !2D Fields
        case ('t0tt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_t0tt,VarID)
        case ('t0mt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_t0mt,VarID)
        case ('m0tt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_m0tt,VarID)
        case ('m0mt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_m0mt,VarID)
        case ('tt0t')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tt0t,VarID)
        case ('tm0t')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tm0t,VarID)
        case ('mt0t')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_mt0t,VarID)
        case ('0ttt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_0ttt,VarID)
        case ('0tmt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_0tmt,VarID)
        case ('0mtt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_0mtt,VarID)
        case ('ttt0')                  !tg3315 for uav,vav,wav etc.
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_ttt0,VarID)
  !3D Fields
        case ('tttt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tttt,VarID)
        case ('mttt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_mttt,VarID)
        case ('tmtt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tmtt,VarID)
        case ('ttmt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_ttmt,VarID)
        case ('mtmt')                                                   ! bss116
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_mtmt,VarID)
        case ('tmmt')                                                   ! bss116
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tmmt,VarID)
        case ('mmtt')                                                   ! bss116
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_mmtt,VarID)
!Soil fields
        case ('tts')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tts ,VarID)
        case ('t0tts')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_t0tts,VarID)
        case ('0ttts')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_0ttts,VarID)
        case ('tttts')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_tttts,VarID)

!Facet information
       case ('ft')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_ft,VarID)
       case ('flt')
          iret=nf90_def_var(ncID,sx(n,1),NF90_FLOAT,dim_flt,VarID)

        case default
        write(0, *) 'nvar', nvar, sx(n,:)
        write(0, *) 'ERROR: Bad dimensional information ',sx(n,:)
        stop 1
        ! call appl_abort(0)
      end select
      if (iret/=0) then
!        write (*,*) 'nvar', nvar, sx(n,:)
!        write (*,*) 'ncID', ncID
        call nchandle_error(iret)
      end if
      iret = nf90_put_att(ncID,VarID,'longname',sx(n,2))
      iret = nf90_put_att(ncID,VarID,'units',sx(n,3))
      iret = nf90_put_att(ncid, VarID, '_FillValue',nc_fillvalue)

    end do
    iret= nf90_enddef(ncID)
  end subroutine define_nc