1 #ifndef USE_NETCDF4_FEATURES
2 #define USE_NETCDF4_FEATURES 1
4 !*------------------------------------------------------------------------------
7 !* Forecast Systems Laboratory
13 !* ADVANCED COMPUTING BRANCH
14 !* SMS/NNT Version: 2.0.0
16 !* This software and its documentation are in the public domain and
17 !* are furnished "as is". The United States government, its
18 !* instrumentalities, officers, employees, and agents make no
19 !* warranty, express or implied, as to the usefulness of the software
20 !* and documentation for any purpose. They assume no
21 !* responsibility (1) for the use of the software and documentation;
22 !* or (2) to provide technical support to users.
24 !* Permission to use, copy, modify, and distribute this software is
25 !* hereby granted, provided that this disclaimer notice appears in
26 !* all copies. All modifications to this software must be clearly
27 !* documented, and are solely the responsibility of the agent making
28 !* the modification. If significant modifications or enhancements
29 !* are made to this software, the SMS Development team
30 !* (sms-info@fsl.noaa.gov) should be notified.
32 !*----------------------------------------------------------------------------
35 ! Author: Jacques Middlecoff jacquesm@fsl.noaa.gov
36 !* Date: October 6, 2000
37 !* NetCDFpar: Adapted by E. Mansell for parallel write/read via netcdf4 interface
39 !*----------------------------------------------------------------------------
43 integer , parameter :: FATAL = 1
44 integer , parameter :: WARN = 1
45 integer , parameter :: WrfDataHandleMax = 99
46 integer , parameter :: MaxDims = 2000 ! = NF_MAX_VARS
48 integer , parameter :: MaxVars = 10000
50 integer , parameter :: MaxVars = 3000
52 integer , parameter :: MaxTimes = 10000
53 integer , parameter :: DateStrLen = 19
54 integer , parameter :: VarNameLen = 31
55 integer , parameter :: NO_DIM = 0
56 integer , parameter :: NVarDims = 4
57 integer , parameter :: NMDVarDims = 2
58 character (8) , parameter :: NO_NAME = 'NULL'
59 character (DateStrLen) , parameter :: ZeroDate = '0000-00-00-00_00_00'
61 #include "wrf_io_flags.h"
63 character (256) :: msg
64 logical :: WrfIOnotInitialized = .true.
66 type :: wrf_data_handle
67 character (255) :: FileName
73 character (5) :: TimesName
75 integer :: CurrentTime !Only used for read
76 integer :: NumberTimes !Only used for read
77 character (DateStrLen), pointer :: Times(:)
79 integer , pointer :: DimLengths(:)
80 integer , pointer :: DimIDs(:)
81 character (31) , pointer :: DimNames(:)
83 character (9) :: DimUnlimName
84 integer , dimension(NVarDims) :: DimID
85 integer , dimension(NVarDims) :: Dimension
86 integer , pointer :: MDVarIDs(:)
87 integer , pointer :: MDVarDimLens(:)
88 character (80) , pointer :: MDVarNames(:)
89 integer , pointer :: VarIDs(:)
90 integer , pointer :: VarDimLens(:,:)
91 character (VarNameLen), pointer :: VarNames(:)
92 integer :: CurrentVariable !Only used for read
94 ! first_operation is set to .TRUE. when a new handle is allocated
95 ! or when open-for-write or open-for-read are committed. It is set
96 ! to .FALSE. when the first field is read or written.
97 logical :: first_operation
100 logical :: use_netcdf_classic
101 logical :: Collective
102 integer :: ind_or_collective
103 end type wrf_data_handle
104 type(wrf_data_handle),target :: WrfDataHandles(WrfDataHandleMax)
105 end module wrf_data_ncpar
107 module ext_ncdpar_support_routines
114 subroutine allocHandle(DataHandle,DH,Comm,Status)
117 include 'wrf_status_codes.h'
118 integer ,intent(out) :: DataHandle
119 type(wrf_data_handle),pointer :: DH
120 integer ,intent(IN) :: Comm
121 integer ,intent(out) :: Status
125 do i=1,WrfDataHandleMax
126 if(WrfDataHandles(i)%Free) then
127 DH => WrfDataHandles(i)
129 allocate(DH%Times(MaxTimes), STAT=stat)
131 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
132 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
133 call wrf_debug ( FATAL , msg)
136 allocate(DH%DimLengths(MaxDims), STAT=stat)
138 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
139 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
140 call wrf_debug ( FATAL , msg)
143 allocate(DH%DimIDs(MaxDims), STAT=stat)
145 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
146 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
147 call wrf_debug ( FATAL , msg)
150 allocate(DH%DimNames(MaxDims), STAT=stat)
152 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
153 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
154 call wrf_debug ( FATAL , msg)
157 allocate(DH%MDVarIDs(MaxVars), STAT=stat)
159 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
160 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
161 call wrf_debug ( FATAL , msg)
164 allocate(DH%MDVarDimLens(MaxVars), STAT=stat)
166 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
167 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
168 call wrf_debug ( FATAL , msg)
171 allocate(DH%MDVarNames(MaxVars), STAT=stat)
173 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
174 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
175 call wrf_debug ( FATAL , msg)
178 allocate(DH%VarIDs(MaxVars), STAT=stat)
180 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
181 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
182 call wrf_debug ( FATAL , msg)
185 allocate(DH%VarDimLens(NVarDims-1,MaxVars), STAT=stat)
187 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
188 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
189 call wrf_debug ( FATAL , msg)
192 allocate(DH%VarNames(MaxVars), STAT=stat)
194 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
195 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
196 call wrf_debug ( FATAL , msg)
201 if(i==WrfDataHandleMax) then
202 Status = WRF_WARN_TOO_MANY_FILES
203 write(msg,*) 'Warning TOO MANY FILES in ',__FILE__,', line', __LINE__
204 call wrf_debug ( WARN , TRIM(msg))
205 write(msg,*) 'Did you call ext_ncdpar_ioinit?'
206 call wrf_debug ( WARN , TRIM(msg))
213 DH%first_operation = .TRUE.
214 DH%R4OnOutput = .false.
216 DH%Collective = .TRUE.
217 DH%ind_or_collective = NF_COLLECTIVE
219 end subroutine allocHandle
221 subroutine deallocHandle(DataHandle, Status)
223 include 'wrf_status_codes.h'
224 integer ,intent(in) :: DataHandle
225 integer ,intent(out) :: Status
226 type(wrf_data_handle),pointer :: DH
230 IF ( DataHandle .GE. 1 .AND. DataHandle .LE. WrfDataHandleMax ) THEN
231 if(.NOT. WrfDataHandles(DataHandle)%Free) then
232 DH => WrfDataHandles(DataHandle)
233 deallocate(DH%Times, STAT=stat)
235 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
236 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
237 call wrf_debug ( FATAL , msg)
240 deallocate(DH%DimLengths, STAT=stat)
242 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
243 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
244 call wrf_debug ( FATAL , msg)
247 deallocate(DH%DimIDs, STAT=stat)
249 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
250 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
251 call wrf_debug ( FATAL , msg)
254 deallocate(DH%DimNames, STAT=stat)
256 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
257 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
258 call wrf_debug ( FATAL , msg)
261 deallocate(DH%MDVarIDs, STAT=stat)
263 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
264 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
265 call wrf_debug ( FATAL , msg)
268 deallocate(DH%MDVarDimLens, STAT=stat)
270 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
271 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
272 call wrf_debug ( FATAL , msg)
275 deallocate(DH%MDVarNames, STAT=stat)
277 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
278 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
279 call wrf_debug ( FATAL , msg)
282 deallocate(DH%VarIDs, STAT=stat)
284 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
285 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
286 call wrf_debug ( FATAL , msg)
289 deallocate(DH%VarDimLens, STAT=stat)
291 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
292 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
293 call wrf_debug ( FATAL , msg)
296 deallocate(DH%VarNames, STAT=stat)
298 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
299 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
300 call wrf_debug ( FATAL , msg)
307 end subroutine deallocHandle
309 subroutine GetDH(DataHandle,DH,Status)
311 include 'wrf_status_codes.h'
312 integer ,intent(in) :: DataHandle
313 type(wrf_data_handle) ,pointer :: DH
314 integer ,intent(out) :: Status
316 if(DataHandle < 1 .or. DataHandle > WrfDataHandleMax) then
317 Status = WRF_WARN_BAD_DATA_HANDLE
320 DH => WrfDataHandles(DataHandle)
322 Status = WRF_WARN_BAD_DATA_HANDLE
329 subroutine DateCheck(Date,Status)
331 include 'wrf_status_codes.h'
332 character*(*) ,intent(in) :: Date
333 integer ,intent(out) :: Status
335 if(len(Date) /= DateStrLen) then
336 Status = WRF_WARN_DATESTR_BAD_LENGTH
341 end subroutine DateCheck
343 subroutine GetName(Element,Var,Name,Status)
345 include 'wrf_status_codes.h'
346 character*(*) ,intent(in) :: Element
347 character*(*) ,intent(in) :: Var
348 character*(*) ,intent(out) :: Name
349 integer ,intent(out) :: Status
350 character (VarNameLen) :: VarName
353 integer, parameter :: upper_to_lower =IACHAR('a')-IACHAR('A')
356 Name = 'MD___'//trim(Element)//VarName
359 if('A'<=c .and. c <='Z') Name(i:i)=achar(iachar(c)+upper_to_lower)
360 if(c=='-'.or.c==':') Name(i:i)='_'
364 end subroutine GetName
366 subroutine GetTimeIndex(IO,DataHandle,DateStr,TimeIndex,Status)
368 include 'wrf_status_codes.h'
370 character (*) ,intent(in) :: IO
371 integer ,intent(in) :: DataHandle
372 character*(*) ,intent(in) :: DateStr
373 integer ,intent(out) :: TimeIndex
374 integer ,intent(out) :: Status
375 type(wrf_data_handle) ,pointer :: DH
381 DH => WrfDataHandles(DataHandle)
382 call DateCheck(DateStr,Status)
383 if(Status /= WRF_NO_ERR) then
384 Status = WRF_WARN_DATESTR_ERROR
385 write(msg,*) 'Warning DATE STRING ERROR in ',__FILE__,', line', __LINE__
386 call wrf_debug ( WARN , TRIM(msg))
389 if(IO == 'write') then
390 TimeIndex = DH%TimeIndex
391 if(TimeIndex <= 0) then
393 elseif(DateStr == DH%Times(TimeIndex)) then
397 TimeIndex = TimeIndex +1
398 if(TimeIndex > MaxTimes) then
399 Status = WRF_WARN_TIME_EOF
400 write(msg,*) 'Warning TIME EOF in ',__FILE__,', line', __LINE__
401 call wrf_debug ( WARN , TRIM(msg))
405 DH%TimeIndex = TimeIndex
406 DH%Times(TimeIndex) = DateStr
408 VStart(2) = TimeIndex
409 VCount(1) = DateStrLen
411 stat = NF_VAR_PAR_ACCESS(DH%NCID,DH%TimesVarID,nf_collective)
412 stat = NF_PUT_VARA_TEXT(DH%NCID,DH%TimesVarID,VStart,VCount,DateStr)
413 call netcdf_err(stat,Status)
414 if(Status /= WRF_NO_ERR) then
415 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
416 call wrf_debug ( WARN , TRIM(msg))
421 if(DH%Times(i)==DateStr) then
427 Status = WRF_WARN_TIME_NF
428 write(msg,*) 'Warning TIME ',DateStr,' NOT FOUND in ',__FILE__,', line', __LINE__
429 call wrf_debug ( WARN , TRIM(msg))
435 end subroutine GetTimeIndex
437 subroutine GetDim(MemoryOrder,NDim,Status)
438 include 'wrf_status_codes.h'
439 character*(*) ,intent(in) :: MemoryOrder
440 integer ,intent(out) :: NDim
441 integer ,intent(out) :: Status
442 character*3 :: MemOrd
444 call LowerCase(MemoryOrder,MemOrd)
446 case ('xyz','xzy','yxz','yzx','zxy','zyx','xsz','xez','ysz','yez')
448 case ('xy','yx','xs','xe','ys','ye','cc')
452 case ('0') ! NDim=0 for scalars. TBH: 20060502
455 Status = WRF_WARN_BAD_MEMORYORDER
460 end subroutine GetDim
462 #ifdef USE_NETCDF4_FEATURES
463 subroutine set_chunking(MemoryOrder,need_chunking)
464 include 'wrf_status_codes.h'
465 character*(*) ,intent(in) :: MemoryOrder
466 integer ,intent(out) :: need_chunking
467 character*3 :: MemOrd
469 call LowerCase(MemoryOrder,MemOrd)
470 if(len(MemOrd) >= 2) then
472 case ('xyz','xzy','yxz','yzx','zxy','zyx','xsz','xez','ysz','yez')
481 end subroutine set_chunking
484 subroutine GetIndices(NDim,Start,End,i1,i2,j1,j2,k1,k2)
485 integer ,intent(in) :: NDim
486 integer ,dimension(*),intent(in) :: Start,End
487 integer ,intent(out) :: i1,i2,j1,j2,k1,k2
495 if(NDim == 0) return ! NDim=0 for scalars. TBH: 20060502
505 end subroutine GetIndices
507 logical function ZeroLengthHorzDim(MemoryOrder,Vector,Status)
509 include 'wrf_status_codes.h'
510 character*(*) ,intent(in) :: MemoryOrder
511 integer,dimension(*) ,intent(in) :: Vector
512 integer ,intent(out) :: Status
514 integer,dimension(NVarDims) :: temp
515 character*3 :: MemOrd
518 call GetDim(MemoryOrder,NDim,Status)
519 temp(1:NDim) = Vector(1:NDim)
520 call LowerCase(MemoryOrder,MemOrd)
521 zero_length = .false.
523 case ('xsz','xez','ysz','yez','xs','xe','ys','ye','z','c')
526 continue ! NDim=0 for scalars. TBH: 20060502
528 zero_length = temp(1) .lt. 1 .or. temp(3) .lt. 1
529 case ('xy','yx','xyz','yxz')
530 zero_length = temp(1) .lt. 1 .or. temp(2) .lt. 1
532 zero_length = temp(2) .lt. 1 .or. temp(3) .lt. 1
534 Status = WRF_WARN_BAD_MEMORYORDER
535 ZeroLengthHorzDim = .true.
539 ZeroLengthHorzDim = zero_length
541 end function ZeroLengthHorzDim
543 subroutine ExtOrder(MemoryOrder,Vector,Status)
545 include 'wrf_status_codes.h'
546 character*(*) ,intent(in) :: MemoryOrder
547 integer,dimension(*) ,intent(inout) :: Vector
548 integer ,intent(out) :: Status
550 integer,dimension(NVarDims) :: temp
551 character*3 :: MemOrd
553 call GetDim(MemoryOrder,NDim,Status)
554 temp(1:NDim) = Vector(1:NDim)
555 call LowerCase(MemoryOrder,MemOrd)
558 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c')
561 continue ! NDim=0 for scalars. TBH: 20060502
583 Status = WRF_WARN_BAD_MEMORYORDER
588 end subroutine ExtOrder
590 subroutine ExtOrderStr(MemoryOrder,Vector,ROVector,Status)
592 include 'wrf_status_codes.h'
593 character*(*) ,intent(in) :: MemoryOrder
594 character*(*),dimension(*) ,intent(in) :: Vector
595 character(80),dimension(NVarDims),intent(out) :: ROVector
596 integer ,intent(out) :: Status
598 character*3 :: MemOrd
600 call GetDim(MemoryOrder,NDim,Status)
601 ROVector(1:NDim) = Vector(1:NDim)
602 call LowerCase(MemoryOrder,MemOrd)
605 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c')
608 continue ! NDim=0 for scalars. TBH: 20060502
610 ROVector(2) = Vector(3)
611 ROVector(3) = Vector(2)
613 ROVector(1) = Vector(2)
614 ROVector(2) = Vector(1)
616 ROVector(1) = Vector(3)
617 ROVector(2) = Vector(1)
618 ROVector(3) = Vector(2)
620 ROVector(1) = Vector(2)
621 ROVector(2) = Vector(3)
622 ROVector(3) = Vector(1)
624 ROVector(1) = Vector(3)
625 ROVector(3) = Vector(1)
627 ROVector(1) = Vector(2)
628 ROVector(2) = Vector(1)
630 Status = WRF_WARN_BAD_MEMORYORDER
635 end subroutine ExtOrderStr
638 subroutine LowerCase(MemoryOrder,MemOrd)
639 character*(*) ,intent(in) :: MemoryOrder
640 character*(*) ,intent(out) :: MemOrd
642 integer ,parameter :: upper_to_lower =IACHAR('a')-IACHAR('A')
647 MemOrd(1:N) = MemoryOrder(1:N)
650 if('A'<=c .and. c <='Z') MemOrd(i:i)=achar(iachar(c)+upper_to_lower)
653 end subroutine LowerCase
655 subroutine UpperCase(MemoryOrder,MemOrd)
656 character*(*) ,intent(in) :: MemoryOrder
657 character*(*) ,intent(out) :: MemOrd
659 integer ,parameter :: lower_to_upper =IACHAR('A')-IACHAR('a')
664 MemOrd(1:N) = MemoryOrder(1:N)
667 if('a'<=c .and. c <='z') MemOrd(i:i)=achar(iachar(c)+lower_to_upper)
670 end subroutine UpperCase
672 subroutine netcdf_err(err,Status)
674 include 'wrf_status_codes.h'
676 integer ,intent(in) :: err
677 integer ,intent(out) :: Status
678 character(len=180) :: errmsg
681 if( err==NF_NOERR )then
684 errmsg = NF_STRERROR(err)
685 write(msg,*) 'NetCDF error: ',errmsg
686 call wrf_debug ( WARN , TRIM(msg))
687 Status = WRF_WARN_NETCDF
690 end subroutine netcdf_err
692 subroutine FieldIO(IO,DataHandle,DateStr,Starts,Length,MemoryOrder &
693 ,FieldType,NCID,VarID,XField,Status)
695 include 'wrf_status_codes.h'
697 character (*) ,intent(in) :: IO
698 integer ,intent(in) :: DataHandle
699 character*(*) ,intent(in) :: DateStr
700 integer,dimension(NVarDims),intent(in) :: Starts
701 integer,dimension(NVarDims),intent(in) :: Length
702 character*(*) ,intent(in) :: MemoryOrder
703 integer ,intent(in) :: FieldType
704 integer ,intent(in) :: NCID
705 integer ,intent(in) :: VarID
706 integer,dimension(*) ,intent(inout) :: XField
707 integer ,intent(out) :: Status
710 integer,dimension(NVarDims) :: VStart
711 integer,dimension(NVarDims) :: VCount
713 call GetTimeIndex(IO,DataHandle,DateStr,TimeIndex,Status)
714 if(Status /= WRF_NO_ERR) then
715 write(msg,*) 'Warning in ',__FILE__,', line', __LINE__
716 call wrf_debug ( WARN , TRIM(msg))
717 write(msg,*) ' Bad time index for DateStr = ',DateStr
718 call wrf_debug ( WARN , TRIM(msg))
721 call GetDim(MemoryOrder,NDim,Status)
724 !jm for parallel netcef VStart(1:NDim) = 1
725 VStart(1:NDim) = Starts(1:NDim)
726 VCount(1:NDim) = Length(1:NDim)
727 VStart(NDim+1) = TimeIndex
730 ! Do not use SELECT statement here as sometimes WRF_REAL=WRF_DOUBLE
731 IF (FieldType == WRF_REAL) THEN
732 call ext_ncdpar_RealFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
733 ELSE IF (FieldType == WRF_DOUBLE) THEN
734 call ext_ncdpar_DoubleFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
735 ELSE IF (FieldType == WRF_INTEGER) THEN
736 call ext_ncdpar_IntFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
737 ELSE IF (FieldType == WRF_LOGICAL) THEN
738 call ext_ncdpar_LogicalFieldIO (IO,NCID,VarID,VStart,VCount,XField,Status)
739 if(Status /= WRF_NO_ERR) return
741 !for wrf_complex, double_complex
742 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
743 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
744 call wrf_debug ( WARN , TRIM(msg))
749 end subroutine FieldIO
751 subroutine Transpose(IO,MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
752 ,XField,x1,x2,y1,y2,z1,z2 &
754 character*(*) ,intent(in) :: IO
755 character*(*) ,intent(in) :: MemoryOrder
756 integer ,intent(in) :: l1,l2,m1,m2,n1,n2
757 integer ,intent(in) :: di
758 integer ,intent(in) :: x1,x2,y1,y2,z1,z2
759 integer ,intent(in) :: i1,i2,j1,j2,k1,k2
760 integer ,intent(inout) :: Field(di,l1:l2,m1:m2,n1:n2)
761 !jm 010827 integer ,intent(inout) :: XField(di,x1:x2,y1:y2,z1:z2)
762 integer ,intent(inout) :: XField(di,(i2-i1+1)*(j2-j1+1)*(k2-k1+1))
763 character*3 :: MemOrd
765 integer ,parameter :: MaxUpperCase=IACHAR('Z')
766 integer :: i,j,k,ix,jx,kx
768 call LowerCase(MemoryOrder,MemOrd)
771 !#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))
772 ! define(`XDEX',($1-``$1''1+1+(``$1''2-``$1''1+1)*(($2-``$2''1)+($3-``$3''1)*(``$2''2-``$2''1+1))))
776 #define DFIELD XField(1:di,XDEX(i,k,j))
777 #include "transpose.code"
778 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c','0')
780 #define DFIELD XField(1:di,XDEX(i,j,k))
781 #include "transpose.code"
784 #define DFIELD XField(1:di,XDEX(j,i,k))
785 #include "transpose.code"
788 #define DFIELD XField(1:di,XDEX(k,i,j))
789 #include "transpose.code"
792 #define DFIELD XField(1:di,XDEX(j,k,i))
793 #include "transpose.code"
796 #define DFIELD XField(1:di,XDEX(k,j,i))
797 #include "transpose.code"
800 #define DFIELD XField(1:di,XDEX(j,i,k))
801 #include "transpose.code"
804 end subroutine Transpose
806 subroutine reorder (MemoryOrder,MemO)
807 character*(*) ,intent(in) :: MemoryOrder
808 character*3 ,intent(out) :: MemO
809 character*3 :: MemOrd
810 integer :: N,i,i1,i2,i3
813 N = len_trim(MemoryOrder)
815 call lowercase(MemoryOrder,MemOrd)
816 ! never invert the boundary codes
817 select case ( MemOrd )
818 case ( 'xsz','xez','ysz','yez' )
826 if(ichar(MemOrd(i:i)) < ichar(MemOrd(i1:i1))) I1 = i
827 if(ichar(MemOrd(i:i)) > ichar(MemOrd(i3:i3))) I3 = i
834 MemO(1:1) = MemoryOrder(i1:i1)
835 MemO(2:2) = MemoryOrder(i2:i2)
836 if(N == 3) MemO(3:3) = MemoryOrder(i3:i3)
837 if(MemOrd(i1:i1) == 's' .or. MemOrd(i1:i1) == 'e') then
838 MemO(1:N-1) = MemO(2:N)
839 MemO(N:N ) = MemoryOrder(i1:i1)
842 end subroutine reorder
844 ! Returns .TRUE. iff it is OK to write time-independent domain metadata to the
845 ! file referenced by DataHandle. If DataHandle is invalid, .FALSE. is
847 LOGICAL FUNCTION ncdpar_ok_to_put_dom_ti( DataHandle )
849 include 'wrf_status_codes.h'
850 INTEGER, INTENT(IN) :: DataHandle
851 CHARACTER*80 :: fname
854 LOGICAL :: dryrun, first_output, retval
855 call ext_ncdpar_inquire_filename( DataHandle, fname, filestate, Status )
856 IF ( Status /= WRF_NO_ERR ) THEN
857 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
859 call wrf_debug ( WARN , TRIM(msg) )
862 dryrun = ( filestate .EQ. WRF_FILE_OPENED_NOT_COMMITTED )
863 first_output = ncdpar_is_first_operation( DataHandle )
864 retval = .NOT. dryrun .AND. first_output
866 ncdpar_ok_to_put_dom_ti = retval
868 END FUNCTION ncdpar_ok_to_put_dom_ti
870 ! Returns .TRUE. iff it is OK to read time-independent domain metadata from the
871 ! file referenced by DataHandle. If DataHandle is invalid, .FALSE. is
873 LOGICAL FUNCTION ncdpar_ok_to_get_dom_ti( DataHandle )
875 include 'wrf_status_codes.h'
876 INTEGER, INTENT(IN) :: DataHandle
877 CHARACTER*80 :: fname
880 LOGICAL :: dryrun, retval
881 call ext_ncdpar_inquire_filename( DataHandle, fname, filestate, Status )
882 IF ( Status /= WRF_NO_ERR ) THEN
883 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
885 call wrf_debug ( WARN , TRIM(msg) )
888 dryrun = ( filestate .EQ. WRF_FILE_OPENED_NOT_COMMITTED )
889 retval = .NOT. dryrun
891 ncdpar_ok_to_get_dom_ti = retval
893 END FUNCTION ncdpar_ok_to_get_dom_ti
895 ! Returns .TRUE. iff nothing has been read from or written to the file
896 ! referenced by DataHandle. If DataHandle is invalid, .FALSE. is returned.
897 LOGICAL FUNCTION ncdpar_is_first_operation( DataHandle )
899 INCLUDE 'wrf_status_codes.h'
900 INTEGER, INTENT(IN) :: DataHandle
901 TYPE(wrf_data_handle) ,POINTER :: DH
904 CALL GetDH( DataHandle, DH, Status )
905 IF ( Status /= WRF_NO_ERR ) THEN
906 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
908 call wrf_debug ( WARN , TRIM(msg) )
911 retval = DH%first_operation
913 ncdpar_is_first_operation = retval
915 END FUNCTION ncdpar_is_first_operation
917 subroutine upgrade_filename(FileName)
920 character*(*), intent(inout) :: FileName
923 do i = 1, len(trim(FileName))
924 if(FileName(i:i) == '-') then
926 else if(FileName(i:i) == ':') then
931 end subroutine upgrade_filename
933 end module ext_ncdpar_support_routines
935 subroutine TransposeToR4a(IO,MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
936 ,XField,x1,x2,y1,y2,z1,z2 &
939 use ext_ncdpar_support_routines
941 character*(*) ,intent(in) :: IO
942 character*(*) ,intent(in) :: MemoryOrder
943 integer ,intent(in) :: l1,l2,m1,m2,n1,n2
944 integer ,intent(in) :: di
945 integer ,intent(in) :: x1,x2,y1,y2,z1,z2
946 integer ,intent(in) :: i1,i2,j1,j2,k1,k2
947 real*8 ,intent(inout) :: Field(di,l1:l2,m1:m2,n1:n2)
948 real*4 ,intent(inout) :: XField(di,(i2-i1+1)*(j2-j1+1)*(k2-k1+1))
949 character*3 :: MemOrd
951 integer ,parameter :: MaxUpperCase=IACHAR('Z')
952 integer :: i,j,k,ix,jx,kx
954 call LowerCase(MemoryOrder,MemOrd)
957 !#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))
958 ! define(`XDEX',($1-``$1''1+1+(``$1''2-``$1''1+1)*(($2-``$2''1)+($3-``$3''1)*(``$2''2-``$2''1+1))))
962 #define DFIELD XField(1:di,XDEX(i,k,j))
963 #include "transpose.code"
964 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c','0')
966 #define DFIELD XField(1:di,XDEX(i,j,k))
967 #include "transpose.code"
970 #define DFIELD XField(1:di,XDEX(j,i,k))
971 #include "transpose.code"
974 #define DFIELD XField(1:di,XDEX(k,i,j))
975 #include "transpose.code"
978 #define DFIELD XField(1:di,XDEX(j,k,i))
979 #include "transpose.code"
982 #define DFIELD XField(1:di,XDEX(k,j,i))
983 #include "transpose.code"
986 #define DFIELD XField(1:di,XDEX(j,i,k))
987 #include "transpose.code"
990 end subroutine TransposeToR4a
992 subroutine ext_ncdpar_open_for_read(DatasetName, Comm1, Comm2, SysDepInfo, DataHandle, Status)
994 use ext_ncdpar_support_routines
996 include 'wrf_status_codes.h'
998 character *(*), INTENT(IN) :: DatasetName
999 integer , INTENT(IN) :: Comm1, Comm2
1000 character *(*), INTENT(IN) :: SysDepInfo
1001 integer , INTENT(OUT) :: DataHandle
1002 integer , INTENT(OUT) :: Status
1003 DataHandle = 0 ! dummy setting to quiet warning message
1004 CALL ext_ncdpar_open_for_read_begin( DatasetName, Comm1, Comm2, SysDepInfo, DataHandle, Status )
1005 IF ( Status .EQ. WRF_NO_ERR ) THEN
1006 CALL ext_ncdpar_open_for_read_commit( DataHandle, Status )
1009 end subroutine ext_ncdpar_open_for_read
1011 !ends training phase; switches internal flag to enable input
1012 !must be paired with call to ext_ncdpar_open_for_read_begin
1013 subroutine ext_ncdpar_open_for_read_commit(DataHandle, Status)
1015 use ext_ncdpar_support_routines
1017 include 'wrf_status_codes.h'
1018 include 'netcdf.inc'
1019 integer, intent(in) :: DataHandle
1020 integer, intent(out) :: Status
1021 type(wrf_data_handle) ,pointer :: DH
1023 if(WrfIOnotInitialized) then
1024 Status = WRF_IO_NOT_INITIALIZED
1025 write(msg,*) 'ext_ncdpar_ioinit was not called ',__FILE__,', line', __LINE__
1026 call wrf_debug ( FATAL , msg)
1029 call GetDH(DataHandle,DH,Status)
1030 if(Status /= WRF_NO_ERR) then
1031 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1032 call wrf_debug ( WARN , TRIM(msg))
1035 DH%FileStatus = WRF_FILE_OPENED_FOR_READ
1036 DH%first_operation = .TRUE.
1039 end subroutine ext_ncdpar_open_for_read_commit
1041 subroutine ext_ncdpar_open_for_read_begin( FileName, Comm, IOComm, SysDepInfo, DataHandle, Status)
1043 use ext_ncdpar_support_routines
1045 include 'wrf_status_codes.h'
1046 include 'netcdf.inc'
1047 character*(*) ,intent(INOUT) :: FileName
1048 integer ,intent(IN) :: Comm
1049 integer ,intent(IN) :: IOComm
1050 character*(*) ,intent(in) :: SysDepInfo
1051 integer ,intent(out) :: DataHandle
1052 integer ,intent(out) :: Status
1053 type(wrf_data_handle) ,pointer :: DH
1056 integer ,allocatable :: Buffer(:)
1058 integer :: StoredDim
1060 integer :: DimIDs(2)
1061 integer :: VStart(2)
1063 integer :: TotalNumVars
1066 character (NF_MAX_NAME) :: Name
1068 #ifdef USE_NETCDF4_FEATURES
1069 integer :: open_mode
1072 !call upgrade_filename(FileName)
1074 if(WrfIOnotInitialized) then
1075 Status = WRF_IO_NOT_INITIALIZED
1076 write(msg,*) 'ext_ncdpar_ioinit was not called ',__FILE__,', line', __LINE__
1077 call wrf_debug ( FATAL , msg)
1080 call allocHandle(DataHandle,DH,Comm,Status)
1081 if(Status /= WRF_NO_ERR) then
1082 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
1083 call wrf_debug ( WARN , TRIM(msg))
1087 ! stat = NF_OPEN_PAR(FileName, NF_NOWRITE, comm, MPI_INFO_NULL, DH%NCID)
1088 stat = NF_OPEN(FileName, NF_NOWRITE, DH%NCID)
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 stat = NF_INQ_VARID(DH%NCID,DH%TimesName,VarID)
1096 call netcdf_err(stat,Status)
1097 if(Status /= WRF_NO_ERR) then
1098 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1099 call wrf_debug ( WARN , TRIM(msg))
1102 stat = NF_INQ_VAR(DH%NCID,VarID,DH%TimesName, XType, StoredDim, DimIDs, NAtts)
1103 call netcdf_err(stat,Status)
1104 if(Status /= WRF_NO_ERR) then
1105 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1106 call wrf_debug ( WARN , TRIM(msg))
1109 if(XType/=NF_CHAR) then
1110 Status = WRF_WARN_TYPE_MISMATCH
1111 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
1112 call wrf_debug ( WARN , TRIM(msg))
1115 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(1),VLen(1))
1116 call netcdf_err(stat,Status)
1117 if(Status /= WRF_NO_ERR) then
1118 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1119 call wrf_debug ( WARN , TRIM(msg))
1122 if(VLen(1) /= DateStrLen) then
1123 Status = WRF_WARN_DATESTR_BAD_LENGTH
1124 write(msg,*) 'Warning DATESTR BAD LENGTH in ',__FILE__,', line', __LINE__
1125 call wrf_debug ( WARN , TRIM(msg))
1128 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(2),VLen(2))
1129 call netcdf_err(stat,Status)
1130 if(Status /= WRF_NO_ERR) then
1131 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1132 call wrf_debug ( WARN , TRIM(msg))
1135 if(VLen(2) > MaxTimes) then
1136 Status = WRF_ERR_FATAL_TOO_MANY_TIMES
1137 write(msg,*) 'Fatal TOO MANY TIME VALUES in ',__FILE__,', line', __LINE__
1138 call wrf_debug ( FATAL , TRIM(msg))
1143 stat = NF_GET_VARA_TEXT(DH%NCID,VarID,VStart,VLen,DH%Times)
1144 call netcdf_err(stat,Status)
1145 if(Status /= WRF_NO_ERR) then
1146 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1147 call wrf_debug ( WARN , TRIM(msg))
1150 stat = NF_INQ_NVARS(DH%NCID,TotalNumVars)
1151 call netcdf_err(stat,Status)
1152 if(Status /= WRF_NO_ERR) then
1153 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1154 call wrf_debug ( WARN , TRIM(msg))
1159 stat = NF_INQ_VARNAME(DH%NCID,i,Name)
1160 call netcdf_err(stat,Status)
1161 if(Status /= WRF_NO_ERR) then
1162 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1163 call wrf_debug ( WARN , TRIM(msg))
1165 elseif(Name(1:5) /= 'md___' .and. Name /= DH%TimesName) then
1167 DH%VarNames(NumVars) = Name
1168 DH%VarIDs(NumVars) = i
1171 DH%NumVars = NumVars
1172 DH%NumberTimes = VLen(2)
1173 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1174 DH%FileName = trim(FileName)
1175 DH%CurrentVariable = 0
1177 DH%TimesVarID = VarID
1180 end subroutine ext_ncdpar_open_for_read_begin
1182 subroutine ext_ncdpar_open_for_update( FileName, Comm, IOComm, SysDepInfo, DataHandle, Status)
1184 use ext_ncdpar_support_routines
1186 include 'wrf_status_codes.h'
1187 include 'netcdf.inc'
1188 character*(*) ,intent(INOUT) :: FileName
1189 integer ,intent(IN) :: Comm
1190 integer ,intent(IN) :: IOComm
1191 character*(*) ,intent(in) :: SysDepInfo
1192 integer ,intent(out) :: DataHandle
1193 integer ,intent(out) :: Status
1194 type(wrf_data_handle) ,pointer :: DH
1197 integer ,allocatable :: Buffer(:)
1199 integer :: StoredDim
1201 integer :: DimIDs(2)
1202 integer :: VStart(2)
1204 integer :: TotalNumVars
1207 character (NF_MAX_NAME) :: Name
1209 !#ifdef USE_NETCDF4_FEATURES
1210 integer :: open_mode
1213 !call upgrade_filename(FileName)
1215 if(WrfIOnotInitialized) then
1216 Status = WRF_IO_NOT_INITIALIZED
1217 write(msg,*) 'ext_ncdpar_ioinit was not called ',__FILE__,', line', __LINE__
1218 call wrf_debug ( FATAL , msg)
1221 call allocHandle(DataHandle,DH,Comm,Status)
1222 if(Status /= WRF_NO_ERR) then
1223 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
1224 call wrf_debug ( WARN , TRIM(msg))
1227 ! stat = NF_OPEN(FileName, NF_WRITE, DH%NCID)
1228 open_mode = ior(NF_MPIIO, NF_WRITE)
1229 stat = NF_OPEN_PAR(FileName, open_mode, comm, MPI_INFO_NULL, DH%NCID)
1230 call netcdf_err(stat,Status)
1231 if(Status /= WRF_NO_ERR) then
1232 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1233 call wrf_debug ( WARN , TRIM(msg))
1236 stat = NF_INQ_VARID(DH%NCID,DH%TimesName,VarID)
1237 call netcdf_err(stat,Status)
1238 if(Status /= WRF_NO_ERR) then
1239 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1240 call wrf_debug ( WARN , TRIM(msg))
1243 stat = NF_INQ_VAR(DH%NCID,VarID,DH%TimesName, XType, StoredDim, DimIDs, NAtts)
1244 call netcdf_err(stat,Status)
1245 if(Status /= WRF_NO_ERR) then
1246 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1247 call wrf_debug ( WARN , TRIM(msg))
1250 if(XType/=NF_CHAR) then
1251 Status = WRF_WARN_TYPE_MISMATCH
1252 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
1253 call wrf_debug ( WARN , TRIM(msg))
1256 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(1),VLen(1))
1257 call netcdf_err(stat,Status)
1258 if(Status /= WRF_NO_ERR) then
1259 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1260 call wrf_debug ( WARN , TRIM(msg))
1263 if(VLen(1) /= DateStrLen) then
1264 Status = WRF_WARN_DATESTR_BAD_LENGTH
1265 write(msg,*) 'Warning DATESTR BAD LENGTH in ',__FILE__,', line', __LINE__
1266 call wrf_debug ( WARN , TRIM(msg))
1269 stat = NF_INQ_DIMLEN(DH%NCID,DimIDs(2),VLen(2))
1270 call netcdf_err(stat,Status)
1271 if(Status /= WRF_NO_ERR) then
1272 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1273 call wrf_debug ( WARN , TRIM(msg))
1276 if(VLen(2) > MaxTimes) then
1277 Status = WRF_ERR_FATAL_TOO_MANY_TIMES
1278 write(msg,*) 'Fatal TOO MANY TIME VALUES in ',__FILE__,', line', __LINE__
1279 call wrf_debug ( FATAL , TRIM(msg))
1284 stat = NF_GET_VARA_TEXT(DH%NCID,VarID,VStart,VLen,DH%Times)
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))
1291 stat = NF_INQ_NVARS(DH%NCID,TotalNumVars)
1292 call netcdf_err(stat,Status)
1293 if(Status /= WRF_NO_ERR) then
1294 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1295 call wrf_debug ( WARN , TRIM(msg))
1300 stat = NF_INQ_VARNAME(DH%NCID,i,Name)
1301 call netcdf_err(stat,Status)
1302 if(Status /= WRF_NO_ERR) then
1303 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1304 call wrf_debug ( WARN , TRIM(msg))
1306 elseif(Name(1:5) /= 'md___' .and. Name /= DH%TimesName) then
1308 DH%VarNames(NumVars) = Name
1309 DH%VarIDs(NumVars) = i
1312 DH%NumVars = NumVars
1313 DH%NumberTimes = VLen(2)
1314 DH%FileStatus = WRF_FILE_OPENED_FOR_UPDATE
1315 DH%FileName = trim(FileName)
1316 DH%CurrentVariable = 0
1318 DH%TimesVarID = VarID
1321 end subroutine ext_ncdpar_open_for_update
1324 SUBROUTINE ext_ncdpar_open_for_write_begin(FileName,Comm,IOComm,SysDepInfo,DataHandle,Status)
1326 use ext_ncdpar_support_routines
1328 include 'wrf_status_codes.h'
1329 include 'netcdf.inc'
1330 character*(*) ,intent(inout) :: FileName
1331 integer ,intent(in) :: Comm
1332 integer ,intent(in) :: IOComm
1333 character*(*) ,intent(in) :: SysDepInfo
1334 integer ,intent(out) :: DataHandle
1335 integer ,intent(out) :: Status
1336 type(wrf_data_handle),pointer :: DH
1339 character (7) :: Buffer
1340 integer :: VDimIDs(2)
1342 #ifdef USE_NETCDF4_FEATURES
1343 integer :: create_mode
1344 integer, parameter :: cache_size = 32, &
1346 cache_preemption = 100
1349 !call upgrade_filename(FileName)
1351 if(WrfIOnotInitialized) then
1352 Status = WRF_IO_NOT_INITIALIZED
1353 write(msg,*) 'ext_ncdpar_open_for_write_begin: ext_ncdpar_ioinit was not called ',__FILE__,', line', __LINE__
1354 call wrf_debug ( FATAL , msg)
1357 call allocHandle(DataHandle,DH,Comm,Status)
1358 if(Status /= WRF_NO_ERR) then
1359 write(msg,*) 'Fatal ALLOCATION ERROR in ext_ncdpar_open_for_write_begin ',__FILE__,', line', __LINE__
1360 call wrf_debug ( FATAL , TRIM(msg))
1365 ! create_mode = ior(NF90_NETCDF4,nf90_iotype)
1366 create_mode = IOR(nf_netcdf4, nf_classic_model)
1367 create_mode = IOR(create_mode, nf_mpiio)
1369 #ifdef USE_NETCDF4_FEATURES
1370 ! create_mode = IOR(nf_netcdf4, nf_classic_model)
1371 if ( DH%use_netcdf_classic ) then
1372 write(msg,*) 'output will be in classic NetCDF format'
1373 call wrf_debug ( WARN , TRIM(msg))
1374 #ifdef WRFIO_ncdpar_NO_LARGE_FILE_SUPPORT
1375 ! stat = NF_CREATE(FileName, NF_CLOBBER, DH%NCID)
1376 create_mode = IOR(create_mode, NF_CLOBBER)
1377 stat = NF_CREATE_PAR(FileName, create_mode, comm, MPI_INFO_NULL, DH%NCID)
1378 ! stat = NF_OPEN_PAR(FileName, NF_NOWRITE, comm, MPI_INFO_NULL, DH%NCID)
1380 ! stat = NF_CREATE(FileName, IOR(NF_CLOBBER,NF_64BIT_OFFSET), DH%NCID)
1381 ! stat = NF_CREATE_PAR(FileName, IOR(NF_CLOBBER,NF_64BIT_OFFSET), comm, MPI_INFO_NULL, DH%NCID)
1382 stat = NF_CREATE_PAR(FileName, create_mode, comm, MPI_INFO_NULL, DH%NCID)
1385 ! create_mode = nf_netcdf4
1386 ! stat = NF_CREATE(FileName, create_mode, DH%NCID)
1387 stat = NF_CREATE_PAR(FileName, create_mode, comm, MPI_INFO_NULL, DH%NCID)
1388 call netcdf_err(stat,Status)
1389 if(Status /= WRF_NO_ERR) then
1390 write(msg,*) 'NetCDF error in NF_CREATE_PAR ext_ncdpar_open_for_write_begin ',__FILE__,', line', __LINE__
1391 call wrf_debug ( WARN , TRIM(msg))
1394 stat = NF_SET_CHUNK_CACHE(cache_size, cache_nelem, cache_preemption)
1397 #ifdef WRFIO_ncdpar_NO_LARGE_FILE_SUPPORT
1398 stat = NF_CREATE(FileName, NF_CLOBBER, DH%NCID)
1400 stat = NF_CREATE(FileName, IOR(NF_CLOBBER,NF_64BIT_OFFSET), DH%NCID)
1403 call netcdf_err(stat,Status)
1404 if(Status /= WRF_NO_ERR) then
1405 write(msg,*) 'NetCDF error in ext_ncdpar_open_for_write_begin ',__FILE__,', line', __LINE__
1406 call wrf_debug ( WARN , TRIM(msg))
1409 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1410 DH%FileName = trim(FileName)
1411 stat = NF_DEF_DIM(DH%NCID,DH%DimUnlimName,NF_UNLIMITED,DH%DimUnlimID)
1412 call netcdf_err(stat,Status)
1413 if(Status /= WRF_NO_ERR) then
1414 write(msg,*) 'NetCDF error in ext_ncdpar_open_for_write_begin ',__FILE__,', line', __LINE__
1415 call wrf_debug ( WARN , TRIM(msg))
1418 DH%VarNames (1:MaxVars) = NO_NAME
1419 DH%MDVarNames(1:MaxVars) = NO_NAME
1421 write(Buffer,FMT="('DIM',i4.4)") i
1422 DH%DimNames (i) = Buffer
1423 DH%DimLengths(i) = NO_DIM
1425 DH%DimNames(1) = 'DateStrLen'
1426 stat = NF_DEF_DIM(DH%NCID,DH%DimNames(1),DateStrLen,DH%DimIDs(1))
1427 call netcdf_err(stat,Status)
1428 if(Status /= WRF_NO_ERR) then
1429 write(msg,*) 'NetCDF error in ext_ncdpar_open_for_write_begin ',__FILE__,', line', __LINE__
1430 call wrf_debug ( WARN , TRIM(msg))
1433 VDimIDs(1) = DH%DimIDs(1)
1434 VDimIDs(2) = DH%DimUnlimID
1435 stat = NF_DEF_VAR(DH%NCID,DH%TimesName,NF_CHAR,2,VDimIDs,DH%TimesVarID)
1436 call netcdf_err(stat,Status)
1437 if(Status /= WRF_NO_ERR) then
1438 write(msg,*) 'NetCDF error in ext_ncdpar_open_for_write_begin ',__FILE__,', line', __LINE__
1439 call wrf_debug ( WARN , TRIM(msg))
1442 DH%DimLengths(1) = DateStrLen
1444 if (index(SysDepInfo,'REAL_OUTPUT_SIZE=4') /= 0) then
1445 DH%R4OnOutput = .true.
1447 !toggle on nofill mode
1448 if (index(SysDepInfo,'NOFILL=.TRUE.') /= 0) then
1453 end subroutine ext_ncdpar_open_for_write_begin
1456 !opens a file for writing or coupler datastream for sending messages.
1457 !no training phase for this version of the open stmt.
1458 subroutine ext_ncdpar_open_for_write (DatasetName, Comm1, Comm2, &
1459 SysDepInfo, DataHandle, Status)
1461 use ext_ncdpar_support_routines
1463 include 'wrf_status_codes.h'
1464 include 'netcdf.inc'
1465 character *(*), intent(in) ::DatasetName
1466 integer , intent(in) ::Comm1, Comm2
1467 character *(*), intent(in) ::SysDepInfo
1468 integer , intent(out) :: DataHandle
1469 integer , intent(out) :: Status
1470 Status=WRF_WARN_NOOP
1471 DataHandle = 0 ! dummy setting to quiet warning message
1473 end subroutine ext_ncdpar_open_for_write
1475 SUBROUTINE ext_ncdpar_open_for_write_commit(DataHandle, Status)
1477 use ext_ncdpar_support_routines
1479 include 'wrf_status_codes.h'
1480 include 'netcdf.inc'
1481 integer ,intent(in) :: DataHandle
1482 integer ,intent(out) :: Status
1483 type(wrf_data_handle),pointer :: DH
1486 integer :: oldmode ! for nf_set_fill, not used
1488 if(WrfIOnotInitialized) then
1489 Status = WRF_IO_NOT_INITIALIZED
1490 write(msg,*) 'ext_ncdpar_open_for_write_commit: ext_ncdpar_ioinit was not called ',__FILE__,', line', __LINE__
1491 call wrf_debug ( FATAL , msg)
1494 call GetDH(DataHandle,DH,Status)
1495 if(Status /= WRF_NO_ERR) then
1496 write(msg,*) 'Warning Status = ',Status,' in ext_ncdpar_open_for_write_commit ',__FILE__,', line', __LINE__
1497 call wrf_debug ( WARN , TRIM(msg))
1500 if ( DH%nofill ) then
1501 Status = NF_SET_FILL(DH%NCID,NF_NOFILL, oldmode )
1502 if(Status /= WRF_NO_ERR) then
1503 write(msg,*) 'Warning Status = ',Status,' from NF_SET_FILL ',__FILE__,', line', __LINE__
1504 call wrf_debug ( WARN , TRIM(msg))
1507 write(msg,*) 'Information: NOFILL being set for writing to ',TRIM(DH%FileName)
1508 call wrf_debug ( WARN , TRIM(msg))
1510 stat = NF_ENDDEF(DH%NCID)
1511 call netcdf_err(stat,Status)
1512 if(Status /= WRF_NO_ERR) then
1513 write(msg,*) 'NetCDF error in ext_ncdpar_open_for_write_commit ',__FILE__,', line', __LINE__
1514 call wrf_debug ( WARN , TRIM(msg))
1517 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1518 DH%first_operation = .TRUE.
1520 end subroutine ext_ncdpar_open_for_write_commit
1522 subroutine ext_ncdpar_ioclose(DataHandle, Status)
1524 use ext_ncdpar_support_routines
1526 include 'wrf_status_codes.h'
1527 include 'netcdf.inc'
1528 integer ,intent(in) :: DataHandle
1529 integer ,intent(out) :: Status
1530 type(wrf_data_handle),pointer :: DH
1533 call GetDH(DataHandle,DH,Status)
1534 if(Status /= WRF_NO_ERR) then
1535 write(msg,*) 'Warning Status = ',Status,' in ext_ncdpar_ioclose ',__FILE__,', line', __LINE__
1536 call wrf_debug ( WARN , TRIM(msg))
1539 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1540 Status = WRF_WARN_FILE_NOT_OPENED
1541 write(msg,*) 'Warning FILE NOT OPENED in ext_ncdpar_ioclose ',__FILE__,', line', __LINE__
1542 call wrf_debug ( WARN , TRIM(msg))
1543 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1544 Status = WRF_WARN_DRYRUN_CLOSE
1545 write(msg,*) 'Warning TRY TO CLOSE DRYRUN in ext_ncdpar_ioclose ',__FILE__,', line', __LINE__
1546 call wrf_debug ( WARN , TRIM(msg))
1547 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1549 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1551 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
1554 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1555 write(msg,*) 'Fatal error BAD FILE STATUS in ext_ncdpar_ioclose ',__FILE__,', line', __LINE__
1556 call wrf_debug ( FATAL , TRIM(msg))
1560 stat = NF_CLOSE(DH%NCID)
1561 call netcdf_err(stat,Status)
1562 if(Status /= WRF_NO_ERR) then
1563 write(msg,*) 'NetCDF error in ext_ncdpar_ioclose ',__FILE__,', line', __LINE__
1564 call wrf_debug ( WARN , TRIM(msg))
1567 CALL deallocHandle( DataHandle, Status )
1570 end subroutine ext_ncdpar_ioclose
1572 subroutine ext_ncdpar_iosync( DataHandle, Status)
1574 use ext_ncdpar_support_routines
1576 include 'wrf_status_codes.h'
1577 include 'netcdf.inc'
1578 integer ,intent(in) :: DataHandle
1579 integer ,intent(out) :: Status
1580 type(wrf_data_handle),pointer :: DH
1583 call GetDH(DataHandle,DH,Status)
1584 if(Status /= WRF_NO_ERR) then
1585 write(msg,*) 'Warning Status = ',Status,' in ext_ncdpar_iosync ',__FILE__,', line', __LINE__
1586 call wrf_debug ( WARN , TRIM(msg))
1589 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1590 Status = WRF_WARN_FILE_NOT_OPENED
1591 write(msg,*) 'Warning FILE NOT OPENED in ext_ncdpar_iosync ',__FILE__,', line', __LINE__
1592 call wrf_debug ( WARN , TRIM(msg))
1593 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1594 Status = WRF_WARN_FILE_NOT_COMMITTED
1595 write(msg,*) 'Warning FILE NOT COMMITTED in ext_ncdpar_iosync ',__FILE__,', line', __LINE__
1596 call wrf_debug ( WARN , TRIM(msg))
1597 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1599 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1602 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1603 write(msg,*) 'Fatal error BAD FILE STATUS in ext_ncdpar_iosync ',__FILE__,', line', __LINE__
1604 call wrf_debug ( FATAL , TRIM(msg))
1607 stat = NF_SYNC(DH%NCID)
1608 call netcdf_err(stat,Status)
1609 if(Status /= WRF_NO_ERR) then
1610 write(msg,*) 'NetCDF error in ext_ncdpar_iosync ',__FILE__,', line', __LINE__
1611 call wrf_debug ( WARN , TRIM(msg))
1615 end subroutine ext_ncdpar_iosync
1619 subroutine ext_ncdpar_redef( DataHandle, Status)
1621 use ext_ncdpar_support_routines
1623 include 'wrf_status_codes.h'
1624 include 'netcdf.inc'
1625 integer ,intent(in) :: DataHandle
1626 integer ,intent(out) :: Status
1627 type(wrf_data_handle),pointer :: DH
1630 call GetDH(DataHandle,DH,Status)
1631 if(Status /= WRF_NO_ERR) then
1632 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1633 call wrf_debug ( WARN , TRIM(msg))
1636 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1637 Status = WRF_WARN_FILE_NOT_OPENED
1638 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1639 call wrf_debug ( WARN , TRIM(msg))
1640 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1641 Status = WRF_WARN_FILE_NOT_COMMITTED
1642 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1643 call wrf_debug ( WARN , TRIM(msg))
1644 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1646 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
1648 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1649 Status = WRF_WARN_FILE_OPEN_FOR_READ
1650 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1651 call wrf_debug ( WARN , TRIM(msg))
1653 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1654 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1655 call wrf_debug ( FATAL , TRIM(msg))
1658 stat = NF_REDEF(DH%NCID)
1659 call netcdf_err(stat,Status)
1660 if(Status /= WRF_NO_ERR) then
1661 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1662 call wrf_debug ( WARN , TRIM(msg))
1665 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1667 end subroutine ext_ncdpar_redef
1669 subroutine ext_ncdpar_enddef( DataHandle, Status)
1671 use ext_ncdpar_support_routines
1673 include 'wrf_status_codes.h'
1674 include 'netcdf.inc'
1675 integer ,intent(in) :: DataHandle
1676 integer ,intent(out) :: Status
1677 type(wrf_data_handle),pointer :: DH
1680 call GetDH(DataHandle,DH,Status)
1681 if(Status /= WRF_NO_ERR) then
1682 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1683 call wrf_debug ( WARN , TRIM(msg))
1686 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1687 Status = WRF_WARN_FILE_NOT_OPENED
1688 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1689 call wrf_debug ( WARN , TRIM(msg))
1690 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1691 Status = WRF_WARN_FILE_NOT_COMMITTED
1692 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1693 call wrf_debug ( WARN , TRIM(msg))
1694 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1696 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1697 Status = WRF_WARN_FILE_OPEN_FOR_READ
1698 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1699 call wrf_debug ( WARN , TRIM(msg))
1701 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1702 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1703 call wrf_debug ( FATAL , TRIM(msg))
1706 stat = NF_ENDDEF(DH%NCID)
1707 call netcdf_err(stat,Status)
1708 if(Status /= WRF_NO_ERR) then
1709 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1710 call wrf_debug ( WARN , TRIM(msg))
1713 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1715 end subroutine ext_ncdpar_enddef
1717 subroutine ext_ncdpar_ioinit(SysDepInfo, Status)
1720 include 'wrf_status_codes.h'
1721 CHARACTER*(*), INTENT(IN) :: SysDepInfo
1722 INTEGER ,INTENT(INOUT) :: Status
1724 WrfIOnotInitialized = .false.
1725 WrfDataHandles(1:WrfDataHandleMax)%Free = .true.
1726 WrfDataHandles(1:WrfDataHandleMax)%TimesName = 'Times'
1727 WrfDataHandles(1:WrfDataHandleMax)%DimUnlimName = 'Time'
1728 WrfDataHandles(1:WrfDataHandleMax)%FileStatus = WRF_FILE_NOT_OPENED
1729 ! if(trim(SysDepInfo) == "use_netcdf_classic" ) then
1730 ! WrfDataHandles(1:WrfDataHandleMax)%use_netcdf_classic = .true.
1732 WrfDataHandles(1:WrfDataHandleMax)%use_netcdf_classic = .false.
1736 end subroutine ext_ncdpar_ioinit
1739 subroutine ext_ncdpar_inquiry (Inquiry, Result, Status)
1742 include 'wrf_status_codes.h'
1743 character *(*), INTENT(IN) :: Inquiry
1744 character *(*), INTENT(OUT) :: Result
1745 integer ,INTENT(INOUT) :: Status
1746 SELECT CASE (Inquiry)
1747 CASE ("RANDOM_WRITE","RANDOM_READ","SEQUENTIAL_WRITE","SEQUENTIAL_READ")
1749 CASE ("OPEN_READ","OPEN_COMMIT_WRITE")
1751 CASE ("OPEN_WRITE","OPEN_COMMIT_READ","PARALLEL_IO")
1753 CASE ("SELF_DESCRIBING","SUPPORT_METADATA","SUPPORT_3D_FIELDS")
1758 Result = 'No Result for that inquiry!'
1762 end subroutine ext_ncdpar_inquiry
1767 subroutine ext_ncdpar_ioexit(Status)
1769 use ext_ncdpar_support_routines
1771 include 'wrf_status_codes.h'
1772 include 'netcdf.inc'
1773 integer , INTENT(INOUT) ::Status
1775 type(wrf_data_handle),pointer :: DH
1778 if(WrfIOnotInitialized) then
1779 Status = WRF_IO_NOT_INITIALIZED
1780 write(msg,*) 'ext_ncdpar_ioinit was not called ',__FILE__,', line', __LINE__
1781 call wrf_debug ( FATAL , msg)
1784 do i=1,WrfDataHandleMax
1785 CALL deallocHandle( i , stat )
1788 end subroutine ext_ncdpar_ioexit
1790 subroutine ext_ncdpar_get_dom_ti_real(DataHandle,Element,Data,Count,OutCount,Status)
1791 #define ROUTINE_TYPE 'REAL'
1792 #define TYPE_DATA real,intent(out) :: Data(*)
1793 #define TYPE_COUNT integer,intent(in) :: Count
1794 #define TYPE_OUTCOUNT integer,intent(out) :: OutCOunt
1795 #define TYPE_BUFFER real,allocatable :: Buffer(:)
1796 #define NF_TYPE NF_FLOAT
1797 #define NF_ROUTINE NF_GET_ATT_REAL
1798 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1799 #include "ext_ncdpar_get_dom_ti.code"
1800 end subroutine ext_ncdpar_get_dom_ti_real
1802 subroutine ext_ncdpar_get_dom_ti_integer(DataHandle,Element,Data,Count,OutCount,Status)
1809 #define ROUTINE_TYPE 'INTEGER'
1810 #define TYPE_DATA integer,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))
1815 #include "ext_ncdpar_get_dom_ti.code"
1816 end subroutine ext_ncdpar_get_dom_ti_integer
1818 subroutine ext_ncdpar_get_dom_ti_double(DataHandle,Element,Data,Count,OutCount,Status)
1825 #define ROUTINE_TYPE 'DOUBLE'
1826 #define TYPE_DATA real*8,intent(out) :: Data(*)
1827 #define TYPE_BUFFER real*8,allocatable :: Buffer(:)
1828 #define NF_TYPE NF_DOUBLE
1829 #define NF_ROUTINE NF_GET_ATT_DOUBLE
1830 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1831 #include "ext_ncdpar_get_dom_ti.code"
1832 end subroutine ext_ncdpar_get_dom_ti_double
1834 subroutine ext_ncdpar_get_dom_ti_logical(DataHandle,Element,Data,Count,OutCount,Status)
1841 #define ROUTINE_TYPE 'LOGICAL'
1842 #define TYPE_DATA logical,intent(out) :: Data(*)
1843 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
1844 #define NF_TYPE NF_INT
1845 #define NF_ROUTINE NF_GET_ATT_INT
1846 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))==1
1847 #include "ext_ncdpar_get_dom_ti.code"
1848 end subroutine ext_ncdpar_get_dom_ti_logical
1850 subroutine ext_ncdpar_get_dom_ti_char(DataHandle,Element,Data,Status)
1854 #undef TYPE_OUTCOUNT
1857 #define ROUTINE_TYPE 'CHAR'
1858 #define TYPE_DATA character*(*),intent(out) :: Data
1860 #define TYPE_OUTCOUNT
1862 #define NF_TYPE NF_CHAR
1864 #include "ext_ncdpar_get_dom_ti.code"
1866 end subroutine ext_ncdpar_get_dom_ti_char
1868 subroutine ext_ncdpar_put_dom_ti_real(DataHandle,Element,Data,Count,Status)
1875 #define ROUTINE_TYPE 'REAL'
1876 #define TYPE_DATA real ,intent(in) :: Data(*)
1877 #define TYPE_COUNT integer,intent(in) :: Count
1878 #define NF_ROUTINE NF_PUT_ATT_REAL
1879 #define ARGS NF_FLOAT,Count,Data
1880 #include "ext_ncdpar_put_dom_ti.code"
1881 end subroutine ext_ncdpar_put_dom_ti_real
1883 subroutine ext_ncdpar_put_dom_ti_integer(DataHandle,Element,Data,Count,Status)
1890 #define ROUTINE_TYPE 'INTEGER'
1891 #define TYPE_DATA integer,intent(in) :: Data(*)
1892 #define TYPE_COUNT integer,intent(in) :: Count
1893 #define NF_ROUTINE NF_PUT_ATT_INT
1894 #define ARGS NF_INT,Count,Data
1895 #include "ext_ncdpar_put_dom_ti.code"
1896 end subroutine ext_ncdpar_put_dom_ti_integer
1898 subroutine ext_ncdpar_put_dom_ti_double(DataHandle,Element,Data,Count,Status)
1905 #define ROUTINE_TYPE 'DOUBLE'
1906 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1907 #define TYPE_COUNT integer,intent(in) :: Count
1908 #define NF_ROUTINE NF_PUT_ATT_DOUBLE
1909 #define ARGS NF_DOUBLE,Count,Data
1910 #include "ext_ncdpar_put_dom_ti.code"
1911 end subroutine ext_ncdpar_put_dom_ti_double
1913 subroutine ext_ncdpar_put_dom_ti_logical(DataHandle,Element,Data,Count,Status)
1919 #define ROUTINE_TYPE 'LOGICAL'
1920 #define TYPE_DATA logical,intent(in) :: Data(*)
1921 #define TYPE_COUNT integer,intent(in) :: Count
1922 #define NF_ROUTINE NF_PUT_ATT_INT
1923 #define ARGS NF_INT,Count,Buffer
1925 #include "ext_ncdpar_put_dom_ti.code"
1926 end subroutine ext_ncdpar_put_dom_ti_logical
1928 subroutine ext_ncdpar_put_dom_ti_char(DataHandle,Element,Data,Status)
1935 #define ROUTINE_TYPE 'CHAR'
1936 #define TYPE_DATA character*(*),intent(in) :: Data
1937 #define TYPE_COUNT integer,parameter :: Count=1
1938 #define NF_ROUTINE NF_PUT_ATT_TEXT
1939 #define ARGS len_trim(Data),Data
1940 #include "ext_ncdpar_put_dom_ti.code"
1941 end subroutine ext_ncdpar_put_dom_ti_char
1943 subroutine ext_ncdpar_put_var_ti_real(DataHandle,Element,Var,Data,Count,Status)
1950 #define ROUTINE_TYPE 'REAL'
1951 #define TYPE_DATA real ,intent(in) :: Data(*)
1952 #define TYPE_COUNT integer ,intent(in) :: Count
1953 #define NF_ROUTINE NF_PUT_ATT_REAL
1954 #define ARGS NF_FLOAT,Count,Data
1955 #include "ext_ncdpar_put_var_ti.code"
1956 end subroutine ext_ncdpar_put_var_ti_real
1958 subroutine ext_ncdpar_put_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,Status)
1967 #define ROUTINE_TYPE 'REAL'
1968 #define TYPE_DATA real ,intent(in) :: Data(*)
1969 #define TYPE_COUNT integer ,intent(in) :: Count
1970 #define NF_ROUTINE NF_PUT_VARA_REAL
1971 #define NF_TYPE NF_FLOAT
1972 #define LENGTH Count
1974 #include "ext_ncdpar_put_var_td.code"
1975 end subroutine ext_ncdpar_put_var_td_real
1977 subroutine ext_ncdpar_put_var_ti_double(DataHandle,Element,Var,Data,Count,Status)
1984 #define ROUTINE_TYPE 'DOUBLE'
1985 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1986 #define TYPE_COUNT integer ,intent(in) :: Count
1987 #define NF_ROUTINE NF_PUT_ATT_DOUBLE
1988 #define ARGS NF_DOUBLE,Count,Data
1989 #include "ext_ncdpar_put_var_ti.code"
1990 end subroutine ext_ncdpar_put_var_ti_double
1992 subroutine ext_ncdpar_put_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,Status)
2001 #define ROUTINE_TYPE 'DOUBLE'
2002 #define TYPE_DATA real*8,intent(in) :: Data(*)
2003 #define TYPE_COUNT integer ,intent(in) :: Count
2004 #define NF_ROUTINE NF_PUT_VARA_DOUBLE
2005 #define NF_TYPE NF_DOUBLE
2006 #define LENGTH Count
2008 #include "ext_ncdpar_put_var_td.code"
2009 end subroutine ext_ncdpar_put_var_td_double
2011 subroutine ext_ncdpar_put_var_ti_integer(DataHandle,Element,Var,Data,Count,Status)
2018 #define ROUTINE_TYPE 'INTEGER'
2019 #define TYPE_DATA integer ,intent(in) :: Data(*)
2020 #define TYPE_COUNT integer ,intent(in) :: Count
2021 #define NF_ROUTINE NF_PUT_ATT_INT
2022 #define ARGS NF_INT,Count,Data
2023 #include "ext_ncdpar_put_var_ti.code"
2024 end subroutine ext_ncdpar_put_var_ti_integer
2026 subroutine ext_ncdpar_put_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,Status)
2035 #define ROUTINE_TYPE 'INTEGER'
2036 #define TYPE_DATA integer ,intent(in) :: Data(*)
2037 #define TYPE_COUNT integer ,intent(in) :: Count
2038 #define NF_ROUTINE NF_PUT_VARA_INT
2039 #define NF_TYPE NF_INT
2040 #define LENGTH Count
2042 #include "ext_ncdpar_put_var_td.code"
2043 end subroutine ext_ncdpar_put_var_td_integer
2045 subroutine ext_ncdpar_put_var_ti_logical(DataHandle,Element,Var,Data,Count,Status)
2051 #define ROUTINE_TYPE 'LOGICAL'
2052 #define TYPE_DATA logical ,intent(in) :: Data(*)
2053 #define TYPE_COUNT integer ,intent(in) :: Count
2054 #define NF_ROUTINE NF_PUT_ATT_INT
2056 #define ARGS NF_INT,Count,Buffer
2057 #include "ext_ncdpar_put_var_ti.code"
2058 end subroutine ext_ncdpar_put_var_ti_logical
2060 subroutine ext_ncdpar_put_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,Status)
2068 #define ROUTINE_TYPE 'LOGICAL'
2069 #define TYPE_DATA logical ,intent(in) :: Data(*)
2070 #define TYPE_COUNT integer ,intent(in) :: Count
2071 #define NF_ROUTINE NF_PUT_VARA_INT
2072 #define NF_TYPE NF_INT
2074 #define LENGTH Count
2076 #include "ext_ncdpar_put_var_td.code"
2077 end subroutine ext_ncdpar_put_var_td_logical
2079 subroutine ext_ncdpar_put_var_ti_char(DataHandle,Element,Var,Data,Status)
2086 #define ROUTINE_TYPE 'CHAR'
2087 #define TYPE_DATA character*(*) ,intent(in) :: Data
2089 #define NF_ROUTINE NF_PUT_ATT_TEXT
2090 #define ARGS len_trim(Data),trim(Data)
2092 #include "ext_ncdpar_put_var_ti.code"
2094 end subroutine ext_ncdpar_put_var_ti_char
2096 subroutine ext_ncdpar_put_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
2105 #define ROUTINE_TYPE 'CHAR'
2106 #define TYPE_DATA character*(*) ,intent(in) :: Data
2108 #define NF_ROUTINE NF_PUT_VARA_TEXT
2109 #define NF_TYPE NF_CHAR
2110 #define LENGTH len(Data)
2111 #include "ext_ncdpar_put_var_td.code"
2112 end subroutine ext_ncdpar_put_var_td_char
2114 subroutine ext_ncdpar_get_var_ti_real(DataHandle,Element,Var,Data,Count,OutCount,Status)
2119 #undef TYPE_OUTCOUNT
2123 #define ROUTINE_TYPE 'REAL'
2124 #define TYPE_DATA real ,intent(out) :: Data(*)
2125 #define TYPE_BUFFER real ,allocatable :: Buffer(:)
2126 #define TYPE_COUNT integer,intent(in) :: Count
2127 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2128 #define NF_TYPE NF_FLOAT
2129 #define NF_ROUTINE NF_GET_ATT_REAL
2130 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2131 #include "ext_ncdpar_get_var_ti.code"
2132 end subroutine ext_ncdpar_get_var_ti_real
2134 subroutine ext_ncdpar_get_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2139 #undef TYPE_OUTCOUNT
2144 #define ROUTINE_TYPE 'REAL'
2145 #define TYPE_DATA real ,intent(out) :: Data(*)
2146 #define TYPE_BUFFER real
2147 #define TYPE_COUNT integer,intent(in) :: Count
2148 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2149 #define NF_TYPE NF_FLOAT
2150 #define NF_ROUTINE NF_GET_VARA_REAL
2151 #define LENGTH min(Count,Len1)
2152 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2153 #include "ext_ncdpar_get_var_td.code"
2154 end subroutine ext_ncdpar_get_var_td_real
2156 subroutine ext_ncdpar_get_var_ti_double(DataHandle,Element,Var,Data,Count,OutCount,Status)
2161 #undef TYPE_OUTCOUNT
2165 #define ROUTINE_TYPE 'DOUBLE'
2166 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2167 #define TYPE_BUFFER real*8 ,allocatable :: Buffer(:)
2168 #define TYPE_COUNT integer,intent(in) :: Count
2169 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2170 #define NF_TYPE NF_DOUBLE
2171 #define NF_ROUTINE NF_GET_ATT_DOUBLE
2172 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2173 #include "ext_ncdpar_get_var_ti.code"
2174 end subroutine ext_ncdpar_get_var_ti_double
2176 subroutine ext_ncdpar_get_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2181 #undef TYPE_OUTCOUNT
2186 #define ROUTINE_TYPE 'DOUBLE'
2187 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2188 #define TYPE_BUFFER real*8
2189 #define TYPE_COUNT integer,intent(in) :: Count
2190 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2191 #define NF_TYPE NF_DOUBLE
2192 #define NF_ROUTINE NF_GET_VARA_DOUBLE
2193 #define LENGTH min(Count,Len1)
2194 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2195 #include "ext_ncdpar_get_var_td.code"
2196 end subroutine ext_ncdpar_get_var_td_double
2198 subroutine ext_ncdpar_get_var_ti_integer(DataHandle,Element,Var,Data,Count,OutCount,Status)
2203 #undef TYPE_OUTCOUNT
2207 #define ROUTINE_TYPE 'INTEGER'
2208 #define TYPE_DATA integer,intent(out) :: Data(*)
2209 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2210 #define TYPE_COUNT integer,intent(in) :: Count
2211 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2212 #define NF_TYPE NF_INT
2213 #define NF_ROUTINE NF_GET_ATT_INT
2214 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2215 #include "ext_ncdpar_get_var_ti.code"
2216 end subroutine ext_ncdpar_get_var_ti_integer
2218 subroutine ext_ncdpar_get_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2223 #undef TYPE_OUTCOUNT
2228 #define ROUTINE_TYPE 'INTEGER'
2229 #define TYPE_DATA integer,intent(out) :: Data(*)
2230 #define TYPE_BUFFER integer
2231 #define TYPE_COUNT integer,intent(in) :: Count
2232 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2233 #define NF_TYPE NF_INT
2234 #define NF_ROUTINE NF_GET_VARA_INT
2235 #define LENGTH min(Count,Len1)
2236 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2237 #include "ext_ncdpar_get_var_td.code"
2238 end subroutine ext_ncdpar_get_var_td_integer
2240 subroutine ext_ncdpar_get_var_ti_logical(DataHandle,Element,Var,Data,Count,OutCount,Status)
2245 #undef TYPE_OUTCOUNT
2249 #define ROUTINE_TYPE 'LOGICAL'
2250 #define TYPE_DATA logical,intent(out) :: Data(*)
2251 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2252 #define TYPE_COUNT integer,intent(in) :: Count
2253 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2254 #define NF_TYPE NF_INT
2255 #define NF_ROUTINE NF_GET_ATT_INT
2256 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))==1
2257 #include "ext_ncdpar_get_var_ti.code"
2258 end subroutine ext_ncdpar_get_var_ti_logical
2260 subroutine ext_ncdpar_get_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2265 #undef TYPE_OUTCOUNT
2270 #define ROUTINE_TYPE 'LOGICAL'
2271 #define TYPE_DATA logical,intent(out) :: Data(*)
2272 #define TYPE_BUFFER integer
2273 #define TYPE_COUNT integer,intent(in) :: Count
2274 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2275 #define NF_TYPE NF_INT
2276 #define NF_ROUTINE NF_GET_VARA_INT
2277 #define LENGTH min(Count,Len1)
2278 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))==1
2279 #include "ext_ncdpar_get_var_td.code"
2280 end subroutine ext_ncdpar_get_var_td_logical
2282 subroutine ext_ncdpar_get_var_ti_char(DataHandle,Element,Var,Data,Status)
2287 #undef TYPE_OUTCOUNT
2291 #define ROUTINE_TYPE 'CHAR'
2292 #define TYPE_DATA character*(*) ,intent(out) :: Data
2294 #define TYPE_COUNT integer :: Count = 1
2295 #define TYPE_OUTCOUNT
2296 #define NF_TYPE NF_CHAR
2297 #define NF_ROUTINE NF_GET_ATT_TEXT
2300 #include "ext_ncdpar_get_var_ti.code"
2302 end subroutine ext_ncdpar_get_var_ti_char
2304 subroutine ext_ncdpar_get_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
2309 #undef TYPE_OUTCOUNT
2313 #define ROUTINE_TYPE 'CHAR'
2314 #define TYPE_DATA character*(*) ,intent(out) :: Data
2315 #define TYPE_BUFFER character (80)
2316 #define TYPE_COUNT integer :: Count = 1
2317 #define TYPE_OUTCOUNT
2318 #define NF_TYPE NF_CHAR
2319 #define NF_ROUTINE NF_GET_VARA_TEXT
2322 #include "ext_ncdpar_get_var_td.code"
2324 end subroutine ext_ncdpar_get_var_td_char
2326 subroutine ext_ncdpar_put_dom_td_real(DataHandle,Element,DateStr,Data,Count,Status)
2327 integer ,intent(in) :: DataHandle
2328 character*(*) ,intent(in) :: Element
2329 character*(*) ,intent(in) :: DateStr
2330 real ,intent(in) :: Data(*)
2331 integer ,intent(in) :: Count
2332 integer ,intent(out) :: Status
2334 call ext_ncdpar_put_var_td_real(DataHandle,Element,DateStr, &
2335 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2337 end subroutine ext_ncdpar_put_dom_td_real
2339 subroutine ext_ncdpar_put_dom_td_integer(DataHandle,Element,DateStr,Data,Count,Status)
2340 integer ,intent(in) :: DataHandle
2341 character*(*) ,intent(in) :: Element
2342 character*(*) ,intent(in) :: DateStr
2343 integer ,intent(in) :: Data(*)
2344 integer ,intent(in) :: Count
2345 integer ,intent(out) :: Status
2347 call ext_ncdpar_put_var_td_integer(DataHandle,Element,DateStr, &
2348 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2350 end subroutine ext_ncdpar_put_dom_td_integer
2352 subroutine ext_ncdpar_put_dom_td_double(DataHandle,Element,DateStr,Data,Count,Status)
2353 integer ,intent(in) :: DataHandle
2354 character*(*) ,intent(in) :: Element
2355 character*(*) ,intent(in) :: DateStr
2356 real*8 ,intent(in) :: Data(*)
2357 integer ,intent(in) :: Count
2358 integer ,intent(out) :: Status
2360 call ext_ncdpar_put_var_td_double(DataHandle,Element,DateStr, &
2361 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2363 end subroutine ext_ncdpar_put_dom_td_double
2365 subroutine ext_ncdpar_put_dom_td_logical(DataHandle,Element,DateStr,Data,Count,Status)
2366 integer ,intent(in) :: DataHandle
2367 character*(*) ,intent(in) :: Element
2368 character*(*) ,intent(in) :: DateStr
2369 logical ,intent(in) :: Data(*)
2370 integer ,intent(in) :: Count
2371 integer ,intent(out) :: Status
2373 call ext_ncdpar_put_var_td_logical(DataHandle,Element,DateStr, &
2374 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2376 end subroutine ext_ncdpar_put_dom_td_logical
2378 subroutine ext_ncdpar_put_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2379 integer ,intent(in) :: DataHandle
2380 character*(*) ,intent(in) :: Element
2381 character*(*) ,intent(in) :: DateStr
2382 character*(*) ,intent(in) :: Data
2383 integer ,intent(out) :: Status
2385 call ext_ncdpar_put_var_td_char(DataHandle,Element,DateStr, &
2386 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2388 end subroutine ext_ncdpar_put_dom_td_char
2390 subroutine ext_ncdpar_get_dom_td_real(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2391 integer ,intent(in) :: DataHandle
2392 character*(*) ,intent(in) :: Element
2393 character*(*) ,intent(in) :: DateStr
2394 real ,intent(out) :: Data(*)
2395 integer ,intent(in) :: Count
2396 integer ,intent(out) :: OutCount
2397 integer ,intent(out) :: Status
2398 call ext_ncdpar_get_var_td_real(DataHandle,Element,DateStr, &
2399 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2401 end subroutine ext_ncdpar_get_dom_td_real
2403 subroutine ext_ncdpar_get_dom_td_integer(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2404 integer ,intent(in) :: DataHandle
2405 character*(*) ,intent(in) :: Element
2406 character*(*) ,intent(in) :: DateStr
2407 integer ,intent(out) :: Data(*)
2408 integer ,intent(in) :: Count
2409 integer ,intent(out) :: OutCount
2410 integer ,intent(out) :: Status
2411 call ext_ncdpar_get_var_td_integer(DataHandle,Element,DateStr, &
2412 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2414 end subroutine ext_ncdpar_get_dom_td_integer
2416 subroutine ext_ncdpar_get_dom_td_double(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2417 integer ,intent(in) :: DataHandle
2418 character*(*) ,intent(in) :: Element
2419 character*(*) ,intent(in) :: DateStr
2420 real*8 ,intent(out) :: Data(*)
2421 integer ,intent(in) :: Count
2422 integer ,intent(out) :: OutCount
2423 integer ,intent(out) :: Status
2424 call ext_ncdpar_get_var_td_double(DataHandle,Element,DateStr, &
2425 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2427 end subroutine ext_ncdpar_get_dom_td_double
2429 subroutine ext_ncdpar_get_dom_td_logical(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2430 integer ,intent(in) :: DataHandle
2431 character*(*) ,intent(in) :: Element
2432 character*(*) ,intent(in) :: DateStr
2433 logical ,intent(out) :: Data(*)
2434 integer ,intent(in) :: Count
2435 integer ,intent(out) :: OutCount
2436 integer ,intent(out) :: Status
2437 call ext_ncdpar_get_var_td_logical(DataHandle,Element,DateStr, &
2438 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2440 end subroutine ext_ncdpar_get_dom_td_logical
2442 subroutine ext_ncdpar_get_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2443 integer ,intent(in) :: DataHandle
2444 character*(*) ,intent(in) :: Element
2445 character*(*) ,intent(in) :: DateStr
2446 character*(*) ,intent(out) :: Data
2447 integer ,intent(out) :: Status
2448 call ext_ncdpar_get_var_td_char(DataHandle,Element,DateStr, &
2449 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2451 end subroutine ext_ncdpar_get_dom_td_char
2453 subroutine ext_ncdpar_write_field(DataHandle,DateStr,Var,Field,FieldTypeIn, &
2454 Comm, IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2455 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2457 use ext_ncdpar_support_routines
2459 include 'wrf_status_codes.h'
2460 include 'netcdf.inc'
2461 integer ,intent(in) :: DataHandle
2462 character*(*) ,intent(in) :: DateStr
2463 character*(*) ,intent(in) :: Var
2464 integer ,intent(inout) :: Field(*)
2465 integer ,intent(in) :: FieldTypeIn
2466 integer ,intent(inout) :: Comm
2467 integer ,intent(inout) :: IOComm
2468 integer ,intent(in) :: DomainDesc
2469 character*(*) ,intent(in) :: MemoryOrdIn
2470 character*(*) ,intent(in) :: Stagger ! Dummy for now
2471 character*(*) ,dimension(*) ,intent(in) :: DimNames
2472 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2473 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2474 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2475 integer ,intent(out) :: Status
2476 integer :: FieldType
2477 character (3) :: MemoryOrder
2478 type(wrf_data_handle) ,pointer :: DH
2481 character (VarNameLen) :: VarName
2482 character (3) :: MemO
2483 character (3) :: UCMemO
2485 integer ,dimension(NVarDims) :: Length_global, Length_native
2486 integer ,dimension(NVarDims) :: Length
2487 integer ,dimension(NVarDims) :: VDimIDs
2488 character(80),dimension(NVarDims) :: RODimNames
2489 integer ,dimension(NVarDims) :: StoredStart
2490 integer ,dimension(:,:,:,:),allocatable :: XField
2494 integer :: i1,i2,j1,j2,k1,k2
2495 integer :: x1,x2,y1,y2,z1,z2
2496 integer :: p1,p2,q1,q2,r1,r2
2497 integer :: l1,l2,m1,m2,n1,n2
2500 character (80) :: NullName
2502 ! Local, possibly adjusted, copies of MemoryStart and MemoryEnd
2503 integer ,dimension(NVarDims) :: lMemoryStart, lMemoryEnd
2505 #ifdef USE_NETCDF4_FEATURES
2506 integer, parameter :: cache_size = 32000000
2507 integer,dimension(NVarDims) :: chunks
2508 integer :: need_chunking
2509 integer :: compression_level
2510 integer :: block_size
2512 integer :: mpi_error_code
2514 MemoryOrder = trim(adjustl(MemoryOrdIn))
2516 call GetDim(MemoryOrder,NDim,Status)
2517 if(Status /= WRF_NO_ERR) then
2518 write(msg,*) 'Warning BAD MEMORY ORDER |',MemoryOrder,'| in ',__FILE__,', line', __LINE__
2519 call wrf_debug ( WARN , TRIM(msg))
2523 call DateCheck(DateStr,Status)
2524 if(Status /= WRF_NO_ERR) then
2525 write(msg,*) 'Warning DATE STRING ERROR |',DateStr,'| in ',__FILE__,', line', __LINE__
2526 call wrf_debug ( WARN , TRIM(msg))
2530 call GetDH(DataHandle,DH,Status)
2531 if(Status /= WRF_NO_ERR) then
2532 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2533 call wrf_debug ( WARN , TRIM(msg))
2538 #ifdef USE_NETCDF4_FEATURES
2539 if ( .not. DH%use_netcdf_classic ) then
2540 call set_chunking(MemoryOrder,need_chunking)
2541 compression_level = 2
2547 if ( DH%R4OnOutput .and. FieldTypeIn == WRF_DOUBLE ) then
2548 FieldType = WRF_REAL
2550 FieldType = FieldTypeIn
2553 write(msg,*)'ext_ncdpar_write_field: called for ',TRIM(Var)
2555 !jm 010827 Length(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
2557 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2559 IF ( ZeroLengthHorzDim(MemoryOrder,Length,Status) ) THEN
2560 write(msg,*)'ext_ncdpar_write_field: zero length dimension in ',TRIM(Var),'. Ignoring'
2561 call wrf_debug ( WARN , TRIM(msg))
2565 Length_native(1:NDim) = Length(1:NDim)
2566 Length_global(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
2567 ! Length_global(1:NDim) = Length(1:NDim)
2569 call ExtOrder(MemoryOrder,Length,Status)
2570 call ExtOrder(MemoryOrder,Length_global,Status)
2572 call ExtOrderStr(MemoryOrder,DimNames,RODimNames,Status)
2574 ! Magic number to identify call from IO server when doing quilting
2575 ! quilting = (MemoryStart(1) == -998899 .AND. MemoryEnd(1) == -998899)
2577 ! lMemoryStart(1:NDim) = 1
2578 ! lMemoryEnd(1:NDim) = Length(1:NDim)
2580 lMemoryStart(1:NDim) = MemoryStart(1:NDim)
2581 lMemoryEnd(1:NDim) = MemoryEnd(1:NDim)
2584 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2585 Status = WRF_WARN_FILE_NOT_OPENED
2586 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2587 call wrf_debug ( WARN , TRIM(msg))
2588 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
2589 Status = WRF_WARN_WRITE_RONLY_FILE
2590 write(msg,*) 'Warning WRITE READ ONLY FILE in ',__FILE__,', line', __LINE__
2591 call wrf_debug ( WARN , TRIM(msg))
2592 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2594 if(DH%VarNames(NVar) == VarName ) then
2595 Status = WRF_WARN_2DRYRUNS_1VARIABLE
2596 write(msg,*) 'Warning 2 DRYRUNS 1 VARIABLE (',TRIM(VarName),') in ',__FILE__,', line', __LINE__
2597 call wrf_debug ( WARN , TRIM(msg))
2599 elseif(DH%VarNames(NVar) == NO_NAME) then
2600 DH%VarNames(NVar) = VarName
2603 elseif(NVar == MaxVars) then
2604 Status = WRF_WARN_TOO_MANY_VARIABLES
2605 write(msg,*) 'Warning TOO MANY VARIABLES in ',__FILE__,', line', __LINE__
2606 call wrf_debug ( WARN , TRIM(msg))
2611 if(RODimNames(j) == NullName .or. RODimNames(j) == '') then
2613 if(DH%DimLengths(i) == Length_global(j)) then
2615 elseif(DH%DimLengths(i) == NO_DIM) then
2616 stat = NF_DEF_DIM(NCID,DH%DimNames(i),Length_global(j),DH%DimIDs(i))
2617 call netcdf_err(stat,Status)
2618 if(Status /= WRF_NO_ERR) then
2619 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2620 call wrf_debug ( WARN , TRIM(msg))
2623 DH%DimLengths(i) = Length_global(j)
2625 elseif(i == MaxDims) then
2626 Status = WRF_WARN_TOO_MANY_DIMS
2627 write(msg,*) 'Warning TOO MANY DIMENSIONS (',i,') in (',TRIM(VarName),') ',__FILE__,', line', __LINE__
2628 call wrf_debug ( WARN , TRIM(msg))
2632 else !look for input name and check if already defined
2635 if (DH%DimNames(i) == RODimNames(j)) then
2636 if (DH%DimLengths(i) == Length_global(j)) then
2640 Status = WRF_WARN_DIMNAME_REDEFINED
2641 write(msg,*) 'Warning DIM ',i,', NAME ',TRIM(DH%DimNames(i)),' REDEFINED by var ', &
2642 TRIM(Var),' ',DH%DimLengths(i),Length_global(j) ,' in ', __FILE__ ,' line', __LINE__
2643 call wrf_debug ( WARN , TRIM(msg))
2650 if (DH%DimLengths(i) == NO_DIM) then
2651 DH%DimNames(i) = RODimNames(j)
2652 stat = NF_DEF_DIM(NCID,DH%DimNames(i),Length_global(j),DH%DimIDs(i))
2653 call netcdf_err(stat,Status)
2654 if(Status /= WRF_NO_ERR) then
2655 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2656 call wrf_debug ( WARN , TRIM(msg))
2659 DH%DimLengths(i) = Length_global(j)
2661 elseif(i == MaxDims) then
2662 Status = WRF_WARN_TOO_MANY_DIMS
2663 write(msg,*) 'Warning TOO MANY DIMENSIONS in ',__FILE__,', line', __LINE__
2664 call wrf_debug ( WARN , TRIM(msg))
2670 VDimIDs(j) = DH%DimIDs(i)
2671 DH%VarDimLens(j,NVar) = Length_global(j)
2673 VDimIDs(NDim+1) = DH%DimUnlimID
2675 ! Do not use SELECT statement here as sometimes WRF_REAL=WRF_DOUBLE
2676 IF (FieldType == WRF_REAL) THEN
2678 ELSE IF (FieldType == WRF_DOUBLE) THEN
2680 ELSE IF (FieldType == WRF_INTEGER) THEN
2682 ELSE IF (FieldType == WRF_LOGICAL) THEN
2685 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
2686 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
2687 call wrf_debug ( WARN , TRIM(msg))
2691 stat = NF_DEF_VAR(NCID,VarName,XType,NDim+1,VDimIDs,VarID)
2692 call netcdf_err(stat,Status)
2693 if(Status /= WRF_NO_ERR) then
2694 write(msg,*) 'ext_ncdpar_write_field: NetCDF error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2695 call wrf_debug ( WARN , TRIM(msg))
2699 #ifdef USE_NETCDF4_FEATURES
2700 if(need_chunking > 0 ) then
2701 chunks(1:NDim) = Length(1:NDim)
2704 IF( need_chunking == 3 ) THEN
2705 IF ( 4*(chunks(3)/4) == chunks(3) ) THEN
2706 chunks(3) = chunks(3)/4
2708 chunks(3) = chunks(3)/4 + 1
2714 chunks(1) = (Length(1) + 1)/2
2715 chunks(2) = (Length(2) + 1)/2
2719 block_size = block_size * chunks(i)
2722 do while (block_size > cache_size)
2723 chunks(1) = (chunks(1) + 1)/2
2724 chunks(2) = (chunks(2) + 1)/2
2728 block_size = block_size * chunks(i)
2734 ! send size from rank 0 so that all use the same value for chunking
2735 CALL MPI_Bcast(chunks, 2, MPI_INTEGER, 0, MPI_COMM_WORLD, mpi_error_code)
2737 ! write(unit=0, fmt='(2x, 3a,i6)') 'file: ', __FILE__, ', line: ', __LINE__
2738 ! write(unit=0, fmt='(2x, 3a)') TRIM(VarName),':'
2739 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length 1 = ', Length(1), ', chunk 1 = ', chunks(1)
2740 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length 2 = ', Length(2), ', chunk 2 = ', chunks(2)
2741 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length NDim+1 = ', Length(NDim+1), ', chunk NDim+1 = ', chunks(NDim+1)
2742 ! write(unit=0, fmt='(10x, a,i6)') 'compression_level = ', compression_level
2745 stat = NF_DEF_VAR_CHUNKING(NCID, VarID, NF_CHUNKED, chunks(1:NDim+1))
2746 call netcdf_err(stat,Status)
2747 if(Status /= WRF_NO_ERR) then
2748 write(msg,*) 'ext_ncdpar_write_field: NetCDF def chunking error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2749 call wrf_debug ( WARN , TRIM(msg))
2754 stat = NF_DEF_VAR_DEFLATE(NCID, VarID, 1, 1, compression_level)
2755 call netcdf_err(stat,Status)
2756 if(Status /= WRF_NO_ERR) then
2757 write(msg,*) 'ext_ncdpar_write_field: NetCDF def compression error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2758 call wrf_debug ( WARN , TRIM(msg))
2764 DH%VarIDs(NVar) = VarID
2765 stat = NF_PUT_ATT_INT(NCID,VarID,'FieldType',NF_INT,1,FieldType)
2766 call netcdf_err(stat,Status)
2767 if(Status /= WRF_NO_ERR) then
2768 write(msg,*) 'ext_ncdpar_write_field: NetCDF error in ',__FILE__,', line', __LINE__
2769 call wrf_debug ( WARN , TRIM(msg))
2772 call reorder(MemoryOrder,MemO)
2773 call uppercase(MemO,UCMemO)
2774 stat = NF_PUT_ATT_TEXT(NCID,VarID,'MemoryOrder',3,UCMemO)
2775 call netcdf_err(stat,Status)
2776 if(Status /= WRF_NO_ERR) then
2777 write(msg,*) 'ext_ncdpar_write_field: NetCDF error in ',__FILE__,', line', __LINE__
2778 call wrf_debug ( WARN , TRIM(msg))
2781 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
2782 do NVar=1,DH%NumVars
2783 if(DH%VarNames(NVar) == VarName) then
2785 elseif(NVar == DH%NumVars) then
2786 Status = WRF_WARN_VAR_NF
2787 write(msg,*) 'Warning VARIABLE NOT FOUND in ',__FILE__,', line', __LINE__
2788 call wrf_debug ( WARN , TRIM(msg))
2792 VarID = DH%VarIDs(NVar)
2794 if(Length_global(j) /= DH%VarDimLens(j,NVar) .AND. DH%FileStatus /= WRF_FILE_OPENED_FOR_UPDATE ) then
2795 Status = WRF_WARN_WRTLEN_NE_DRRUNLEN
2796 write(msg,*) 'Warning LENGTH != DRY RUN LENGTH for |', &
2797 VarName,'| dim ',j,' in ',__FILE__,', line', __LINE__
2798 call wrf_debug ( WARN , TRIM(msg))
2799 write(msg,*) ' LENGTH ',Length_global(j),' DRY RUN LENGTH ',DH%VarDimLens(j,NVar)
2800 call wrf_debug ( WARN , TRIM(msg))
2802 !jm 010825 elseif(DomainStart(j) < MemoryStart(j)) then
2803 elseif(PatchStart(j) < lMemoryStart(j)) then
2804 Status = WRF_WARN_DIMENSION_ERROR
2805 write(msg,*) 'Warning DIMENSION ERROR for |',VarName, &
2806 '| in ',__FILE__,', line', __LINE__
2807 call wrf_debug ( WARN , TRIM(msg))
2812 call GetIndices(NDim,lMemoryStart,lMemoryEnd,l1,l2,m1,m2,n1,n2)
2813 call GetIndices(NDim,StoredStart,Length ,x1,x2,y1,y2,z1,z2)
2814 call GetIndices(NDim,StoredStart,Length_native ,p1,p2,q1,q2,r1,r2)
2815 call GetIndices(NDim,PatchStart, PatchEnd ,i1,i2,j1,j2,k1,k2)
2817 if(FieldType == WRF_DOUBLE) di=2
2818 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
2820 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
2821 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
2822 call wrf_debug ( FATAL , TRIM(msg))
2825 if (DH%R4OnOutput .and. FieldTypeIn == WRF_DOUBLE) then
2826 call TransposeToR4a('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2827 ,XField,x1,x2,y1,y2,z1,z2 &
2828 ,i1,i2,j1,j2,k1,k2 )
2830 call Transpose('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2831 ,XField,x1,x2,y1,y2,z1,z2 &
2832 ,i1,i2,j1,j2,k1,k2 )
2834 StoredStart(1:NDim) = PatchStart(1:NDim)
2835 call ExtOrder(MemoryOrder,StoredStart,Status)
2836 call FieldIO('write',DataHandle,DateStr,StoredStart,Length,MemoryOrder, &
2837 FieldType,NCID,VarID,XField,Status)
2838 if(Status /= WRF_NO_ERR) then
2839 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2840 call wrf_debug ( WARN , TRIM(msg))
2843 deallocate(XField, STAT=stat)
2845 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
2846 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
2847 call wrf_debug ( FATAL , TRIM(msg))
2851 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2852 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2853 call wrf_debug ( FATAL , TRIM(msg))
2855 DH%first_operation = .FALSE.
2857 end subroutine ext_ncdpar_write_field
2860 subroutine ext_ncdpar_write_field_orig(DataHandle,DateStr,Var,Field,FieldTypeIn, &
2861 Comm, IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2862 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2864 use ext_ncdpar_support_routines
2866 include 'wrf_status_codes.h'
2867 include 'netcdf.inc'
2868 integer ,intent(in) :: DataHandle
2869 character*(*) ,intent(in) :: DateStr
2870 character*(*) ,intent(in) :: Var
2871 integer ,intent(inout) :: Field(*)
2872 integer ,intent(in) :: FieldTypeIn
2873 integer ,intent(inout) :: Comm
2874 integer ,intent(inout) :: IOComm
2875 integer ,intent(in) :: DomainDesc
2876 character*(*) ,intent(in) :: MemoryOrdIn
2877 character*(*) ,intent(in) :: Stagger ! Dummy for now
2878 character*(*) ,dimension(*) ,intent(in) :: DimNames
2879 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2880 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2881 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2882 integer ,intent(out) :: Status
2883 integer :: FieldType
2884 character (3) :: MemoryOrder
2885 type(wrf_data_handle) ,pointer :: DH
2888 character (VarNameLen) :: VarName
2889 character (3) :: MemO
2890 character (3) :: UCMemO
2892 integer ,dimension(NVarDims) :: Length
2893 integer ,dimension(NVarDims) :: VDimIDs
2894 character(80),dimension(NVarDims) :: RODimNames
2895 integer ,dimension(NVarDims) :: StoredStart
2896 integer ,dimension(:,:,:,:),allocatable :: XField
2900 integer :: i1,i2,j1,j2,k1,k2
2901 integer :: x1,x2,y1,y2,z1,z2
2902 integer :: l1,l2,m1,m2,n1,n2
2905 character (80) :: NullName
2908 #ifdef USE_NETCDF4_FEATURES
2909 integer, parameter :: cache_size = 32000000
2910 integer,dimension(NVarDims) :: chunks
2911 integer :: need_chunking
2912 integer :: compression_level
2913 integer :: block_size
2916 MemoryOrder = trim(adjustl(MemoryOrdIn))
2918 call GetDim(MemoryOrder,NDim,Status)
2919 if(Status /= WRF_NO_ERR) then
2920 write(msg,*) 'Warning BAD MEMORY ORDER |',MemoryOrder,'| in ',__FILE__,', line', __LINE__
2921 call wrf_debug ( WARN , TRIM(msg))
2925 call DateCheck(DateStr,Status)
2926 if(Status /= WRF_NO_ERR) then
2927 write(msg,*) 'Warning DATE STRING ERROR |',DateStr,'| in ',__FILE__,', line', __LINE__
2928 call wrf_debug ( WARN , TRIM(msg))
2932 call GetDH(DataHandle,DH,Status)
2933 if(Status /= WRF_NO_ERR) then
2934 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2935 call wrf_debug ( WARN , TRIM(msg))
2940 #ifdef USE_NETCDF4_FEATURES
2941 if ( .not. DH%use_netcdf_classic ) then
2942 call set_chunking(MemoryOrder,need_chunking)
2943 compression_level = 2
2949 if ( DH%R4OnOutput .and. FieldTypeIn == WRF_DOUBLE ) then
2950 FieldType = WRF_REAL
2952 FieldType = FieldTypeIn
2955 write(msg,*)'ext_ncdpar_write_field: called for ',TRIM(Var)
2957 !jm 010827 Length(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
2959 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2961 IF ( ZeroLengthHorzDim(MemoryOrder,Length,Status) ) THEN
2962 write(msg,*)'ext_ncdpar_write_field: zero length dimension in ',TRIM(Var),'. Ignoring'
2963 call wrf_debug ( WARN , TRIM(msg))
2967 call ExtOrder(MemoryOrder,Length,Status)
2968 call ExtOrderStr(MemoryOrder,DimNames,RODimNames,Status)
2969 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2970 Status = WRF_WARN_FILE_NOT_OPENED
2971 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2972 call wrf_debug ( WARN , TRIM(msg))
2973 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
2974 Status = WRF_WARN_WRITE_RONLY_FILE
2975 write(msg,*) 'Warning WRITE READ ONLY FILE in ',__FILE__,', line', __LINE__
2976 call wrf_debug ( WARN , TRIM(msg))
2977 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2979 if(DH%VarNames(NVar) == VarName ) then
2980 Status = WRF_WARN_2DRYRUNS_1VARIABLE
2981 write(msg,*) 'Warning 2 DRYRUNS 1 VARIABLE in ',__FILE__,', line', __LINE__
2982 call wrf_debug ( WARN , TRIM(msg))
2984 elseif(DH%VarNames(NVar) == NO_NAME) then
2985 DH%VarNames(NVar) = VarName
2988 elseif(NVar == MaxVars) then
2989 Status = WRF_WARN_TOO_MANY_VARIABLES
2990 write(msg,*) 'Warning TOO MANY VARIABLES in ',__FILE__,', line', __LINE__
2991 call wrf_debug ( WARN , TRIM(msg))
2996 if(RODimNames(j) == NullName .or. RODimNames(j) == '') then
2998 if(DH%DimLengths(i) == Length(j)) then
3000 elseif(DH%DimLengths(i) == NO_DIM) then
3001 stat = NF_DEF_DIM(NCID,DH%DimNames(i),Length(j),DH%DimIDs(i))
3002 call netcdf_err(stat,Status)
3003 if(Status /= WRF_NO_ERR) then
3004 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3005 call wrf_debug ( WARN , TRIM(msg))
3008 DH%DimLengths(i) = Length(j)
3010 elseif(i == MaxDims) then
3011 Status = WRF_WARN_TOO_MANY_DIMS
3012 write(msg,*) 'Warning TOO MANY DIMENSIONS in ',__FILE__,', line', __LINE__
3013 call wrf_debug ( WARN , TRIM(msg))
3017 else !look for input name and check if already defined
3020 if (DH%DimNames(i) == RODimNames(j)) then
3021 if (DH%DimLengths(i) == Length(j)) then
3025 Status = WRF_WARN_DIMNAME_REDEFINED
3026 write(msg,*) 'Warning DIM ',i,', NAME ',TRIM(DH%DimNames(i)),' REDEFINED by var ', &
3027 TRIM(Var),' ',DH%DimLengths(i),Length(j) ,' in ', __FILE__ ,' line', __LINE__
3028 call wrf_debug ( WARN , TRIM(msg))
3035 if (DH%DimLengths(i) == NO_DIM) then
3036 DH%DimNames(i) = RODimNames(j)
3037 stat = NF_DEF_DIM(NCID,DH%DimNames(i),Length(j),DH%DimIDs(i))
3038 call netcdf_err(stat,Status)
3039 if(Status /= WRF_NO_ERR) then
3040 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3041 call wrf_debug ( WARN , TRIM(msg))
3044 DH%DimLengths(i) = Length(j)
3046 elseif(i == MaxDims) then
3047 Status = WRF_WARN_TOO_MANY_DIMS
3048 write(msg,*) 'Warning TOO MANY DIMENSIONS in ',__FILE__,', line', __LINE__
3049 call wrf_debug ( WARN , TRIM(msg))
3055 VDimIDs(j) = DH%DimIDs(i)
3056 DH%VarDimLens(j,NVar) = Length(j)
3058 VDimIDs(NDim+1) = DH%DimUnlimID
3060 ! Do not use SELECT statement here as sometimes WRF_REAL=WRF_DOUBLE
3061 IF (FieldType == WRF_REAL) THEN
3063 ELSE IF (FieldType == WRF_DOUBLE) THEN
3065 ELSE IF (FieldType == WRF_INTEGER) THEN
3067 ELSE IF (FieldType == WRF_LOGICAL) THEN
3070 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
3071 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
3072 call wrf_debug ( WARN , TRIM(msg))
3076 stat = NF_DEF_VAR(NCID,VarName,XType,NDim+1,VDimIDs,VarID)
3077 call netcdf_err(stat,Status)
3078 if(Status /= WRF_NO_ERR) then
3079 write(msg,*) 'ext_ncdpar_write_field: NetCDF error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
3080 call wrf_debug ( WARN , TRIM(msg))
3084 #ifdef USE_NETCDF4_FEATURES
3085 if(need_chunking > 0 ) then
3086 chunks(1:NDim) = Length(1:NDim)
3088 chunks(1) = (Length(1) + 1)/2
3089 chunks(2) = (Length(2) + 1)/2
3093 block_size = block_size * chunks(i)
3096 do while (block_size > cache_size)
3097 chunks(1) = (chunks(1) + 1)/2
3098 chunks(2) = (chunks(2) + 1)/2
3102 block_size = block_size * chunks(i)
3106 ! write(unit=0, fmt='(2x, 3a,i6)') 'file: ', __FILE__, ', line: ', __LINE__
3107 ! write(unit=0, fmt='(2x, 3a)') TRIM(VarName),':'
3108 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length 1 = ', Length(1), ', chunk 1 = ', chunks(1)
3109 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length 2 = ', Length(2), ', chunk 2 = ', chunks(2)
3110 ! write(unit=0, fmt='(10x, 2(a,i6))') 'length NDim+1 = ', Length(NDim+1), ', chunk NDim+1 = ', chunks(NDim+1)
3111 ! write(unit=0, fmt='(10x, a,i6)') 'compression_level = ', compression_level
3113 stat = NF_DEF_VAR_CHUNKING(NCID, VarID, NF_CHUNKED, chunks(1:NDim+1))
3114 call netcdf_err(stat,Status)
3115 if(Status /= WRF_NO_ERR) then
3116 write(msg,*) 'ext_ncdpar_write_field: NetCDF def chunking error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
3117 call wrf_debug ( WARN , TRIM(msg))
3122 stat = NF_DEF_VAR_DEFLATE(NCID, VarID, 1, 1, compression_level)
3123 call netcdf_err(stat,Status)
3124 if(Status /= WRF_NO_ERR) then
3125 write(msg,*) 'ext_ncdpar_write_field: NetCDF def compression error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
3126 call wrf_debug ( WARN , TRIM(msg))
3133 DH%VarIDs(NVar) = VarID
3134 stat = NF_PUT_ATT_INT(NCID,VarID,'FieldType',NF_INT,1,FieldType)
3135 call netcdf_err(stat,Status)
3136 if(Status /= WRF_NO_ERR) then
3137 write(msg,*) 'ext_ncdpar_write_field: NetCDF error in ',__FILE__,', line', __LINE__
3138 call wrf_debug ( WARN , TRIM(msg))
3141 call reorder(MemoryOrder,MemO)
3142 call uppercase(MemO,UCMemO)
3143 stat = NF_PUT_ATT_TEXT(NCID,VarID,'MemoryOrder',3,UCMemO)
3144 call netcdf_err(stat,Status)
3145 if(Status /= WRF_NO_ERR) then
3146 write(msg,*) 'ext_ncdpar_write_field: NetCDF error in ',__FILE__,', line', __LINE__
3147 call wrf_debug ( WARN , TRIM(msg))
3150 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3151 do NVar=1,DH%NumVars
3152 if(DH%VarNames(NVar) == VarName) then
3154 elseif(NVar == DH%NumVars) then
3155 Status = WRF_WARN_VAR_NF
3156 write(msg,*) 'Warning VARIABLE NOT FOUND in ',__FILE__,', line', __LINE__
3157 call wrf_debug ( WARN , TRIM(msg))
3161 VarID = DH%VarIDs(NVar)
3163 if(Length(j) /= DH%VarDimLens(j,NVar) .AND. DH%FileStatus /= WRF_FILE_OPENED_FOR_UPDATE ) then
3164 Status = WRF_WARN_WRTLEN_NE_DRRUNLEN
3165 write(msg,*) 'Warning LENGTH != DRY RUN LENGTH for |', &
3166 VarName,'| dim ',j,' in ',__FILE__,', line', __LINE__
3167 call wrf_debug ( WARN , TRIM(msg))
3168 write(msg,*) ' LENGTH ',Length(j),' DRY RUN LENGTH ',DH%VarDimLens(j,NVar)
3169 call wrf_debug ( WARN , TRIM(msg))
3171 !jm 010825 elseif(DomainStart(j) < MemoryStart(j)) then
3172 elseif(PatchStart(j) < MemoryStart(j)) then
3173 Status = WRF_WARN_DIMENSION_ERROR
3174 write(msg,*) 'Warning DIMENSION ERROR for |',VarName, &
3175 '| in ',__FILE__,', line', __LINE__
3176 call wrf_debug ( WARN , TRIM(msg))
3181 call GetIndices(NDim,MemoryStart,MemoryEnd,l1,l2,m1,m2,n1,n2)
3182 call GetIndices(NDim,StoredStart,Length ,x1,x2,y1,y2,z1,z2)
3183 call GetIndices(NDim,PatchStart, PatchEnd ,i1,i2,j1,j2,k1,k2)
3185 if(FieldType == WRF_DOUBLE) di=2
3186 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
3188 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
3189 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
3190 call wrf_debug ( FATAL , TRIM(msg))
3193 if (DH%R4OnOutput .and. FieldTypeIn == WRF_DOUBLE) then
3194 call TransposeToR4a('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
3195 ,XField,x1,x2,y1,y2,z1,z2 &
3196 ,i1,i2,j1,j2,k1,k2 )
3198 call Transpose('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
3199 ,XField,x1,x2,y1,y2,z1,z2 &
3200 ,i1,i2,j1,j2,k1,k2 )
3202 call FieldIO('write',DataHandle,DateStr,StoredStart,Length,MemoryOrder, &
3203 FieldType,NCID,VarID,XField,Status)
3204 ! call FieldIO('write',DataHandle,DateStr,Length,MemoryOrder, &
3205 ! FieldType,NCID,VarID,XField,Status)
3206 if(Status /= WRF_NO_ERR) then
3207 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3208 call wrf_debug ( WARN , TRIM(msg))
3211 deallocate(XField, STAT=stat)
3213 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
3214 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
3215 call wrf_debug ( FATAL , TRIM(msg))
3219 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3220 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3221 call wrf_debug ( FATAL , TRIM(msg))
3223 DH%first_operation = .FALSE.
3225 end subroutine ext_ncdpar_write_field_orig
3227 subroutine ext_ncdpar_read_field(DataHandle,DateStr,Var,Field,FieldType,Comm, &
3228 IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
3229 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
3231 use ext_ncdpar_support_routines
3233 include 'wrf_status_codes.h'
3234 include 'netcdf.inc'
3235 integer ,intent(in) :: DataHandle
3236 character*(*) ,intent(in) :: DateStr
3237 character*(*) ,intent(in) :: Var
3238 integer ,intent(out) :: Field(*)
3239 integer ,intent(in) :: FieldType
3240 integer ,intent(inout) :: Comm
3241 integer ,intent(inout) :: IOComm
3242 integer ,intent(in) :: DomainDesc
3243 character*(*) ,intent(in) :: MemoryOrdIn
3244 character*(*) ,intent(in) :: Stagger ! Dummy for now
3245 character*(*) , dimension (*) ,intent(in) :: DimNames
3246 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
3247 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
3248 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
3249 integer ,intent(out) :: Status
3250 character (3) :: MemoryOrder
3251 character (NF_MAX_NAME) :: dimname
3252 type(wrf_data_handle) ,pointer :: DH
3255 character (VarNameLen) :: VarName
3257 integer ,dimension(NVarDims) :: VCount
3258 integer ,dimension(NVarDims) :: VStart
3259 integer ,dimension(NVarDims) :: Length
3260 integer ,dimension(NVarDims) :: VDimIDs
3261 integer ,dimension(NVarDims) :: MemS
3262 integer ,dimension(NVarDims) :: MemE
3263 integer ,dimension(NVarDims) :: StoredStart
3264 integer ,dimension(NVarDims) :: StoredLen
3265 integer(KIND=MPI_OFFSET_KIND) ,dimension(NVarDims) :: StoredLen_okind
3266 integer ,dimension(:,:,:,:) ,allocatable :: XField
3269 integer :: i1,i2,j1,j2,k1,k2
3270 integer :: x1,x2,y1,y2,z1,z2
3271 integer :: l1,l2,m1,m2,n1,n2
3272 character (VarNameLen) :: Name
3274 integer :: StoredDim
3281 MemoryOrder = trim(adjustl(MemoryOrdIn))
3282 call GetDim(MemoryOrder,NDim,Status)
3283 if(Status /= WRF_NO_ERR) then
3284 write(msg,*) 'Warning BAD MEMORY ORDER |',TRIM(MemoryOrder),'| for |', &
3285 TRIM(Var),'| in ext_ncdpar_read_field ',__FILE__,', line', __LINE__
3286 call wrf_debug ( WARN , TRIM(msg))
3289 call DateCheck(DateStr,Status)
3290 if(Status /= WRF_NO_ERR) then
3291 write(msg,*) 'Warning DATE STRING ERROR |',TRIM(DateStr),'| for |',TRIM(Var), &
3292 '| in ext_ncdpar_read_field ',__FILE__,', line', __LINE__
3293 call wrf_debug ( WARN , TRIM(msg))
3297 call GetDH(DataHandle,DH,Status)
3298 if(Status /= WRF_NO_ERR) then
3299 write(msg,*) 'Warning Status = ',Status,' in ext_ncdpar_read_field ',__FILE__,', line', __LINE__
3300 call wrf_debug ( WARN , TRIM(msg))
3303 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3304 Status = WRF_WARN_FILE_NOT_OPENED
3305 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3306 call wrf_debug ( WARN , TRIM(msg))
3307 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3308 ! jm it is okay to have a dry run read. means read is called between ofrb and ofrc. Just return.
3309 ! Status = WRF_WARN_DRYRUN_READ
3310 ! write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3311 ! call wrf_debug ( WARN , TRIM(msg))
3314 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3315 Status = WRF_WARN_READ_WONLY_FILE
3316 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3317 call wrf_debug ( WARN , TRIM(msg))
3318 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
3321 !jm Length(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
3322 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
3323 StoredStart(1:NDim) = PatchStart(1:NDim)
3325 call ExtOrder(MemoryOrder,Length,Status)
3326 stat = NF_INQ_VARID(NCID,VarName,VarID)
3327 call netcdf_err(stat,Status)
3328 if(Status /= WRF_NO_ERR) then
3329 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__,' Varname ',Varname
3330 call wrf_debug ( WARN , TRIM(msg))
3333 stat = NF_INQ_VAR(NCID,VarID,Name,XType,StoredDim,VDimIDs,NAtts)
3334 call netcdf_err(stat,Status)
3335 if(Status /= WRF_NO_ERR) then
3336 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3337 call wrf_debug ( WARN , TRIM(msg))
3340 stat = NF_GET_ATT_INT(NCID,VarID,'FieldType',FType)
3341 call netcdf_err(stat,Status)
3342 if(Status /= WRF_NO_ERR) then
3343 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3344 call wrf_debug ( WARN , TRIM(msg))
3347 ! allow coercion between double and single prec real
3348 !jm if(FieldType /= Ftype) then
3349 if( (FieldType == WRF_REAL .OR. FieldType == WRF_DOUBLE) ) then
3350 if ( .NOT. (Ftype == WRF_REAL .OR. Ftype == WRF_DOUBLE )) then
3351 Status = WRF_WARN_TYPE_MISMATCH
3352 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
3353 call wrf_debug ( WARN , TRIM(msg))
3356 else if(FieldType /= Ftype) then
3357 Status = WRF_WARN_TYPE_MISMATCH
3358 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
3359 call wrf_debug ( WARN , TRIM(msg))
3363 ! Do not use SELECT statement here as sometimes WRF_REAL=WRF_DOUBLE
3364 IF (FieldType == WRF_REAL) THEN
3365 ! allow coercion between double and single prec real
3366 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
3367 Status = WRF_WARN_TYPE_MISMATCH
3368 write(msg,*) 'Warning REAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
3370 ELSE IF (FieldType == WRF_DOUBLE) THEN
3371 ! allow coercion between double and single prec real
3372 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
3373 Status = WRF_WARN_TYPE_MISMATCH
3374 write(msg,*) 'Warning DOUBLE TYPE MISMATCH in ',__FILE__,', line', __LINE__
3376 ELSE IF (FieldType == WRF_INTEGER) THEN
3377 if(XType /= NF_INT) then
3378 Status = WRF_WARN_TYPE_MISMATCH
3379 write(msg,*) 'Warning INTEGER TYPE MISMATCH in ',__FILE__,', line', __LINE__
3381 ELSE IF (FieldType == WRF_LOGICAL) THEN
3382 if(XType /= NF_INT) then
3383 Status = WRF_WARN_TYPE_MISMATCH
3384 write(msg,*) 'Warning LOGICAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
3387 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
3388 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
3391 if(Status /= WRF_NO_ERR) then
3392 call wrf_debug ( WARN , TRIM(msg))
3395 ! NDim=0 for scalars. Handle read of old NDim=1 files. TBH: 20060502
3396 IF ( ( NDim == 0 ) .AND. ( StoredDim == 2 ) ) THEN
3397 stat = NF_INQ_DIMNAME(NCID,VDimIDs(1),dimname)
3398 call netcdf_err(stat,Status)
3399 if(Status /= WRF_NO_ERR) then
3400 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3401 call wrf_debug ( WARN , TRIM(msg))
3404 IF ( dimname(1:10) == 'ext_scalar' ) THEN
3409 if(StoredDim /= NDim+1) then
3410 Status = WRF_ERR_FATAL_BAD_VARIABLE_DIM
3411 write(msg,*) 'Fatal error BAD VARIABLE DIMENSION in ext_ncdpar_read_field ',TRIM(Var),TRIM(DateStr)
3412 call wrf_debug ( FATAL , msg)
3413 write(msg,*) ' StoredDim ', StoredDim, ' .NE. NDim+1 ', NDim+1
3414 call wrf_debug ( FATAL , msg)
3418 stat = NF_INQ_DIMLEN(NCID,VDimIDs(j),StoredLen(j))
3419 call netcdf_err(stat,Status)
3420 if(Status /= WRF_NO_ERR) then
3421 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3422 call wrf_debug ( WARN , TRIM(msg))
3425 if(Length(j) > StoredLen(j)) then
3426 Status = WRF_WARN_READ_PAST_EOF
3427 write(msg,*) 'Warning READ PAST EOF in ext_ncdpar_read_field of ',TRIM(Var),Length(j),'>',StoredLen(j)
3428 call wrf_debug ( WARN , TRIM(msg))
3430 elseif(Length(j) <= 0) then
3431 Status = WRF_WARN_ZERO_LENGTH_READ
3432 write(msg,*) 'Warning ZERO LENGTH READ in ',__FILE__,', line', __LINE__
3433 call wrf_debug ( WARN , TRIM(msg))
3435 elseif(DomainStart(j) < MemoryStart(j)) then
3436 Status = WRF_WARN_DIMENSION_ERROR
3437 write(msg,*) 'Warning dim ',j,' DomainStart (',DomainStart(j), &
3438 ') < MemoryStart (',MemoryStart(j),') in ',__FILE__,', line', __LINE__
3439 call wrf_debug ( WARN , TRIM(msg))
3445 call GetIndices(NDim,MemoryStart,MemoryEnd,l1,l2,m1,m2,n1,n2)
3446 call GetIndices(NDim,StoredStart,StoredLen,x1,x2,y1,y2,z1,z2)
3447 !jm call GetIndices(NDim,DomainStart,DomainEnd,i1,i2,j1,j2,k1,k2)
3448 call GetIndices(NDim,PatchStart,PatchEnd,i1,i2,j1,j2,k1,k2)
3449 StoredStart(1:NDim) = PatchStart(1:NDim)
3450 call ExtOrder(MemoryOrder,StoredStart,Status)
3453 if(FieldType == WRF_DOUBLE) di=2
3454 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
3456 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
3457 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
3458 call wrf_debug ( FATAL , msg)
3461 call FieldIO('read',DataHandle,DateStr,StoredStart,Length,MemoryOrder, &
3462 FieldType,NCID,VarID,XField,Status)
3463 if(Status /= WRF_NO_ERR) then
3464 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3465 call wrf_debug ( WARN , TRIM(msg))
3468 call Transpose('read',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
3469 ,XField,x1,x2,y1,y2,z1,z2 &
3470 ,i1,i2,j1,j2,k1,k2 )
3471 deallocate(XField, STAT=stat)
3473 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
3474 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
3475 call wrf_debug ( FATAL , msg)
3479 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3480 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3481 call wrf_debug ( FATAL , msg)
3483 DH%first_operation = .FALSE.
3485 end subroutine ext_ncdpar_read_field
3487 subroutine ext_ncdpar_inquire_opened( DataHandle, FileName , FileStatus, Status )
3489 use ext_ncdpar_support_routines
3491 include 'wrf_status_codes.h'
3492 integer ,intent(in) :: DataHandle
3493 character*(*) ,intent(inout) :: FileName
3494 integer ,intent(out) :: FileStatus
3495 integer ,intent(out) :: Status
3496 type(wrf_data_handle) ,pointer :: DH
3498 !call upgrade_filename(FileName)
3500 call GetDH(DataHandle,DH,Status)
3501 if(Status /= WRF_NO_ERR) then
3502 FileStatus = WRF_FILE_NOT_OPENED
3505 if(trim(FileName) /= trim(DH%FileName)) then
3506 FileStatus = WRF_FILE_NOT_OPENED
3508 FileStatus = DH%FileStatus
3512 end subroutine ext_ncdpar_inquire_opened
3514 subroutine ext_ncdpar_inquire_filename( Datahandle, FileName, FileStatus, Status )
3516 use ext_ncdpar_support_routines
3518 include 'wrf_status_codes.h'
3519 integer ,intent(in) :: DataHandle
3520 character*(*) ,intent(out) :: FileName
3521 integer ,intent(out) :: FileStatus
3522 integer ,intent(out) :: Status
3523 type(wrf_data_handle) ,pointer :: DH
3524 FileStatus = WRF_FILE_NOT_OPENED
3525 call GetDH(DataHandle,DH,Status)
3526 if(Status /= WRF_NO_ERR) then
3527 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3528 call wrf_debug ( WARN , TRIM(msg))
3531 FileName = trim(DH%FileName)
3532 FileStatus = DH%FileStatus
3535 end subroutine ext_ncdpar_inquire_filename
3537 subroutine ext_ncdpar_set_time(DataHandle, DateStr, Status)
3539 use ext_ncdpar_support_routines
3541 include 'wrf_status_codes.h'
3542 integer ,intent(in) :: DataHandle
3543 character*(*) ,intent(in) :: DateStr
3544 integer ,intent(out) :: Status
3545 type(wrf_data_handle) ,pointer :: DH
3548 call DateCheck(DateStr,Status)
3549 if(Status /= WRF_NO_ERR) then
3550 write(msg,*) 'Warning DATE STRING ERROR in ',__FILE__,', line', __LINE__
3551 call wrf_debug ( WARN , TRIM(msg))
3554 call GetDH(DataHandle,DH,Status)
3555 if(Status /= WRF_NO_ERR) then
3556 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3557 call wrf_debug ( WARN , TRIM(msg))
3560 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3561 Status = WRF_WARN_FILE_NOT_OPENED
3562 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3563 call wrf_debug ( WARN , TRIM(msg))
3564 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3565 Status = WRF_WARN_FILE_NOT_COMMITTED
3566 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
3567 call wrf_debug ( WARN , TRIM(msg))
3568 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3569 Status = WRF_WARN_READ_WONLY_FILE
3570 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3571 call wrf_debug ( WARN , TRIM(msg))
3572 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
3574 if(DH%Times(i)==DateStr) then
3578 if(i==MaxTimes) then
3579 Status = WRF_WARN_TIME_NF
3583 DH%CurrentVariable = 0
3586 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3587 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3588 call wrf_debug ( FATAL , msg)
3591 end subroutine ext_ncdpar_set_time
3593 subroutine ext_ncdpar_get_next_time(DataHandle, DateStr, Status)
3595 use ext_ncdpar_support_routines
3597 include 'wrf_status_codes.h'
3598 integer ,intent(in) :: DataHandle
3599 character*(*) ,intent(out) :: DateStr
3600 integer ,intent(out) :: Status
3601 type(wrf_data_handle) ,pointer :: DH
3603 call GetDH(DataHandle,DH,Status)
3604 if(Status /= WRF_NO_ERR) then
3605 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3606 call wrf_debug ( WARN , TRIM(msg))
3609 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3610 Status = WRF_WARN_FILE_NOT_OPENED
3611 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3612 call wrf_debug ( WARN , TRIM(msg))
3613 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3614 Status = WRF_WARN_DRYRUN_READ
3615 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3616 call wrf_debug ( WARN , TRIM(msg))
3617 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3618 Status = WRF_WARN_READ_WONLY_FILE
3619 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3620 call wrf_debug ( WARN , TRIM(msg))
3621 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
3622 if(DH%CurrentTime >= DH%NumberTimes) then
3623 Status = WRF_WARN_TIME_EOF
3626 DH%CurrentTime = DH%CurrentTime +1
3627 DateStr = DH%Times(DH%CurrentTime)
3628 DH%CurrentVariable = 0
3631 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3632 write(msg,*) 'DH%FileStatus ',DH%FileStatus
3633 call wrf_debug ( FATAL , msg)
3634 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3635 call wrf_debug ( FATAL , msg)
3638 end subroutine ext_ncdpar_get_next_time
3640 subroutine ext_ncdpar_get_previous_time(DataHandle, DateStr, Status)
3642 use ext_ncdpar_support_routines
3644 include 'wrf_status_codes.h'
3645 integer ,intent(in) :: DataHandle
3646 character*(*) ,intent(out) :: DateStr
3647 integer ,intent(out) :: Status
3648 type(wrf_data_handle) ,pointer :: DH
3650 call GetDH(DataHandle,DH,Status)
3651 if(Status /= WRF_NO_ERR) then
3652 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3653 call wrf_debug ( WARN , TRIM(msg))
3656 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3657 Status = WRF_WARN_FILE_NOT_OPENED
3658 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3659 call wrf_debug ( WARN , TRIM(msg))
3660 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3661 Status = WRF_WARN_DRYRUN_READ
3662 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3663 call wrf_debug ( WARN , TRIM(msg))
3664 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3665 Status = WRF_WARN_READ_WONLY_FILE
3666 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3667 call wrf_debug ( WARN , TRIM(msg))
3668 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
3669 if(DH%CurrentTime.GT.0) then
3670 DH%CurrentTime = DH%CurrentTime -1
3672 DateStr = DH%Times(DH%CurrentTime)
3673 DH%CurrentVariable = 0
3676 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3677 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3678 call wrf_debug ( FATAL , msg)
3681 end subroutine ext_ncdpar_get_previous_time
3683 subroutine ext_ncdpar_get_next_var(DataHandle, VarName, Status)
3685 use ext_ncdpar_support_routines
3687 include 'wrf_status_codes.h'
3688 include 'netcdf.inc'
3689 integer ,intent(in) :: DataHandle
3690 character*(*) ,intent(out) :: VarName
3691 integer ,intent(out) :: Status
3692 type(wrf_data_handle) ,pointer :: DH
3694 character (80) :: Name
3696 call GetDH(DataHandle,DH,Status)
3697 if(Status /= WRF_NO_ERR) then
3698 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3699 call wrf_debug ( WARN , TRIM(msg))
3702 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3703 Status = WRF_WARN_FILE_NOT_OPENED
3704 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3705 call wrf_debug ( WARN , TRIM(msg))
3706 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3707 Status = WRF_WARN_DRYRUN_READ
3708 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3709 call wrf_debug ( WARN , TRIM(msg))
3710 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3711 Status = WRF_WARN_READ_WONLY_FILE
3712 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3713 call wrf_debug ( WARN , TRIM(msg))
3714 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3716 DH%CurrentVariable = DH%CurrentVariable +1
3717 if(DH%CurrentVariable > DH%NumVars) then
3718 Status = WRF_WARN_VAR_EOF
3721 VarName = DH%VarNames(DH%CurrentVariable)
3724 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3725 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3726 call wrf_debug ( FATAL , msg)
3729 end subroutine ext_ncdpar_get_next_var
3731 subroutine ext_ncdpar_end_of_frame(DataHandle, Status)
3733 use ext_ncdpar_support_routines
3735 include 'netcdf.inc'
3736 include 'wrf_status_codes.h'
3737 integer ,intent(in) :: DataHandle
3738 integer ,intent(out) :: Status
3739 type(wrf_data_handle) ,pointer :: DH
3741 call GetDH(DataHandle,DH,Status)
3743 end subroutine ext_ncdpar_end_of_frame
3745 ! NOTE: For scalar variables NDim is set to zero and DomainStart and
3746 ! NOTE: DomainEnd are left unmodified.
3747 subroutine ext_ncdpar_get_var_info(DataHandle,Name,NDim,MemoryOrder,Stagger,DomainStart,DomainEnd,WrfType,Status)
3749 use ext_ncdpar_support_routines
3751 include 'netcdf.inc'
3752 include 'wrf_status_codes.h'
3753 integer ,intent(in) :: DataHandle
3754 character*(*) ,intent(in) :: Name
3755 integer ,intent(out) :: NDim
3756 character*(*) ,intent(out) :: MemoryOrder
3757 character*(*) :: Stagger ! Dummy for now
3758 integer ,dimension(*) ,intent(out) :: DomainStart, DomainEnd
3759 integer ,intent(out) :: WrfType
3760 integer ,intent(out) :: Status
3761 type(wrf_data_handle) ,pointer :: DH
3763 integer ,dimension(NVarDims) :: VDimIDs
3768 call GetDH(DataHandle,DH,Status)
3769 if(Status /= WRF_NO_ERR) then
3770 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3771 call wrf_debug ( WARN , TRIM(msg))
3774 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3775 Status = WRF_WARN_FILE_NOT_OPENED
3776 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3777 call wrf_debug ( WARN , TRIM(msg))
3779 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3780 Status = WRF_WARN_DRYRUN_READ
3781 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3782 call wrf_debug ( WARN , TRIM(msg))
3784 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3785 Status = WRF_WARN_READ_WONLY_FILE
3786 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3787 call wrf_debug ( WARN , TRIM(msg))
3789 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3790 stat = NF_INQ_VARID(DH%NCID,Name,VarID)
3791 call netcdf_err(stat,Status)
3792 if(Status /= WRF_NO_ERR) then
3793 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3794 call wrf_debug ( WARN , TRIM(msg))
3797 stat = NF_INQ_VARTYPE(DH%NCID,VarID,XType)
3798 call netcdf_err(stat,Status)
3799 if(Status /= WRF_NO_ERR) then
3800 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3801 call wrf_debug ( WARN , TRIM(msg))
3804 stat = NF_GET_ATT_INT(DH%NCID,VarID,'FieldType',WrfType)
3805 call netcdf_err(stat,Status)
3806 if(Status /= WRF_NO_ERR) then
3807 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3808 call wrf_debug ( WARN , TRIM(msg))
3813 Status = WRF_WARN_BAD_DATA_TYPE
3814 write(msg,*) 'Warning BYTE IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3815 call wrf_debug ( WARN , TRIM(msg))
3818 Status = WRF_WARN_BAD_DATA_TYPE
3819 write(msg,*) 'Warning CHAR IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3820 call wrf_debug ( WARN , TRIM(msg))
3823 Status = WRF_WARN_BAD_DATA_TYPE
3824 write(msg,*) 'Warning SHORT IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3825 call wrf_debug ( WARN , TRIM(msg))
3828 if(WrfType /= WRF_INTEGER .and. WrfType /= WRF_LOGICAL) then
3829 Status = WRF_WARN_BAD_DATA_TYPE
3830 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3831 call wrf_debug ( WARN , TRIM(msg))
3835 if(WrfType /= WRF_REAL) then
3836 Status = WRF_WARN_BAD_DATA_TYPE
3837 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3838 call wrf_debug ( WARN , TRIM(msg))
3842 if(WrfType /= WRF_DOUBLE) then
3843 Status = WRF_WARN_BAD_DATA_TYPE
3844 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3845 call wrf_debug ( WARN , TRIM(msg))
3849 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
3850 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
3851 call wrf_debug ( WARN , TRIM(msg))
3855 stat = NF_GET_ATT_TEXT(DH%NCID,VarID,'MemoryOrder',MemoryOrder)
3856 call netcdf_err(stat,Status)
3857 if(Status /= WRF_NO_ERR) then
3858 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3859 call wrf_debug ( WARN , TRIM(msg))
3862 call GetDim(MemoryOrder,NDim,Status)
3863 if(Status /= WRF_NO_ERR) then
3864 write(msg,*) 'Warning BAD MEMORY ORDER ',TRIM(MemoryOrder),' in ',__FILE__,', line', __LINE__
3865 call wrf_debug ( WARN , TRIM(msg))
3868 stat = NF_INQ_VARDIMID(DH%NCID,VarID,VDimIDs)
3869 call netcdf_err(stat,Status)
3870 if(Status /= WRF_NO_ERR) then
3871 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3872 call wrf_debug ( WARN , TRIM(msg))
3877 stat = NF_INQ_DIMLEN(DH%NCID,VDimIDs(j),DomainEnd(j))
3878 call netcdf_err(stat,Status)
3879 if(Status /= WRF_NO_ERR) then
3880 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3881 call wrf_debug ( WARN , TRIM(msg))
3886 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3887 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3888 call wrf_debug ( FATAL , msg)
3891 end subroutine ext_ncdpar_get_var_info
3893 subroutine ext_ncdpar_warning_str( Code, ReturnString, Status)
3895 use ext_ncdpar_support_routines
3897 include 'netcdf.inc'
3898 include 'wrf_status_codes.h'
3900 integer , intent(in) ::Code
3901 character *(*), intent(out) :: ReturnString
3902 integer, intent(out) ::Status
3906 ReturnString='No error'
3910 ReturnString= 'File not found (or file is incomplete)'
3914 ReturnString='Metadata not found'
3918 ReturnString= 'Timestamp not found'
3922 ReturnString= 'No more timestamps'
3926 ReturnString= 'Variable not found'
3930 ReturnString= 'No more variables for the current time'
3934 ReturnString= 'Too many open files'
3938 ReturnString= 'Data type mismatch'
3942 ReturnString= 'Attempt to write read-only file'
3946 ReturnString= 'Attempt to read write-only file'
3950 ReturnString= 'Attempt to access unopened file'
3954 ReturnString= 'Attempt to do 2 trainings for 1 variable'
3958 ReturnString= 'Attempt to read past EOF'
3962 ReturnString= 'Bad data handle'
3966 ReturnString= 'Write length not equal to training length'
3970 ReturnString= 'More dimensions requested than training'
3974 ReturnString= 'Attempt to read more data than exists'
3978 ReturnString= 'Input dimensions inconsistent'
3982 ReturnString= 'Input MemoryOrder not recognized'
3986 ReturnString= 'A dimension name with 2 different lengths'
3990 ReturnString= 'String longer than provided storage'
3994 ReturnString= 'Function not supportable'
3998 ReturnString= 'Package implements this routine as NOOP'
4002 !netcdf-specific warning messages
4004 ReturnString= 'Bad data type'
4008 ReturnString= 'File not committed'
4012 ReturnString= 'File is opened for reading'
4016 ReturnString= 'Attempt to write metadata after open commit'
4020 ReturnString= 'I/O not initialized'
4024 ReturnString= 'Too many variables requested'
4028 ReturnString= 'Attempt to close file during a dry run'
4032 ReturnString= 'Date string not 19 characters in length'
4036 ReturnString= 'Attempt to read zero length words'
4040 ReturnString= 'Data type not found'
4044 ReturnString= 'Badly formatted date string'
4048 ReturnString= 'Attempt at read during a dry run'
4052 ReturnString= 'Attempt to get zero words'
4056 ReturnString= 'Attempt to put zero length words'
4060 ReturnString= 'NetCDF error'
4064 ReturnString= 'Requested length <= 1'
4068 ReturnString= 'More data available than requested'
4072 ReturnString= 'New date less than previous date'
4077 ReturnString= 'This warning code is not supported or handled directly by WRF and NetCDF. &
4078 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
4079 & to be calling a package-specific routine to return a message for this warning code.'
4084 end subroutine ext_ncdpar_warning_str
4086 !returns message string for all WRF and netCDF warning/error status codes
4087 !Other i/o packages must provide their own routines to return their own status messages
4088 subroutine ext_ncdpar_error_str( Code, ReturnString, Status)
4090 use ext_ncdpar_support_routines
4092 include 'netcdf.inc'
4093 include 'wrf_status_codes.h'
4095 integer , intent(in) ::Code
4096 character *(*), intent(out) :: ReturnString
4097 integer, intent(out) ::Status
4101 ReturnString= 'Allocation Error'
4105 ReturnString= 'Deallocation Error'
4109 ReturnString= 'Bad File Status'
4113 ReturnString= 'Variable on disk is not 3D'
4117 ReturnString= 'Metadata on disk is not 1D'
4121 ReturnString= 'Time dimension too small'
4125 ReturnString= 'This error code is not supported or handled directly by WRF and NetCDF. &
4126 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
4127 & to be calling a package-specific routine to return a message for this error code.'
4132 end subroutine ext_ncdpar_error_str