1 !*------------------------------------------------------------------------------
4 !* Forecast Systems Laboratory
10 !* ADVANCED COMPUTING BRANCH
11 !* SMS/NNT Version: 2.0.0
13 !* This software and its documentation are in the public domain and
14 !* are furnished "as is". The United States government, its
15 !* instrumentalities, officers, employees, and agents make no
16 !* warranty, express or implied, as to the usefulness of the software
17 !* and documentation for any purpose. They assume no
18 !* responsibility (1) for the use of the software and documentation;
19 !* or (2) to provide technical support to users.
21 !* Permission to use, copy, modify, and distribute this software is
22 !* hereby granted, provided that this disclaimer notice appears in
23 !* all copies. All modifications to this software must be clearly
24 !* documented, and are solely the responsibility of the agent making
25 !* the modification. If significant modifications or enhancements
26 !* are made to this software, the SMS Development team
27 !* (sms-info@fsl.noaa.gov) should be notified.
29 !*----------------------------------------------------------------------------
32 ! Author: Jacques Middlecoff jacquesm@fsl.noaa.gov
33 !* Date: October 6, 2000
35 !*----------------------------------------------------------------------------
39 integer , parameter :: FATAL = 1
40 integer , parameter :: WARN = 1
41 integer , parameter :: WrfDataHandleMax = 99
42 integer , parameter :: MaxDims = 2000 ! = NF_MAX_VARS
44 integer , parameter :: MaxVars = 10000
46 integer , parameter :: MaxVars = 3000
48 integer , parameter :: MaxTimes = 60000
49 integer , parameter :: DateStrLen = 19
50 integer , parameter :: VarNameLen = 31
51 integer , parameter :: NO_DIM = 0
52 integer , parameter :: NVarDims = 4
53 integer , parameter :: NMDVarDims = 2
54 character (8) , parameter :: NO_NAME = 'NULL'
55 character (DateStrLen) , parameter :: ZeroDate = '0000-00-00-00:00:00'
57 #include "wrf_io_flags.h"
59 character (256) :: msg
60 logical :: WrfIOnotInitialized = .true.
62 type :: wrf_data_handle
63 character (255) :: FileName
69 character (5) :: TimesName
71 integer :: CurrentTime !Only used for read
72 integer :: NumberTimes !Only used for read
73 character (DateStrLen), pointer :: Times(:)
75 integer , pointer :: DimLengths(:)
76 integer , pointer :: DimIDs(:)
77 character (31) , pointer :: DimNames(:)
79 character (9) :: DimUnlimName
80 integer , dimension(NVarDims) :: DimID
81 integer , dimension(NVarDims) :: Dimension
82 integer , pointer :: MDVarIDs(:)
83 integer , pointer :: MDVarDimLens(:)
84 character (80) , pointer :: MDVarNames(:)
85 integer , pointer :: VarIDs(:)
86 integer , pointer :: VarDimLens(:,:)
87 character (VarNameLen), pointer :: VarNames(:)
88 integer :: CurrentVariable !Only used for read
90 ! first_operation is set to .TRUE. when a new handle is allocated
91 ! or when open-for-write or open-for-read are committed. It is set
92 ! to .FALSE. when the first field is read or written.
93 logical :: first_operation
96 logical :: use_netcdf_classic
97 end type wrf_data_handle
98 type(wrf_data_handle),target :: WrfDataHandles(WrfDataHandleMax)
101 module ext_ncd_support_routines
107 subroutine allocHandle(DataHandle,DH,Comm,Status)
109 include 'wrf_status_codes.h'
110 integer ,intent(out) :: DataHandle
111 type(wrf_data_handle),pointer :: DH
112 integer ,intent(IN) :: Comm
113 integer ,intent(out) :: Status
117 do i=1,WrfDataHandleMax
118 if(WrfDataHandles(i)%Free) then
119 DH => WrfDataHandles(i)
121 allocate(DH%Times(MaxTimes), STAT=stat)
123 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
124 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
125 call wrf_debug ( FATAL , msg)
128 allocate(DH%DimLengths(MaxDims), STAT=stat)
130 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
131 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
132 call wrf_debug ( FATAL , msg)
135 allocate(DH%DimIDs(MaxDims), STAT=stat)
137 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
138 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
139 call wrf_debug ( FATAL , msg)
142 allocate(DH%DimNames(MaxDims), STAT=stat)
144 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
145 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
146 call wrf_debug ( FATAL , msg)
149 allocate(DH%MDVarIDs(MaxVars), STAT=stat)
151 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
152 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
153 call wrf_debug ( FATAL , msg)
156 allocate(DH%MDVarDimLens(MaxVars), STAT=stat)
158 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
159 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
160 call wrf_debug ( FATAL , msg)
163 allocate(DH%MDVarNames(MaxVars), STAT=stat)
165 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
166 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
167 call wrf_debug ( FATAL , msg)
170 allocate(DH%VarIDs(MaxVars), STAT=stat)
172 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
173 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
174 call wrf_debug ( FATAL , msg)
177 allocate(DH%VarDimLens(NVarDims-1,MaxVars), STAT=stat)
179 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
180 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
181 call wrf_debug ( FATAL , msg)
184 allocate(DH%VarNames(MaxVars), STAT=stat)
186 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
187 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
188 call wrf_debug ( FATAL , msg)
193 if(i==WrfDataHandleMax) then
194 Status = WRF_WARN_TOO_MANY_FILES
195 write(msg,*) 'Warning TOO MANY FILES in ',__FILE__,', line', __LINE__
196 call wrf_debug ( WARN , TRIM(msg))
197 write(msg,*) 'Did you call ext_ncd_ioinit?'
198 call wrf_debug ( WARN , TRIM(msg))
205 DH%first_operation = .TRUE.
206 DH%R4OnOutput = .false.
209 end subroutine allocHandle
211 subroutine deallocHandle(DataHandle, Status)
213 include 'wrf_status_codes.h'
214 integer ,intent(in) :: DataHandle
215 integer ,intent(out) :: Status
216 type(wrf_data_handle),pointer :: DH
220 IF ( DataHandle .GE. 1 .AND. DataHandle .LE. WrfDataHandleMax ) THEN
221 if(.NOT. WrfDataHandles(DataHandle)%Free) then
222 DH => WrfDataHandles(DataHandle)
223 deallocate(DH%Times, STAT=stat)
225 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
226 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
227 call wrf_debug ( FATAL , msg)
230 deallocate(DH%DimLengths, STAT=stat)
232 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
233 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
234 call wrf_debug ( FATAL , msg)
237 deallocate(DH%DimIDs, STAT=stat)
239 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
240 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
241 call wrf_debug ( FATAL , msg)
244 deallocate(DH%DimNames, STAT=stat)
246 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
247 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
248 call wrf_debug ( FATAL , msg)
251 deallocate(DH%MDVarIDs, STAT=stat)
253 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
254 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
255 call wrf_debug ( FATAL , msg)
258 deallocate(DH%MDVarDimLens, STAT=stat)
260 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
261 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
262 call wrf_debug ( FATAL , msg)
265 deallocate(DH%MDVarNames, STAT=stat)
267 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
268 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
269 call wrf_debug ( FATAL , msg)
272 deallocate(DH%VarIDs, STAT=stat)
274 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
275 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
276 call wrf_debug ( FATAL , msg)
279 deallocate(DH%VarDimLens, STAT=stat)
281 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
282 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
283 call wrf_debug ( FATAL , msg)
286 deallocate(DH%VarNames, STAT=stat)
288 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
289 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
290 call wrf_debug ( FATAL , msg)
297 end subroutine deallocHandle
299 subroutine GetDH(DataHandle,DH,Status)
301 include 'wrf_status_codes.h'
302 integer ,intent(in) :: DataHandle
303 type(wrf_data_handle) ,pointer :: DH
304 integer ,intent(out) :: Status
306 if(DataHandle < 1 .or. DataHandle > WrfDataHandleMax) then
307 Status = WRF_WARN_BAD_DATA_HANDLE
310 DH => WrfDataHandles(DataHandle)
312 Status = WRF_WARN_BAD_DATA_HANDLE
319 subroutine DateCheck(Date,Status)
321 include 'wrf_status_codes.h'
322 character*(*) ,intent(in) :: Date
323 integer ,intent(out) :: Status
325 if(len(Date) /= DateStrLen) then
326 Status = WRF_WARN_DATESTR_BAD_LENGTH
331 end subroutine DateCheck
333 subroutine GetName(Element,Var,Name,Status)
335 include 'wrf_status_codes.h'
336 character*(*) ,intent(in) :: Element
337 character*(*) ,intent(in) :: Var
338 character*(*) ,intent(out) :: Name
339 integer ,intent(out) :: Status
340 character (VarNameLen) :: VarName
343 integer, parameter :: upper_to_lower =IACHAR('a')-IACHAR('A')
346 Name = 'MD___'//trim(Element)//VarName
349 if('A'<=c .and. c <='Z') Name(i:i)=achar(iachar(c)+upper_to_lower)
350 if(c=='-'.or.c==':') Name(i:i)='_'
354 end subroutine GetName
356 subroutine GetTimeIndex(IO,DataHandle,DateStr,TimeIndex,Status)
358 include 'wrf_status_codes.h'
360 character (*) ,intent(in) :: IO
361 integer ,intent(in) :: DataHandle
362 character*(*) ,intent(in) :: DateStr
363 integer ,intent(out) :: TimeIndex
364 integer ,intent(out) :: Status
365 type(wrf_data_handle) ,pointer :: DH
371 DH => WrfDataHandles(DataHandle)
372 call DateCheck(DateStr,Status)
373 if(Status /= WRF_NO_ERR) then
374 Status = WRF_WARN_DATESTR_ERROR
375 write(msg,*) 'Warning DATE STRING ERROR in ',__FILE__,', line', __LINE__
376 call wrf_debug ( WARN , TRIM(msg))
379 if(IO == 'write') then
380 TimeIndex = DH%TimeIndex
381 if(TimeIndex <= 0) then
383 elseif(DateStr == DH%Times(TimeIndex)) then
387 TimeIndex = TimeIndex +1
388 if(TimeIndex > MaxTimes) then
389 Status = WRF_WARN_TIME_EOF
390 write(msg,*) 'Warning TIME EOF in ',__FILE__,', line', __LINE__
391 call wrf_debug ( WARN , TRIM(msg))
395 DH%TimeIndex = TimeIndex
396 DH%Times(TimeIndex) = DateStr
398 VStart(2) = TimeIndex
399 VCount(1) = DateStrLen
401 stat = NF_PUT_VARA_TEXT(DH%NCID,DH%TimesVarID,VStart,VCount,DateStr)
402 call netcdf_err(stat,Status)
403 if(Status /= WRF_NO_ERR) then
404 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
405 call wrf_debug ( WARN , TRIM(msg))
410 if(DH%Times(i)==DateStr) then
416 Status = WRF_WARN_TIME_NF
417 write(msg,*) 'Warning TIME ',DateStr,' NOT FOUND in ',__FILE__,', line', __LINE__
418 call wrf_debug ( WARN , TRIM(msg))
424 end subroutine GetTimeIndex
426 subroutine GetDim(MemoryOrder,NDim,Status)
427 include 'wrf_status_codes.h'
428 character*(*) ,intent(in) :: MemoryOrder
429 integer ,intent(out) :: NDim
430 integer ,intent(out) :: Status
431 character*3 :: MemOrd
433 call LowerCase(MemoryOrder,MemOrd)
435 case ('xyz','xzy','yxz','yzx','zxy','zyx','xsz','xez','ysz','yez')
437 case ('xy','yx','xs','xe','ys','ye','cc')
441 case ('0') ! NDim=0 for scalars. TBH: 20060502
444 Status = WRF_WARN_BAD_MEMORYORDER
449 end subroutine GetDim
451 #ifdef USE_NETCDF4_FEATURES
452 subroutine set_chunking(MemoryOrder,need_chunking)
453 include 'wrf_status_codes.h'
454 character*(*) ,intent(in) :: MemoryOrder
455 logical ,intent(out) :: need_chunking
456 character*3 :: MemOrd
458 call LowerCase(MemoryOrder,MemOrd)
459 if(len(MemOrd) >= 2) then
461 case ('xyz','xzy','yxz','yzx','zxy','zyx','xsz','xez','ysz','yez')
462 need_chunking = .true.
464 need_chunking = .true.
466 need_chunking = .false.
470 end subroutine set_chunking
473 subroutine GetIndices(NDim,Start,End,i1,i2,j1,j2,k1,k2)
474 integer ,intent(in) :: NDim
475 integer ,dimension(*),intent(in) :: Start,End
476 integer ,intent(out) :: i1,i2,j1,j2,k1,k2
484 if(NDim == 0) return ! NDim=0 for scalars. TBH: 20060502
494 end subroutine GetIndices
496 logical function ZeroLengthHorzDim(MemoryOrder,Vector,Status)
498 include 'wrf_status_codes.h'
499 character*(*) ,intent(in) :: MemoryOrder
500 integer,dimension(*) ,intent(in) :: Vector
501 integer ,intent(out) :: Status
503 integer,dimension(NVarDims) :: temp
504 character*3 :: MemOrd
507 call GetDim(MemoryOrder,NDim,Status)
508 temp(1:NDim) = Vector(1:NDim)
509 call LowerCase(MemoryOrder,MemOrd)
510 zero_length = .false.
512 case ('xsz','xez','ysz','yez','xs','xe','ys','ye','z','c')
515 continue ! NDim=0 for scalars. TBH: 20060502
517 zero_length = temp(1) .lt. 1 .or. temp(3) .lt. 1
518 case ('xy','yx','xyz','yxz')
519 zero_length = temp(1) .lt. 1 .or. temp(2) .lt. 1
521 zero_length = temp(2) .lt. 1 .or. temp(3) .lt. 1
523 Status = WRF_WARN_BAD_MEMORYORDER
524 ZeroLengthHorzDim = .true.
528 ZeroLengthHorzDim = zero_length
530 end function ZeroLengthHorzDim
532 subroutine ExtOrder(MemoryOrder,Vector,Status)
534 include 'wrf_status_codes.h'
535 character*(*) ,intent(in) :: MemoryOrder
536 integer,dimension(*) ,intent(inout) :: Vector
537 integer ,intent(out) :: Status
539 integer,dimension(NVarDims) :: temp
540 character*3 :: MemOrd
542 call GetDim(MemoryOrder,NDim,Status)
543 temp(1:NDim) = Vector(1:NDim)
544 call LowerCase(MemoryOrder,MemOrd)
547 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c')
550 continue ! NDim=0 for scalars. TBH: 20060502
572 Status = WRF_WARN_BAD_MEMORYORDER
577 end subroutine ExtOrder
579 subroutine ExtOrderStr(MemoryOrder,Vector,ROVector,Status)
581 include 'wrf_status_codes.h'
582 character*(*) ,intent(in) :: MemoryOrder
583 character*(*),dimension(*) ,intent(in) :: Vector
584 character(80),dimension(NVarDims),intent(out) :: ROVector
585 integer ,intent(out) :: Status
587 character*3 :: MemOrd
589 call GetDim(MemoryOrder,NDim,Status)
590 ROVector(1:NDim) = Vector(1:NDim)
591 call LowerCase(MemoryOrder,MemOrd)
594 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c')
597 continue ! NDim=0 for scalars. TBH: 20060502
599 ROVector(2) = Vector(3)
600 ROVector(3) = Vector(2)
602 ROVector(1) = Vector(2)
603 ROVector(2) = Vector(1)
605 ROVector(1) = Vector(3)
606 ROVector(2) = Vector(1)
607 ROVector(3) = Vector(2)
609 ROVector(1) = Vector(2)
610 ROVector(2) = Vector(3)
611 ROVector(3) = Vector(1)
613 ROVector(1) = Vector(3)
614 ROVector(3) = Vector(1)
616 ROVector(1) = Vector(2)
617 ROVector(2) = Vector(1)
619 Status = WRF_WARN_BAD_MEMORYORDER
624 end subroutine ExtOrderStr
627 subroutine LowerCase(MemoryOrder,MemOrd)
628 character*(*) ,intent(in) :: MemoryOrder
629 character*(*) ,intent(out) :: MemOrd
631 integer ,parameter :: upper_to_lower =IACHAR('a')-IACHAR('A')
636 MemOrd(1:N) = MemoryOrder(1:N)
639 if('A'<=c .and. c <='Z') MemOrd(i:i)=achar(iachar(c)+upper_to_lower)
642 end subroutine LowerCase
644 subroutine UpperCase(MemoryOrder,MemOrd)
645 character*(*) ,intent(in) :: MemoryOrder
646 character*(*) ,intent(out) :: MemOrd
648 integer ,parameter :: lower_to_upper =IACHAR('A')-IACHAR('a')
653 MemOrd(1:N) = MemoryOrder(1:N)
656 if('a'<=c .and. c <='z') MemOrd(i:i)=achar(iachar(c)+lower_to_upper)
659 end subroutine UpperCase
661 subroutine netcdf_err(err,Status)
663 include 'wrf_status_codes.h'
665 integer ,intent(in) :: err
666 integer ,intent(out) :: Status
667 character(len=80) :: errmsg
670 if( err==NF_NOERR )then
673 errmsg = NF_STRERROR(err)
674 write(msg,*) 'NetCDF error: ',errmsg
675 call wrf_debug ( WARN , TRIM(msg))
676 Status = WRF_WARN_NETCDF
679 end subroutine netcdf_err
681 subroutine FieldIO(IO,DataHandle,DateStr,Length,MemoryOrder &
682 ,FieldType,NCID,VarID,XField,Status)
684 include 'wrf_status_codes.h'
686 character (*) ,intent(in) :: IO
687 integer ,intent(in) :: DataHandle
688 character*(*) ,intent(in) :: DateStr
689 integer,dimension(NVarDims),intent(in) :: Length
690 character*(*) ,intent(in) :: MemoryOrder
691 integer ,intent(in) :: FieldType
692 integer ,intent(in) :: NCID
693 integer ,intent(in) :: VarID
694 integer,dimension(*) ,intent(inout) :: XField
695 integer ,intent(out) :: Status
698 integer,dimension(NVarDims) :: VStart
699 integer,dimension(NVarDims) :: VCount
701 call GetTimeIndex(IO,DataHandle,DateStr,TimeIndex,Status)
702 if(Status /= WRF_NO_ERR) then
703 write(msg,*) 'Warning in ',__FILE__,', line', __LINE__
704 call wrf_debug ( WARN , TRIM(msg))
705 write(msg,*) ' Bad time index for DateStr = ',DateStr
706 call wrf_debug ( WARN , TRIM(msg))
709 call GetDim(MemoryOrder,NDim,Status)
713 VCount(1:NDim) = Length(1:NDim)
714 VStart(NDim+1) = TimeIndex
717 ! Do not use SELECT statement here as sometimes WRF_REAL=WRF_DOUBLE
718 IF (FieldType == WRF_REAL) THEN
719 call ext_ncd_RealFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
720 ELSE IF (FieldType == WRF_DOUBLE) THEN
721 call ext_ncd_DoubleFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
722 ELSE IF (FieldType == WRF_INTEGER) THEN
723 call ext_ncd_IntFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
724 ELSE IF (FieldType == WRF_LOGICAL) THEN
725 call ext_ncd_LogicalFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
726 if(Status /= WRF_NO_ERR) return
728 !for wrf_complex, double_complex
729 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
730 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
731 call wrf_debug ( WARN , TRIM(msg))
736 end subroutine FieldIO
738 subroutine Transpose(IO,MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
739 ,XField,x1,x2,y1,y2,z1,z2 &
741 character*(*) ,intent(in) :: IO
742 character*(*) ,intent(in) :: MemoryOrder
743 integer ,intent(in) :: l1,l2,m1,m2,n1,n2
744 integer ,intent(in) :: di
745 integer ,intent(in) :: x1,x2,y1,y2,z1,z2
746 integer ,intent(in) :: i1,i2,j1,j2,k1,k2
747 integer ,intent(inout) :: Field(di,l1:l2,m1:m2,n1:n2)
748 !jm 010827 integer ,intent(inout) :: XField(di,x1:x2,y1:y2,z1:z2)
749 integer ,intent(inout) :: XField(di,(i2-i1+1)*(j2-j1+1)*(k2-k1+1))
750 character*3 :: MemOrd
752 integer ,parameter :: MaxUpperCase=IACHAR('Z')
753 integer :: i,j,k,ix,jx,kx
755 call LowerCase(MemoryOrder,MemOrd)
758 !#define XDEX(A,B,C) A-A ## 1+1+(A ## 2-A ## 1+1)*((B-B ## 1)+(C-C ## 1)*(B ## 2-B ## 1+1))
759 ! define(`XDEX',($1-``$1''1+1+(``$1''2-``$1''1+1)*(($2-``$2''1)+($3-``$3''1)*(``$2''2-``$2''1+1))))
763 #define DFIELD XField(1:di,XDEX(i,k,j))
764 #include "transpose.code"
765 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c','0')
767 #define DFIELD XField(1:di,XDEX(i,j,k))
768 #include "transpose.code"
771 #define DFIELD XField(1:di,XDEX(j,i,k))
772 #include "transpose.code"
775 #define DFIELD XField(1:di,XDEX(k,i,j))
776 #include "transpose.code"
779 #define DFIELD XField(1:di,XDEX(j,k,i))
780 #include "transpose.code"
783 #define DFIELD XField(1:di,XDEX(k,j,i))
784 #include "transpose.code"
787 #define DFIELD XField(1:di,XDEX(j,i,k))
788 #include "transpose.code"
791 end subroutine Transpose
793 subroutine reorder (MemoryOrder,MemO)
794 character*(*) ,intent(in) :: MemoryOrder
795 character*3 ,intent(out) :: MemO
796 character*3 :: MemOrd
797 integer :: N,i,i1,i2,i3
800 N = len_trim(MemoryOrder)
802 call lowercase(MemoryOrder,MemOrd)
803 ! never invert the boundary codes
804 select case ( MemOrd )
805 case ( 'xsz','xez','ysz','yez' )
813 if(ichar(MemOrd(i:i)) < ichar(MemOrd(i1:i1))) I1 = i
814 if(ichar(MemOrd(i:i)) > ichar(MemOrd(i3:i3))) I3 = i
821 MemO(1:1) = MemoryOrder(i1:i1)
822 MemO(2:2) = MemoryOrder(i2:i2)
823 if(N == 3) MemO(3:3) = MemoryOrder(i3:i3)
824 if(MemOrd(i1:i1) == 's' .or. MemOrd(i1:i1) == 'e') then
825 MemO(1:N-1) = MemO(2:N)
826 MemO(N:N ) = MemoryOrder(i1:i1)
829 end subroutine reorder
831 ! Returns .TRUE. iff it is OK to write time-independent domain metadata to the
832 ! file referenced by DataHandle. If DataHandle is invalid, .FALSE. is
834 LOGICAL FUNCTION ncd_ok_to_put_dom_ti( DataHandle )
836 include 'wrf_status_codes.h'
837 INTEGER, INTENT(IN) :: DataHandle
838 CHARACTER*80 :: fname
841 LOGICAL :: dryrun, first_output, retval
842 call ext_ncd_inquire_filename( DataHandle, fname, filestate, Status )
843 IF ( Status /= WRF_NO_ERR ) THEN
844 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
846 call wrf_debug ( WARN , TRIM(msg) )
849 dryrun = ( filestate .EQ. WRF_FILE_OPENED_NOT_COMMITTED )
850 first_output = ncd_is_first_operation( DataHandle )
851 retval = .NOT. dryrun .AND. first_output
853 ncd_ok_to_put_dom_ti = retval
855 END FUNCTION ncd_ok_to_put_dom_ti
857 ! Returns .TRUE. iff it is OK to read time-independent domain metadata from the
858 ! file referenced by DataHandle. If DataHandle is invalid, .FALSE. is
860 LOGICAL FUNCTION ncd_ok_to_get_dom_ti( DataHandle )
862 include 'wrf_status_codes.h'
863 INTEGER, INTENT(IN) :: DataHandle
864 CHARACTER*80 :: fname
867 LOGICAL :: dryrun, retval
868 call ext_ncd_inquire_filename( DataHandle, fname, filestate, Status )
869 IF ( Status /= WRF_NO_ERR ) THEN
870 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
872 call wrf_debug ( WARN , TRIM(msg) )
875 dryrun = ( filestate .EQ. WRF_FILE_OPENED_NOT_COMMITTED )
876 retval = .NOT. dryrun
878 ncd_ok_to_get_dom_ti = retval
880 END FUNCTION ncd_ok_to_get_dom_ti
882 ! Returns .TRUE. iff nothing has been read from or written to the file
883 ! referenced by DataHandle. If DataHandle is invalid, .FALSE. is returned.
884 LOGICAL FUNCTION ncd_is_first_operation( DataHandle )
886 INCLUDE 'wrf_status_codes.h'
887 INTEGER, INTENT(IN) :: DataHandle
888 TYPE(wrf_data_handle) ,POINTER :: DH
891 CALL GetDH( DataHandle, DH, Status )
892 IF ( Status /= WRF_NO_ERR ) THEN
893 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
895 call wrf_debug ( WARN , TRIM(msg) )
898 retval = DH%first_operation
900 ncd_is_first_operation = retval
902 END FUNCTION ncd_is_first_operation
904 subroutine upgrade_filename(FileName)
907 character*(*), intent(inout) :: FileName
910 do i = 1, len(trim(FileName))
911 if(FileName(i:i) == '-') then
913 else if(FileName(i:i) == ':') then
918 end subroutine upgrade_filename
920 end module ext_ncd_support_routines
922 subroutine TransposeToR4(IO,MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
923 ,XField,x1,x2,y1,y2,z1,z2 &
926 use ext_ncd_support_routines
928 character*(*) ,intent(in) :: IO
929 character*(*) ,intent(in) :: MemoryOrder
930 integer ,intent(in) :: l1,l2,m1,m2,n1,n2
931 integer ,intent(in) :: di
932 integer ,intent(in) :: x1,x2,y1,y2,z1,z2
933 integer ,intent(in) :: i1,i2,j1,j2,k1,k2
934 real*8 ,intent(inout) :: Field(di,l1:l2,m1:m2,n1:n2)
935 real*4 ,intent(inout) :: XField(di,(i2-i1+1)*(j2-j1+1)*(k2-k1+1))
936 character*3 :: MemOrd
938 integer ,parameter :: MaxUpperCase=IACHAR('Z')
939 integer :: i,j,k,ix,jx,kx
941 call LowerCase(MemoryOrder,MemOrd)
944 !#define XDEX(A,B,C) A-A ## 1+1+(A ## 2-A ## 1+1)*((B-B ## 1)+(C-C ## 1)*(B ## 2-B ## 1+1))
945 ! define(`XDEX',($1-``$1''1+1+(``$1''2-``$1''1+1)*(($2-``$2''1)+($3-``$3''1)*(``$2''2-``$2''1+1))))
949 #define DFIELD XField(1:di,XDEX(i,k,j))
950 #include "transpose.code"
951 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c','0')
953 #define DFIELD XField(1:di,XDEX(i,j,k))
954 #include "transpose.code"
957 #define DFIELD XField(1:di,XDEX(j,i,k))
958 #include "transpose.code"
961 #define DFIELD XField(1:di,XDEX(k,i,j))
962 #include "transpose.code"
965 #define DFIELD XField(1:di,XDEX(j,k,i))
966 #include "transpose.code"
969 #define DFIELD XField(1:di,XDEX(k,j,i))
970 #include "transpose.code"
973 #define DFIELD XField(1:di,XDEX(j,i,k))
974 #include "transpose.code"
977 end subroutine TransposeToR4
979 subroutine ext_ncd_open_for_read(DatasetName, Comm1, Comm2, SysDepInfo, DataHandle, Status)
981 use ext_ncd_support_routines
983 include 'wrf_status_codes.h'
985 character *(*), INTENT(IN) :: DatasetName
986 integer , INTENT(IN) :: Comm1, Comm2
987 character *(*), INTENT(IN) :: SysDepInfo
988 integer , INTENT(OUT) :: DataHandle
989 integer , INTENT(OUT) :: Status
990 DataHandle = 0 ! dummy setting to quiet warning message
991 CALL ext_ncd_open_for_read_begin( DatasetName, Comm1, Comm2, SysDepInfo, DataHandle, Status )
992 IF ( Status .EQ. WRF_NO_ERR ) THEN
993 CALL ext_ncd_open_for_read_commit( DataHandle, Status )
996 end subroutine ext_ncd_open_for_read
998 !ends training phase; switches internal flag to enable input
999 !must be paired with call to ext_ncd_open_for_read_begin
1000 subroutine ext_ncd_open_for_read_commit(DataHandle, Status)
1002 use ext_ncd_support_routines
1004 include 'wrf_status_codes.h'
1005 include 'netcdf.inc'
1006 integer, intent(in) :: DataHandle
1007 integer, intent(out) :: Status
1008 type(wrf_data_handle) ,pointer :: DH
1010 if(WrfIOnotInitialized) then
1011 Status = WRF_IO_NOT_INITIALIZED
1012 write(msg,*) 'ext_ncd_ioinit was not called ',__FILE__,', line', __LINE__
1013 call wrf_debug ( FATAL , msg)
1016 call GetDH(DataHandle,DH,Status)
1017 if(Status /= WRF_NO_ERR) then
1018 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1019 call wrf_debug ( WARN , TRIM(msg))
1022 DH%FileStatus = WRF_FILE_OPENED_FOR_READ
1023 DH%first_operation = .TRUE.
1026 end subroutine ext_ncd_open_for_read_commit
1028 subroutine ext_ncd_open_for_read_begin( FileName, Comm, IOComm, SysDepInfo, DataHandle, Status)
1030 use ext_ncd_support_routines
1032 include 'wrf_status_codes.h'
1033 include 'netcdf.inc'
1034 character*(*) ,intent(INOUT) :: FileName
1035 integer ,intent(IN) :: Comm
1036 integer ,intent(IN) :: IOComm
1037 character*(*) ,intent(in) :: SysDepInfo
1038 integer ,intent(out) :: DataHandle
1039 integer ,intent(out) :: Status
1040 type(wrf_data_handle) ,pointer :: DH
1043 integer ,allocatable :: Buffer(:)
1045 integer :: StoredDim
1047 integer :: DimIDs(2)
1048 integer :: VStart(2)
1050 integer :: TotalNumVars
1053 character (NF_MAX_NAME) :: Name
1055 #ifdef USE_NETCDF4_FEATURES
1056 integer :: open_mode
1059 !call upgrade_filename(FileName)
1061 if(WrfIOnotInitialized) then
1062 Status = WRF_IO_NOT_INITIALIZED
1063 write(msg,*) 'ext_ncd_ioinit was not called ',__FILE__,', line', __LINE__
1064 call wrf_debug ( FATAL , msg)
1067 call allocHandle(DataHandle,DH,Comm,Status)
1068 if(Status /= WRF_NO_ERR) then
1069 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
1070 call wrf_debug ( WARN , TRIM(msg))
1074 stat = NF_OPEN(FileName, NF_NOWRITE, DH%NCID)
1075 call netcdf_err(stat,Status)
1076 if(Status /= WRF_NO_ERR) then
1077 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1078 call wrf_debug ( WARN , TRIM(msg))
1081 stat = NF_INQ_VARID(DH%NCID,DH%TimesName,VarID)
1082 call netcdf_err(stat,Status)
1083 if(Status /= WRF_NO_ERR) then
1084 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1085 call wrf_debug ( WARN , TRIM(msg))
1088 stat = NF_INQ_VAR(DH%NCID,VarID,DH%TimesName, XType, StoredDim, DimIDs, NAtts)
1089 call netcdf_err(stat,Status)
1090 if(Status /= WRF_NO_ERR) then
1091 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1092 call wrf_debug ( WARN , TRIM(msg))
1095 if(XType/=NF_CHAR) then
1096 Status = WRF_WARN_TYPE_MISMATCH
1097 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
1098 call wrf_debug ( WARN , TRIM(msg))
1101 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(1),VLen(1))
1102 call netcdf_err(stat,Status)
1103 if(Status /= WRF_NO_ERR) then
1104 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1105 call wrf_debug ( WARN , TRIM(msg))
1108 if(VLen(1) /= DateStrLen) then
1109 Status = WRF_WARN_DATESTR_BAD_LENGTH
1110 write(msg,*) 'Warning DATESTR BAD LENGTH in ',__FILE__,', line', __LINE__
1111 call wrf_debug ( WARN , TRIM(msg))
1114 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(2),VLen(2))
1115 call netcdf_err(stat,Status)
1116 if(Status /= WRF_NO_ERR) then
1117 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1118 call wrf_debug ( WARN , TRIM(msg))
1121 if(VLen(2) > MaxTimes) then
1122 Status = WRF_ERR_FATAL_TOO_MANY_TIMES
1123 write(msg,*) 'Fatal TOO MANY TIME VALUES in ',__FILE__,', line', __LINE__
1124 call wrf_debug ( FATAL , TRIM(msg))
1129 stat = NF_GET_VARA_TEXT(DH%NCID,VarID,VStart,VLen,DH%Times)
1130 call netcdf_err(stat,Status)
1131 if(Status /= WRF_NO_ERR) then
1132 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1133 call wrf_debug ( WARN , TRIM(msg))
1136 stat = NF_INQ_NVARS(DH%NCID,TotalNumVars)
1137 call netcdf_err(stat,Status)
1138 if(Status /= WRF_NO_ERR) then
1139 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1140 call wrf_debug ( WARN , TRIM(msg))
1145 stat = NF_INQ_VARNAME(DH%NCID,i,Name)
1146 call netcdf_err(stat,Status)
1147 if(Status /= WRF_NO_ERR) then
1148 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1149 call wrf_debug ( WARN , TRIM(msg))
1151 elseif(Name(1:5) /= 'md___' .and. Name /= DH%TimesName) then
1153 DH%VarNames(NumVars) = Name
1154 DH%VarIDs(NumVars) = i
1157 DH%NumVars = NumVars
1158 DH%NumberTimes = VLen(2)
1159 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1160 DH%FileName = trim(FileName)
1161 DH%CurrentVariable = 0
1163 DH%TimesVarID = VarID
1166 end subroutine ext_ncd_open_for_read_begin
1168 subroutine ext_ncd_open_for_update( FileName, Comm, IOComm, SysDepInfo, DataHandle, Status)
1170 use ext_ncd_support_routines
1172 include 'wrf_status_codes.h'
1173 include 'netcdf.inc'
1174 character*(*) ,intent(INOUT) :: FileName
1175 integer ,intent(IN) :: Comm
1176 integer ,intent(IN) :: IOComm
1177 character*(*) ,intent(in) :: SysDepInfo
1178 integer ,intent(out) :: DataHandle
1179 integer ,intent(out) :: Status
1180 type(wrf_data_handle) ,pointer :: DH
1183 integer ,allocatable :: Buffer(:)
1185 integer :: StoredDim
1187 integer :: DimIDs(2)
1188 integer :: VStart(2)
1190 integer :: TotalNumVars
1193 character (NF_MAX_NAME) :: Name
1195 #ifdef USE_NETCDF4_FEATURES
1196 integer :: open_mode
1199 !call upgrade_filename(FileName)
1201 if(WrfIOnotInitialized) then
1202 Status = WRF_IO_NOT_INITIALIZED
1203 write(msg,*) 'ext_ncd_ioinit was not called ',__FILE__,', line', __LINE__
1204 call wrf_debug ( FATAL , msg)
1207 call allocHandle(DataHandle,DH,Comm,Status)
1208 if(Status /= WRF_NO_ERR) then
1209 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
1210 call wrf_debug ( WARN , TRIM(msg))
1213 stat = NF_OPEN(FileName, NF_WRITE, DH%NCID)
1214 call netcdf_err(stat,Status)
1215 if(Status /= WRF_NO_ERR) then
1216 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1217 call wrf_debug ( WARN , TRIM(msg))
1220 stat = NF_INQ_VARID(DH%NCID,DH%TimesName,VarID)
1221 call netcdf_err(stat,Status)
1222 if(Status /= WRF_NO_ERR) then
1223 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1224 call wrf_debug ( WARN , TRIM(msg))
1227 stat = NF_INQ_VAR(DH%NCID,VarID,DH%TimesName, XType, StoredDim, DimIDs, NAtts)
1228 call netcdf_err(stat,Status)
1229 if(Status /= WRF_NO_ERR) then
1230 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1231 call wrf_debug ( WARN , TRIM(msg))
1234 if(XType/=NF_CHAR) then
1235 Status = WRF_WARN_TYPE_MISMATCH
1236 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
1237 call wrf_debug ( WARN , TRIM(msg))
1240 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(1),VLen(1))
1241 call netcdf_err(stat,Status)
1242 if(Status /= WRF_NO_ERR) then
1243 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1244 call wrf_debug ( WARN , TRIM(msg))
1247 if(VLen(1) /= DateStrLen) then
1248 Status = WRF_WARN_DATESTR_BAD_LENGTH
1249 write(msg,*) 'Warning DATESTR BAD LENGTH in ',__FILE__,', line', __LINE__
1250 call wrf_debug ( WARN , TRIM(msg))
1253 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(2),VLen(2))
1254 call netcdf_err(stat,Status)
1255 if(Status /= WRF_NO_ERR) then
1256 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1257 call wrf_debug ( WARN , TRIM(msg))
1260 if(VLen(2) > MaxTimes) then
1261 Status = WRF_ERR_FATAL_TOO_MANY_TIMES
1262 write(msg,*) 'Fatal TOO MANY TIME VALUES in ',__FILE__,', line', __LINE__
1263 call wrf_debug ( FATAL , TRIM(msg))
1268 stat = NF_GET_VARA_TEXT(DH%NCID,VarID,VStart,VLen,DH%Times)
1269 call netcdf_err(stat,Status)
1270 if(Status /= WRF_NO_ERR) then
1271 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1272 call wrf_debug ( WARN , TRIM(msg))
1275 stat = NF_INQ_NVARS(DH%NCID,TotalNumVars)
1276 call netcdf_err(stat,Status)
1277 if(Status /= WRF_NO_ERR) then
1278 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1279 call wrf_debug ( WARN , TRIM(msg))
1284 stat = NF_INQ_VARNAME(DH%NCID,i,Name)
1285 call netcdf_err(stat,Status)
1286 if(Status /= WRF_NO_ERR) then
1287 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1288 call wrf_debug ( WARN , TRIM(msg))
1290 elseif(Name(1:5) /= 'md___' .and. Name /= DH%TimesName) then
1292 DH%VarNames(NumVars) = Name
1293 DH%VarIDs(NumVars) = i
1296 DH%NumVars = NumVars
1297 DH%NumberTimes = VLen(2)
1298 DH%FileStatus = WRF_FILE_OPENED_FOR_UPDATE
1299 DH%FileName = trim(FileName)
1300 DH%CurrentVariable = 0
1302 DH%TimesVarID = VarID
1305 end subroutine ext_ncd_open_for_update
1308 SUBROUTINE ext_ncd_open_for_write_begin(FileName,Comm,IOComm,SysDepInfo,DataHandle,Status)
1310 use ext_ncd_support_routines
1312 include 'wrf_status_codes.h'
1313 include 'netcdf.inc'
1314 character*(*) ,intent(inout) :: FileName
1315 integer ,intent(in) :: Comm
1316 integer ,intent(in) :: IOComm
1317 character*(*) ,intent(in) :: SysDepInfo
1318 integer ,intent(out) :: DataHandle
1319 integer ,intent(out) :: Status
1320 type(wrf_data_handle),pointer :: DH
1323 character (7) :: Buffer
1324 integer :: VDimIDs(2)
1326 #ifdef USE_NETCDF4_FEATURES
1327 integer :: create_mode
1328 integer, parameter :: cache_size = 32, &
1330 cache_preemption = 100
1333 !call upgrade_filename(FileName)
1335 if(WrfIOnotInitialized) then
1336 Status = WRF_IO_NOT_INITIALIZED
1337 write(msg,*) 'ext_ncd_open_for_write_begin: ext_ncd_ioinit was not called ',__FILE__,', line', __LINE__
1338 call wrf_debug ( FATAL , msg)
1341 call allocHandle(DataHandle,DH,Comm,Status)
1342 if(Status /= WRF_NO_ERR) then
1343 write(msg,*) 'Fatal ALLOCATION ERROR in ext_ncd_open_for_write_begin ',__FILE__,', line', __LINE__
1344 call wrf_debug ( FATAL , TRIM(msg))
1349 #ifdef USE_NETCDF4_FEATURES
1350 ! create_mode = IOR(nf_netcdf4, nf_classic_model)
1351 if ( DH%use_netcdf_classic ) then
1352 write(msg,*) 'output will be in classic NetCDF format'
1353 call wrf_debug ( WARN , TRIM(msg))
1354 #ifdef WRFIO_NCD_NO_LARGE_FILE_SUPPORT
1355 stat = NF_CREATE(FileName, NF_CLOBBER, DH%NCID)
1357 stat = NF_CREATE(FileName, IOR(NF_CLOBBER,NF_64BIT_OFFSET), DH%NCID)
1360 create_mode = nf_netcdf4
1361 stat = NF_CREATE(FileName, create_mode, DH%NCID)
1362 stat = NF_SET_CHUNK_CACHE(cache_size, cache_nelem, cache_preemption)
1365 #ifdef WRFIO_NCD_NO_LARGE_FILE_SUPPORT
1366 stat = NF_CREATE(FileName, NF_CLOBBER, DH%NCID)
1368 stat = NF_CREATE(FileName, IOR(NF_CLOBBER,NF_64BIT_OFFSET), DH%NCID)
1371 call netcdf_err(stat,Status)
1372 if(Status /= WRF_NO_ERR) then
1373 write(msg,*) 'NetCDF error in ext_ncd_open_for_write_begin ',__FILE__,', line', __LINE__
1374 call wrf_debug ( WARN , TRIM(msg))
1377 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1378 DH%FileName = trim(FileName)
1379 stat = NF_DEF_DIM(DH%NCID,DH%DimUnlimName,NF_UNLIMITED,DH%DimUnlimID)
1380 call netcdf_err(stat,Status)
1381 if(Status /= WRF_NO_ERR) then
1382 write(msg,*) 'NetCDF error in ext_ncd_open_for_write_begin ',__FILE__,', line', __LINE__
1383 call wrf_debug ( WARN , TRIM(msg))
1386 DH%VarNames (1:MaxVars) = NO_NAME
1387 DH%MDVarNames(1:MaxVars) = NO_NAME
1389 write(Buffer,FMT="('DIM',i4.4)") i
1390 DH%DimNames (i) = Buffer
1391 DH%DimLengths(i) = NO_DIM
1393 DH%DimNames(1) = 'DateStrLen'
1394 stat = NF_DEF_DIM(DH%NCID,DH%DimNames(1),DateStrLen,DH%DimIDs(1))
1395 call netcdf_err(stat,Status)
1396 if(Status /= WRF_NO_ERR) then
1397 write(msg,*) 'NetCDF error in ext_ncd_open_for_write_begin ',__FILE__,', line', __LINE__
1398 call wrf_debug ( WARN , TRIM(msg))
1401 VDimIDs(1) = DH%DimIDs(1)
1402 VDimIDs(2) = DH%DimUnlimID
1403 stat = NF_DEF_VAR(DH%NCID,DH%TimesName,NF_CHAR,2,VDimIDs,DH%TimesVarID)
1404 call netcdf_err(stat,Status)
1405 if(Status /= WRF_NO_ERR) then
1406 write(msg,*) 'NetCDF error in ext_ncd_open_for_write_begin ',__FILE__,', line', __LINE__
1407 call wrf_debug ( WARN , TRIM(msg))
1410 DH%DimLengths(1) = DateStrLen
1412 if (index(SysDepInfo,'REAL_OUTPUT_SIZE=4') /= 0) then
1413 DH%R4OnOutput = .true.
1415 !toggle on nofill mode
1416 if (index(SysDepInfo,'NOFILL=.TRUE.') /= 0) then
1421 end subroutine ext_ncd_open_for_write_begin
1424 !opens a file for writing or coupler datastream for sending messages.
1425 !no training phase for this version of the open stmt.
1426 subroutine ext_ncd_open_for_write (DatasetName, Comm1, Comm2, &
1427 SysDepInfo, DataHandle, Status)
1429 use ext_ncd_support_routines
1431 include 'wrf_status_codes.h'
1432 include 'netcdf.inc'
1433 character *(*), intent(in) ::DatasetName
1434 integer , intent(in) ::Comm1, Comm2
1435 character *(*), intent(in) ::SysDepInfo
1436 integer , intent(out) :: DataHandle
1437 integer , intent(out) :: Status
1438 Status=WRF_WARN_NOOP
1439 DataHandle = 0 ! dummy setting to quiet warning message
1441 end subroutine ext_ncd_open_for_write
1443 SUBROUTINE ext_ncd_open_for_write_commit(DataHandle, Status)
1445 use ext_ncd_support_routines
1447 include 'wrf_status_codes.h'
1448 include 'netcdf.inc'
1449 integer ,intent(in) :: DataHandle
1450 integer ,intent(out) :: Status
1451 type(wrf_data_handle),pointer :: DH
1454 integer :: oldmode ! for nf_set_fill, not used
1456 if(WrfIOnotInitialized) then
1457 Status = WRF_IO_NOT_INITIALIZED
1458 write(msg,*) 'ext_ncd_open_for_write_commit: ext_ncd_ioinit was not called ',__FILE__,', line', __LINE__
1459 call wrf_debug ( FATAL , msg)
1462 call GetDH(DataHandle,DH,Status)
1463 if(Status /= WRF_NO_ERR) then
1464 write(msg,*) 'Warning Status = ',Status,' in ext_ncd_open_for_write_commit ',__FILE__,', line', __LINE__
1465 call wrf_debug ( WARN , TRIM(msg))
1468 if ( DH%nofill ) then
1469 Status = NF_SET_FILL(DH%NCID,NF_NOFILL, oldmode )
1470 if(Status /= WRF_NO_ERR) then
1471 write(msg,*) 'Warning Status = ',Status,' from NF_SET_FILL ',__FILE__,', line', __LINE__
1472 call wrf_debug ( WARN , TRIM(msg))
1475 write(msg,*) 'Information: NOFILL being set for writing to ',TRIM(DH%FileName)
1476 call wrf_debug ( WARN , TRIM(msg))
1478 stat = NF_ENDDEF(DH%NCID)
1479 call netcdf_err(stat,Status)
1480 if(Status /= WRF_NO_ERR) then
1481 write(msg,*) 'NetCDF error in ext_ncd_open_for_write_commit ',__FILE__,', line', __LINE__
1482 call wrf_debug ( WARN , TRIM(msg))
1485 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1486 DH%first_operation = .TRUE.
1488 end subroutine ext_ncd_open_for_write_commit
1490 subroutine ext_ncd_ioclose(DataHandle, Status)
1492 use ext_ncd_support_routines
1494 include 'wrf_status_codes.h'
1495 include 'netcdf.inc'
1496 integer ,intent(in) :: DataHandle
1497 integer ,intent(out) :: Status
1498 type(wrf_data_handle),pointer :: DH
1501 call GetDH(DataHandle,DH,Status)
1502 if(Status /= WRF_NO_ERR) then
1503 write(msg,*) 'Warning Status = ',Status,' in ext_ncd_ioclose ',__FILE__,', line', __LINE__
1504 call wrf_debug ( WARN , TRIM(msg))
1507 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1508 Status = WRF_WARN_FILE_NOT_OPENED
1509 write(msg,*) 'Warning FILE NOT OPENED in ext_ncd_ioclose ',__FILE__,', line', __LINE__
1510 call wrf_debug ( WARN , TRIM(msg))
1511 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1512 Status = WRF_WARN_DRYRUN_CLOSE
1513 write(msg,*) 'Warning TRY TO CLOSE DRYRUN in ext_ncd_ioclose ',__FILE__,', line', __LINE__
1514 call wrf_debug ( WARN , TRIM(msg))
1515 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1517 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1519 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
1522 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1523 write(msg,*) 'Fatal error BAD FILE STATUS in ext_ncd_ioclose ',__FILE__,', line', __LINE__
1524 call wrf_debug ( FATAL , TRIM(msg))
1528 stat = NF_CLOSE(DH%NCID)
1529 call netcdf_err(stat,Status)
1530 if(Status /= WRF_NO_ERR) then
1531 write(msg,*) 'NetCDF error in ext_ncd_ioclose ',__FILE__,', line', __LINE__
1532 call wrf_debug ( WARN , TRIM(msg))
1535 CALL deallocHandle( DataHandle, Status )
1538 end subroutine ext_ncd_ioclose
1540 subroutine ext_ncd_iosync( DataHandle, Status)
1542 use ext_ncd_support_routines
1544 include 'wrf_status_codes.h'
1545 include 'netcdf.inc'
1546 integer ,intent(in) :: DataHandle
1547 integer ,intent(out) :: Status
1548 type(wrf_data_handle),pointer :: DH
1551 call GetDH(DataHandle,DH,Status)
1552 if(Status /= WRF_NO_ERR) then
1553 write(msg,*) 'Warning Status = ',Status,' in ext_ncd_iosync ',__FILE__,', line', __LINE__
1554 call wrf_debug ( WARN , TRIM(msg))
1557 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1558 Status = WRF_WARN_FILE_NOT_OPENED
1559 write(msg,*) 'Warning FILE NOT OPENED in ext_ncd_iosync ',__FILE__,', line', __LINE__
1560 call wrf_debug ( WARN , TRIM(msg))
1561 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1562 Status = WRF_WARN_FILE_NOT_COMMITTED
1563 write(msg,*) 'Warning FILE NOT COMMITTED in ext_ncd_iosync ',__FILE__,', line', __LINE__
1564 call wrf_debug ( WARN , TRIM(msg))
1565 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1567 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1570 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1571 write(msg,*) 'Fatal error BAD FILE STATUS in ext_ncd_iosync ',__FILE__,', line', __LINE__
1572 call wrf_debug ( FATAL , TRIM(msg))
1575 stat = NF_SYNC(DH%NCID)
1576 call netcdf_err(stat,Status)
1577 if(Status /= WRF_NO_ERR) then
1578 write(msg,*) 'NetCDF error in ext_ncd_iosync ',__FILE__,', line', __LINE__
1579 call wrf_debug ( WARN , TRIM(msg))
1583 end subroutine ext_ncd_iosync
1587 subroutine ext_ncd_redef( DataHandle, Status)
1589 use ext_ncd_support_routines
1591 include 'wrf_status_codes.h'
1592 include 'netcdf.inc'
1593 integer ,intent(in) :: DataHandle
1594 integer ,intent(out) :: Status
1595 type(wrf_data_handle),pointer :: DH
1598 call GetDH(DataHandle,DH,Status)
1599 if(Status /= WRF_NO_ERR) then
1600 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1601 call wrf_debug ( WARN , TRIM(msg))
1604 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1605 Status = WRF_WARN_FILE_NOT_OPENED
1606 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1607 call wrf_debug ( WARN , TRIM(msg))
1608 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1609 Status = WRF_WARN_FILE_NOT_COMMITTED
1610 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1611 call wrf_debug ( WARN , TRIM(msg))
1612 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1614 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
1616 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1617 Status = WRF_WARN_FILE_OPEN_FOR_READ
1618 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1619 call wrf_debug ( WARN , TRIM(msg))
1621 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1622 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1623 call wrf_debug ( FATAL , TRIM(msg))
1626 stat = NF_REDEF(DH%NCID)
1627 call netcdf_err(stat,Status)
1628 if(Status /= WRF_NO_ERR) then
1629 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1630 call wrf_debug ( WARN , TRIM(msg))
1633 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1635 end subroutine ext_ncd_redef
1637 subroutine ext_ncd_enddef( DataHandle, Status)
1639 use ext_ncd_support_routines
1641 include 'wrf_status_codes.h'
1642 include 'netcdf.inc'
1643 integer ,intent(in) :: DataHandle
1644 integer ,intent(out) :: Status
1645 type(wrf_data_handle),pointer :: DH
1648 call GetDH(DataHandle,DH,Status)
1649 if(Status /= WRF_NO_ERR) then
1650 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1651 call wrf_debug ( WARN , TRIM(msg))
1654 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1655 Status = WRF_WARN_FILE_NOT_OPENED
1656 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1657 call wrf_debug ( WARN , TRIM(msg))
1658 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1659 Status = WRF_WARN_FILE_NOT_COMMITTED
1660 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1661 call wrf_debug ( WARN , TRIM(msg))
1662 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1664 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1665 Status = WRF_WARN_FILE_OPEN_FOR_READ
1666 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1667 call wrf_debug ( WARN , TRIM(msg))
1669 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1670 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1671 call wrf_debug ( FATAL , TRIM(msg))
1674 stat = NF_ENDDEF(DH%NCID)
1675 call netcdf_err(stat,Status)
1676 if(Status /= WRF_NO_ERR) then
1677 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1678 call wrf_debug ( WARN , TRIM(msg))
1681 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1683 end subroutine ext_ncd_enddef
1685 subroutine ext_ncd_ioinit(SysDepInfo, Status)
1688 include 'wrf_status_codes.h'
1689 CHARACTER*(*), INTENT(IN) :: SysDepInfo
1690 INTEGER ,INTENT(INOUT) :: Status
1692 WrfIOnotInitialized = .false.
1693 WrfDataHandles(1:WrfDataHandleMax)%Free = .true.
1694 WrfDataHandles(1:WrfDataHandleMax)%TimesName = 'Times'
1695 WrfDataHandles(1:WrfDataHandleMax)%DimUnlimName = 'Time'
1696 WrfDataHandles(1:WrfDataHandleMax)%FileStatus = WRF_FILE_NOT_OPENED
1697 if(trim(SysDepInfo) == "use_netcdf_classic" ) then
1698 WrfDataHandles(1:WrfDataHandleMax)%use_netcdf_classic = .true.
1700 WrfDataHandles(1:WrfDataHandleMax)%use_netcdf_classic = .false.
1704 end subroutine ext_ncd_ioinit
1707 subroutine ext_ncd_inquiry (Inquiry, Result, Status)
1710 include 'wrf_status_codes.h'
1711 character *(*), INTENT(IN) :: Inquiry
1712 character *(*), INTENT(OUT) :: Result
1713 integer ,INTENT(INOUT) :: Status
1714 SELECT CASE (Inquiry)
1715 CASE ("RANDOM_WRITE","RANDOM_READ","SEQUENTIAL_WRITE","SEQUENTIAL_READ")
1717 CASE ("OPEN_READ","OPEN_COMMIT_WRITE")
1719 CASE ("OPEN_WRITE","OPEN_COMMIT_READ","PARALLEL_IO")
1721 CASE ("SELF_DESCRIBING","SUPPORT_METADATA","SUPPORT_3D_FIELDS")
1726 Result = 'No Result for that inquiry!'
1730 end subroutine ext_ncd_inquiry
1735 subroutine ext_ncd_ioexit(Status)
1737 use ext_ncd_support_routines
1739 include 'wrf_status_codes.h'
1740 include 'netcdf.inc'
1741 integer , INTENT(INOUT) ::Status
1743 type(wrf_data_handle),pointer :: DH
1746 if(WrfIOnotInitialized) then
1747 Status = WRF_IO_NOT_INITIALIZED
1748 write(msg,*) 'ext_ncd_ioinit was not called ',__FILE__,', line', __LINE__
1749 call wrf_debug ( FATAL , msg)
1752 do i=1,WrfDataHandleMax
1753 CALL deallocHandle( i , stat )
1756 end subroutine ext_ncd_ioexit
1758 subroutine ext_ncd_get_dom_ti_real(DataHandle,Element,Data,Count,OutCount,Status)
1759 #define ROUTINE_TYPE 'REAL'
1760 #define TYPE_DATA real,intent(out) :: Data(*)
1761 #define TYPE_COUNT integer,intent(in) :: Count
1762 #define TYPE_OUTCOUNT integer,intent(out) :: OutCOunt
1763 #define TYPE_BUFFER real,allocatable :: Buffer(:)
1764 #define NF_TYPE NF_FLOAT
1765 #define NF_ROUTINE NF_GET_ATT_REAL
1766 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1767 #include "ext_ncd_get_dom_ti.code"
1768 end subroutine ext_ncd_get_dom_ti_real
1770 subroutine ext_ncd_get_dom_ti_integer(DataHandle,Element,Data,Count,OutCount,Status)
1777 #define ROUTINE_TYPE 'INTEGER'
1778 #define TYPE_DATA integer,intent(out) :: Data(*)
1779 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
1780 #define NF_TYPE NF_INT
1781 #define NF_ROUTINE NF_GET_ATT_INT
1782 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1783 #include "ext_ncd_get_dom_ti.code"
1784 end subroutine ext_ncd_get_dom_ti_integer
1786 subroutine ext_ncd_get_dom_ti_double(DataHandle,Element,Data,Count,OutCount,Status)
1793 #define ROUTINE_TYPE 'DOUBLE'
1794 #define TYPE_DATA real*8,intent(out) :: Data(*)
1795 #define TYPE_BUFFER real*8,allocatable :: Buffer(:)
1796 #define NF_TYPE NF_DOUBLE
1797 #define NF_ROUTINE NF_GET_ATT_DOUBLE
1798 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1799 #include "ext_ncd_get_dom_ti.code"
1800 end subroutine ext_ncd_get_dom_ti_double
1802 subroutine ext_ncd_get_dom_ti_logical(DataHandle,Element,Data,Count,OutCount,Status)
1809 #define ROUTINE_TYPE 'LOGICAL'
1810 #define TYPE_DATA logical,intent(out) :: Data(*)
1811 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
1812 #define NF_TYPE NF_INT
1813 #define NF_ROUTINE NF_GET_ATT_INT
1814 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))==1
1815 #include "ext_ncd_get_dom_ti.code"
1816 end subroutine ext_ncd_get_dom_ti_logical
1818 subroutine ext_ncd_get_dom_ti_char(DataHandle,Element,Data,Status)
1822 #undef TYPE_OUTCOUNT
1825 #define ROUTINE_TYPE 'CHAR'
1826 #define TYPE_DATA character*(*),intent(out) :: Data
1828 #define TYPE_OUTCOUNT
1830 #define NF_TYPE NF_CHAR
1832 #include "ext_ncd_get_dom_ti.code"
1834 end subroutine ext_ncd_get_dom_ti_char
1836 subroutine ext_ncd_put_dom_ti_real(DataHandle,Element,Data,Count,Status)
1843 #define ROUTINE_TYPE 'REAL'
1844 #define TYPE_DATA real ,intent(in) :: Data(*)
1845 #define TYPE_COUNT integer,intent(in) :: Count
1846 #define NF_ROUTINE NF_PUT_ATT_REAL
1847 #define ARGS NF_FLOAT,Count,Data
1848 #include "ext_ncd_put_dom_ti.code"
1849 end subroutine ext_ncd_put_dom_ti_real
1851 subroutine ext_ncd_put_dom_ti_integer(DataHandle,Element,Data,Count,Status)
1858 #define ROUTINE_TYPE 'INTEGER'
1859 #define TYPE_DATA integer,intent(in) :: Data(*)
1860 #define TYPE_COUNT integer,intent(in) :: Count
1861 #define NF_ROUTINE NF_PUT_ATT_INT
1862 #define ARGS NF_INT,Count,Data
1863 #include "ext_ncd_put_dom_ti.code"
1864 end subroutine ext_ncd_put_dom_ti_integer
1866 subroutine ext_ncd_put_dom_ti_double(DataHandle,Element,Data,Count,Status)
1873 #define ROUTINE_TYPE 'DOUBLE'
1874 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1875 #define TYPE_COUNT integer,intent(in) :: Count
1876 #define NF_ROUTINE NF_PUT_ATT_DOUBLE
1877 #define ARGS NF_DOUBLE,Count,Data
1878 #include "ext_ncd_put_dom_ti.code"
1879 end subroutine ext_ncd_put_dom_ti_double
1881 subroutine ext_ncd_put_dom_ti_logical(DataHandle,Element,Data,Count,Status)
1887 #define ROUTINE_TYPE 'LOGICAL'
1888 #define TYPE_DATA logical,intent(in) :: Data(*)
1889 #define TYPE_COUNT integer,intent(in) :: Count
1890 #define NF_ROUTINE NF_PUT_ATT_INT
1891 #define ARGS NF_INT,Count,Buffer
1893 #include "ext_ncd_put_dom_ti.code"
1894 end subroutine ext_ncd_put_dom_ti_logical
1896 subroutine ext_ncd_put_dom_ti_char(DataHandle,Element,Data,Status)
1903 #define ROUTINE_TYPE 'CHAR'
1904 #define TYPE_DATA character*(*),intent(in) :: Data
1905 #define TYPE_COUNT integer,parameter :: Count=1
1906 #define NF_ROUTINE NF_PUT_ATT_TEXT
1907 #define ARGS len_trim(Data),Data
1908 #include "ext_ncd_put_dom_ti.code"
1909 end subroutine ext_ncd_put_dom_ti_char
1911 subroutine ext_ncd_put_var_ti_real(DataHandle,Element,Var,Data,Count,Status)
1918 #define ROUTINE_TYPE 'REAL'
1919 #define TYPE_DATA real ,intent(in) :: Data(*)
1920 #define TYPE_COUNT integer ,intent(in) :: Count
1921 #define NF_ROUTINE NF_PUT_ATT_REAL
1922 #define ARGS NF_FLOAT,Count,Data
1923 #include "ext_ncd_put_var_ti.code"
1924 end subroutine ext_ncd_put_var_ti_real
1926 subroutine ext_ncd_put_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,Status)
1935 #define ROUTINE_TYPE 'REAL'
1936 #define TYPE_DATA real ,intent(in) :: Data(*)
1937 #define TYPE_COUNT integer ,intent(in) :: Count
1938 #define NF_ROUTINE NF_PUT_VARA_REAL
1939 #define NF_TYPE NF_FLOAT
1940 #define LENGTH Count
1942 #include "ext_ncd_put_var_td.code"
1943 end subroutine ext_ncd_put_var_td_real
1945 subroutine ext_ncd_put_var_ti_double(DataHandle,Element,Var,Data,Count,Status)
1952 #define ROUTINE_TYPE 'DOUBLE'
1953 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1954 #define TYPE_COUNT integer ,intent(in) :: Count
1955 #define NF_ROUTINE NF_PUT_ATT_DOUBLE
1956 #define ARGS NF_DOUBLE,Count,Data
1957 #include "ext_ncd_put_var_ti.code"
1958 end subroutine ext_ncd_put_var_ti_double
1960 subroutine ext_ncd_put_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,Status)
1969 #define ROUTINE_TYPE 'DOUBLE'
1970 #define TYPE_DATA real*8,intent(in) :: Data(*)
1971 #define TYPE_COUNT integer ,intent(in) :: Count
1972 #define NF_ROUTINE NF_PUT_VARA_DOUBLE
1973 #define NF_TYPE NF_DOUBLE
1974 #define LENGTH Count
1976 #include "ext_ncd_put_var_td.code"
1977 end subroutine ext_ncd_put_var_td_double
1979 subroutine ext_ncd_put_var_ti_integer(DataHandle,Element,Var,Data,Count,Status)
1986 #define ROUTINE_TYPE 'INTEGER'
1987 #define TYPE_DATA integer ,intent(in) :: Data(*)
1988 #define TYPE_COUNT integer ,intent(in) :: Count
1989 #define NF_ROUTINE NF_PUT_ATT_INT
1990 #define ARGS NF_INT,Count,Data
1991 #include "ext_ncd_put_var_ti.code"
1992 end subroutine ext_ncd_put_var_ti_integer
1994 subroutine ext_ncd_put_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,Status)
2003 #define ROUTINE_TYPE 'INTEGER'
2004 #define TYPE_DATA integer ,intent(in) :: Data(*)
2005 #define TYPE_COUNT integer ,intent(in) :: Count
2006 #define NF_ROUTINE NF_PUT_VARA_INT
2007 #define NF_TYPE NF_INT
2008 #define LENGTH Count
2010 #include "ext_ncd_put_var_td.code"
2011 end subroutine ext_ncd_put_var_td_integer
2013 subroutine ext_ncd_put_var_ti_logical(DataHandle,Element,Var,Data,Count,Status)
2019 #define ROUTINE_TYPE 'LOGICAL'
2020 #define TYPE_DATA logical ,intent(in) :: Data(*)
2021 #define TYPE_COUNT integer ,intent(in) :: Count
2022 #define NF_ROUTINE NF_PUT_ATT_INT
2024 #define ARGS NF_INT,Count,Buffer
2025 #include "ext_ncd_put_var_ti.code"
2026 end subroutine ext_ncd_put_var_ti_logical
2028 subroutine ext_ncd_put_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,Status)
2036 #define ROUTINE_TYPE 'LOGICAL'
2037 #define TYPE_DATA logical ,intent(in) :: Data(*)
2038 #define TYPE_COUNT integer ,intent(in) :: Count
2039 #define NF_ROUTINE NF_PUT_VARA_INT
2040 #define NF_TYPE NF_INT
2042 #define LENGTH Count
2044 #include "ext_ncd_put_var_td.code"
2045 end subroutine ext_ncd_put_var_td_logical
2047 subroutine ext_ncd_put_var_ti_char(DataHandle,Element,Var,Data,Status)
2054 #define ROUTINE_TYPE 'CHAR'
2055 #define TYPE_DATA character*(*) ,intent(in) :: Data
2057 #define NF_ROUTINE NF_PUT_ATT_TEXT
2058 #define ARGS len_trim(Data),trim(Data)
2060 #include "ext_ncd_put_var_ti.code"
2062 end subroutine ext_ncd_put_var_ti_char
2064 subroutine ext_ncd_put_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
2073 #define ROUTINE_TYPE 'CHAR'
2074 #define TYPE_DATA character*(*) ,intent(in) :: Data
2076 #define NF_ROUTINE NF_PUT_VARA_TEXT
2077 #define NF_TYPE NF_CHAR
2078 #define LENGTH len(Data)
2079 #include "ext_ncd_put_var_td.code"
2080 end subroutine ext_ncd_put_var_td_char
2082 subroutine ext_ncd_get_var_ti_real(DataHandle,Element,Var,Data,Count,OutCount,Status)
2087 #undef TYPE_OUTCOUNT
2091 #define ROUTINE_TYPE 'REAL'
2092 #define TYPE_DATA real ,intent(out) :: Data(*)
2093 #define TYPE_BUFFER real ,allocatable :: Buffer(:)
2094 #define TYPE_COUNT integer,intent(in) :: Count
2095 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2096 #define NF_TYPE NF_FLOAT
2097 #define NF_ROUTINE NF_GET_ATT_REAL
2098 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2099 #include "ext_ncd_get_var_ti.code"
2100 end subroutine ext_ncd_get_var_ti_real
2102 subroutine ext_ncd_get_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2107 #undef TYPE_OUTCOUNT
2112 #define ROUTINE_TYPE 'REAL'
2113 #define TYPE_DATA real ,intent(out) :: Data(*)
2114 #define TYPE_BUFFER real
2115 #define TYPE_COUNT integer,intent(in) :: Count
2116 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2117 #define NF_TYPE NF_FLOAT
2118 #define NF_ROUTINE NF_GET_VARA_REAL
2119 #define LENGTH min(Count,Len1)
2120 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2121 #include "ext_ncd_get_var_td.code"
2122 end subroutine ext_ncd_get_var_td_real
2124 subroutine ext_ncd_get_var_ti_double(DataHandle,Element,Var,Data,Count,OutCount,Status)
2129 #undef TYPE_OUTCOUNT
2133 #define ROUTINE_TYPE 'DOUBLE'
2134 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2135 #define TYPE_BUFFER real*8 ,allocatable :: Buffer(:)
2136 #define TYPE_COUNT integer,intent(in) :: Count
2137 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2138 #define NF_TYPE NF_DOUBLE
2139 #define NF_ROUTINE NF_GET_ATT_DOUBLE
2140 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2141 #include "ext_ncd_get_var_ti.code"
2142 end subroutine ext_ncd_get_var_ti_double
2144 subroutine ext_ncd_get_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2149 #undef TYPE_OUTCOUNT
2154 #define ROUTINE_TYPE 'DOUBLE'
2155 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2156 #define TYPE_BUFFER real*8
2157 #define TYPE_COUNT integer,intent(in) :: Count
2158 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2159 #define NF_TYPE NF_DOUBLE
2160 #define NF_ROUTINE NF_GET_VARA_DOUBLE
2161 #define LENGTH min(Count,Len1)
2162 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2163 #include "ext_ncd_get_var_td.code"
2164 end subroutine ext_ncd_get_var_td_double
2166 subroutine ext_ncd_get_var_ti_integer(DataHandle,Element,Var,Data,Count,OutCount,Status)
2171 #undef TYPE_OUTCOUNT
2175 #define ROUTINE_TYPE 'INTEGER'
2176 #define TYPE_DATA integer,intent(out) :: Data(*)
2177 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2178 #define TYPE_COUNT integer,intent(in) :: Count
2179 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2180 #define NF_TYPE NF_INT
2181 #define NF_ROUTINE NF_GET_ATT_INT
2182 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2183 #include "ext_ncd_get_var_ti.code"
2184 end subroutine ext_ncd_get_var_ti_integer
2186 subroutine ext_ncd_get_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2191 #undef TYPE_OUTCOUNT
2196 #define ROUTINE_TYPE 'INTEGER'
2197 #define TYPE_DATA integer,intent(out) :: Data(*)
2198 #define TYPE_BUFFER integer
2199 #define TYPE_COUNT integer,intent(in) :: Count
2200 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2201 #define NF_TYPE NF_INT
2202 #define NF_ROUTINE NF_GET_VARA_INT
2203 #define LENGTH min(Count,Len1)
2204 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2205 #include "ext_ncd_get_var_td.code"
2206 end subroutine ext_ncd_get_var_td_integer
2208 subroutine ext_ncd_get_var_ti_logical(DataHandle,Element,Var,Data,Count,OutCount,Status)
2213 #undef TYPE_OUTCOUNT
2217 #define ROUTINE_TYPE 'LOGICAL'
2218 #define TYPE_DATA logical,intent(out) :: Data(*)
2219 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2220 #define TYPE_COUNT integer,intent(in) :: Count
2221 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2222 #define NF_TYPE NF_INT
2223 #define NF_ROUTINE NF_GET_ATT_INT
2224 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))==1
2225 #include "ext_ncd_get_var_ti.code"
2226 end subroutine ext_ncd_get_var_ti_logical
2228 subroutine ext_ncd_get_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2233 #undef TYPE_OUTCOUNT
2238 #define ROUTINE_TYPE 'LOGICAL'
2239 #define TYPE_DATA logical,intent(out) :: Data(*)
2240 #define TYPE_BUFFER integer
2241 #define TYPE_COUNT integer,intent(in) :: Count
2242 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2243 #define NF_TYPE NF_INT
2244 #define NF_ROUTINE NF_GET_VARA_INT
2245 #define LENGTH min(Count,Len1)
2246 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))==1
2247 #include "ext_ncd_get_var_td.code"
2248 end subroutine ext_ncd_get_var_td_logical
2250 subroutine ext_ncd_get_var_ti_char(DataHandle,Element,Var,Data,Status)
2255 #undef TYPE_OUTCOUNT
2259 #define ROUTINE_TYPE 'CHAR'
2260 #define TYPE_DATA character*(*) ,intent(out) :: Data
2262 #define TYPE_COUNT integer :: Count = 1
2263 #define TYPE_OUTCOUNT
2264 #define NF_TYPE NF_CHAR
2265 #define NF_ROUTINE NF_GET_ATT_TEXT
2268 #include "ext_ncd_get_var_ti.code"
2270 end subroutine ext_ncd_get_var_ti_char
2272 subroutine ext_ncd_get_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
2277 #undef TYPE_OUTCOUNT
2281 #define ROUTINE_TYPE 'CHAR'
2282 #define TYPE_DATA character*(*) ,intent(out) :: Data
2283 #define TYPE_BUFFER character (80)
2284 #define TYPE_COUNT integer :: Count = 1
2285 #define TYPE_OUTCOUNT
2286 #define NF_TYPE NF_CHAR
2287 #define NF_ROUTINE NF_GET_VARA_TEXT
2290 #include "ext_ncd_get_var_td.code"
2292 end subroutine ext_ncd_get_var_td_char
2294 subroutine ext_ncd_put_dom_td_real(DataHandle,Element,DateStr,Data,Count,Status)
2295 integer ,intent(in) :: DataHandle
2296 character*(*) ,intent(in) :: Element
2297 character*(*) ,intent(in) :: DateStr
2298 real ,intent(in) :: Data(*)
2299 integer ,intent(in) :: Count
2300 integer ,intent(out) :: Status
2302 call ext_ncd_put_var_td_real(DataHandle,Element,DateStr, &
2303 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2305 end subroutine ext_ncd_put_dom_td_real
2307 subroutine ext_ncd_put_dom_td_integer(DataHandle,Element,DateStr,Data,Count,Status)
2308 integer ,intent(in) :: DataHandle
2309 character*(*) ,intent(in) :: Element
2310 character*(*) ,intent(in) :: DateStr
2311 integer ,intent(in) :: Data(*)
2312 integer ,intent(in) :: Count
2313 integer ,intent(out) :: Status
2315 call ext_ncd_put_var_td_integer(DataHandle,Element,DateStr, &
2316 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2318 end subroutine ext_ncd_put_dom_td_integer
2320 subroutine ext_ncd_put_dom_td_double(DataHandle,Element,DateStr,Data,Count,Status)
2321 integer ,intent(in) :: DataHandle
2322 character*(*) ,intent(in) :: Element
2323 character*(*) ,intent(in) :: DateStr
2324 real*8 ,intent(in) :: Data(*)
2325 integer ,intent(in) :: Count
2326 integer ,intent(out) :: Status
2328 call ext_ncd_put_var_td_double(DataHandle,Element,DateStr, &
2329 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2331 end subroutine ext_ncd_put_dom_td_double
2333 subroutine ext_ncd_put_dom_td_logical(DataHandle,Element,DateStr,Data,Count,Status)
2334 integer ,intent(in) :: DataHandle
2335 character*(*) ,intent(in) :: Element
2336 character*(*) ,intent(in) :: DateStr
2337 logical ,intent(in) :: Data(*)
2338 integer ,intent(in) :: Count
2339 integer ,intent(out) :: Status
2341 call ext_ncd_put_var_td_logical(DataHandle,Element,DateStr, &
2342 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2344 end subroutine ext_ncd_put_dom_td_logical
2346 subroutine ext_ncd_put_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2347 integer ,intent(in) :: DataHandle
2348 character*(*) ,intent(in) :: Element
2349 character*(*) ,intent(in) :: DateStr
2350 character*(*) ,intent(in) :: Data
2351 integer ,intent(out) :: Status
2353 call ext_ncd_put_var_td_char(DataHandle,Element,DateStr, &
2354 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2356 end subroutine ext_ncd_put_dom_td_char
2358 subroutine ext_ncd_get_dom_td_real(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2359 integer ,intent(in) :: DataHandle
2360 character*(*) ,intent(in) :: Element
2361 character*(*) ,intent(in) :: DateStr
2362 real ,intent(out) :: Data(*)
2363 integer ,intent(in) :: Count
2364 integer ,intent(out) :: OutCount
2365 integer ,intent(out) :: Status
2366 call ext_ncd_get_var_td_real(DataHandle,Element,DateStr, &
2367 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2369 end subroutine ext_ncd_get_dom_td_real
2371 subroutine ext_ncd_get_dom_td_integer(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2372 integer ,intent(in) :: DataHandle
2373 character*(*) ,intent(in) :: Element
2374 character*(*) ,intent(in) :: DateStr
2375 integer ,intent(out) :: Data(*)
2376 integer ,intent(in) :: Count
2377 integer ,intent(out) :: OutCount
2378 integer ,intent(out) :: Status
2379 call ext_ncd_get_var_td_integer(DataHandle,Element,DateStr, &
2380 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2382 end subroutine ext_ncd_get_dom_td_integer
2384 subroutine ext_ncd_get_dom_td_double(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2385 integer ,intent(in) :: DataHandle
2386 character*(*) ,intent(in) :: Element
2387 character*(*) ,intent(in) :: DateStr
2388 real*8 ,intent(out) :: Data(*)
2389 integer ,intent(in) :: Count
2390 integer ,intent(out) :: OutCount
2391 integer ,intent(out) :: Status
2392 call ext_ncd_get_var_td_double(DataHandle,Element,DateStr, &
2393 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2395 end subroutine ext_ncd_get_dom_td_double
2397 subroutine ext_ncd_get_dom_td_logical(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2398 integer ,intent(in) :: DataHandle
2399 character*(*) ,intent(in) :: Element
2400 character*(*) ,intent(in) :: DateStr
2401 logical ,intent(out) :: Data(*)
2402 integer ,intent(in) :: Count
2403 integer ,intent(out) :: OutCount
2404 integer ,intent(out) :: Status
2405 call ext_ncd_get_var_td_logical(DataHandle,Element,DateStr, &
2406 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2408 end subroutine ext_ncd_get_dom_td_logical
2410 subroutine ext_ncd_get_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2411 integer ,intent(in) :: DataHandle
2412 character*(*) ,intent(in) :: Element
2413 character*(*) ,intent(in) :: DateStr
2414 character*(*) ,intent(out) :: Data
2415 integer ,intent(out) :: Status
2416 call ext_ncd_get_var_td_char(DataHandle,Element,DateStr, &
2417 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2419 end subroutine ext_ncd_get_dom_td_char
2421 subroutine ext_ncd_write_field(DataHandle,DateStr,Var,Field,FieldTypeIn, &
2422 Comm, IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2423 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2425 use ext_ncd_support_routines
2427 include 'wrf_status_codes.h'
2428 include 'netcdf.inc'
2429 integer ,intent(in) :: DataHandle
2430 character*(*) ,intent(in) :: DateStr
2431 character*(*) ,intent(in) :: Var
2432 integer ,intent(inout) :: Field(*)
2433 integer ,intent(in) :: FieldTypeIn
2434 integer ,intent(inout) :: Comm
2435 integer ,intent(inout) :: IOComm
2436 integer ,intent(in) :: DomainDesc
2437 character*(*) ,intent(in) :: MemoryOrdIn
2438 character*(*) ,intent(in) :: Stagger ! Dummy for now
2439 character*(*) ,dimension(*) ,intent(in) :: DimNames
2440 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2441 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2442 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2443 integer ,intent(out) :: Status
2444 integer :: FieldType
2445 character (3) :: MemoryOrder
2446 type(wrf_data_handle) ,pointer :: DH
2449 character (VarNameLen) :: VarName
2450 character (3) :: MemO
2451 character (3) :: UCMemO
2453 integer ,dimension(NVarDims) :: Length
2454 integer ,dimension(NVarDims) :: VDimIDs
2455 character(80),dimension(NVarDims) :: RODimNames
2456 integer ,dimension(NVarDims) :: StoredStart
2457 integer ,dimension(:,:,:,:),allocatable :: XField
2461 integer :: i1,i2,j1,j2,k1,k2
2462 integer :: x1,x2,y1,y2,z1,z2
2463 integer :: l1,l2,m1,m2,n1,n2
2466 character (80) :: NullName
2469 #ifdef USE_NETCDF4_FEATURES
2470 integer, parameter :: cache_size = 32000000
2471 integer,dimension(NVarDims) :: chunks
2472 logical :: need_chunking
2473 integer :: compression_level
2474 integer :: block_size
2477 MemoryOrder = trim(adjustl(MemoryOrdIn))
2479 call GetDim(MemoryOrder,NDim,Status)
2480 if(Status /= WRF_NO_ERR) then
2481 write(msg,*) 'Warning BAD MEMORY ORDER |',MemoryOrder,'| in ',__FILE__,', line', __LINE__
2482 call wrf_debug ( WARN , TRIM(msg))
2486 call DateCheck(DateStr,Status)
2487 if(Status /= WRF_NO_ERR) then
2488 write(msg,*) 'Warning DATE STRING ERROR |',DateStr,'| in ',__FILE__,', line', __LINE__
2489 call wrf_debug ( WARN , TRIM(msg))
2493 call GetDH(DataHandle,DH,Status)
2494 if(Status /= WRF_NO_ERR) then
2495 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2496 call wrf_debug ( WARN , TRIM(msg))
2501 #ifdef USE_NETCDF4_FEATURES
2502 if ( .not. DH%use_netcdf_classic ) then
2503 call set_chunking(MemoryOrder,need_chunking)
2504 compression_level = 2
2506 need_chunking = .false.
2510 if ( DH%R4OnOutput .and. FieldTypeIn == WRF_DOUBLE ) then
2511 FieldType = WRF_REAL
2513 FieldType = FieldTypeIn
2516 write(msg,*)'ext_ncd_write_field: called for ',TRIM(Var)
2518 !jm 010827 Length(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
2520 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2522 IF ( ZeroLengthHorzDim(MemoryOrder,Length,Status) ) THEN
2523 write(msg,*)'ext_ncd_write_field: zero length dimension in ',TRIM(Var),'. Ignoring'
2524 call wrf_debug ( WARN , TRIM(msg))
2528 call ExtOrder(MemoryOrder,Length,Status)
2529 call ExtOrderStr(MemoryOrder,DimNames,RODimNames,Status)
2530 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2531 Status = WRF_WARN_FILE_NOT_OPENED
2532 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2533 call wrf_debug ( WARN , TRIM(msg))
2534 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
2535 Status = WRF_WARN_WRITE_RONLY_FILE
2536 write(msg,*) 'Warning WRITE READ ONLY FILE in ',__FILE__,', line', __LINE__
2537 call wrf_debug ( WARN , TRIM(msg))
2538 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2540 if(DH%VarNames(NVar) == VarName ) then
2541 Status = WRF_WARN_2DRYRUNS_1VARIABLE
2542 write(msg,*) 'Warning 2 DRYRUNS 1 VARIABLE in ',__FILE__,', line', __LINE__
2543 call wrf_debug ( WARN , TRIM(msg))
2545 elseif(DH%VarNames(NVar) == NO_NAME) then
2546 DH%VarNames(NVar) = VarName
2549 elseif(NVar == MaxVars) then
2550 Status = WRF_WARN_TOO_MANY_VARIABLES
2551 write(msg,*) 'Warning TOO MANY VARIABLES in ',__FILE__,', line', __LINE__
2552 call wrf_debug ( WARN , TRIM(msg))
2557 if(RODimNames(j) == NullName .or. RODimNames(j) == '') then
2559 if(DH%DimLengths(i) == Length(j)) then
2561 elseif(DH%DimLengths(i) == NO_DIM) then
2562 stat = NF_DEF_DIM(NCID,DH%DimNames(i),Length(j),DH%DimIDs(i))
2563 call netcdf_err(stat,Status)
2564 if(Status /= WRF_NO_ERR) then
2565 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2566 call wrf_debug ( WARN , TRIM(msg))
2569 DH%DimLengths(i) = Length(j)
2571 elseif(i == MaxDims) then
2572 Status = WRF_WARN_TOO_MANY_DIMS
2573 write(msg,*) 'Warning TOO MANY DIMENSIONS in ',__FILE__,', line', __LINE__
2574 call wrf_debug ( WARN , TRIM(msg))
2578 else !look for input name and check if already defined
2581 if (DH%DimNames(i) == RODimNames(j)) then
2582 if (DH%DimLengths(i) == Length(j)) then
2586 Status = WRF_WARN_DIMNAME_REDEFINED
2587 write(msg,*) 'Warning DIM ',i,', NAME ',TRIM(DH%DimNames(i)),' REDEFINED by var ', &
2588 TRIM(Var),' ',DH%DimLengths(i),Length(j) ,' in ', __FILE__ ,' line', __LINE__
2589 call wrf_debug ( WARN , TRIM(msg))
2596 if (DH%DimLengths(i) == NO_DIM) then
2597 DH%DimNames(i) = RODimNames(j)
2598 stat = NF_DEF_DIM(NCID,DH%DimNames(i),Length(j),DH%DimIDs(i))
2599 call netcdf_err(stat,Status)
2600 if(Status /= WRF_NO_ERR) then
2601 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2602 call wrf_debug ( WARN , TRIM(msg))
2605 DH%DimLengths(i) = Length(j)
2607 elseif(i == MaxDims) then
2608 Status = WRF_WARN_TOO_MANY_DIMS
2609 write(msg,*) 'Warning TOO MANY DIMENSIONS in ',__FILE__,', line', __LINE__
2610 call wrf_debug ( WARN , TRIM(msg))
2616 VDimIDs(j) = DH%DimIDs(i)
2617 DH%VarDimLens(j,NVar) = Length(j)
2619 VDimIDs(NDim+1) = DH%DimUnlimID
2621 ! Do not use SELECT statement here as sometimes WRF_REAL=WRF_DOUBLE
2622 IF (FieldType == WRF_REAL) THEN
2624 ELSE IF (FieldType == WRF_DOUBLE) THEN
2626 ELSE IF (FieldType == WRF_INTEGER) THEN
2628 ELSE IF (FieldType == WRF_LOGICAL) THEN
2631 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
2632 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
2633 call wrf_debug ( WARN , TRIM(msg))
2637 stat = NF_DEF_VAR(NCID,VarName,XType,NDim+1,VDimIDs,VarID)
2638 call netcdf_err(stat,Status)
2639 if(Status /= WRF_NO_ERR) then
2640 write(msg,*) 'ext_ncd_write_field: NetCDF error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2641 call wrf_debug ( WARN , TRIM(msg))
2645 #ifdef USE_NETCDF4_FEATURES
2646 if(need_chunking) then
2647 chunks(1:NDim) = Length(1:NDim)
2649 chunks(1) = (Length(1) + 1)/2
2650 chunks(2) = (Length(2) + 1)/2
2654 block_size = block_size * chunks(i)
2657 do while (block_size > cache_size)
2658 chunks(1) = (chunks(1) + 1)/2
2659 chunks(2) = (chunks(2) + 1)/2
2663 block_size = block_size * chunks(i)
2667 ! write(unit=0, fmt='(2x, 3a,i6)') 'file: ', __FILE__, ', line: ', __LINE__
2668 ! write(unit=0, fmt='(2x, 3a)') TRIM(VarName),':'
2669 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length 1 = ', Length(1), ', chunk 1 = ', chunks(1)
2670 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length 2 = ', Length(2), ', chunk 2 = ', chunks(2)
2671 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length NDim+1 = ', Length(NDim+1), ', chunk NDim+1 = ', chunks(NDim+1)
2672 ! write(unit=0, fmt='(10x, a,i6)') 'compression_level = ', compression_level
2674 stat = NF_DEF_VAR_CHUNKING(NCID, VarID, NF_CHUNKED, chunks(1:NDim+1))
2675 call netcdf_err(stat,Status)
2676 if(Status /= WRF_NO_ERR) then
2677 write(msg,*) 'ext_ncd_write_field: NetCDF def chunking error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2678 call wrf_debug ( WARN , TRIM(msg))
2682 stat = NF_DEF_VAR_DEFLATE(NCID, VarID, 1, 1, compression_level)
2683 call netcdf_err(stat,Status)
2684 if(Status /= WRF_NO_ERR) then
2685 write(msg,*) 'ext_ncd_write_field: NetCDF def compression error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2686 call wrf_debug ( WARN , TRIM(msg))
2692 DH%VarIDs(NVar) = VarID
2693 stat = NF_PUT_ATT_INT(NCID,VarID,'FieldType',NF_INT,1,FieldType)
2694 call netcdf_err(stat,Status)
2695 if(Status /= WRF_NO_ERR) then
2696 write(msg,*) 'ext_ncd_write_field: NetCDF error in ',__FILE__,', line', __LINE__
2697 call wrf_debug ( WARN , TRIM(msg))
2700 call reorder(MemoryOrder,MemO)
2701 call uppercase(MemO,UCMemO)
2702 stat = NF_PUT_ATT_TEXT(NCID,VarID,'MemoryOrder',3,UCMemO)
2703 call netcdf_err(stat,Status)
2704 if(Status /= WRF_NO_ERR) then
2705 write(msg,*) 'ext_ncd_write_field: NetCDF error in ',__FILE__,', line', __LINE__
2706 call wrf_debug ( WARN , TRIM(msg))
2709 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
2710 do NVar=1,DH%NumVars
2711 if(DH%VarNames(NVar) == VarName) then
2713 elseif(NVar == DH%NumVars) then
2714 Status = WRF_WARN_VAR_NF
2715 write(msg,*) 'Warning VARIABLE NOT FOUND in ',__FILE__,', line', __LINE__
2716 call wrf_debug ( WARN , TRIM(msg))
2720 VarID = DH%VarIDs(NVar)
2722 if(Length(j) /= DH%VarDimLens(j,NVar) .AND. DH%FileStatus /= WRF_FILE_OPENED_FOR_UPDATE ) then
2723 Status = WRF_WARN_WRTLEN_NE_DRRUNLEN
2724 write(msg,*) 'Warning LENGTH != DRY RUN LENGTH for |', &
2725 VarName,'| dim ',j,' in ',__FILE__,', line', __LINE__
2726 call wrf_debug ( WARN , TRIM(msg))
2727 write(msg,*) ' LENGTH ',Length(j),' DRY RUN LENGTH ',DH%VarDimLens(j,NVar)
2728 call wrf_debug ( WARN , TRIM(msg))
2730 !jm 010825 elseif(DomainStart(j) < MemoryStart(j)) then
2731 elseif(PatchStart(j) < MemoryStart(j)) then
2732 Status = WRF_WARN_DIMENSION_ERROR
2733 write(msg,*) 'Warning DIMENSION ERROR for |',VarName, &
2734 '| in ',__FILE__,', line', __LINE__
2735 call wrf_debug ( WARN , TRIM(msg))
2740 call GetIndices(NDim,MemoryStart,MemoryEnd,l1,l2,m1,m2,n1,n2)
2741 call GetIndices(NDim,StoredStart,Length ,x1,x2,y1,y2,z1,z2)
2742 call GetIndices(NDim,PatchStart, PatchEnd ,i1,i2,j1,j2,k1,k2)
2744 if(FieldType == WRF_DOUBLE) di=2
2745 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
2747 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
2748 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
2749 call wrf_debug ( FATAL , TRIM(msg))
2752 if (DH%R4OnOutput .and. FieldTypeIn == WRF_DOUBLE) then
2753 call TransposeToR4('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2754 ,XField,x1,x2,y1,y2,z1,z2 &
2755 ,i1,i2,j1,j2,k1,k2 )
2757 call Transpose('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2758 ,XField,x1,x2,y1,y2,z1,z2 &
2759 ,i1,i2,j1,j2,k1,k2 )
2761 call FieldIO('write',DataHandle,DateStr,Length,MemoryOrder, &
2762 FieldType,NCID,VarID,XField,Status)
2763 if(Status /= WRF_NO_ERR) then
2764 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2765 call wrf_debug ( WARN , TRIM(msg))
2768 deallocate(XField, STAT=stat)
2770 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
2771 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
2772 call wrf_debug ( FATAL , TRIM(msg))
2776 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2777 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2778 call wrf_debug ( FATAL , TRIM(msg))
2780 DH%first_operation = .FALSE.
2782 end subroutine ext_ncd_write_field
2784 subroutine ext_ncd_read_field(DataHandle,DateStr,Var,Field,FieldType,Comm, &
2785 IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2786 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2788 use ext_ncd_support_routines
2790 include 'wrf_status_codes.h'
2791 include 'netcdf.inc'
2792 integer ,intent(in) :: DataHandle
2793 character*(*) ,intent(in) :: DateStr
2794 character*(*) ,intent(in) :: Var
2795 integer ,intent(out) :: Field(*)
2796 integer ,intent(in) :: FieldType
2797 integer ,intent(inout) :: Comm
2798 integer ,intent(inout) :: IOComm
2799 integer ,intent(in) :: DomainDesc
2800 character*(*) ,intent(in) :: MemoryOrdIn
2801 character*(*) ,intent(in) :: Stagger ! Dummy for now
2802 character*(*) , dimension (*) ,intent(in) :: DimNames
2803 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2804 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2805 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2806 integer ,intent(out) :: Status
2807 character (3) :: MemoryOrder
2808 character (NF_MAX_NAME) :: dimname
2809 type(wrf_data_handle) ,pointer :: DH
2812 character (VarNameLen) :: VarName
2814 integer ,dimension(NVarDims) :: VCount
2815 integer ,dimension(NVarDims) :: VStart
2816 integer ,dimension(NVarDims) :: Length
2817 integer ,dimension(NVarDims) :: VDimIDs
2818 integer ,dimension(NVarDims) :: MemS
2819 integer ,dimension(NVarDims) :: MemE
2820 integer ,dimension(NVarDims) :: StoredStart
2821 integer ,dimension(NVarDims) :: StoredLen
2822 integer ,dimension(:,:,:,:) ,allocatable :: XField
2825 integer :: i1,i2,j1,j2,k1,k2
2826 integer :: x1,x2,y1,y2,z1,z2
2827 integer :: l1,l2,m1,m2,n1,n2
2828 character (VarNameLen) :: Name
2830 integer :: StoredDim
2837 MemoryOrder = trim(adjustl(MemoryOrdIn))
2838 call GetDim(MemoryOrder,NDim,Status)
2839 if(Status /= WRF_NO_ERR) then
2840 write(msg,*) 'Warning BAD MEMORY ORDER |',TRIM(MemoryOrder),'| for |', &
2841 TRIM(Var),'| in ext_ncd_read_field ',__FILE__,', line', __LINE__
2842 call wrf_debug ( WARN , TRIM(msg))
2845 call DateCheck(DateStr,Status)
2846 if(Status /= WRF_NO_ERR) then
2847 write(msg,*) 'Warning DATE STRING ERROR |',TRIM(DateStr),'| for |',TRIM(Var), &
2848 '| in ext_ncd_read_field ',__FILE__,', line', __LINE__
2849 call wrf_debug ( WARN , TRIM(msg))
2853 call GetDH(DataHandle,DH,Status)
2854 if(Status /= WRF_NO_ERR) then
2855 write(msg,*) 'Warning Status = ',Status,' in ext_ncd_read_field ',__FILE__,', line', __LINE__
2856 call wrf_debug ( WARN , TRIM(msg))
2859 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2860 Status = WRF_WARN_FILE_NOT_OPENED
2861 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2862 call wrf_debug ( WARN , TRIM(msg))
2863 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2864 ! jm it is okay to have a dry run read. means read is called between ofrb and ofrc. Just return.
2865 ! Status = WRF_WARN_DRYRUN_READ
2866 ! write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
2867 ! call wrf_debug ( WARN , TRIM(msg))
2870 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
2871 Status = WRF_WARN_READ_WONLY_FILE
2872 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
2873 call wrf_debug ( WARN , TRIM(msg))
2874 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
2877 !jm Length(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
2878 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2879 call ExtOrder(MemoryOrder,Length,Status)
2880 stat = NF_INQ_VARID(NCID,VarName,VarID)
2881 call netcdf_err(stat,Status)
2882 if(Status /= WRF_NO_ERR) then
2883 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__,' Varname ',Varname
2884 call wrf_debug ( WARN , TRIM(msg))
2887 stat = NF_INQ_VAR(NCID,VarID,Name,XType,StoredDim,VDimIDs,NAtts)
2888 call netcdf_err(stat,Status)
2889 if(Status /= WRF_NO_ERR) then
2890 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2891 call wrf_debug ( WARN , TRIM(msg))
2894 stat = NF_GET_ATT_INT(NCID,VarID,'FieldType',FType)
2895 call netcdf_err(stat,Status)
2896 if(Status /= WRF_NO_ERR) then
2897 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2898 call wrf_debug ( WARN , TRIM(msg))
2901 ! allow coercion between double and single prec real
2902 !jm if(FieldType /= Ftype) then
2903 if( (FieldType == WRF_REAL .OR. FieldType == WRF_DOUBLE) ) then
2904 if ( .NOT. (Ftype == WRF_REAL .OR. Ftype == WRF_DOUBLE )) then
2905 Status = WRF_WARN_TYPE_MISMATCH
2906 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
2907 call wrf_debug ( WARN , TRIM(msg))
2910 else if(FieldType /= Ftype) then
2911 Status = WRF_WARN_TYPE_MISMATCH
2912 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
2913 call wrf_debug ( WARN , TRIM(msg))
2917 ! Do not use SELECT statement here as sometimes WRF_REAL=WRF_DOUBLE
2918 IF (FieldType == WRF_REAL) THEN
2919 ! allow coercion between double and single prec real
2920 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
2921 Status = WRF_WARN_TYPE_MISMATCH
2922 write(msg,*) 'Warning REAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
2924 ELSE IF (FieldType == WRF_DOUBLE) THEN
2925 ! allow coercion between double and single prec real
2926 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
2927 Status = WRF_WARN_TYPE_MISMATCH
2928 write(msg,*) 'Warning DOUBLE TYPE MISMATCH in ',__FILE__,', line', __LINE__
2930 ELSE IF (FieldType == WRF_INTEGER) THEN
2931 if(XType /= NF_INT) then
2932 Status = WRF_WARN_TYPE_MISMATCH
2933 write(msg,*) 'Warning INTEGER TYPE MISMATCH in ',__FILE__,', line', __LINE__
2935 ELSE IF (FieldType == WRF_LOGICAL) THEN
2936 if(XType /= NF_INT) then
2937 Status = WRF_WARN_TYPE_MISMATCH
2938 write(msg,*) 'Warning LOGICAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
2941 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
2942 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
2945 if(Status /= WRF_NO_ERR) then
2946 call wrf_debug ( WARN , TRIM(msg))
2949 ! NDim=0 for scalars. Handle read of old NDim=1 files. TBH: 20060502
2950 IF ( ( NDim == 0 ) .AND. ( StoredDim == 2 ) ) THEN
2951 stat = NF_INQ_DIMNAME(NCID,VDimIDs(1),dimname)
2952 call netcdf_err(stat,Status)
2953 if(Status /= WRF_NO_ERR) then
2954 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2955 call wrf_debug ( WARN , TRIM(msg))
2958 IF ( dimname(1:10) == 'ext_scalar' ) THEN
2963 if(StoredDim /= NDim+1) then
2964 Status = WRF_ERR_FATAL_BAD_VARIABLE_DIM
2965 write(msg,*) 'Fatal error BAD VARIABLE DIMENSION in ext_ncd_read_field ',TRIM(Var),TRIM(DateStr)
2966 call wrf_debug ( FATAL , msg)
2967 write(msg,*) ' StoredDim ', StoredDim, ' .NE. NDim+1 ', NDim+1
2968 call wrf_debug ( FATAL , msg)
2972 stat = NF_INQ_DIMLEN(NCID,VDimIDs(j),StoredLen(j))
2973 call netcdf_err(stat,Status)
2974 if(Status /= WRF_NO_ERR) then
2975 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2976 call wrf_debug ( WARN , TRIM(msg))
2979 if(Length(j) > StoredLen(j)) then
2980 Status = WRF_WARN_READ_PAST_EOF
2981 write(msg,*) 'Warning READ PAST EOF in ext_ncd_read_field of ',TRIM(Var),Length(j),'>',StoredLen(j)
2982 call wrf_debug ( WARN , TRIM(msg))
2984 elseif(Length(j) <= 0) then
2985 Status = WRF_WARN_ZERO_LENGTH_READ
2986 write(msg,*) 'Warning ZERO LENGTH READ in ',__FILE__,', line', __LINE__
2987 call wrf_debug ( WARN , TRIM(msg))
2989 elseif(DomainStart(j) < MemoryStart(j)) then
2990 Status = WRF_WARN_DIMENSION_ERROR
2991 write(msg,*) 'Warning dim ',j,' DomainStart (',DomainStart(j), &
2992 ') < MemoryStart (',MemoryStart(j),') in ',__FILE__,', line', __LINE__
2993 call wrf_debug ( WARN , TRIM(msg))
2999 call GetIndices(NDim,MemoryStart,MemoryEnd,l1,l2,m1,m2,n1,n2)
3000 call GetIndices(NDim,StoredStart,StoredLen,x1,x2,y1,y2,z1,z2)
3001 !jm call GetIndices(NDim,DomainStart,DomainEnd,i1,i2,j1,j2,k1,k2)
3002 call GetIndices(NDim,PatchStart,PatchEnd,i1,i2,j1,j2,k1,k2)
3005 if(FieldType == WRF_DOUBLE) di=2
3006 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
3008 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
3009 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
3010 call wrf_debug ( FATAL , msg)
3013 call FieldIO('read',DataHandle,DateStr,Length,MemoryOrder, &
3014 FieldType,NCID,VarID,XField,Status)
3015 if(Status /= WRF_NO_ERR) then
3016 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3017 call wrf_debug ( WARN , TRIM(msg))
3020 call Transpose('read',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
3021 ,XField,x1,x2,y1,y2,z1,z2 &
3022 ,i1,i2,j1,j2,k1,k2 )
3023 deallocate(XField, STAT=stat)
3025 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
3026 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
3027 call wrf_debug ( FATAL , msg)
3031 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3032 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3033 call wrf_debug ( FATAL , msg)
3035 DH%first_operation = .FALSE.
3037 end subroutine ext_ncd_read_field
3039 subroutine ext_ncd_inquire_opened( DataHandle, FileName , FileStatus, Status )
3041 use ext_ncd_support_routines
3043 include 'wrf_status_codes.h'
3044 integer ,intent(in) :: DataHandle
3045 character*(*) ,intent(inout) :: FileName
3046 integer ,intent(out) :: FileStatus
3047 integer ,intent(out) :: Status
3048 type(wrf_data_handle) ,pointer :: DH
3050 !call upgrade_filename(FileName)
3052 call GetDH(DataHandle,DH,Status)
3053 if(Status /= WRF_NO_ERR) then
3054 FileStatus = WRF_FILE_NOT_OPENED
3057 if(trim(FileName) /= trim(DH%FileName)) then
3058 FileStatus = WRF_FILE_NOT_OPENED
3060 FileStatus = DH%FileStatus
3064 end subroutine ext_ncd_inquire_opened
3066 subroutine ext_ncd_inquire_filename( Datahandle, FileName, FileStatus, Status )
3068 use ext_ncd_support_routines
3070 include 'wrf_status_codes.h'
3071 integer ,intent(in) :: DataHandle
3072 character*(*) ,intent(out) :: FileName
3073 integer ,intent(out) :: FileStatus
3074 integer ,intent(out) :: Status
3075 type(wrf_data_handle) ,pointer :: DH
3076 FileStatus = WRF_FILE_NOT_OPENED
3077 call GetDH(DataHandle,DH,Status)
3078 if(Status /= WRF_NO_ERR) then
3079 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3080 call wrf_debug ( WARN , TRIM(msg))
3083 FileName = trim(DH%FileName)
3084 FileStatus = DH%FileStatus
3087 end subroutine ext_ncd_inquire_filename
3089 subroutine ext_ncd_set_time(DataHandle, DateStr, Status)
3091 use ext_ncd_support_routines
3093 include 'wrf_status_codes.h'
3094 integer ,intent(in) :: DataHandle
3095 character*(*) ,intent(in) :: DateStr
3096 integer ,intent(out) :: Status
3097 type(wrf_data_handle) ,pointer :: DH
3100 call DateCheck(DateStr,Status)
3101 if(Status /= WRF_NO_ERR) then
3102 write(msg,*) 'Warning DATE STRING ERROR in ',__FILE__,', line', __LINE__
3103 call wrf_debug ( WARN , TRIM(msg))
3106 call GetDH(DataHandle,DH,Status)
3107 if(Status /= WRF_NO_ERR) then
3108 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3109 call wrf_debug ( WARN , TRIM(msg))
3112 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3113 Status = WRF_WARN_FILE_NOT_OPENED
3114 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3115 call wrf_debug ( WARN , TRIM(msg))
3116 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3117 Status = WRF_WARN_FILE_NOT_COMMITTED
3118 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
3119 call wrf_debug ( WARN , TRIM(msg))
3120 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3121 Status = WRF_WARN_READ_WONLY_FILE
3122 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3123 call wrf_debug ( WARN , TRIM(msg))
3124 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
3126 if(DH%Times(i)==DateStr) then
3130 if(i==MaxTimes) then
3131 Status = WRF_WARN_TIME_NF
3135 DH%CurrentVariable = 0
3138 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3139 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3140 call wrf_debug ( FATAL , msg)
3143 end subroutine ext_ncd_set_time
3145 subroutine ext_ncd_get_next_time(DataHandle, DateStr, Status)
3147 use ext_ncd_support_routines
3149 include 'wrf_status_codes.h'
3150 integer ,intent(in) :: DataHandle
3151 character*(*) ,intent(out) :: DateStr
3152 integer ,intent(out) :: Status
3153 type(wrf_data_handle) ,pointer :: DH
3155 call GetDH(DataHandle,DH,Status)
3156 if(Status /= WRF_NO_ERR) then
3157 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3158 call wrf_debug ( WARN , TRIM(msg))
3161 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3162 Status = WRF_WARN_FILE_NOT_OPENED
3163 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3164 call wrf_debug ( WARN , TRIM(msg))
3165 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3166 Status = WRF_WARN_DRYRUN_READ
3167 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3168 call wrf_debug ( WARN , TRIM(msg))
3169 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3170 Status = WRF_WARN_READ_WONLY_FILE
3171 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3172 call wrf_debug ( WARN , TRIM(msg))
3173 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
3174 if(DH%CurrentTime >= DH%NumberTimes) then
3175 Status = WRF_WARN_TIME_EOF
3178 DH%CurrentTime = DH%CurrentTime +1
3179 DateStr = DH%Times(DH%CurrentTime)
3180 DH%CurrentVariable = 0
3183 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3184 write(msg,*) 'DH%FileStatus ',DH%FileStatus
3185 call wrf_debug ( FATAL , msg)
3186 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3187 call wrf_debug ( FATAL , msg)
3190 end subroutine ext_ncd_get_next_time
3192 subroutine ext_ncd_get_previous_time(DataHandle, DateStr, Status)
3194 use ext_ncd_support_routines
3196 include 'wrf_status_codes.h'
3197 integer ,intent(in) :: DataHandle
3198 character*(*) ,intent(out) :: DateStr
3199 integer ,intent(out) :: Status
3200 type(wrf_data_handle) ,pointer :: DH
3202 call GetDH(DataHandle,DH,Status)
3203 if(Status /= WRF_NO_ERR) then
3204 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3205 call wrf_debug ( WARN , TRIM(msg))
3208 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3209 Status = WRF_WARN_FILE_NOT_OPENED
3210 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3211 call wrf_debug ( WARN , TRIM(msg))
3212 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3213 Status = WRF_WARN_DRYRUN_READ
3214 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3215 call wrf_debug ( WARN , TRIM(msg))
3216 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3217 Status = WRF_WARN_READ_WONLY_FILE
3218 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3219 call wrf_debug ( WARN , TRIM(msg))
3220 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
3221 if(DH%CurrentTime.GT.0) then
3222 DH%CurrentTime = DH%CurrentTime -1
3224 DateStr = DH%Times(MAX(1,DH%CurrentTime))
3225 DH%CurrentVariable = 0
3228 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3229 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3230 call wrf_debug ( FATAL , msg)
3233 end subroutine ext_ncd_get_previous_time
3235 subroutine ext_ncd_get_next_var(DataHandle, VarName, Status)
3237 use ext_ncd_support_routines
3239 include 'wrf_status_codes.h'
3240 include 'netcdf.inc'
3241 integer ,intent(in) :: DataHandle
3242 character*(*) ,intent(out) :: VarName
3243 integer ,intent(out) :: Status
3244 type(wrf_data_handle) ,pointer :: DH
3246 character (80) :: Name
3248 call GetDH(DataHandle,DH,Status)
3249 if(Status /= WRF_NO_ERR) then
3250 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3251 call wrf_debug ( WARN , TRIM(msg))
3254 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3255 Status = WRF_WARN_FILE_NOT_OPENED
3256 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3257 call wrf_debug ( WARN , TRIM(msg))
3258 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3259 Status = WRF_WARN_DRYRUN_READ
3260 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3261 call wrf_debug ( WARN , TRIM(msg))
3262 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3263 Status = WRF_WARN_READ_WONLY_FILE
3264 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3265 call wrf_debug ( WARN , TRIM(msg))
3266 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3268 DH%CurrentVariable = DH%CurrentVariable +1
3269 if(DH%CurrentVariable > DH%NumVars) then
3270 Status = WRF_WARN_VAR_EOF
3273 VarName = DH%VarNames(DH%CurrentVariable)
3276 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3277 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3278 call wrf_debug ( FATAL , msg)
3281 end subroutine ext_ncd_get_next_var
3283 subroutine ext_ncd_end_of_frame(DataHandle, Status)
3285 use ext_ncd_support_routines
3287 include 'netcdf.inc'
3288 include 'wrf_status_codes.h'
3289 integer ,intent(in) :: DataHandle
3290 integer ,intent(out) :: Status
3291 type(wrf_data_handle) ,pointer :: DH
3293 call GetDH(DataHandle,DH,Status)
3295 end subroutine ext_ncd_end_of_frame
3297 ! NOTE: For scalar variables NDim is set to zero and DomainStart and
3298 ! NOTE: DomainEnd are left unmodified.
3299 subroutine ext_ncd_get_var_info(DataHandle,Name,NDim,MemoryOrder,Stagger,DomainStart,DomainEnd,WrfType,Status)
3301 use ext_ncd_support_routines
3303 include 'netcdf.inc'
3304 include 'wrf_status_codes.h'
3305 integer ,intent(in) :: DataHandle
3306 character*(*) ,intent(in) :: Name
3307 integer ,intent(out) :: NDim
3308 character*(*) ,intent(out) :: MemoryOrder
3309 character*(*) :: Stagger ! Dummy for now
3310 integer ,dimension(*) ,intent(out) :: DomainStart, DomainEnd
3311 integer ,intent(out) :: WrfType
3312 integer ,intent(out) :: Status
3313 type(wrf_data_handle) ,pointer :: DH
3315 integer ,dimension(NVarDims) :: VDimIDs
3320 call GetDH(DataHandle,DH,Status)
3321 if(Status /= WRF_NO_ERR) then
3322 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3323 call wrf_debug ( WARN , TRIM(msg))
3326 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3327 Status = WRF_WARN_FILE_NOT_OPENED
3328 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3329 call wrf_debug ( WARN , TRIM(msg))
3331 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3332 Status = WRF_WARN_DRYRUN_READ
3333 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3334 call wrf_debug ( WARN , TRIM(msg))
3336 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3337 Status = WRF_WARN_READ_WONLY_FILE
3338 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3339 call wrf_debug ( WARN , TRIM(msg))
3341 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3342 stat = NF_INQ_VARID(DH%NCID,Name,VarID)
3343 call netcdf_err(stat,Status)
3344 if(Status /= WRF_NO_ERR) then
3345 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3346 call wrf_debug ( WARN , TRIM(msg))
3349 stat = NF_INQ_VARTYPE(DH%NCID,VarID,XType)
3350 call netcdf_err(stat,Status)
3351 if(Status /= WRF_NO_ERR) then
3352 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3353 call wrf_debug ( WARN , TRIM(msg))
3356 stat = NF_GET_ATT_INT(DH%NCID,VarID,'FieldType',WrfType)
3357 call netcdf_err(stat,Status)
3358 if(Status /= WRF_NO_ERR) then
3359 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3360 call wrf_debug ( WARN , TRIM(msg))
3365 Status = WRF_WARN_BAD_DATA_TYPE
3366 write(msg,*) 'Warning BYTE IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3367 call wrf_debug ( WARN , TRIM(msg))
3370 Status = WRF_WARN_BAD_DATA_TYPE
3371 write(msg,*) 'Warning CHAR IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3372 call wrf_debug ( WARN , TRIM(msg))
3375 Status = WRF_WARN_BAD_DATA_TYPE
3376 write(msg,*) 'Warning SHORT IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3377 call wrf_debug ( WARN , TRIM(msg))
3380 if(WrfType /= WRF_INTEGER .and. WrfType /= WRF_LOGICAL) then
3381 Status = WRF_WARN_BAD_DATA_TYPE
3382 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3383 call wrf_debug ( WARN , TRIM(msg))
3387 if(WrfType /= WRF_REAL) then
3388 Status = WRF_WARN_BAD_DATA_TYPE
3389 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3390 call wrf_debug ( WARN , TRIM(msg))
3394 if(WrfType /= WRF_DOUBLE) then
3395 Status = WRF_WARN_BAD_DATA_TYPE
3396 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3397 call wrf_debug ( WARN , TRIM(msg))
3401 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
3402 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
3403 call wrf_debug ( WARN , TRIM(msg))
3407 stat = NF_GET_ATT_TEXT(DH%NCID,VarID,'MemoryOrder',MemoryOrder)
3408 call netcdf_err(stat,Status)
3409 if(Status /= WRF_NO_ERR) then
3410 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3411 call wrf_debug ( WARN , TRIM(msg))
3414 call GetDim(MemoryOrder,NDim,Status)
3415 if(Status /= WRF_NO_ERR) then
3416 write(msg,*) 'Warning BAD MEMORY ORDER ',TRIM(MemoryOrder),' in ',__FILE__,', line', __LINE__
3417 call wrf_debug ( WARN , TRIM(msg))
3420 stat = NF_INQ_VARDIMID(DH%NCID,VarID,VDimIDs)
3421 call netcdf_err(stat,Status)
3422 if(Status /= WRF_NO_ERR) then
3423 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3424 call wrf_debug ( WARN , TRIM(msg))
3429 stat = NF_INQ_DIMLEN(DH%NCID,VDimIDs(j),DomainEnd(j))
3430 call netcdf_err(stat,Status)
3431 if(Status /= WRF_NO_ERR) then
3432 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3433 call wrf_debug ( WARN , TRIM(msg))
3438 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3439 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3440 call wrf_debug ( FATAL , msg)
3443 end subroutine ext_ncd_get_var_info
3445 subroutine ext_ncd_warning_str( Code, ReturnString, Status)
3447 use ext_ncd_support_routines
3449 include 'netcdf.inc'
3450 include 'wrf_status_codes.h'
3452 integer , intent(in) ::Code
3453 character *(*), intent(out) :: ReturnString
3454 integer, intent(out) ::Status
3458 ReturnString='No error'
3462 ReturnString= 'File not found (or file is incomplete)'
3466 ReturnString='Metadata not found'
3470 ReturnString= 'Timestamp not found'
3474 ReturnString= 'No more timestamps'
3478 ReturnString= 'Variable not found'
3482 ReturnString= 'No more variables for the current time'
3486 ReturnString= 'Too many open files'
3490 ReturnString= 'Data type mismatch'
3494 ReturnString= 'Attempt to write read-only file'
3498 ReturnString= 'Attempt to read write-only file'
3502 ReturnString= 'Attempt to access unopened file'
3506 ReturnString= 'Attempt to do 2 trainings for 1 variable'
3510 ReturnString= 'Attempt to read past EOF'
3514 ReturnString= 'Bad data handle'
3518 ReturnString= 'Write length not equal to training length'
3522 ReturnString= 'More dimensions requested than training'
3526 ReturnString= 'Attempt to read more data than exists'
3530 ReturnString= 'Input dimensions inconsistent'
3534 ReturnString= 'Input MemoryOrder not recognized'
3538 ReturnString= 'A dimension name with 2 different lengths'
3542 ReturnString= 'String longer than provided storage'
3546 ReturnString= 'Function not supportable'
3550 ReturnString= 'Package implements this routine as NOOP'
3554 !netcdf-specific warning messages
3556 ReturnString= 'Bad data type'
3560 ReturnString= 'File not committed'
3564 ReturnString= 'File is opened for reading'
3568 ReturnString= 'Attempt to write metadata after open commit'
3572 ReturnString= 'I/O not initialized'
3576 ReturnString= 'Too many variables requested'
3580 ReturnString= 'Attempt to close file during a dry run'
3584 ReturnString= 'Date string not 19 characters in length'
3588 ReturnString= 'Attempt to read zero length words'
3592 ReturnString= 'Data type not found'
3596 ReturnString= 'Badly formatted date string'
3600 ReturnString= 'Attempt at read during a dry run'
3604 ReturnString= 'Attempt to get zero words'
3608 ReturnString= 'Attempt to put zero length words'
3612 ReturnString= 'NetCDF error'
3616 ReturnString= 'Requested length <= 1'
3620 ReturnString= 'More data available than requested'
3624 ReturnString= 'New date less than previous date'
3629 ReturnString= 'This warning code is not supported or handled directly by WRF and NetCDF. &
3630 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
3631 & to be calling a package-specific routine to return a message for this warning code.'
3636 end subroutine ext_ncd_warning_str
3638 !returns message string for all WRF and netCDF warning/error status codes
3639 !Other i/o packages must provide their own routines to return their own status messages
3640 subroutine ext_ncd_error_str( Code, ReturnString, Status)
3642 use ext_ncd_support_routines
3644 include 'netcdf.inc'
3645 include 'wrf_status_codes.h'
3647 integer , intent(in) ::Code
3648 character *(*), intent(out) :: ReturnString
3649 integer, intent(out) ::Status
3653 ReturnString= 'Allocation Error'
3657 ReturnString= 'Deallocation Error'
3661 ReturnString= 'Bad File Status'
3665 ReturnString= 'Variable on disk is not 3D'
3669 ReturnString= 'Metadata on disk is not 1D'
3673 ReturnString= 'Time dimension too small'
3677 ReturnString= 'This error code is not supported or handled directly by WRF and NetCDF. &
3678 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
3679 & to be calling a package-specific routine to return a message for this error code.'
3684 end subroutine ext_ncd_error_str