1 !/*------------------------------------------------------------------------------
4 !* Forecast Systems Laboratory
10 !* ADVANCED COMPUTING BRANCH
11 !* SMS/NNT Version: 2.0.0
13 !* This software and its documentation are in the public domain and
14 !* are furnished "as is". The United States government, its
15 !* instrumentalities, officers, employees, and agents make no
16 !* warranty, express or implied, as to the usefulness of the software
17 !* and documentation for any purpose. They assume no
18 !* responsibility (1) for the use of the software and documentation;
19 !* or (2) to provide technical support to users.
21 !* Permission to use, copy, modify, and distribute this software is
22 !* hereby granted, provided that this disclaimer notice appears in
23 !* all copies. All modifications to this software must be clearly
24 !* documented, and are solely the responsibility of the agent making
25 !* the modification. If significant modifications or enhancements
26 !* are made to this software, the SMS Development team
27 !* (sms-info@fsl.noaa.gov) should be notified.
29 !*----------------------------------------------------------------------------
32 ! Author: Jacques Middlecoff jacquesm@fsl.noaa.gov
33 !* Date: October 6, 2000
35 !*----------------------------------------------------------------------------
40 integer , parameter :: FATAL = 0
41 integer , parameter :: WARN = 0
42 integer , parameter :: WrfDataHandleMax = 99
43 integer , parameter :: MaxDims = 2000 ! = NF_MAX_VARS
44 integer , parameter :: MaxVars = 3000
45 integer , parameter :: MaxTimes = 60000
46 integer , parameter :: DateStrLen = 19
47 integer , parameter :: VarNameLen = 31
48 integer , parameter :: NO_DIM = 0
49 integer , parameter :: NVarDims = 4
50 integer , parameter :: NMDVarDims = 2
51 character (8) , parameter :: NO_NAME = 'NULL'
52 character (DateStrLen) , parameter :: ZeroDate = '0000-00-00-00:00:00'
54 #include "wrf_io_flags.h"
56 character (256) :: msg
57 logical :: WrfIOnotInitialized = .true.
59 type :: wrf_data_handle
60 character (255) :: FileName
66 character (5) :: TimesName
68 integer :: CurrentTime !Only used for read
69 integer :: NumberTimes !Only used for read
70 character (DateStrLen), pointer :: Times(:)
72 integer , pointer :: DimLengths(:)
73 integer , pointer :: DimIDs(:)
74 character (31) , pointer :: DimNames(:)
76 character (9) :: DimUnlimName
77 integer , dimension(NVarDims) :: DimID
78 integer , dimension(NVarDims) :: Dimension
79 integer , pointer :: MDVarIDs(:)
80 integer , pointer :: MDVarDimLens(:)
81 character (80) , pointer :: MDVarNames(:)
82 integer , pointer :: VarIDs(:)
83 integer , pointer :: VarDimLens(:,:)
84 character (VarNameLen), pointer :: VarNames(:)
85 integer :: CurrentVariable !Only used for read
87 ! first_operation is set to .TRUE. when a new handle is allocated
88 ! or when open-for-write or open-for-read are committed. It is set
89 ! to .FALSE. when the first field is read or written.
90 logical :: first_operation
91 ! Whether pnetcdf file is in collective (.true.) or independent mode
92 ! Collective mode is the default.
94 end type wrf_data_handle
95 type(wrf_data_handle),target :: WrfDataHandles(WrfDataHandleMax)
96 end module wrf_data_pnc
98 module ext_pnc_support_routines
105 integer(KIND=MPI_OFFSET_KIND) function i2offset(i)
109 end function i2offset
111 subroutine allocHandle(DataHandle,DH,Comm,Status)
113 include 'wrf_status_codes.h'
114 integer ,intent(out) :: DataHandle
115 type(wrf_data_handle),pointer :: DH
116 integer ,intent(IN) :: Comm
117 integer ,intent(out) :: Status
121 do i=1,WrfDataHandleMax
122 if(WrfDataHandles(i)%Free) then
123 DH => WrfDataHandles(i)
125 allocate(DH%Times(MaxTimes), STAT=stat)
127 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
128 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
129 call wrf_debug ( FATAL , msg)
132 allocate(DH%DimLengths(MaxDims), STAT=stat)
134 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
135 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
136 call wrf_debug ( FATAL , msg)
139 allocate(DH%DimIDs(MaxDims), STAT=stat)
141 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
142 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
143 call wrf_debug ( FATAL , msg)
146 allocate(DH%DimNames(MaxDims), STAT=stat)
148 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
149 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
150 call wrf_debug ( FATAL , msg)
153 allocate(DH%MDVarIDs(MaxVars), STAT=stat)
155 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
156 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
157 call wrf_debug ( FATAL , msg)
160 allocate(DH%MDVarDimLens(MaxVars), STAT=stat)
162 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
163 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
164 call wrf_debug ( FATAL , msg)
167 allocate(DH%MDVarNames(MaxVars), STAT=stat)
169 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
170 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
171 call wrf_debug ( FATAL , msg)
174 allocate(DH%VarIDs(MaxVars), STAT=stat)
176 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
177 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
178 call wrf_debug ( FATAL , msg)
181 allocate(DH%VarDimLens(NVarDims-1,MaxVars), STAT=stat)
183 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
184 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
185 call wrf_debug ( FATAL , msg)
188 allocate(DH%VarNames(MaxVars), STAT=stat)
190 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
191 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
192 call wrf_debug ( FATAL , msg)
197 if(i==WrfDataHandleMax) then
198 Status = WRF_WARN_TOO_MANY_FILES
199 write(msg,*) 'Warning TOO MANY FILES in ',__FILE__,', line', __LINE__
200 call wrf_debug ( WARN , TRIM(msg))
201 write(msg,*) 'Did you call ext_pnc_ioinit?'
202 call wrf_debug ( WARN , TRIM(msg))
209 DH%first_operation = .TRUE.
210 DH%Collective = .TRUE.
212 end subroutine allocHandle
214 subroutine deallocHandle(DataHandle, Status)
216 include 'wrf_status_codes.h'
217 integer ,intent(in) :: DataHandle
218 integer ,intent(out) :: Status
219 type(wrf_data_handle),pointer :: DH
223 IF ( DataHandle .GE. 1 .AND. DataHandle .LE. WrfDataHandleMax ) THEN
224 if(.NOT. WrfDataHandles(DataHandle)%Free) then
225 DH => WrfDataHandles(DataHandle)
226 deallocate(DH%Times, STAT=stat)
228 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
229 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
230 call wrf_debug ( FATAL , msg)
233 deallocate(DH%DimLengths, 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%DimIDs, 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%DimNames, STAT=stat)
249 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
250 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
251 call wrf_debug ( FATAL , msg)
254 deallocate(DH%MDVarIDs, STAT=stat)
256 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
257 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
258 call wrf_debug ( FATAL , msg)
261 deallocate(DH%MDVarDimLens, 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%MDVarNames, 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%VarIDs, 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%VarDimLens, 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%VarNames, 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)
300 end subroutine deallocHandle
302 subroutine GetDH(DataHandle,DH,Status)
304 include 'wrf_status_codes.h'
305 integer ,intent(in) :: DataHandle
306 type(wrf_data_handle) ,pointer :: DH
307 integer ,intent(out) :: Status
309 if(DataHandle < 1 .or. DataHandle > WrfDataHandleMax) then
310 Status = WRF_WARN_BAD_DATA_HANDLE
313 DH => WrfDataHandles(DataHandle)
315 Status = WRF_WARN_BAD_DATA_HANDLE
322 subroutine DateCheck(Date,Status)
324 include 'wrf_status_codes.h'
325 character*(*) ,intent(in) :: Date
326 integer ,intent(out) :: Status
328 if(len(Date) /= DateStrLen) then
329 Status = WRF_WARN_DATESTR_BAD_LENGTH
334 end subroutine DateCheck
336 subroutine GetName(Element,Var,Name,Status)
338 include 'wrf_status_codes.h'
339 character*(*) ,intent(in) :: Element
340 character*(*) ,intent(in) :: Var
341 character*(*) ,intent(out) :: Name
342 integer ,intent(out) :: Status
343 character (VarNameLen) :: VarName
346 integer, parameter :: upper_to_lower =IACHAR('a')-IACHAR('A')
349 Name = 'MD___'//trim(Element)//VarName
352 if('A'<=c .and. c <='Z') Name(i:i)=achar(iachar(c)+upper_to_lower)
353 if(c=='-'.or.c==':') Name(i:i)='_'
357 end subroutine GetName
359 subroutine GetTimeIndex(IO,DataHandle,DateStr,TimeIndex,Status)
361 include 'wrf_status_codes.h'
362 # include "pnetcdf.inc"
363 character (*) ,intent(in) :: IO
364 integer ,intent(in) :: DataHandle
365 character*(*) ,intent(in) :: DateStr
366 integer ,intent(out) :: TimeIndex
367 integer ,intent(out) :: Status
368 type(wrf_data_handle) ,pointer :: DH
369 integer(KIND=MPI_OFFSET_KIND) :: VStart(2)
370 integer(KIND=MPI_OFFSET_KIND) :: VCount(2)
374 DH => WrfDataHandles(DataHandle)
375 call DateCheck(DateStr,Status)
376 if(Status /= WRF_NO_ERR) then
377 Status = WRF_WARN_DATESTR_ERROR
378 write(msg,*) 'Warning DATE STRING ERROR in ',__FILE__,', line', __LINE__
379 call wrf_debug ( WARN , TRIM(msg))
382 if(IO == 'write') then
383 TimeIndex = DH%TimeIndex
384 if(TimeIndex <= 0) then
386 elseif(DateStr == DH%Times(TimeIndex)) then
390 TimeIndex = TimeIndex +1
391 if(TimeIndex > MaxTimes) then
392 Status = WRF_WARN_TIME_EOF
393 write(msg,*) 'Warning TIME EOF in ',__FILE__,', line', __LINE__
394 call wrf_debug ( WARN , TRIM(msg))
398 DH%TimeIndex = TimeIndex
399 DH%Times(TimeIndex) = DateStr
401 VStart(2) = TimeIndex
402 VCount(1) = DateStrLen
404 stat = NFMPI_PUT_VARA_TEXT_ALL(DH%NCID,DH%TimesVarID,VStart,VCount,DateStr)
405 call netcdf_err(stat,Status)
406 if(Status /= WRF_NO_ERR) then
407 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
408 call wrf_debug ( WARN , TRIM(msg))
413 if(DH%Times(i)==DateStr) then
419 Status = WRF_WARN_TIME_NF
420 write(msg,*) 'Warning TIME ',DateStr,' NOT FOUND in ',__FILE__,', line', __LINE__
421 call wrf_debug ( WARN , TRIM(msg))
427 end subroutine GetTimeIndex
429 subroutine GetDim(MemoryOrder,NDim,Status)
430 include 'wrf_status_codes.h'
431 character*(*) ,intent(in) :: MemoryOrder
432 integer ,intent(out) :: NDim
433 integer ,intent(out) :: Status
434 character*3 :: MemOrd
436 call LowerCase(MemoryOrder,MemOrd)
438 case ('xyz','xzy','yxz','yzx','zxy','zyx','xsz','xez','ysz','yez')
440 case ('xy','yx','xs','xe','ys','ye')
444 case ('0') ! NDim=0 for scalars. TBH: 20060502
447 print *, 'memory order = ',MemOrd,' ',MemoryOrder
448 Status = WRF_WARN_BAD_MEMORYORDER
453 end subroutine GetDim
455 subroutine GetIndices(NDim,Start,End,i1,i2,j1,j2,k1,k2)
456 integer ,intent(in) :: NDim
457 integer ,dimension(*),intent(in) :: Start,End
458 integer ,intent(out) :: i1,i2,j1,j2,k1,k2
466 if(NDim == 0) return ! NDim=0 for scalars. TBH: 20060502
476 end subroutine GetIndices
478 logical function ZeroLengthHorzDim(MemoryOrder,Vector,Status)
480 include 'wrf_status_codes.h'
481 character*(*) ,intent(in) :: MemoryOrder
482 integer,dimension(*) ,intent(in) :: Vector
483 integer ,intent(out) :: Status
485 integer,dimension(NVarDims) :: temp
486 character*3 :: MemOrd
489 call GetDim(MemoryOrder,NDim,Status)
490 temp(1:NDim) = Vector(1:NDim)
491 call LowerCase(MemoryOrder,MemOrd)
492 zero_length = .false.
494 case ('xsz','xez','ysz','yez','xs','xe','ys','ye','z','c')
497 continue ! NDim=0 for scalars. TBH: 20060502
499 zero_length = temp(1) .lt. 1 .or. temp(3) .lt. 1
500 case ('xy','yx','xyz','yxz')
501 zero_length = temp(1) .lt. 1 .or. temp(2) .lt. 1
503 zero_length = temp(2) .lt. 1 .or. temp(3) .lt. 1
505 Status = WRF_WARN_BAD_MEMORYORDER
506 ZeroLengthHorzDim = .true.
510 ZeroLengthHorzDim = zero_length
512 end function ZeroLengthHorzDim
514 subroutine ExtOrder(MemoryOrder,Vector,Status)
516 include 'wrf_status_codes.h'
517 character*(*) ,intent(in) :: MemoryOrder
518 integer,dimension(*) ,intent(inout) :: Vector
519 integer ,intent(out) :: Status
521 integer,dimension(NVarDims) :: temp
522 character*3 :: MemOrd
524 call GetDim(MemoryOrder,NDim,Status)
525 temp(1:NDim) = Vector(1:NDim)
526 call LowerCase(MemoryOrder,MemOrd)
529 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c')
532 continue ! NDim=0 for scalars. TBH: 20060502
554 Status = WRF_WARN_BAD_MEMORYORDER
559 end subroutine ExtOrder
561 subroutine ExtOrderStr(MemoryOrder,Vector,ROVector,Status)
563 include 'wrf_status_codes.h'
564 character*(*) ,intent(in) :: MemoryOrder
565 character*(*),dimension(*) ,intent(in) :: Vector
566 character(80),dimension(NVarDims),intent(out) :: ROVector
567 integer ,intent(out) :: Status
569 character*3 :: MemOrd
571 call GetDim(MemoryOrder,NDim,Status)
572 ROVector(1:NDim) = Vector(1:NDim)
573 call LowerCase(MemoryOrder,MemOrd)
576 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c')
579 continue ! NDim=0 for scalars. TBH: 20060502
581 ROVector(2) = Vector(3)
582 ROVector(3) = Vector(2)
584 ROVector(1) = Vector(2)
585 ROVector(2) = Vector(1)
587 ROVector(1) = Vector(3)
588 ROVector(2) = Vector(1)
589 ROVector(3) = Vector(2)
591 ROVector(1) = Vector(2)
592 ROVector(2) = Vector(3)
593 ROVector(3) = Vector(1)
595 ROVector(1) = Vector(3)
596 ROVector(3) = Vector(1)
598 ROVector(1) = Vector(2)
599 ROVector(2) = Vector(1)
601 Status = WRF_WARN_BAD_MEMORYORDER
606 end subroutine ExtOrderStr
609 subroutine LowerCase(MemoryOrder,MemOrd)
610 character*(*) ,intent(in) :: MemoryOrder
611 character*(*) ,intent(out) :: MemOrd
613 integer ,parameter :: upper_to_lower =IACHAR('a')-IACHAR('A')
618 MemOrd(1:N) = MemoryOrder(1:N)
621 if('A'<=c .and. c <='Z') MemOrd(i:i)=achar(iachar(c)+upper_to_lower)
624 end subroutine LowerCase
626 subroutine UpperCase(MemoryOrder,MemOrd)
627 character*(*) ,intent(in) :: MemoryOrder
628 character*(*) ,intent(out) :: MemOrd
630 integer ,parameter :: lower_to_upper =IACHAR('A')-IACHAR('a')
635 MemOrd(1:N) = MemoryOrder(1:N)
638 if('a'<=c .and. c <='z') MemOrd(i:i)=achar(iachar(c)+lower_to_upper)
641 end subroutine UpperCase
643 subroutine netcdf_err(err,Status)
645 include 'wrf_status_codes.h'
646 # include "pnetcdf.inc"
647 integer ,intent(in) :: err
648 integer ,intent(out) :: Status
649 character(len=80) :: errmsg
652 if( err==NF_NOERR )then
655 errmsg = NFMPI_STRERROR(err)
656 write(msg,*) 'NetCDF error: ',errmsg
657 call wrf_debug ( WARN , TRIM(msg))
658 Status = WRF_WARN_NETCDF
661 end subroutine netcdf_err
663 subroutine FieldIO(IO,DataHandle,DateStr,Starts,Length,MemoryOrder &
664 ,FieldType,NCID,VarID,XField,Status)
666 include 'wrf_status_codes.h'
667 # include "pnetcdf.inc"
668 character (*) ,intent(in) :: IO
669 integer ,intent(in) :: DataHandle
670 character*(*) ,intent(in) :: DateStr
671 integer,dimension(NVarDims),intent(in) :: Starts
672 integer,dimension(NVarDims),intent(in) :: Length
673 character*(*) ,intent(in) :: MemoryOrder
674 integer ,intent(in) :: FieldType
675 integer ,intent(in) :: NCID
676 integer ,intent(in) :: VarID
677 integer,dimension(*) ,intent(inout) :: XField
678 integer ,intent(out) :: Status
681 integer,dimension(NVarDims) :: VStart
682 integer,dimension(NVarDims) :: VCount
684 call GetTimeIndex(IO,DataHandle,DateStr,TimeIndex,Status)
685 if(Status /= WRF_NO_ERR) then
686 write(msg,*) 'Warning in ',__FILE__,', line', __LINE__
687 call wrf_debug ( WARN , TRIM(msg))
688 write(msg,*) ' Bad time index for DateStr = ',DateStr
689 call wrf_debug ( WARN , TRIM(msg))
692 call GetDim(MemoryOrder,NDim,Status)
695 !jm for parallel netcef VStart(1:NDim) = 1
696 VStart(1:NDim) = Starts(1:NDim)
697 VCount(1:NDim) = Length(1:NDim)
698 VStart(NDim+1) = TimeIndex
700 select case (FieldType)
702 call ext_pnc_RealFieldIO (WrfDataHandles(DataHandle)%Collective, &
703 IO,NCID,VarID,VStart,VCount,XField,Status)
705 call ext_pnc_DoubleFieldIO (WrfDataHandles(DataHandle)%Collective, &
706 IO,NCID,VarID,VStart,VCount,XField,Status)
708 call ext_pnc_IntFieldIO (WrfDataHandles(DataHandle)%Collective, &
709 IO,NCID,VarID,VStart,VCount,XField,Status)
711 call ext_pnc_LogicalFieldIO (WrfDataHandles(DataHandle)%Collective, &
712 IO,NCID,VarID,VStart,VCount,XField,Status)
713 if(Status /= WRF_NO_ERR) return
715 !for wrf_complex, double_complex
716 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
717 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
718 call wrf_debug ( WARN , TRIM(msg))
722 end subroutine FieldIO
724 subroutine Transpose(IO,MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
725 ,XField,x1,x2,y1,y2,z1,z2 &
727 character*(*) ,intent(in) :: IO
728 character*(*) ,intent(in) :: MemoryOrder
729 integer ,intent(in) :: l1,l2,m1,m2,n1,n2
730 integer ,intent(in) :: di
731 integer ,intent(in) :: x1,x2,y1,y2,z1,z2
732 integer ,intent(in) :: i1,i2,j1,j2,k1,k2
733 integer ,intent(inout) :: Field(di,l1:l2,m1:m2,n1:n2)
734 !jm 010827 integer ,intent(inout) :: XField(di,x1:x2,y1:y2,z1:z2)
735 integer ,intent(inout) :: XField(di,(i2-i1+1)*(j2-j1+1)*(k2-k1+1))
736 character*3 :: MemOrd
738 integer ,parameter :: MaxUpperCase=IACHAR('Z')
739 integer :: i,j,k,ix,jx,kx
741 call LowerCase(MemoryOrder,MemOrd)
744 !#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))
745 ! define(`XDEX',($1-``$1''1+1+(``$1''2-``$1''1+1)*(($2-``$2''1)+($3-``$3''1)*(``$2''2-``$2''1+1))))
749 #define DFIELD XField(1:di,XDEX(i,k,j))
750 #include "transpose.code"
751 case ('xyz','xsz','xez','ysz','yez','xy','xs','xe','ys','ye','z','c','0')
753 #define DFIELD XField(1:di,XDEX(i,j,k))
754 #include "transpose.code"
757 #define DFIELD XField(1:di,XDEX(j,i,k))
758 #include "transpose.code"
761 #define DFIELD XField(1:di,XDEX(k,i,j))
762 #include "transpose.code"
765 #define DFIELD XField(1:di,XDEX(j,k,i))
766 #include "transpose.code"
769 #define DFIELD XField(1:di,XDEX(k,j,i))
770 #include "transpose.code"
773 #define DFIELD XField(1:di,XDEX(j,i,k))
774 #include "transpose.code"
777 end subroutine Transpose
779 subroutine reorder (MemoryOrder,MemO)
780 character*(*) ,intent(in) :: MemoryOrder
781 character*3 ,intent(out) :: MemO
782 character*3 :: MemOrd
783 integer :: N,i,i1,i2,i3
786 N = len_trim(MemoryOrder)
788 call lowercase(MemoryOrder,MemOrd)
789 ! never invert the boundary codes
790 select case ( MemOrd )
791 case ( 'xsz','xez','ysz','yez' )
799 if(ichar(MemOrd(i:i)) < ichar(MemOrd(i1:i1))) I1 = i
800 if(ichar(MemOrd(i:i)) > ichar(MemOrd(i3:i3))) I3 = i
807 MemO(1:1) = MemoryOrder(i1:i1)
808 MemO(2:2) = MemoryOrder(i2:i2)
809 if(N == 3) MemO(3:3) = MemoryOrder(i3:i3)
810 if(MemOrd(i1:i1) == 's' .or. MemOrd(i1:i1) == 'e') then
811 MemO(1:N-1) = MemO(2:N)
812 MemO(N:N ) = MemoryOrder(i1:i1)
815 end subroutine reorder
817 ! Returns .TRUE. iff it is OK to write time-independent domain metadata to the
818 ! file referenced by DataHandle. If DataHandle is invalid, .FALSE. is
820 LOGICAL FUNCTION ncd_ok_to_put_dom_ti( DataHandle )
822 include 'wrf_status_codes.h'
823 INTEGER, INTENT(IN) :: DataHandle
824 CHARACTER*80 :: fname
827 LOGICAL :: dryrun, first_output, retval
828 call ext_pnc_inquire_filename( DataHandle, fname, filestate, Status )
829 IF ( Status /= WRF_NO_ERR ) THEN
830 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
832 call wrf_debug ( WARN , TRIM(msg) )
835 dryrun = ( filestate .EQ. WRF_FILE_OPENED_NOT_COMMITTED )
836 first_output = ncd_is_first_operation( DataHandle )
837 ! retval = .NOT. dryrun .AND. first_output
840 ncd_ok_to_put_dom_ti = retval
842 END FUNCTION ncd_ok_to_put_dom_ti
844 ! Returns .TRUE. iff it is OK to read time-independent domain metadata from the
845 ! file referenced by DataHandle. If DataHandle is invalid, .FALSE. is
847 LOGICAL FUNCTION ncd_ok_to_get_dom_ti( DataHandle )
849 include 'wrf_status_codes.h'
850 INTEGER, INTENT(IN) :: DataHandle
851 CHARACTER*80 :: fname
854 LOGICAL :: dryrun, retval
855 call ext_pnc_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 retval = .NOT. dryrun
865 ncd_ok_to_get_dom_ti = retval
867 END FUNCTION ncd_ok_to_get_dom_ti
869 ! Returns .TRUE. iff nothing has been read from or written to the file
870 ! referenced by DataHandle. If DataHandle is invalid, .FALSE. is returned.
871 LOGICAL FUNCTION ncd_is_first_operation( DataHandle )
873 INCLUDE 'wrf_status_codes.h'
874 INTEGER, INTENT(IN) :: DataHandle
875 TYPE(wrf_data_handle) ,POINTER :: DH
878 CALL GetDH( DataHandle, DH, Status )
879 IF ( Status /= WRF_NO_ERR ) THEN
880 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__, &
882 call wrf_debug ( WARN , TRIM(msg) )
885 retval = DH%first_operation
887 ncd_is_first_operation = retval
889 END FUNCTION ncd_is_first_operation
891 end module ext_pnc_support_routines
893 subroutine ext_pnc_open_for_read(DatasetName, Comm1, Comm2, SysDepInfo, DataHandle, Status)
895 use ext_pnc_support_routines
897 include 'wrf_status_codes.h'
898 # include "pnetcdf.inc"
899 character *(*), INTENT(IN) :: DatasetName
900 integer , INTENT(IN) :: Comm1, Comm2
901 character *(*), INTENT(IN) :: SysDepInfo
902 integer , INTENT(OUT) :: DataHandle
903 integer , INTENT(OUT) :: Status
904 DataHandle = 0 ! dummy setting to quiet warning message
905 CALL ext_pnc_open_for_read_begin( DatasetName, Comm1, Comm2, SysDepInfo, DataHandle, Status )
906 IF ( Status .EQ. WRF_NO_ERR ) THEN
907 CALL ext_pnc_open_for_read_commit( DataHandle, Status )
910 end subroutine ext_pnc_open_for_read
912 !ends training phase; switches internal flag to enable input
913 !must be paired with call to ext_pnc_open_for_read_begin
914 subroutine ext_pnc_open_for_read_commit(DataHandle, Status)
916 use ext_pnc_support_routines
918 include 'wrf_status_codes.h'
919 # include "pnetcdf.inc"
920 integer, intent(in) :: DataHandle
921 integer, intent(out) :: Status
922 type(wrf_data_handle) ,pointer :: DH
924 if(WrfIOnotInitialized) then
925 Status = WRF_IO_NOT_INITIALIZED
926 write(msg,*) 'ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
927 call wrf_debug ( FATAL , msg)
930 call GetDH(DataHandle,DH,Status)
931 if(Status /= WRF_NO_ERR) then
932 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
933 call wrf_debug ( WARN , TRIM(msg))
936 DH%FileStatus = WRF_FILE_OPENED_FOR_READ
937 DH%first_operation = .TRUE.
940 end subroutine ext_pnc_open_for_read_commit
942 subroutine ext_pnc_open_for_read_begin( FileName, Comm, IOComm, SysDepInfo, DataHandle, Status)
944 use ext_pnc_support_routines
946 include 'wrf_status_codes.h'
947 # include "pnetcdf.inc"
948 character*(*) ,intent(IN) :: FileName
949 integer ,intent(IN) :: Comm
950 integer ,intent(IN) :: IOComm
951 character*(*) ,intent(in) :: SysDepInfo
952 integer ,intent(out) :: DataHandle
953 integer ,intent(out) :: Status
954 type(wrf_data_handle) ,pointer :: DH
957 integer ,allocatable :: Buffer(:)
962 integer(KIND=MPI_OFFSET_KIND) :: VStart(2)
963 integer(KIND=MPI_OFFSET_KIND) :: VLen(2)
964 integer :: TotalNumVars
967 character (NF_MAX_NAME) :: Name
969 if(WrfIOnotInitialized) then
970 Status = WRF_IO_NOT_INITIALIZED
971 write(msg,*) 'ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
972 call wrf_debug ( FATAL , msg)
975 call allocHandle(DataHandle,DH,Comm,Status)
976 if(Status /= WRF_NO_ERR) then
977 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
978 call wrf_debug ( WARN , TRIM(msg))
981 stat = NFMPI_OPEN(Comm, FileName, NF_NOWRITE, MPI_INFO_NULL, DH%NCID)
982 call netcdf_err(stat,Status)
983 if(Status /= WRF_NO_ERR) then
984 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
985 call wrf_debug ( WARN , TRIM(msg))
988 stat = NFMPI_INQ_VARID(DH%NCID,DH%TimesName,VarID)
989 call netcdf_err(stat,Status)
990 if(Status /= WRF_NO_ERR) then
991 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
992 call wrf_debug ( WARN , TRIM(msg))
995 stat = NFMPI_INQ_VAR(DH%NCID,VarID,DH%TimesName, XType, StoredDim, DimIDs, NAtts)
996 call netcdf_err(stat,Status)
997 if(Status /= WRF_NO_ERR) then
998 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
999 call wrf_debug ( WARN , TRIM(msg))
1002 if(XType/=NF_CHAR) then
1003 Status = WRF_WARN_TYPE_MISMATCH
1004 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
1005 call wrf_debug ( WARN , TRIM(msg))
1008 stat = NFMPI_INQ_DIMLEN(DH%NCID,DimIDs(1),VLen(1))
1009 call netcdf_err(stat,Status)
1010 if(Status /= WRF_NO_ERR) then
1011 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1012 call wrf_debug ( WARN , TRIM(msg))
1015 if(VLen(1) /= DateStrLen) then
1016 Status = WRF_WARN_DATESTR_BAD_LENGTH
1017 write(msg,*) 'Warning DATESTR BAD LENGTH in ',__FILE__,', line', __LINE__
1018 call wrf_debug ( WARN , TRIM(msg))
1021 stat = NFMPI_INQ_DIMLEN(DH%NCID,DimIDs(2),VLen(2))
1022 call netcdf_err(stat,Status)
1023 if(Status /= WRF_NO_ERR) then
1024 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1025 call wrf_debug ( WARN , TRIM(msg))
1028 if(VLen(2) > MaxTimes) then
1029 Status = WRF_ERR_FATAL_TOO_MANY_TIMES
1030 write(msg,*) 'Fatal TOO MANY TIME VALUES in ',__FILE__,', line', __LINE__
1031 call wrf_debug ( FATAL , TRIM(msg))
1036 stat = NFMPI_GET_VARA_TEXT_ALL(DH%NCID,VarID,VStart,VLen,DH%Times)
1037 call netcdf_err(stat,Status)
1038 if(Status /= WRF_NO_ERR) then
1039 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1040 call wrf_debug ( WARN , TRIM(msg))
1043 stat = NFMPI_INQ_NVARS(DH%NCID,TotalNumVars)
1044 call netcdf_err(stat,Status)
1045 if(Status /= WRF_NO_ERR) then
1046 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1047 call wrf_debug ( WARN , TRIM(msg))
1052 stat = NFMPI_INQ_VARNAME(DH%NCID,i,Name)
1053 call netcdf_err(stat,Status)
1054 if(Status /= WRF_NO_ERR) then
1055 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1056 call wrf_debug ( WARN , TRIM(msg))
1058 elseif(Name(1:5) /= 'md___' .and. Name /= DH%TimesName) then
1060 DH%VarNames(NumVars) = Name
1061 DH%VarIDs(NumVars) = i
1064 DH%NumVars = NumVars
1065 DH%NumberTimes = VLen(2)
1066 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1067 DH%FileName = FileName
1068 DH%CurrentVariable = 0
1070 DH%TimesVarID = VarID
1073 end subroutine ext_pnc_open_for_read_begin
1075 subroutine ext_pnc_open_for_update( FileName, Comm, IOComm, SysDepInfo, DataHandle, Status)
1077 use ext_pnc_support_routines
1079 include 'wrf_status_codes.h'
1080 # include "pnetcdf.inc"
1081 character*(*) ,intent(IN) :: FileName
1082 integer ,intent(IN) :: Comm
1083 integer ,intent(IN) :: IOComm
1084 character*(*) ,intent(in) :: SysDepInfo
1085 integer ,intent(out) :: DataHandle
1086 integer ,intent(out) :: Status
1087 type(wrf_data_handle) ,pointer :: DH
1090 integer ,allocatable :: Buffer(:)
1092 integer :: StoredDim
1094 integer :: DimIDs(2)
1095 integer :: VStart(2)
1097 integer :: TotalNumVars
1100 character (NF_MAX_NAME) :: Name
1102 if(WrfIOnotInitialized) then
1103 Status = WRF_IO_NOT_INITIALIZED
1104 write(msg,*) 'ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
1105 call wrf_debug ( FATAL , msg)
1108 call allocHandle(DataHandle,DH,Comm,Status)
1109 if(Status /= WRF_NO_ERR) then
1110 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
1111 call wrf_debug ( WARN , TRIM(msg))
1114 stat = NFMPI_OPEN(Comm, FileName, NF_WRITE, MPI_INFO_NULL, DH%NCID)
1115 call netcdf_err(stat,Status)
1116 if(Status /= WRF_NO_ERR) then
1117 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1118 call wrf_debug ( WARN , TRIM(msg))
1121 stat = NFMPI_INQ_VARID(DH%NCID,DH%TimesName,VarID)
1122 call netcdf_err(stat,Status)
1123 if(Status /= WRF_NO_ERR) then
1124 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1125 call wrf_debug ( WARN , TRIM(msg))
1128 stat = NFMPI_INQ_VAR(DH%NCID,VarID,DH%TimesName, XType, StoredDim, DimIDs, NAtts)
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(XType/=NF_CHAR) then
1136 Status = WRF_WARN_TYPE_MISMATCH
1137 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
1138 call wrf_debug ( WARN , TRIM(msg))
1141 stat = NFMPI_INQ_DIMLEN(DH%NCID,DimIDs(1),VLen(1))
1142 call netcdf_err(stat,Status)
1143 if(Status /= WRF_NO_ERR) then
1144 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1145 call wrf_debug ( WARN , TRIM(msg))
1148 if(VLen(1) /= DateStrLen) then
1149 Status = WRF_WARN_DATESTR_BAD_LENGTH
1150 write(msg,*) 'Warning DATESTR BAD LENGTH in ',__FILE__,', line', __LINE__
1151 call wrf_debug ( WARN , TRIM(msg))
1154 stat = NFMPI_INQ_DIMLEN(DH%NCID,DimIDs(2),VLen(2))
1155 call netcdf_err(stat,Status)
1156 if(Status /= WRF_NO_ERR) then
1157 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1158 call wrf_debug ( WARN , TRIM(msg))
1161 if(VLen(2) > MaxTimes) then
1162 Status = WRF_ERR_FATAL_TOO_MANY_TIMES
1163 write(msg,*) 'Fatal TOO MANY TIME VALUES in ',__FILE__,', line', __LINE__
1164 call wrf_debug ( FATAL , TRIM(msg))
1169 stat = NFMPI_GET_VARA_TEXT_ALL(DH%NCID,VarID,VStart,VLen,DH%Times)
1170 call netcdf_err(stat,Status)
1171 if(Status /= WRF_NO_ERR) then
1172 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1173 call wrf_debug ( WARN , TRIM(msg))
1176 stat = NFMPI_INQ_NVARS(DH%NCID,TotalNumVars)
1177 call netcdf_err(stat,Status)
1178 if(Status /= WRF_NO_ERR) then
1179 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1180 call wrf_debug ( WARN , TRIM(msg))
1185 stat = NFMPI_INQ_VARNAME(DH%NCID,i,Name)
1186 call netcdf_err(stat,Status)
1187 if(Status /= WRF_NO_ERR) then
1188 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1189 call wrf_debug ( WARN , TRIM(msg))
1191 elseif(Name(1:5) /= 'md___' .and. Name /= DH%TimesName) then
1193 DH%VarNames(NumVars) = Name
1194 DH%VarIDs(NumVars) = i
1197 DH%NumVars = NumVars
1198 DH%NumberTimes = VLen(2)
1199 DH%FileStatus = WRF_FILE_OPENED_FOR_UPDATE
1200 DH%FileName = FileName
1201 DH%CurrentVariable = 0
1203 DH%TimesVarID = VarID
1206 end subroutine ext_pnc_open_for_update
1209 SUBROUTINE ext_pnc_open_for_write_begin(FileName,Comm,IOComm,SysDepInfo,DataHandle,Status)
1211 use ext_pnc_support_routines
1213 include 'wrf_status_codes.h'
1214 # include "pnetcdf.inc"
1215 character*(*) ,intent(in) :: FileName
1216 integer ,intent(in) :: Comm
1217 integer ,intent(in) :: IOComm
1218 character*(*) ,intent(in) :: SysDepInfo
1219 integer ,intent(out) :: DataHandle
1220 integer ,intent(out) :: Status
1221 type(wrf_data_handle),pointer :: DH
1224 character (7) :: Buffer
1225 integer :: VDimIDs(2)
1226 integer :: info, ierr ! added for Blue Gene (see NF_CREAT below)
1227 character*1024 :: newFileName
1229 integer local_communicator_x, ntasks_x
1231 if(WrfIOnotInitialized) then
1232 Status = WRF_IO_NOT_INITIALIZED
1233 write(msg,*) 'ext_pnc_open_for_write_begin: ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
1234 call wrf_debug ( FATAL , msg)
1237 call allocHandle(DataHandle,DH,Comm,Status)
1238 if(Status /= WRF_NO_ERR) then
1239 write(msg,*) 'Fatal ALLOCATION ERROR in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1240 call wrf_debug ( FATAL , TRIM(msg))
1247 call mpi_info_create( info, ierr )
1249 CALL mpi_info_set(info,"romio_ds_write","disable", ierr) ; write(0,*)'mpi_info_set write returns ',ierr
1250 CALL mpi_info_set(info,"romio_ds_read","disable", ierr) ; write(0,*)'mpi_info_set read returns ',ierr
1254 ! Remove the dash/underscore change to filenames for pnetcdf...
1255 write(newFileName, fmt="(a)") TRIM(ADJUSTL(FileName))
1256 do i = 1, len_trim(newFileName)
1257 ! if(newFileName(i:i) == '-') newFileName(i:i) = '_'
1258 if(newFileName(i:i) == ':') newFileName(i:i) = '_'
1260 stat = NFMPI_CREATE(Comm, newFileName, IOR(NF_CLOBBER, NF_64BIT_OFFSET), info, DH%NCID)
1261 ! stat = NFMPI_CREATE(Comm, newFileName, NF_64BIT_OFFSET, info, DH%NCID)
1262 call mpi_info_free( info, ierr)
1265 ! rob latham suggested hint
1267 call mpi_info_create( info, ierr )
1268 ! call mpi_info_set(info,'cd_buffer_size','4194304',ierr)
1269 call mpi_info_set(info,'cd_buffer_size','8388608',ierr)
1270 stat = NFMPI_CREATE(Comm, FileName, IOR(NF_CLOBBER, NF_64BIT_OFFSET), info, DH%NCID)
1271 call mpi_info_free( info, ierr)
1276 call netcdf_err(stat,Status)
1277 if(Status /= WRF_NO_ERR) then
1278 write(msg,*) 'NetCDF error in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1279 call wrf_debug ( WARN , TRIM(msg))
1282 ! JPE added for performance
1283 stat = NFMPI_SET_FILL(DH%NCID, NF_NOFILL, i)
1285 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1286 DH%FileName = FileName
1287 stat = NFMPI_DEF_DIM(DH%NCID,DH%DimUnlimName,i2offset(NF_UNLIMITED),DH%DimUnlimID)
1288 call netcdf_err(stat,Status)
1289 if(Status /= WRF_NO_ERR) then
1290 write(msg,*) 'NetCDF error in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1291 call wrf_debug ( WARN , TRIM(msg))
1294 DH%VarNames (1:MaxVars) = NO_NAME
1295 DH%MDVarNames(1:MaxVars) = NO_NAME
1297 write(Buffer,FMT="('DIM',i4.4)") i
1298 DH%DimNames (i) = Buffer
1299 DH%DimLengths(i) = NO_DIM
1301 DH%DimNames(1) = 'DateStrLen'
1302 stat = NFMPI_DEF_DIM(DH%NCID,DH%DimNames(1),i2offset(DateStrLen),DH%DimIDs(1))
1303 call netcdf_err(stat,Status)
1304 if(Status /= WRF_NO_ERR) then
1305 write(msg,*) 'NetCDF error in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1306 call wrf_debug ( WARN , TRIM(msg))
1309 VDimIDs(1) = DH%DimIDs(1)
1310 VDimIDs(2) = DH%DimUnlimID
1311 stat = NFMPI_DEF_VAR(DH%NCID,DH%TimesName,NF_CHAR,2,VDimIDs,DH%TimesVarID)
1312 call netcdf_err(stat,Status)
1313 if(Status /= WRF_NO_ERR) then
1314 write(msg,*) 'NetCDF error in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1315 call wrf_debug ( WARN , TRIM(msg))
1318 DH%DimLengths(1) = DateStrLen
1320 end subroutine ext_pnc_open_for_write_begin
1323 !opens a file for writing or coupler datastream for sending messages.
1324 !no training phase for this version of the open stmt.
1325 subroutine ext_pnc_open_for_write (DatasetName, Comm1, Comm2, &
1326 SysDepInfo, DataHandle, Status)
1328 use ext_pnc_support_routines
1330 include 'wrf_status_codes.h'
1331 # include "pnetcdf.inc"
1332 character *(*), intent(in) ::DatasetName
1333 integer , intent(in) ::Comm1, Comm2
1334 character *(*), intent(in) ::SysDepInfo
1335 integer , intent(out) :: DataHandle
1336 integer , intent(out) :: Status
1337 Status=WRF_WARN_NOOP
1338 DataHandle = 0 ! dummy setting to quiet warning message
1340 end subroutine ext_pnc_open_for_write
1342 SUBROUTINE ext_pnc_open_for_write_commit(DataHandle, Status)
1344 use ext_pnc_support_routines
1346 include 'wrf_status_codes.h'
1347 # include "pnetcdf.inc"
1348 integer ,intent(in) :: DataHandle
1349 integer ,intent(out) :: Status
1350 type(wrf_data_handle),pointer :: DH
1354 if(WrfIOnotInitialized) then
1355 Status = WRF_IO_NOT_INITIALIZED
1356 write(msg,*) 'ext_pnc_open_for_write_commit: ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
1357 call wrf_debug ( FATAL , msg)
1360 call GetDH(DataHandle,DH,Status)
1361 if(Status /= WRF_NO_ERR) then
1362 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_open_for_write_commit ',__FILE__,', line', __LINE__
1363 call wrf_debug ( WARN , TRIM(msg))
1366 stat = NFMPI_ENDDEF(DH%NCID)
1367 call netcdf_err(stat,Status)
1368 if(Status /= WRF_NO_ERR) then
1369 write(msg,*) 'NetCDF error (',stat,') from NFMPI_ENDDEF in ext_pnc_open_for_write_commit ',__FILE__,', line', __LINE__
1370 call wrf_debug ( WARN , TRIM(msg))
1373 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1374 DH%first_operation = .TRUE.
1376 end subroutine ext_pnc_open_for_write_commit
1378 subroutine ext_pnc_ioclose(DataHandle, Status)
1380 use ext_pnc_support_routines
1382 include 'wrf_status_codes.h'
1383 # include "pnetcdf.inc"
1384 integer ,intent(in) :: DataHandle
1385 integer ,intent(out) :: Status
1386 type(wrf_data_handle),pointer :: DH
1389 call GetDH(DataHandle,DH,Status)
1390 if(Status /= WRF_NO_ERR) then
1391 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1392 call wrf_debug ( WARN , TRIM(msg))
1395 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1396 Status = WRF_WARN_FILE_NOT_OPENED
1397 write(msg,*) 'Warning FILE NOT OPENED in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1398 call wrf_debug ( WARN , TRIM(msg))
1399 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1400 Status = WRF_WARN_DRYRUN_CLOSE
1401 write(msg,*) 'Warning TRY TO CLOSE DRYRUN in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1402 call wrf_debug ( WARN , TRIM(msg))
1403 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1405 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1407 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
1410 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1411 write(msg,*) 'Fatal error BAD FILE STATUS in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1412 call wrf_debug ( FATAL , TRIM(msg))
1416 stat = NFMPI_CLOSE(DH%NCID)
1417 call netcdf_err(stat,Status)
1418 if(Status /= WRF_NO_ERR) then
1419 write(msg,*) 'NetCDF error in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1420 call wrf_debug ( WARN , TRIM(msg))
1423 CALL deallocHandle( DataHandle, Status )
1426 end subroutine ext_pnc_ioclose
1428 subroutine ext_pnc_iosync( DataHandle, Status)
1430 use ext_pnc_support_routines
1432 include 'wrf_status_codes.h'
1433 # include "pnetcdf.inc"
1434 integer ,intent(in) :: DataHandle
1435 integer ,intent(out) :: Status
1436 type(wrf_data_handle),pointer :: DH
1439 call GetDH(DataHandle,DH,Status)
1440 if(Status /= WRF_NO_ERR) then
1441 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_iosync ',__FILE__,', line', __LINE__
1442 call wrf_debug ( WARN , TRIM(msg))
1445 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1446 Status = WRF_WARN_FILE_NOT_OPENED
1447 write(msg,*) 'Warning FILE NOT OPENED in ext_pnc_iosync ',__FILE__,', line', __LINE__
1448 call wrf_debug ( WARN , TRIM(msg))
1449 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1450 Status = WRF_WARN_FILE_NOT_COMMITTED
1451 write(msg,*) 'Warning FILE NOT COMMITTED in ext_pnc_iosync ',__FILE__,', line', __LINE__
1452 call wrf_debug ( WARN , TRIM(msg))
1453 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1455 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1458 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1459 write(msg,*) 'Fatal error BAD FILE STATUS in ext_pnc_iosync ',__FILE__,', line', __LINE__
1460 call wrf_debug ( FATAL , TRIM(msg))
1463 stat = NFMPI_SYNC(DH%NCID)
1464 call netcdf_err(stat,Status)
1465 if(Status /= WRF_NO_ERR) then
1466 write(msg,*) 'NetCDF error in ext_pnc_iosync ',__FILE__,', line', __LINE__
1467 call wrf_debug ( WARN , TRIM(msg))
1471 end subroutine ext_pnc_iosync
1475 subroutine ext_pnc_redef( DataHandle, Status)
1477 use ext_pnc_support_routines
1479 include 'wrf_status_codes.h'
1480 # include "pnetcdf.inc"
1481 integer ,intent(in) :: DataHandle
1482 integer ,intent(out) :: Status
1483 type(wrf_data_handle),pointer :: DH
1486 call GetDH(DataHandle,DH,Status)
1487 if(Status /= WRF_NO_ERR) then
1488 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1489 call wrf_debug ( WARN , TRIM(msg))
1492 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1493 Status = WRF_WARN_FILE_NOT_OPENED
1494 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1495 call wrf_debug ( WARN , TRIM(msg))
1496 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1497 Status = WRF_WARN_FILE_NOT_COMMITTED
1498 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1499 call wrf_debug ( WARN , TRIM(msg))
1500 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1502 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1503 Status = WRF_WARN_FILE_OPEN_FOR_READ
1504 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1505 call wrf_debug ( WARN , TRIM(msg))
1507 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1508 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1509 call wrf_debug ( FATAL , TRIM(msg))
1512 stat = NFMPI_REDEF(DH%NCID)
1513 call netcdf_err(stat,Status)
1514 if(Status /= WRF_NO_ERR) then
1515 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1516 call wrf_debug ( WARN , TRIM(msg))
1519 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1521 end subroutine ext_pnc_redef
1523 subroutine ext_pnc_enddef( DataHandle, Status)
1525 use ext_pnc_support_routines
1527 include 'wrf_status_codes.h'
1528 # include "pnetcdf.inc"
1529 integer ,intent(in) :: DataHandle
1530 integer ,intent(out) :: Status
1531 type(wrf_data_handle),pointer :: DH
1534 call GetDH(DataHandle,DH,Status)
1535 if(Status /= WRF_NO_ERR) then
1536 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1537 call wrf_debug ( WARN , TRIM(msg))
1540 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1541 Status = WRF_WARN_FILE_NOT_OPENED
1542 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1543 call wrf_debug ( WARN , TRIM(msg))
1544 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1545 Status = WRF_WARN_FILE_NOT_COMMITTED
1546 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1547 call wrf_debug ( WARN , TRIM(msg))
1548 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1550 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1551 Status = WRF_WARN_FILE_OPEN_FOR_READ
1552 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1553 call wrf_debug ( WARN , TRIM(msg))
1555 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1556 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1557 call wrf_debug ( FATAL , TRIM(msg))
1560 stat = NFMPI_ENDDEF(DH%NCID)
1561 call netcdf_err(stat,Status)
1562 if(Status /= WRF_NO_ERR) then
1563 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1564 call wrf_debug ( WARN , TRIM(msg))
1567 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1569 end subroutine ext_pnc_enddef
1571 subroutine ext_pnc_ioinit(SysDepInfo, Status)
1574 include 'wrf_status_codes.h'
1575 CHARACTER*(*), INTENT(IN) :: SysDepInfo
1576 INTEGER ,INTENT(INOUT) :: Status
1578 WrfIOnotInitialized = .false.
1579 WrfDataHandles(1:WrfDataHandleMax)%Free = .true.
1580 WrfDataHandles(1:WrfDataHandleMax)%TimesName = 'Times'
1581 WrfDataHandles(1:WrfDataHandleMax)%DimUnlimName = 'Time'
1582 WrfDataHandles(1:WrfDataHandleMax)%FileStatus = WRF_FILE_NOT_OPENED
1585 end subroutine ext_pnc_ioinit
1588 subroutine ext_pnc_inquiry (Inquiry, Result, Status)
1591 include 'wrf_status_codes.h'
1592 character *(*), INTENT(IN) :: Inquiry
1593 character *(*), INTENT(OUT) :: Result
1594 integer ,INTENT(INOUT) :: Status
1595 SELECT CASE (Inquiry)
1596 CASE ("RANDOM_WRITE","RANDOM_READ","SEQUENTIAL_WRITE","SEQUENTIAL_READ")
1598 CASE ("OPEN_READ","OPEN_COMMIT_WRITE")
1600 CASE ("OPEN_WRITE","OPEN_COMMIT_READ","PARALLEL_IO")
1602 CASE ("SELF_DESCRIBING","SUPPORT_METADATA","SUPPORT_3D_FIELDS")
1607 Result = 'No Result for that inquiry!'
1611 end subroutine ext_pnc_inquiry
1616 subroutine ext_pnc_ioexit(Status)
1618 use ext_pnc_support_routines
1620 include 'wrf_status_codes.h'
1621 # include "pnetcdf.inc"
1622 integer , INTENT(INOUT) ::Status
1624 type(wrf_data_handle),pointer :: DH
1627 if(WrfIOnotInitialized) then
1628 Status = WRF_IO_NOT_INITIALIZED
1629 write(msg,*) 'ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
1630 call wrf_debug ( FATAL , msg)
1633 do i=1,WrfDataHandleMax
1634 CALL deallocHandle( i , stat )
1637 end subroutine ext_pnc_ioexit
1639 subroutine ext_pnc_get_dom_ti_real(DataHandle,Element,Data,Count,OutCount,Status)
1640 #define ROUTINE_TYPE 'REAL'
1641 #define TYPE_DATA real,intent(out) :: Data(*)
1642 #define TYPE_COUNT integer,intent(in) :: Count
1643 #define TYPE_OUTCOUNT integer,intent(out) :: OutCOunt
1644 #define TYPE_BUFFER real,allocatable :: Buffer(:)
1645 #define NF_TYPE NF_FLOAT
1646 #define NF_ROUTINE NFMPI_GET_ATT_REAL
1647 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1648 #include "ext_pnc_get_dom_ti.code"
1649 end subroutine ext_pnc_get_dom_ti_real
1651 subroutine ext_pnc_get_dom_ti_integer(DataHandle,Element,Data,Count,OutCount,Status)
1658 #define ROUTINE_TYPE 'INTEGER'
1659 #define TYPE_DATA integer,intent(out) :: Data(*)
1660 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
1661 #define NF_TYPE NF_INT
1662 #define NF_ROUTINE NFMPI_GET_ATT_INT
1663 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1664 #include "ext_pnc_get_dom_ti.code"
1665 end subroutine ext_pnc_get_dom_ti_integer
1667 subroutine ext_pnc_get_dom_ti_double(DataHandle,Element,Data,Count,OutCount,Status)
1674 #define ROUTINE_TYPE 'DOUBLE'
1675 #define TYPE_DATA real*8,intent(out) :: Data(*)
1676 #define TYPE_BUFFER real*8,allocatable :: Buffer(:)
1677 #define NF_TYPE NF_DOUBLE
1678 #define NF_ROUTINE NFMPI_GET_ATT_DOUBLE
1679 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1680 #include "ext_pnc_get_dom_ti.code"
1681 end subroutine ext_pnc_get_dom_ti_double
1683 subroutine ext_pnc_get_dom_ti_logical(DataHandle,Element,Data,Count,OutCount,Status)
1690 #define ROUTINE_TYPE 'LOGICAL'
1691 #define TYPE_DATA logical,intent(out) :: Data(*)
1692 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
1693 #define NF_TYPE NF_INT
1694 #define NF_ROUTINE NFMPI_GET_ATT_INT
1695 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))==1
1696 #include "ext_pnc_get_dom_ti.code"
1697 end subroutine ext_pnc_get_dom_ti_logical
1699 subroutine ext_pnc_get_dom_ti_char(DataHandle,Element,Data,Status)
1703 #undef TYPE_OUTCOUNT
1706 #define ROUTINE_TYPE 'CHAR'
1707 #define TYPE_DATA character*(*),intent(out) :: Data
1709 #define TYPE_OUTCOUNT
1711 #define NF_TYPE NF_CHAR
1713 #include "ext_pnc_get_dom_ti.code"
1715 end subroutine ext_pnc_get_dom_ti_char
1717 subroutine ext_pnc_put_dom_ti_real(DataHandle,Element,Data,Count,Status)
1724 #define ROUTINE_TYPE 'REAL'
1725 #define TYPE_DATA real ,intent(in) :: Data(*)
1726 #define TYPE_COUNT integer,intent(in) :: Count
1727 #define NF_ROUTINE NFMPI_PUT_ATT_REAL
1728 #define ARGS NF_FLOAT,i2offset(Count),Data
1729 #include "ext_pnc_put_dom_ti.code"
1730 end subroutine ext_pnc_put_dom_ti_real
1732 subroutine ext_pnc_put_dom_ti_integer(DataHandle,Element,Data,Count,Status)
1739 #define ROUTINE_TYPE 'INTEGER'
1740 #define TYPE_DATA integer,intent(in) :: Data(*)
1741 #define TYPE_COUNT integer,intent(in) :: Count
1742 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1743 #define ARGS NF_INT,i2offset(Count),Data
1744 #include "ext_pnc_put_dom_ti.code"
1745 end subroutine ext_pnc_put_dom_ti_integer
1747 subroutine ext_pnc_put_dom_ti_double(DataHandle,Element,Data,Count,Status)
1754 #define ROUTINE_TYPE 'DOUBLE'
1755 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1756 #define TYPE_COUNT integer,intent(in) :: Count
1757 #define NF_ROUTINE NFMPI_PUT_ATT_DOUBLE
1758 #define ARGS NF_DOUBLE,i2offset(Count),Data
1759 #include "ext_pnc_put_dom_ti.code"
1760 end subroutine ext_pnc_put_dom_ti_double
1762 subroutine ext_pnc_put_dom_ti_logical(DataHandle,Element,Data,Count,Status)
1768 #define ROUTINE_TYPE 'LOGICAL'
1769 #define TYPE_DATA logical,intent(in) :: Data(*)
1770 #define TYPE_COUNT integer,intent(in) :: Count
1771 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1772 #define ARGS NF_INT,i2offset(Count),Buffer
1774 #include "ext_pnc_put_dom_ti.code"
1775 end subroutine ext_pnc_put_dom_ti_logical
1777 subroutine ext_pnc_put_dom_ti_char(DataHandle,Element,Data,Status)
1784 #define ROUTINE_TYPE 'CHAR'
1785 #define TYPE_DATA character*(*),intent(in) :: Data
1786 #define TYPE_COUNT integer,parameter :: Count=1
1787 #define NF_ROUTINE NFMPI_PUT_ATT_TEXT
1788 #define ARGS i2offset(len_trim(Data)),Data
1789 #include "ext_pnc_put_dom_ti.code"
1790 end subroutine ext_pnc_put_dom_ti_char
1792 subroutine ext_pnc_put_var_ti_real(DataHandle,Element,Var,Data,Count,Status)
1799 #define ROUTINE_TYPE 'REAL'
1800 #define TYPE_DATA real ,intent(in) :: Data(*)
1801 #define TYPE_COUNT integer ,intent(in) :: Count
1802 #define NF_ROUTINE NFMPI_PUT_ATT_REAL
1803 #define ARGS NF_FLOAT,i2offset(Count),Data
1804 #include "ext_pnc_put_var_ti.code"
1805 end subroutine ext_pnc_put_var_ti_real
1807 subroutine ext_pnc_put_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,Status)
1816 #define ROUTINE_TYPE 'REAL'
1817 #define TYPE_DATA real ,intent(in) :: Data(*)
1818 #define TYPE_COUNT integer ,intent(in) :: Count
1819 #define NF_ROUTINE NFMPI_PUT_VARA_REAL_ALL
1820 #define NF_TYPE NF_FLOAT
1821 #define LENGTH Count
1823 #include "ext_pnc_put_var_td.code"
1824 end subroutine ext_pnc_put_var_td_real
1826 subroutine ext_pnc_put_var_ti_double(DataHandle,Element,Var,Data,Count,Status)
1833 #define ROUTINE_TYPE 'DOUBLE'
1834 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1835 #define TYPE_COUNT integer ,intent(in) :: Count
1836 #define NF_ROUTINE NFMPI_PUT_ATT_DOUBLE
1837 #define ARGS NF_DOUBLE,i2offset(Count),Data
1838 #include "ext_pnc_put_var_ti.code"
1839 end subroutine ext_pnc_put_var_ti_double
1841 subroutine ext_pnc_put_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,Status)
1850 #define ROUTINE_TYPE 'DOUBLE'
1851 #define TYPE_DATA real*8,intent(in) :: Data(*)
1852 #define TYPE_COUNT integer ,intent(in) :: Count
1853 #define NF_ROUTINE NFMPI_PUT_VARA_DOUBLE_ALL
1854 #define NF_TYPE NF_DOUBLE
1855 #define LENGTH Count
1857 #include "ext_pnc_put_var_td.code"
1858 end subroutine ext_pnc_put_var_td_double
1860 subroutine ext_pnc_put_var_ti_integer(DataHandle,Element,Var,Data,Count,Status)
1867 #define ROUTINE_TYPE 'INTEGER'
1868 #define TYPE_DATA integer ,intent(in) :: Data(*)
1869 #define TYPE_COUNT integer ,intent(in) :: Count
1870 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1871 #define ARGS NF_INT,i2offset(Count),Data
1872 #include "ext_pnc_put_var_ti.code"
1873 end subroutine ext_pnc_put_var_ti_integer
1875 subroutine ext_pnc_put_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,Status)
1884 #define ROUTINE_TYPE 'INTEGER'
1885 #define TYPE_DATA integer ,intent(in) :: Data(*)
1886 #define TYPE_COUNT integer ,intent(in) :: Count
1887 #define NF_ROUTINE NFMPI_PUT_VARA_INT_ALL
1888 #define NF_TYPE NF_INT
1889 #define LENGTH Count
1891 #include "ext_pnc_put_var_td.code"
1892 end subroutine ext_pnc_put_var_td_integer
1894 subroutine ext_pnc_put_var_ti_logical(DataHandle,Element,Var,Data,Count,Status)
1900 #define ROUTINE_TYPE 'LOGICAL'
1901 #define TYPE_DATA logical ,intent(in) :: Data(*)
1902 #define TYPE_COUNT integer ,intent(in) :: Count
1903 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1905 #define ARGS NF_INT,i2offset(Count),Buffer
1906 #include "ext_pnc_put_var_ti.code"
1907 end subroutine ext_pnc_put_var_ti_logical
1909 subroutine ext_pnc_put_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,Status)
1917 #define ROUTINE_TYPE 'LOGICAL'
1918 #define TYPE_DATA logical ,intent(in) :: Data(*)
1919 #define TYPE_COUNT integer ,intent(in) :: Count
1920 #define NF_ROUTINE NFMPI_PUT_VARA_INT_ALL
1921 #define NF_TYPE NF_INT
1923 #define LENGTH Count
1925 #include "ext_pnc_put_var_td.code"
1926 end subroutine ext_pnc_put_var_td_logical
1928 subroutine ext_pnc_put_var_ti_char(DataHandle,Element,Var,Data,Status)
1935 #define ROUTINE_TYPE 'CHAR'
1936 #define TYPE_DATA character*(*) ,intent(in) :: Data
1938 #define NF_ROUTINE NFMPI_PUT_ATT_TEXT
1939 #define ARGS i2offset(len_trim(Data)),trim(Data)
1941 #include "ext_pnc_put_var_ti.code"
1943 end subroutine ext_pnc_put_var_ti_char
1945 subroutine ext_pnc_put_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
1954 #define ROUTINE_TYPE 'CHAR'
1955 #define TYPE_DATA character*(*) ,intent(in) :: Data
1957 #define NF_ROUTINE NFMPI_PUT_VARA_TEXT_ALL
1958 #define NF_TYPE NF_CHAR
1959 #define LENGTH len(Data)
1960 #include "ext_pnc_put_var_td.code"
1961 end subroutine ext_pnc_put_var_td_char
1963 subroutine ext_pnc_get_var_ti_real(DataHandle,Element,Var,Data,Count,OutCount,Status)
1968 #undef TYPE_OUTCOUNT
1972 #define ROUTINE_TYPE 'REAL'
1973 #define TYPE_DATA real ,intent(out) :: Data(*)
1974 #define TYPE_BUFFER real ,allocatable :: Buffer(:)
1975 #define TYPE_COUNT integer,intent(in) :: Count
1976 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
1977 #define NF_TYPE NF_FLOAT
1978 #define NF_ROUTINE NFMPI_GET_ATT_REAL
1979 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
1980 #include "ext_pnc_get_var_ti.code"
1981 end subroutine ext_pnc_get_var_ti_real
1983 subroutine ext_pnc_get_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
1988 #undef TYPE_OUTCOUNT
1993 #define ROUTINE_TYPE 'REAL'
1994 #define TYPE_DATA real ,intent(out) :: Data(*)
1995 #define TYPE_BUFFER real
1996 #define TYPE_COUNT integer,intent(in) :: Count
1997 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
1998 #define NF_TYPE NF_FLOAT
1999 #define NF_ROUTINE NFMPI_GET_VARA_REAL_ALL
2000 #define LENGTH min(Count,Len1)
2001 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2002 #include "ext_pnc_get_var_td.code"
2003 end subroutine ext_pnc_get_var_td_real
2005 subroutine ext_pnc_get_var_ti_double(DataHandle,Element,Var,Data,Count,OutCount,Status)
2010 #undef TYPE_OUTCOUNT
2014 #define ROUTINE_TYPE 'DOUBLE'
2015 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2016 #define TYPE_BUFFER real*8 ,allocatable :: Buffer(:)
2017 #define TYPE_COUNT integer,intent(in) :: Count
2018 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2019 #define NF_TYPE NF_DOUBLE
2020 #define NF_ROUTINE NFMPI_GET_ATT_DOUBLE
2021 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2022 #include "ext_pnc_get_var_ti.code"
2023 end subroutine ext_pnc_get_var_ti_double
2025 subroutine ext_pnc_get_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2030 #undef TYPE_OUTCOUNT
2035 #define ROUTINE_TYPE 'DOUBLE'
2036 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2037 #define TYPE_BUFFER real*8
2038 #define TYPE_COUNT integer,intent(in) :: Count
2039 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2040 #define NF_TYPE NF_DOUBLE
2041 #define NF_ROUTINE NFMPI_GET_VARA_DOUBLE_ALL
2042 #define LENGTH min(Count,Len1)
2043 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2044 #include "ext_pnc_get_var_td.code"
2045 end subroutine ext_pnc_get_var_td_double
2047 subroutine ext_pnc_get_var_ti_integer(DataHandle,Element,Var,Data,Count,OutCount,Status)
2052 #undef TYPE_OUTCOUNT
2056 #define ROUTINE_TYPE 'INTEGER'
2057 #define TYPE_DATA integer,intent(out) :: Data(*)
2058 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2059 #define TYPE_COUNT integer,intent(in) :: Count
2060 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2061 #define NF_TYPE NF_INT
2062 #define NF_ROUTINE NFMPI_GET_ATT_INT
2063 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2064 #include "ext_pnc_get_var_ti.code"
2065 end subroutine ext_pnc_get_var_ti_integer
2067 subroutine ext_pnc_get_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2072 #undef TYPE_OUTCOUNT
2077 #define ROUTINE_TYPE 'INTEGER'
2078 #define TYPE_DATA integer,intent(out) :: Data(*)
2079 #define TYPE_BUFFER integer
2080 #define TYPE_COUNT integer,intent(in) :: Count
2081 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2082 #define NF_TYPE NF_INT
2083 #define NF_ROUTINE NFMPI_GET_VARA_INT_ALL
2084 #define LENGTH min(Count,Len1)
2085 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2086 #include "ext_pnc_get_var_td.code"
2087 end subroutine ext_pnc_get_var_td_integer
2089 subroutine ext_pnc_get_var_ti_logical(DataHandle,Element,Var,Data,Count,OutCount,Status)
2094 #undef TYPE_OUTCOUNT
2098 #define ROUTINE_TYPE 'LOGICAL'
2099 #define TYPE_DATA logical,intent(out) :: Data(*)
2100 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2101 #define TYPE_COUNT integer,intent(in) :: Count
2102 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2103 #define NF_TYPE NF_INT
2104 #define NF_ROUTINE NFMPI_GET_ATT_INT
2105 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))==1
2106 #include "ext_pnc_get_var_ti.code"
2107 end subroutine ext_pnc_get_var_ti_logical
2109 subroutine ext_pnc_get_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2114 #undef TYPE_OUTCOUNT
2119 #define ROUTINE_TYPE 'LOGICAL'
2120 #define TYPE_DATA logical,intent(out) :: Data(*)
2121 #define TYPE_BUFFER integer
2122 #define TYPE_COUNT integer,intent(in) :: Count
2123 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2124 #define NF_TYPE NF_INT
2125 #define NF_ROUTINE NFMPI_GET_VARA_INT_ALL
2126 #define LENGTH min(Count,Len1)
2127 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))==1
2128 #include "ext_pnc_get_var_td.code"
2129 end subroutine ext_pnc_get_var_td_logical
2131 subroutine ext_pnc_get_var_ti_char(DataHandle,Element,Var,Data,Status)
2136 #undef TYPE_OUTCOUNT
2140 #define ROUTINE_TYPE 'CHAR'
2141 #define TYPE_DATA character*(*) ,intent(out) :: Data
2143 #define TYPE_COUNT integer :: Count = 1
2144 #define TYPE_OUTCOUNT
2145 #define NF_TYPE NF_CHAR
2146 #define NF_ROUTINE NFMPI_GET_ATT_TEXT
2149 #include "ext_pnc_get_var_ti.code"
2151 end subroutine ext_pnc_get_var_ti_char
2153 subroutine ext_pnc_get_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
2158 #undef TYPE_OUTCOUNT
2162 #define ROUTINE_TYPE 'CHAR'
2163 #define TYPE_DATA character*(*) ,intent(out) :: Data
2164 #define TYPE_BUFFER character (80)
2165 #define TYPE_COUNT integer :: Count = 1
2166 #define TYPE_OUTCOUNT
2167 #define NF_TYPE NF_CHAR
2168 #define NF_ROUTINE NFMPI_GET_VARA_TEXT_ALL
2171 #include "ext_pnc_get_var_td.code"
2173 end subroutine ext_pnc_get_var_td_char
2175 subroutine ext_pnc_put_dom_td_real(DataHandle,Element,DateStr,Data,Count,Status)
2176 integer ,intent(in) :: DataHandle
2177 character*(*) ,intent(in) :: Element
2178 character*(*) ,intent(in) :: DateStr
2179 real ,intent(in) :: Data(*)
2180 integer ,intent(in) :: Count
2181 integer ,intent(out) :: Status
2183 call ext_pnc_put_var_td_real(DataHandle,Element,DateStr, &
2184 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2186 end subroutine ext_pnc_put_dom_td_real
2188 subroutine ext_pnc_put_dom_td_integer(DataHandle,Element,DateStr,Data,Count,Status)
2189 integer ,intent(in) :: DataHandle
2190 character*(*) ,intent(in) :: Element
2191 character*(*) ,intent(in) :: DateStr
2192 integer ,intent(in) :: Data(*)
2193 integer ,intent(in) :: Count
2194 integer ,intent(out) :: Status
2196 call ext_pnc_put_var_td_integer(DataHandle,Element,DateStr, &
2197 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2199 end subroutine ext_pnc_put_dom_td_integer
2201 subroutine ext_pnc_put_dom_td_double(DataHandle,Element,DateStr,Data,Count,Status)
2202 integer ,intent(in) :: DataHandle
2203 character*(*) ,intent(in) :: Element
2204 character*(*) ,intent(in) :: DateStr
2205 real*8 ,intent(in) :: Data(*)
2206 integer ,intent(in) :: Count
2207 integer ,intent(out) :: Status
2209 call ext_pnc_put_var_td_double(DataHandle,Element,DateStr, &
2210 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2212 end subroutine ext_pnc_put_dom_td_double
2214 subroutine ext_pnc_put_dom_td_logical(DataHandle,Element,DateStr,Data,Count,Status)
2215 integer ,intent(in) :: DataHandle
2216 character*(*) ,intent(in) :: Element
2217 character*(*) ,intent(in) :: DateStr
2218 logical ,intent(in) :: Data(*)
2219 integer ,intent(in) :: Count
2220 integer ,intent(out) :: Status
2222 call ext_pnc_put_var_td_logical(DataHandle,Element,DateStr, &
2223 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2225 end subroutine ext_pnc_put_dom_td_logical
2227 subroutine ext_pnc_put_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2228 integer ,intent(in) :: DataHandle
2229 character*(*) ,intent(in) :: Element
2230 character*(*) ,intent(in) :: DateStr
2231 character*(*) ,intent(in) :: Data
2232 integer ,intent(out) :: Status
2234 call ext_pnc_put_var_td_char(DataHandle,Element,DateStr, &
2235 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2237 end subroutine ext_pnc_put_dom_td_char
2239 subroutine ext_pnc_get_dom_td_real(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2240 integer ,intent(in) :: DataHandle
2241 character*(*) ,intent(in) :: Element
2242 character*(*) ,intent(in) :: DateStr
2243 real ,intent(out) :: Data(*)
2244 integer ,intent(in) :: Count
2245 integer ,intent(out) :: OutCount
2246 integer ,intent(out) :: Status
2247 call ext_pnc_get_var_td_real(DataHandle,Element,DateStr, &
2248 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2250 end subroutine ext_pnc_get_dom_td_real
2252 subroutine ext_pnc_get_dom_td_integer(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2253 integer ,intent(in) :: DataHandle
2254 character*(*) ,intent(in) :: Element
2255 character*(*) ,intent(in) :: DateStr
2256 integer ,intent(out) :: Data(*)
2257 integer ,intent(in) :: Count
2258 integer ,intent(out) :: OutCount
2259 integer ,intent(out) :: Status
2260 call ext_pnc_get_var_td_integer(DataHandle,Element,DateStr, &
2261 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2263 end subroutine ext_pnc_get_dom_td_integer
2265 subroutine ext_pnc_get_dom_td_double(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2266 integer ,intent(in) :: DataHandle
2267 character*(*) ,intent(in) :: Element
2268 character*(*) ,intent(in) :: DateStr
2269 real*8 ,intent(out) :: Data(*)
2270 integer ,intent(in) :: Count
2271 integer ,intent(out) :: OutCount
2272 integer ,intent(out) :: Status
2273 call ext_pnc_get_var_td_double(DataHandle,Element,DateStr, &
2274 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2276 end subroutine ext_pnc_get_dom_td_double
2278 subroutine ext_pnc_get_dom_td_logical(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2279 integer ,intent(in) :: DataHandle
2280 character*(*) ,intent(in) :: Element
2281 character*(*) ,intent(in) :: DateStr
2282 logical ,intent(out) :: Data(*)
2283 integer ,intent(in) :: Count
2284 integer ,intent(out) :: OutCount
2285 integer ,intent(out) :: Status
2286 call ext_pnc_get_var_td_logical(DataHandle,Element,DateStr, &
2287 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2289 end subroutine ext_pnc_get_dom_td_logical
2291 subroutine ext_pnc_get_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2292 integer ,intent(in) :: DataHandle
2293 character*(*) ,intent(in) :: Element
2294 character*(*) ,intent(in) :: DateStr
2295 character*(*) ,intent(out) :: Data
2296 integer ,intent(out) :: Status
2297 call ext_pnc_get_var_td_char(DataHandle,Element,DateStr, &
2298 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2300 end subroutine ext_pnc_get_dom_td_char
2303 subroutine ext_pnc_write_field(DataHandle,DateStr,Var,Field,FieldType,Comm, &
2304 IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2305 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2307 use ext_pnc_support_routines
2309 include 'wrf_status_codes.h'
2310 # include "pnetcdf.inc"
2311 integer ,intent(in) :: DataHandle
2312 character*(*) ,intent(in) :: DateStr
2313 character*(*) ,intent(in) :: Var
2314 integer ,intent(inout) :: Field(*)
2315 integer ,intent(in) :: FieldType
2316 integer ,intent(inout) :: Comm
2317 integer ,intent(inout) :: IOComm
2318 integer ,intent(in) :: DomainDesc
2319 character*(*) ,intent(in) :: MemoryOrdIn
2320 character*(*) ,intent(in) :: Stagger ! Dummy for now
2321 character*(*) ,dimension(*) ,intent(in) :: DimNames
2322 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2323 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2324 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2325 integer ,intent(out) :: Status
2326 character (3) :: MemoryOrder
2327 type(wrf_data_handle) ,pointer :: DH
2330 character (VarNameLen) :: VarName
2331 character (3) :: MemO
2332 character (3) :: UCMemO
2334 integer ,dimension(NVarDims) :: Length_global, Length_native
2335 integer ,dimension(NVarDims) :: Length
2336 integer ,dimension(NVarDims) :: VDimIDs
2337 character(80),dimension(NVarDims) :: RODimNames
2338 integer ,dimension(NVarDims) :: StoredStart
2339 integer ,dimension(:,:,:,:),allocatable :: XField
2343 integer :: i1,i2,j1,j2,k1,k2
2344 integer :: x1,x2,y1,y2,z1,z2
2345 integer :: p1,p2,q1,q2,r1,r2
2346 integer :: l1,l2,m1,m2,n1,n2
2349 character (80) :: NullName
2352 ! Local, possibly adjusted, copies of MemoryStart and MemoryEnd
2353 integer ,dimension(NVarDims) :: lMemoryStart, lMemoryEnd
2354 MemoryOrder = trim(adjustl(MemoryOrdIn))
2356 call GetDim(MemoryOrder,NDim,Status)
2357 if(Status /= WRF_NO_ERR) then
2358 write(msg,*) 'Warning BAD MEMORY ORDER |',MemoryOrder,'| in ',__FILE__,', line', __LINE__
2359 call wrf_debug ( WARN , TRIM(msg))
2362 call DateCheck(DateStr,Status)
2363 if(Status /= WRF_NO_ERR) then
2364 write(msg,*) 'Warning DATE STRING ERROR |',DateStr,'| in ',__FILE__,', line', __LINE__
2365 call wrf_debug ( WARN , TRIM(msg))
2369 call GetDH(DataHandle,DH,Status)
2370 if(Status /= WRF_NO_ERR) then
2371 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2372 call wrf_debug ( WARN , TRIM(msg))
2377 write(msg,*)'ext_pnc_write_field: called for ',TRIM(Var)
2378 CALL wrf_debug( 100, msg )
2381 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2382 Length_native(1:NDim) = Length(1:NDim)
2383 Length_global(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
2385 call ExtOrder(MemoryOrder,Length,Status)
2386 call ExtOrder(MemoryOrder,Length_global,Status)
2388 call ExtOrderStr(MemoryOrder,DimNames,RODimNames,Status)
2390 ! Magic number to identify call from IO server when doing quilting
2391 quilting = (MemoryStart(1) == -998899 .AND. MemoryEnd(1) == -998899)
2393 lMemoryStart(1:NDim) = 1
2394 lMemoryEnd(1:NDim) = Length(1:NDim)
2396 lMemoryStart(1:NDim) = MemoryStart(1:NDim)
2397 lMemoryEnd(1:NDim) = MemoryEnd(1:NDim)
2400 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2401 Status = WRF_WARN_FILE_NOT_OPENED
2402 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2403 call wrf_debug ( WARN , TRIM(msg))
2404 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
2405 Status = WRF_WARN_WRITE_RONLY_FILE
2406 write(msg,*) 'Warning WRITE READ ONLY FILE in ',__FILE__,', line', __LINE__
2407 call wrf_debug ( WARN , TRIM(msg))
2408 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2410 if(DH%VarNames(NVar) == VarName ) then
2411 Status = WRF_WARN_2DRYRUNS_1VARIABLE
2412 write(msg,*) 'Warning 2 DRYRUNS 1 VARIABLE (',TRIM(VarName),') in ',__FILE__,', line', __LINE__
2413 call wrf_debug ( WARN , TRIM(msg))
2415 elseif(DH%VarNames(NVar) == NO_NAME) then
2416 DH%VarNames(NVar) = VarName
2419 elseif(NVar == MaxVars) then
2420 Status = WRF_WARN_TOO_MANY_VARIABLES
2421 write(msg,*) 'Warning TOO MANY VARIABLES in ',__FILE__,', line', __LINE__
2422 call wrf_debug ( WARN , TRIM(msg))
2427 if(RODimNames(j) == NullName .or. RODimNames(j) == '') then
2429 if(DH%DimLengths(i) == Length_global(j)) then
2431 elseif(DH%DimLengths(i) == NO_DIM) then
2432 stat = NFMPI_DEF_DIM(NCID,DH%DimNames(i),i2offset(Length_global(j)),DH%DimIDs(i))
2433 call netcdf_err(stat,Status)
2434 if(Status /= WRF_NO_ERR) then
2435 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2436 call wrf_debug ( WARN , TRIM(msg))
2439 DH%DimLengths(i) = Length_global(j)
2441 elseif(i == MaxDims) then
2442 Status = WRF_WARN_TOO_MANY_DIMS
2443 write(msg,*) 'Warning TOO MANY DIMENSIONS (',i,') in (',TRIM(VarName),') ',__FILE__,', line', __LINE__
2444 call wrf_debug ( WARN , TRIM(msg))
2448 else !look for input name and check if already defined
2451 if (DH%DimNames(i) == RODimNames(j)) then
2452 if (DH%DimLengths(i) == Length_global(j)) then
2456 Status = WRF_WARN_DIMNAME_REDEFINED
2457 write(msg,*) 'Warning DIM ',i,', NAME ',TRIM(DH%DimNames(i)),' REDEFINED by var ', &
2458 TRIM(Var),' ',DH%DimLengths(i),Length_global(j) ,' in ', __FILE__ ,' line', __LINE__
2459 call wrf_debug ( WARN , TRIM(msg))
2466 if (DH%DimLengths(i) == NO_DIM) then
2467 DH%DimNames(i) = RODimNames(j)
2468 stat = NFMPI_DEF_DIM(NCID,DH%DimNames(i),i2offset(Length_global(j)),DH%DimIDs(i))
2469 call netcdf_err(stat,Status)
2470 if(Status /= WRF_NO_ERR) then
2471 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2472 call wrf_debug ( WARN , TRIM(msg))
2475 DH%DimLengths(i) = Length_global(j)
2477 elseif(i == MaxDims) then
2478 Status = WRF_WARN_TOO_MANY_DIMS
2479 write(msg,*) 'Warning TOO MANY DIMENSIONS in ',__FILE__,', line', __LINE__
2480 call wrf_debug ( WARN , TRIM(msg))
2486 VDimIDs(j) = DH%DimIDs(i)
2487 DH%VarDimLens(j,NVar) = Length_global(j)
2489 VDimIDs(NDim+1) = DH%DimUnlimID
2490 select case (FieldType)
2500 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
2501 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
2502 call wrf_debug ( WARN , TRIM(msg))
2507 stat = NFMPI_DEF_VAR(NCID,VarName,XType,NDim+1,VDimIDs,VarID)
2508 call netcdf_err(stat,Status)
2509 if(Status /= WRF_NO_ERR) then
2510 write(msg,*) 'ext_pnc_write_field: NetCDF error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2511 call wrf_debug ( WARN , TRIM(msg))
2514 DH%VarIDs(NVar) = VarID
2515 stat = NFMPI_PUT_ATT_INT(NCID,VarID,'FieldType',NF_INT,i2offset(1),FieldType)
2516 call netcdf_err(stat,Status)
2517 if(Status /= WRF_NO_ERR) then
2518 write(msg,*) 'ext_pnc_write_field: NetCDF error in ',__FILE__,', line', __LINE__
2519 call wrf_debug ( WARN , TRIM(msg))
2522 call reorder(MemoryOrder,MemO)
2523 call uppercase(MemO,UCMemO)
2524 stat = NFMPI_PUT_ATT_TEXT(NCID,VarID,'MemoryOrder',i2offset(3),UCMemO)
2525 call netcdf_err(stat,Status)
2526 if(Status /= WRF_NO_ERR) then
2527 write(msg,*) 'ext_pnc_write_field: NetCDF error in ',__FILE__,', line', __LINE__
2528 call wrf_debug ( WARN , TRIM(msg))
2531 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
2532 do NVar=1,DH%NumVars
2533 if(DH%VarNames(NVar) == VarName) then
2535 elseif(NVar == DH%NumVars) then
2536 Status = WRF_WARN_VAR_NF
2537 write(msg,*) 'Warning VARIABLE NOT FOUND in ',__FILE__,', line', __LINE__
2538 call wrf_debug ( WARN , TRIM(msg))
2542 VarID = DH%VarIDs(NVar)
2544 if(Length_global(j) /= DH%VarDimLens(j,NVar) .AND. DH%FileStatus /= WRF_FILE_OPENED_FOR_UPDATE ) then
2545 Status = WRF_WARN_WRTLEN_NE_DRRUNLEN
2546 write(msg,*) 'Warning LENGTH != DRY RUN LENGTH for |', &
2547 VarName,'| dim ',j,' in ',__FILE__,', line', __LINE__
2548 call wrf_debug ( WARN , TRIM(msg))
2549 write(msg,*) ' LENGTH ',Length_global(j),' DRY RUN LENGTH ',DH%VarDimLens(j,NVar)
2550 call wrf_debug ( WARN , TRIM(msg))
2552 !jm 061024 elseif(PatchStart(j) < MemoryStart(j)) then
2553 !jm elseif(DomainStart(j) < MemoryStart(j)) then
2554 elseif(PatchStart(j) < lMemoryStart(j)) then
2555 Status = WRF_WARN_DIMENSION_ERROR
2556 write(msg,*) 'Warning DIMENSION ERROR for |',VarName, &
2557 '| in ',__FILE__,', line', __LINE__
2558 call wrf_debug ( WARN , TRIM(msg))
2563 call GetIndices(NDim,lMemoryStart,lMemoryEnd,l1,l2,m1,m2,n1,n2)
2564 call GetIndices(NDim,StoredStart,Length ,x1,x2,y1,y2,z1,z2)
2565 call GetIndices(NDim,StoredStart,Length_native ,p1,p2,q1,q2,r1,r2)
2566 call GetIndices(NDim,PatchStart, PatchEnd ,i1,i2,j1,j2,k1,k2)
2568 if(FieldType == WRF_DOUBLE) di=2
2569 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
2571 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
2572 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
2573 call wrf_debug ( FATAL , TRIM(msg))
2578 WRITE(msg,*) 'ARPDBG: MemoryStart = ',lMemoryStart(1:NDim)
2579 CALL wrf_message(msg)
2580 WRITE(msg,*) 'ARPDBG: lMemoryEnd = ',lMemoryEnd(1:NDim)
2581 CALL wrf_message(msg)
2582 WRITE(msg,*) 'ARPDBG: Length = ',Length(1:NDim)
2583 CALL wrf_message(msg)
2587 ! Don't pass in PatchStart and PatchEnd here since we want to
2588 ! take transpose of whole patch of data which has been sent to
2589 ! the IO server and passed down to us.
2590 ! JM: the field and patch dimensions must be reordered or xpose is a noop
2591 call Transpose('write',MemoryOrder,di, Field,p1,p2,q1,q2,r1,r2 &
2592 ,XField,x1,x2,y1,y2,z1,z2 &
2593 ,p1,p2,q1,q2,r1,r2 )
2595 call Transpose('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2596 ,XField,x1,x2,y1,y2,z1,z2 &
2597 ,i1,i2,j1,j2,k1,k2 )
2600 StoredStart(1:NDim) = PatchStart(1:NDim)
2601 call ExtOrder(MemoryOrder,StoredStart,Status)
2602 call FieldIO('write',DataHandle,DateStr,StoredStart,Length,MemoryOrder, &
2603 FieldType,NCID,VarID,XField,Status)
2604 if(Status /= WRF_NO_ERR) then
2605 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2606 call wrf_debug ( WARN , TRIM(msg))
2609 deallocate(XField, STAT=stat)
2611 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
2612 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
2613 call wrf_debug ( FATAL , TRIM(msg))
2617 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2618 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2619 call wrf_debug ( FATAL , TRIM(msg))
2621 DH%first_operation = .FALSE.
2623 end subroutine ext_pnc_write_field
2625 subroutine ext_pnc_read_field(DataHandle,DateStr,Var,Field,FieldType,Comm, &
2626 IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2627 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2629 use ext_pnc_support_routines
2631 include 'wrf_status_codes.h'
2632 # include "pnetcdf.inc"
2633 integer ,intent(in) :: DataHandle
2634 character*(*) ,intent(in) :: DateStr
2635 character*(*) ,intent(in) :: Var
2636 integer ,intent(out) :: Field(*)
2637 integer ,intent(in) :: FieldType
2638 integer ,intent(inout) :: Comm
2639 integer ,intent(inout) :: IOComm
2640 integer ,intent(in) :: DomainDesc
2641 character*(*) ,intent(in) :: MemoryOrdIn
2642 character*(*) ,intent(in) :: Stagger ! Dummy for now
2643 character*(*) , dimension (*) ,intent(in) :: DimNames
2644 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2645 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2646 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2647 integer ,intent(out) :: Status
2648 character (3) :: MemoryOrder
2649 character (NF_MAX_NAME) :: dimname
2650 type(wrf_data_handle) ,pointer :: DH
2653 character (VarNameLen) :: VarName
2655 integer ,dimension(NVarDims) :: VCount
2656 integer ,dimension(NVarDims) :: VStart
2657 integer ,dimension(NVarDims) :: Length
2658 integer ,dimension(NVarDims) :: VDimIDs
2659 integer ,dimension(NVarDims) :: MemS
2660 integer ,dimension(NVarDims) :: MemE
2661 integer ,dimension(NVarDims) :: StoredStart
2662 integer ,dimension(NVarDims) :: StoredLen
2663 integer(KIND=MPI_OFFSET_KIND) ,dimension(NVarDims) :: StoredLen_okind
2664 integer ,dimension(:,:,:,:) ,allocatable :: XField
2667 integer :: i1,i2,j1,j2,k1,k2
2668 integer :: x1,x2,y1,y2,z1,z2
2669 integer :: l1,l2,m1,m2,n1,n2
2670 character (VarNameLen) :: Name
2672 integer :: StoredDim
2674 integer(KIND=MPI_OFFSET_KIND) :: Len
2679 MemoryOrder = trim(adjustl(MemoryOrdIn))
2680 call GetDim(MemoryOrder,NDim,Status)
2681 if(Status /= WRF_NO_ERR) then
2682 write(msg,*) 'Warning BAD MEMORY ORDER |',TRIM(MemoryOrder),'| for |', &
2683 TRIM(Var),'| in ext_pnc_read_field ',__FILE__,', line', __LINE__
2684 call wrf_debug ( WARN , TRIM(msg))
2687 call DateCheck(DateStr,Status)
2688 if(Status /= WRF_NO_ERR) then
2689 write(msg,*) 'Warning DATE STRING ERROR |',TRIM(DateStr),'| for |',TRIM(Var), &
2690 '| in ext_pnc_read_field ',__FILE__,', line', __LINE__
2691 call wrf_debug ( WARN , TRIM(msg))
2695 call GetDH(DataHandle,DH,Status)
2696 if(Status /= WRF_NO_ERR) then
2697 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_read_field ',__FILE__,', line', __LINE__
2698 call wrf_debug ( WARN , TRIM(msg))
2701 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2702 Status = WRF_WARN_FILE_NOT_OPENED
2703 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2704 call wrf_debug ( WARN , TRIM(msg))
2705 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2706 ! jm it is okay to have a dry run read. means read is called between ofrb and ofrc. Just return.
2707 ! Status = WRF_WARN_DRYRUN_READ
2708 ! write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
2709 ! call wrf_debug ( WARN , TRIM(msg))
2712 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
2713 Status = WRF_WARN_READ_WONLY_FILE
2714 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
2715 call wrf_debug ( WARN , TRIM(msg))
2716 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
2719 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2720 StoredStart(1:NDim) = PatchStart(1:NDim)
2722 call ExtOrder(MemoryOrder,Length,Status)
2724 stat = NFMPI_INQ_VARID(NCID,VarName,VarID)
2725 call netcdf_err(stat,Status)
2726 if(Status /= WRF_NO_ERR) then
2727 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__,' Varname ',Varname
2728 call wrf_debug ( WARN , TRIM(msg))
2731 stat = NFMPI_INQ_VAR(NCID,VarID,Name,XType,StoredDim,VDimIDs,NAtts)
2732 call netcdf_err(stat,Status)
2733 if(Status /= WRF_NO_ERR) then
2734 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2735 call wrf_debug ( WARN , TRIM(msg))
2738 stat = NFMPI_GET_ATT_INT(NCID,VarID,'FieldType',FType)
2739 call netcdf_err(stat,Status)
2740 if(Status /= WRF_NO_ERR) then
2741 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2742 call wrf_debug ( WARN , TRIM(msg))
2745 ! allow coercion between double and single prec real
2746 !jm if(FieldType /= Ftype) then
2747 if( (FieldType == WRF_REAL .OR. FieldType == WRF_DOUBLE) ) then
2748 if ( .NOT. (Ftype == WRF_REAL .OR. Ftype == WRF_DOUBLE )) then
2749 Status = WRF_WARN_TYPE_MISMATCH
2750 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
2751 call wrf_debug ( WARN , TRIM(msg))
2754 else if(FieldType /= Ftype) then
2755 Status = WRF_WARN_TYPE_MISMATCH
2756 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
2757 call wrf_debug ( WARN , TRIM(msg))
2760 select case (FieldType)
2762 ! allow coercion between double and single prec real
2763 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
2764 Status = WRF_WARN_TYPE_MISMATCH
2765 write(msg,*) 'Warning REAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
2768 ! allow coercion between double and single prec real
2769 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
2770 Status = WRF_WARN_TYPE_MISMATCH
2771 write(msg,*) 'Warning DOUBLE TYPE MISMATCH in ',__FILE__,', line', __LINE__
2774 if(XType /= NF_INT) then
2775 Status = WRF_WARN_TYPE_MISMATCH
2776 write(msg,*) 'Warning INTEGER TYPE MISMATCH in ',__FILE__,', line', __LINE__
2779 if(XType /= NF_INT) then
2780 Status = WRF_WARN_TYPE_MISMATCH
2781 write(msg,*) 'Warning LOGICAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
2784 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
2785 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
2787 if(Status /= WRF_NO_ERR) then
2788 call wrf_debug ( WARN , TRIM(msg))
2791 ! NDim=0 for scalars. Handle read of old NDim=1 files. TBH: 20060502
2792 IF ( ( NDim == 0 ) .AND. ( StoredDim == 2 ) ) THEN
2793 stat = NFMPI_INQ_DIMNAME(NCID,VDimIDs(1),dimname)
2794 call netcdf_err(stat,Status)
2795 if(Status /= WRF_NO_ERR) then
2796 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2797 call wrf_debug ( WARN , TRIM(msg))
2800 IF ( dimname(1:10) == 'ext_scalar' ) THEN
2805 if(StoredDim /= NDim+1) then
2806 Status = WRF_ERR_FATAL_BAD_VARIABLE_DIM
2807 write(msg,*) 'Fatal error BAD VARIABLE DIMENSION in ext_pnc_read_field ',TRIM(Var),TRIM(DateStr)
2808 call wrf_debug ( FATAL , msg)
2809 write(msg,*) ' StoredDim ', StoredDim, ' .NE. NDim+1 ', NDim+1
2810 call wrf_debug ( FATAL , msg)
2814 stat = NFMPI_INQ_DIMLEN(NCID,VDimIDs(j),StoredLen_okind(j))
2815 call netcdf_err(stat,Status)
2816 if(Status /= WRF_NO_ERR) then
2817 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2818 call wrf_debug ( WARN , TRIM(msg))
2821 StoredLen(j) = StoredLen_okind(j)
2822 if(Length(j) > StoredLen(j)) then
2823 Status = WRF_WARN_READ_PAST_EOF
2824 write(msg,*) 'Warning READ PAST EOF in ext_pnc_read_field of ',TRIM(Var),Length(j),'>',StoredLen(j)
2825 call wrf_debug ( WARN , TRIM(msg))
2827 elseif(Length(j) <= 0) then
2828 Status = WRF_WARN_ZERO_LENGTH_READ
2829 write(msg,*) 'Warning ZERO LENGTH READ in ',__FILE__,', line', __LINE__
2830 call wrf_debug ( WARN , TRIM(msg))
2836 call GetIndices(NDim,MemoryStart,MemoryEnd,l1,l2,m1,m2,n1,n2)
2837 call GetIndices(NDim,StoredStart,Length,x1,x2,y1,y2,z1,z2)
2838 !jm call GetIndices(NDim,DomainStart,DomainEnd,i1,i2,j1,j2,k1,k2)
2839 call GetIndices(NDim,PatchStart,PatchEnd,i1,i2,j1,j2,k1,k2)
2841 StoredStart(1:NDim) = PatchStart(1:NDim)
2842 call ExtOrder(MemoryOrder,StoredStart,Status)
2845 if(FieldType == WRF_DOUBLE) di=2
2846 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
2848 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
2849 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
2850 call wrf_debug ( FATAL , msg)
2853 call FieldIO('read',DataHandle,DateStr,StoredStart,Length,MemoryOrder, &
2854 FieldType,NCID,VarID,XField,Status)
2855 if(Status /= WRF_NO_ERR) then
2856 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2857 call wrf_debug ( WARN , TRIM(msg))
2860 call Transpose('read',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2861 ,XField,x1,x2,y1,y2,z1,z2 &
2862 ,i1,i2,j1,j2,k1,k2 )
2863 deallocate(XField, STAT=stat)
2865 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
2866 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
2867 call wrf_debug ( FATAL , msg)
2871 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2872 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2873 call wrf_debug ( FATAL , msg)
2875 DH%first_operation = .FALSE.
2877 end subroutine ext_pnc_read_field
2879 subroutine ext_pnc_inquire_opened( DataHandle, FileName , FileStatus, Status )
2881 use ext_pnc_support_routines
2883 include 'wrf_status_codes.h'
2884 integer ,intent(in) :: DataHandle
2885 character*(*) ,intent(in) :: FileName
2886 integer ,intent(out) :: FileStatus
2887 integer ,intent(out) :: Status
2888 type(wrf_data_handle) ,pointer :: DH
2890 call GetDH(DataHandle,DH,Status)
2891 if(Status /= WRF_NO_ERR) then
2892 FileStatus = WRF_FILE_NOT_OPENED
2895 if(FileName /= DH%FileName) then
2896 FileStatus = WRF_FILE_NOT_OPENED
2898 FileStatus = DH%FileStatus
2902 end subroutine ext_pnc_inquire_opened
2904 subroutine ext_pnc_inquire_filename( Datahandle, FileName, FileStatus, Status )
2906 use ext_pnc_support_routines
2908 include 'wrf_status_codes.h'
2909 integer ,intent(in) :: DataHandle
2910 character*(*) ,intent(out) :: FileName
2911 integer ,intent(out) :: FileStatus
2912 integer ,intent(out) :: Status
2913 type(wrf_data_handle) ,pointer :: DH
2914 FileStatus = WRF_FILE_NOT_OPENED
2915 call GetDH(DataHandle,DH,Status)
2916 if(Status /= WRF_NO_ERR) then
2917 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2918 call wrf_debug ( WARN , TRIM(msg))
2921 FileName = DH%FileName
2922 FileStatus = DH%FileStatus
2925 end subroutine ext_pnc_inquire_filename
2927 subroutine ext_pnc_set_time(DataHandle, DateStr, Status)
2929 use ext_pnc_support_routines
2931 include 'wrf_status_codes.h'
2932 integer ,intent(in) :: DataHandle
2933 character*(*) ,intent(in) :: DateStr
2934 integer ,intent(out) :: Status
2935 type(wrf_data_handle) ,pointer :: DH
2938 call DateCheck(DateStr,Status)
2939 if(Status /= WRF_NO_ERR) then
2940 write(msg,*) 'Warning DATE STRING ERROR in ',__FILE__,', line', __LINE__
2941 call wrf_debug ( WARN , TRIM(msg))
2944 call GetDH(DataHandle,DH,Status)
2945 if(Status /= WRF_NO_ERR) then
2946 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2947 call wrf_debug ( WARN , TRIM(msg))
2950 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2951 Status = WRF_WARN_FILE_NOT_OPENED
2952 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2953 call wrf_debug ( WARN , TRIM(msg))
2954 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2955 Status = WRF_WARN_FILE_NOT_COMMITTED
2956 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
2957 call wrf_debug ( WARN , TRIM(msg))
2958 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
2959 Status = WRF_WARN_READ_WONLY_FILE
2960 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
2961 call wrf_debug ( WARN , TRIM(msg))
2962 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
2964 if(DH%Times(i)==DateStr) then
2968 if(i==MaxTimes) then
2969 Status = WRF_WARN_TIME_NF
2973 DH%CurrentVariable = 0
2976 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2977 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2978 call wrf_debug ( FATAL , msg)
2981 end subroutine ext_pnc_set_time
2983 subroutine ext_pnc_get_next_time(DataHandle, DateStr, Status)
2985 use ext_pnc_support_routines
2987 include 'wrf_status_codes.h'
2988 integer ,intent(in) :: DataHandle
2989 character*(*) ,intent(out) :: DateStr
2990 integer ,intent(out) :: Status
2991 type(wrf_data_handle) ,pointer :: DH
2993 call GetDH(DataHandle,DH,Status)
2994 if(Status /= WRF_NO_ERR) then
2995 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2996 call wrf_debug ( WARN , TRIM(msg))
2999 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3000 Status = WRF_WARN_FILE_NOT_OPENED
3001 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3002 call wrf_debug ( WARN , TRIM(msg))
3003 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3004 Status = WRF_WARN_DRYRUN_READ
3005 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3006 call wrf_debug ( WARN , TRIM(msg))
3007 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3008 Status = WRF_WARN_READ_WONLY_FILE
3009 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3010 call wrf_debug ( WARN , TRIM(msg))
3011 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
3012 if(DH%CurrentTime >= DH%NumberTimes) then
3013 write(msg,*) 'Warning ext_pnc_get_next_time: DH%CurrentTime >= DH%NumberTimes ',DH%CurrentTime,DH%NumberTimes
3014 call wrf_debug ( WARN , TRIM(msg))
3015 Status = WRF_WARN_TIME_EOF
3018 DH%CurrentTime = DH%CurrentTime +1
3019 DateStr = DH%Times(DH%CurrentTime)
3020 DH%CurrentVariable = 0
3023 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3024 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3025 call wrf_debug ( FATAL , msg)
3028 end subroutine ext_pnc_get_next_time
3030 subroutine ext_pnc_get_previous_time(DataHandle, DateStr, Status)
3032 use ext_pnc_support_routines
3034 include 'wrf_status_codes.h'
3035 integer ,intent(in) :: DataHandle
3036 character*(*) ,intent(out) :: DateStr
3037 integer ,intent(out) :: Status
3038 type(wrf_data_handle) ,pointer :: DH
3040 call GetDH(DataHandle,DH,Status)
3041 if(Status /= WRF_NO_ERR) then
3042 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3043 call wrf_debug ( WARN , TRIM(msg))
3046 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3047 Status = WRF_WARN_FILE_NOT_OPENED
3048 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3049 call wrf_debug ( WARN , TRIM(msg))
3050 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3051 Status = WRF_WARN_DRYRUN_READ
3052 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3053 call wrf_debug ( WARN , TRIM(msg))
3054 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3055 Status = WRF_WARN_READ_WONLY_FILE
3056 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3057 call wrf_debug ( WARN , TRIM(msg))
3058 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
3059 if(DH%CurrentTime.GT.0) then
3060 DH%CurrentTime = DH%CurrentTime -1
3062 DateStr = DH%Times(DH%CurrentTime)
3063 DH%CurrentVariable = 0
3066 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3067 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3068 call wrf_debug ( FATAL , msg)
3071 end subroutine ext_pnc_get_previous_time
3073 subroutine ext_pnc_get_next_var(DataHandle, VarName, Status)
3075 use ext_pnc_support_routines
3077 include 'wrf_status_codes.h'
3078 # include "pnetcdf.inc"
3079 integer ,intent(in) :: DataHandle
3080 character*(*) ,intent(out) :: VarName
3081 integer ,intent(out) :: Status
3082 type(wrf_data_handle) ,pointer :: DH
3084 character (80) :: Name
3086 call GetDH(DataHandle,DH,Status)
3087 if(Status /= WRF_NO_ERR) then
3088 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3089 call wrf_debug ( WARN , TRIM(msg))
3092 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3093 Status = WRF_WARN_FILE_NOT_OPENED
3094 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3095 call wrf_debug ( WARN , TRIM(msg))
3096 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3097 Status = WRF_WARN_DRYRUN_READ
3098 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3099 call wrf_debug ( WARN , TRIM(msg))
3100 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3101 Status = WRF_WARN_READ_WONLY_FILE
3102 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3103 call wrf_debug ( WARN , TRIM(msg))
3104 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3106 DH%CurrentVariable = DH%CurrentVariable +1
3107 if(DH%CurrentVariable > DH%NumVars) then
3108 Status = WRF_WARN_VAR_EOF
3111 VarName = DH%VarNames(DH%CurrentVariable)
3114 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3115 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3116 call wrf_debug ( FATAL , msg)
3119 end subroutine ext_pnc_get_next_var
3121 subroutine ext_pnc_end_of_frame(DataHandle, Status)
3123 use ext_pnc_support_routines
3125 # include "pnetcdf.inc"
3126 include 'wrf_status_codes.h'
3127 integer ,intent(in) :: DataHandle
3128 integer ,intent(out) :: Status
3129 type(wrf_data_handle) ,pointer :: DH
3131 call GetDH(DataHandle,DH,Status)
3133 end subroutine ext_pnc_end_of_frame
3135 ! NOTE: For scalar variables NDim is set to zero and DomainStart and
3136 ! NOTE: DomainEnd are left unmodified.
3137 subroutine ext_pnc_get_var_info(DataHandle,Name,NDim,MemoryOrder,Stagger,DomainStart,DomainEnd,WrfType,Status)
3139 use ext_pnc_support_routines
3141 # include "pnetcdf.inc"
3142 include 'wrf_status_codes.h'
3143 integer ,intent(in) :: DataHandle
3144 character*(*) ,intent(in) :: Name
3145 integer ,intent(out) :: NDim
3146 character*(*) ,intent(out) :: MemoryOrder
3147 character*(*) :: Stagger ! Dummy for now
3148 integer ,dimension(*) ,intent(out) :: DomainStart, DomainEnd
3149 integer ,intent(out) :: WrfType
3150 integer ,intent(out) :: Status
3151 type(wrf_data_handle) ,pointer :: DH
3153 integer ,dimension(NVarDims) :: VDimIDs
3158 call GetDH(DataHandle,DH,Status)
3159 if(Status /= WRF_NO_ERR) then
3160 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3161 call wrf_debug ( WARN , TRIM(msg))
3164 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3165 Status = WRF_WARN_FILE_NOT_OPENED
3166 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3167 call wrf_debug ( WARN , TRIM(msg))
3169 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3170 Status = WRF_WARN_DRYRUN_READ
3171 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3172 call wrf_debug ( WARN , TRIM(msg))
3174 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3175 Status = WRF_WARN_READ_WONLY_FILE
3176 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3177 call wrf_debug ( WARN , TRIM(msg))
3179 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3180 stat = NFMPI_INQ_VARID(DH%NCID,Name,VarID)
3181 call netcdf_err(stat,Status)
3182 if(Status /= WRF_NO_ERR) then
3183 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3184 call wrf_debug ( WARN , TRIM(msg))
3187 stat = NFMPI_INQ_VARTYPE(DH%NCID,VarID,XType)
3188 call netcdf_err(stat,Status)
3189 if(Status /= WRF_NO_ERR) then
3190 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3191 call wrf_debug ( WARN , TRIM(msg))
3194 stat = NFMPI_GET_ATT_INT(DH%NCID,VarID,'FieldType',WrfType)
3195 call netcdf_err(stat,Status)
3196 if(Status /= WRF_NO_ERR) then
3197 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3198 call wrf_debug ( WARN , TRIM(msg))
3203 Status = WRF_WARN_BAD_DATA_TYPE
3204 write(msg,*) 'Warning BYTE IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3205 call wrf_debug ( WARN , TRIM(msg))
3208 Status = WRF_WARN_BAD_DATA_TYPE
3209 write(msg,*) 'Warning CHAR IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3210 call wrf_debug ( WARN , TRIM(msg))
3213 Status = WRF_WARN_BAD_DATA_TYPE
3214 write(msg,*) 'Warning SHORT IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3215 call wrf_debug ( WARN , TRIM(msg))
3218 if(WrfType /= WRF_INTEGER .and. WrfType /= WRF_LOGICAL) then
3219 Status = WRF_WARN_BAD_DATA_TYPE
3220 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3221 call wrf_debug ( WARN , TRIM(msg))
3225 if(WrfType /= WRF_REAL) then
3226 Status = WRF_WARN_BAD_DATA_TYPE
3227 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3228 call wrf_debug ( WARN , TRIM(msg))
3232 if(WrfType /= WRF_DOUBLE) then
3233 Status = WRF_WARN_BAD_DATA_TYPE
3234 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3235 call wrf_debug ( WARN , TRIM(msg))
3239 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
3240 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
3241 call wrf_debug ( WARN , TRIM(msg))
3245 stat = NFMPI_GET_ATT_TEXT(DH%NCID,VarID,'MemoryOrder',MemoryOrder)
3246 call netcdf_err(stat,Status)
3247 if(Status /= WRF_NO_ERR) then
3248 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3249 call wrf_debug ( WARN , TRIM(msg))
3252 call GetDim(MemoryOrder,NDim,Status)
3253 if(Status /= WRF_NO_ERR) then
3254 write(msg,*) 'Warning BAD MEMORY ORDER ',TRIM(MemoryOrder),' in ',__FILE__,', line', __LINE__
3255 call wrf_debug ( WARN , TRIM(msg))
3258 stat = NFMPI_INQ_VARDIMID(DH%NCID,VarID,VDimIDs)
3259 call netcdf_err(stat,Status)
3260 if(Status /= WRF_NO_ERR) then
3261 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3262 call wrf_debug ( WARN , TRIM(msg))
3267 stat = NFMPI_INQ_DIMLEN(DH%NCID,VDimIDs(j),DomainEnd(j))
3268 call netcdf_err(stat,Status)
3269 if(Status /= WRF_NO_ERR) then
3270 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3271 call wrf_debug ( WARN , TRIM(msg))
3276 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3277 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3278 call wrf_debug ( FATAL , msg)
3281 end subroutine ext_pnc_get_var_info
3283 subroutine ext_pnc_warning_str( Code, ReturnString, Status)
3285 use ext_pnc_support_routines
3287 # include "pnetcdf.inc"
3288 include 'wrf_status_codes.h'
3290 integer , intent(in) ::Code
3291 character *(*), intent(out) :: ReturnString
3292 integer, intent(out) ::Status
3296 ReturnString='No error'
3300 ReturnString= 'File not found (or file is incomplete)'
3304 ReturnString='Metadata not found'
3308 ReturnString= 'Timestamp not found'
3312 ReturnString= 'No more timestamps'
3316 ReturnString= 'Variable not found'
3320 ReturnString= 'No more variables for the current time'
3324 ReturnString= 'Too many open files'
3328 ReturnString= 'Data type mismatch'
3332 ReturnString= 'Attempt to write read-only file'
3336 ReturnString= 'Attempt to read write-only file'
3340 ReturnString= 'Attempt to access unopened file'
3344 ReturnString= 'Attempt to do 2 trainings for 1 variable'
3348 ReturnString= 'Attempt to read past EOF'
3352 ReturnString= 'Bad data handle'
3356 ReturnString= 'Write length not equal to training length'
3360 ReturnString= 'More dimensions requested than training'
3364 ReturnString= 'Attempt to read more data than exists'
3368 ReturnString= 'Input dimensions inconsistent'
3372 ReturnString= 'Input MemoryOrder not recognized'
3376 ReturnString= 'A dimension name with 2 different lengths'
3380 ReturnString= 'String longer than provided storage'
3384 ReturnString= 'Function not supportable'
3388 ReturnString= 'Package implements this routine as NOOP'
3392 !netcdf-specific warning messages
3394 ReturnString= 'Bad data type'
3398 ReturnString= 'File not committed'
3402 ReturnString= 'File is opened for reading'
3406 ReturnString= 'Attempt to write metadata after open commit'
3410 ReturnString= 'I/O not initialized'
3414 ReturnString= 'Too many variables requested'
3418 ReturnString= 'Attempt to close file during a dry run'
3422 ReturnString= 'Date string not 19 characters in length'
3426 ReturnString= 'Attempt to read zero length words'
3430 ReturnString= 'Data type not found'
3434 ReturnString= 'Badly formatted date string'
3438 ReturnString= 'Attempt at read during a dry run'
3442 ReturnString= 'Attempt to get zero words'
3446 ReturnString= 'Attempt to put zero length words'
3450 ReturnString= 'NetCDF error'
3454 ReturnString= 'Requested length <= 1'
3458 ReturnString= 'More data available than requested'
3462 ReturnString= 'New date less than previous date'
3467 ReturnString= 'This warning code is not supported or handled directly by WRF and NetCDF. &
3468 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
3469 & to be calling a package-specific routine to return a message for this warning code.'
3474 end subroutine ext_pnc_warning_str
3477 !returns message string for all WRF and netCDF warning/error status codes
3478 !Other i/o packages must provide their own routines to return their own status messages
3479 subroutine ext_pnc_error_str( Code, ReturnString, Status)
3481 use ext_pnc_support_routines
3483 # include "pnetcdf.inc"
3484 include 'wrf_status_codes.h'
3486 integer , intent(in) ::Code
3487 character *(*), intent(out) :: ReturnString
3488 integer, intent(out) ::Status
3492 ReturnString= 'Allocation Error'
3496 ReturnString= 'Deallocation Error'
3500 ReturnString= 'Bad File Status'
3504 ReturnString= 'Variable on disk is not 3D'
3508 ReturnString= 'Metadata on disk is not 1D'
3512 ReturnString= 'Time dimension too small'
3516 ReturnString= 'This error code is not supported or handled directly by WRF and NetCDF. &
3517 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
3518 & to be calling a package-specific routine to return a message for this error code.'
3523 end subroutine ext_pnc_error_str
3527 subroutine ext_pnc_end_independent_mode(DataHandle, Status)
3529 use ext_pnc_support_routines
3530 include 'wrf_status_codes.h'
3531 # include "pnetcdf.inc"
3532 integer ,intent(in) :: DataHandle
3533 integer ,intent(out) :: Status
3534 type(wrf_data_handle) ,pointer :: DH
3537 DH => WrfDataHandles(DataHandle)
3538 DH%Collective = .TRUE.
3539 stat = NFMPI_END_INDEP_DATA(DH%NCID)
3541 call netcdf_err(stat,Status)
3542 if(Status /= WRF_NO_ERR) then
3543 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3544 call wrf_debug ( WARN , TRIM(msg))
3548 end subroutine ext_pnc_end_independent_mode
3550 subroutine ext_pnc_start_independent_mode(DataHandle, Status)
3552 use ext_pnc_support_routines
3553 include 'wrf_status_codes.h'
3554 # include "pnetcdf.inc"
3555 integer ,intent(in) :: DataHandle
3556 integer ,intent(out) :: Status
3557 type(wrf_data_handle) ,pointer :: DH
3560 DH => WrfDataHandles(DataHandle)
3561 DH%Collective = .FALSE.
3562 stat = NFMPI_BEGIN_INDEP_DATA(DH%NCID)
3563 call netcdf_err(stat,Status)
3564 if(Status /= WRF_NO_ERR) then
3565 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3566 call wrf_debug ( WARN , TRIM(msg))
3571 end subroutine ext_pnc_start_independent_mode