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 = 10000
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*128 :: idstr,ntasks_x_str,loccomm_str
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
1252 stat = NFMPI_CREATE(Comm, FileName, IOR(NF_CLOBBER, NF_64BIT_OFFSET), info, DH%NCID)
1253 call mpi_info_free( info, ierr)
1256 ! rob latham suggested hint
1258 call mpi_info_create( info, ierr )
1259 ! call mpi_info_set(info,'cd_buffer_size','4194304',ierr)
1260 call mpi_info_set(info,'cd_buffer_size','8388608',ierr)
1261 stat = NFMPI_CREATE(Comm, FileName, IOR(NF_CLOBBER, NF_64BIT_OFFSET), info, DH%NCID)
1262 call mpi_info_free( info, ierr)
1267 call netcdf_err(stat,Status)
1268 if(Status /= WRF_NO_ERR) then
1269 write(msg,*) 'NetCDF error in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1270 call wrf_debug ( WARN , TRIM(msg))
1273 ! JPE added for performance
1274 stat = NFMPI_SET_FILL(DH%NCID, NF_NOFILL, i)
1276 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1277 DH%FileName = FileName
1278 stat = NFMPI_DEF_DIM(DH%NCID,DH%DimUnlimName,i2offset(NF_UNLIMITED),DH%DimUnlimID)
1279 call netcdf_err(stat,Status)
1280 if(Status /= WRF_NO_ERR) then
1281 write(msg,*) 'NetCDF error in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1282 call wrf_debug ( WARN , TRIM(msg))
1285 DH%VarNames (1:MaxVars) = NO_NAME
1286 DH%MDVarNames(1:MaxVars) = NO_NAME
1288 write(Buffer,FMT="('DIM',i4.4)") i
1289 DH%DimNames (i) = Buffer
1290 DH%DimLengths(i) = NO_DIM
1292 DH%DimNames(1) = 'DateStrLen'
1293 stat = NFMPI_DEF_DIM(DH%NCID,DH%DimNames(1),i2offset(DateStrLen),DH%DimIDs(1))
1294 call netcdf_err(stat,Status)
1295 if(Status /= WRF_NO_ERR) then
1296 write(msg,*) 'NetCDF error in ext_pnc_open_for_write_begin ',__FILE__,', line', __LINE__
1297 call wrf_debug ( WARN , TRIM(msg))
1300 VDimIDs(1) = DH%DimIDs(1)
1301 VDimIDs(2) = DH%DimUnlimID
1302 stat = NFMPI_DEF_VAR(DH%NCID,DH%TimesName,NF_CHAR,2,VDimIDs,DH%TimesVarID)
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 DH%DimLengths(1) = DateStrLen
1311 end subroutine ext_pnc_open_for_write_begin
1314 !opens a file for writing or coupler datastream for sending messages.
1315 !no training phase for this version of the open stmt.
1316 subroutine ext_pnc_open_for_write (DatasetName, Comm1, Comm2, &
1317 SysDepInfo, DataHandle, Status)
1319 use ext_pnc_support_routines
1321 include 'wrf_status_codes.h'
1322 # include "pnetcdf.inc"
1323 character *(*), intent(in) ::DatasetName
1324 integer , intent(in) ::Comm1, Comm2
1325 character *(*), intent(in) ::SysDepInfo
1326 integer , intent(out) :: DataHandle
1327 integer , intent(out) :: Status
1328 Status=WRF_WARN_NOOP
1329 DataHandle = 0 ! dummy setting to quiet warning message
1331 end subroutine ext_pnc_open_for_write
1333 SUBROUTINE ext_pnc_open_for_write_commit(DataHandle, Status)
1335 use ext_pnc_support_routines
1337 include 'wrf_status_codes.h'
1338 # include "pnetcdf.inc"
1339 integer ,intent(in) :: DataHandle
1340 integer ,intent(out) :: Status
1341 type(wrf_data_handle),pointer :: DH
1345 if(WrfIOnotInitialized) then
1346 Status = WRF_IO_NOT_INITIALIZED
1347 write(msg,*) 'ext_pnc_open_for_write_commit: ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
1348 call wrf_debug ( FATAL , msg)
1351 call GetDH(DataHandle,DH,Status)
1352 if(Status /= WRF_NO_ERR) then
1353 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_open_for_write_commit ',__FILE__,', line', __LINE__
1354 call wrf_debug ( WARN , TRIM(msg))
1357 stat = NFMPI_ENDDEF(DH%NCID)
1358 call netcdf_err(stat,Status)
1359 if(Status /= WRF_NO_ERR) then
1360 write(msg,*) 'NetCDF error (',stat,') from NFMPI_ENDDEF in ext_pnc_open_for_write_commit ',__FILE__,', line', __LINE__
1361 call wrf_debug ( WARN , TRIM(msg))
1364 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1365 DH%first_operation = .TRUE.
1367 end subroutine ext_pnc_open_for_write_commit
1369 subroutine ext_pnc_ioclose(DataHandle, Status)
1371 use ext_pnc_support_routines
1373 include 'wrf_status_codes.h'
1374 # include "pnetcdf.inc"
1375 integer ,intent(in) :: DataHandle
1376 integer ,intent(out) :: Status
1377 type(wrf_data_handle),pointer :: DH
1380 call GetDH(DataHandle,DH,Status)
1381 if(Status /= WRF_NO_ERR) then
1382 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1383 call wrf_debug ( WARN , TRIM(msg))
1386 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1387 Status = WRF_WARN_FILE_NOT_OPENED
1388 write(msg,*) 'Warning FILE NOT OPENED in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1389 call wrf_debug ( WARN , TRIM(msg))
1390 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1391 Status = WRF_WARN_DRYRUN_CLOSE
1392 write(msg,*) 'Warning TRY TO CLOSE DRYRUN in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1393 call wrf_debug ( WARN , TRIM(msg))
1394 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1396 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1398 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
1401 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1402 write(msg,*) 'Fatal error BAD FILE STATUS in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1403 call wrf_debug ( FATAL , TRIM(msg))
1407 stat = NFMPI_CLOSE(DH%NCID)
1408 call netcdf_err(stat,Status)
1409 if(Status /= WRF_NO_ERR) then
1410 write(msg,*) 'NetCDF error in ext_pnc_ioclose ',__FILE__,', line', __LINE__
1411 call wrf_debug ( WARN , TRIM(msg))
1414 CALL deallocHandle( DataHandle, Status )
1417 end subroutine ext_pnc_ioclose
1419 subroutine ext_pnc_iosync( DataHandle, Status)
1421 use ext_pnc_support_routines
1423 include 'wrf_status_codes.h'
1424 # include "pnetcdf.inc"
1425 integer ,intent(in) :: DataHandle
1426 integer ,intent(out) :: Status
1427 type(wrf_data_handle),pointer :: DH
1430 call GetDH(DataHandle,DH,Status)
1431 if(Status /= WRF_NO_ERR) then
1432 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_iosync ',__FILE__,', line', __LINE__
1433 call wrf_debug ( WARN , TRIM(msg))
1436 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1437 Status = WRF_WARN_FILE_NOT_OPENED
1438 write(msg,*) 'Warning FILE NOT OPENED in ext_pnc_iosync ',__FILE__,', line', __LINE__
1439 call wrf_debug ( WARN , TRIM(msg))
1440 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1441 Status = WRF_WARN_FILE_NOT_COMMITTED
1442 write(msg,*) 'Warning FILE NOT COMMITTED in ext_pnc_iosync ',__FILE__,', line', __LINE__
1443 call wrf_debug ( WARN , TRIM(msg))
1444 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1446 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1449 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1450 write(msg,*) 'Fatal error BAD FILE STATUS in ext_pnc_iosync ',__FILE__,', line', __LINE__
1451 call wrf_debug ( FATAL , TRIM(msg))
1454 stat = NFMPI_SYNC(DH%NCID)
1455 call netcdf_err(stat,Status)
1456 if(Status /= WRF_NO_ERR) then
1457 write(msg,*) 'NetCDF error in ext_pnc_iosync ',__FILE__,', line', __LINE__
1458 call wrf_debug ( WARN , TRIM(msg))
1462 end subroutine ext_pnc_iosync
1466 subroutine ext_pnc_redef( DataHandle, Status)
1468 use ext_pnc_support_routines
1470 include 'wrf_status_codes.h'
1471 # include "pnetcdf.inc"
1472 integer ,intent(in) :: DataHandle
1473 integer ,intent(out) :: Status
1474 type(wrf_data_handle),pointer :: DH
1477 call GetDH(DataHandle,DH,Status)
1478 if(Status /= WRF_NO_ERR) then
1479 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1480 call wrf_debug ( WARN , TRIM(msg))
1483 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1484 Status = WRF_WARN_FILE_NOT_OPENED
1485 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1486 call wrf_debug ( WARN , TRIM(msg))
1487 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1488 Status = WRF_WARN_FILE_NOT_COMMITTED
1489 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1490 call wrf_debug ( WARN , TRIM(msg))
1491 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1493 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1494 Status = WRF_WARN_FILE_OPEN_FOR_READ
1495 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1496 call wrf_debug ( WARN , TRIM(msg))
1498 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1499 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1500 call wrf_debug ( FATAL , TRIM(msg))
1503 stat = NFMPI_REDEF(DH%NCID)
1504 call netcdf_err(stat,Status)
1505 if(Status /= WRF_NO_ERR) then
1506 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1507 call wrf_debug ( WARN , TRIM(msg))
1510 DH%FileStatus = WRF_FILE_OPENED_NOT_COMMITTED
1512 end subroutine ext_pnc_redef
1514 subroutine ext_pnc_enddef( DataHandle, Status)
1516 use ext_pnc_support_routines
1518 include 'wrf_status_codes.h'
1519 # include "pnetcdf.inc"
1520 integer ,intent(in) :: DataHandle
1521 integer ,intent(out) :: Status
1522 type(wrf_data_handle),pointer :: DH
1525 call GetDH(DataHandle,DH,Status)
1526 if(Status /= WRF_NO_ERR) then
1527 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1528 call wrf_debug ( WARN , TRIM(msg))
1531 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
1532 Status = WRF_WARN_FILE_NOT_OPENED
1533 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
1534 call wrf_debug ( WARN , TRIM(msg))
1535 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
1536 Status = WRF_WARN_FILE_NOT_COMMITTED
1537 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
1538 call wrf_debug ( WARN , TRIM(msg))
1539 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
1541 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
1542 Status = WRF_WARN_FILE_OPEN_FOR_READ
1543 write(msg,*) 'Warning FILE OPEN FOR READ in ',__FILE__,', line', __LINE__
1544 call wrf_debug ( WARN , TRIM(msg))
1546 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
1547 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
1548 call wrf_debug ( FATAL , TRIM(msg))
1551 stat = NFMPI_ENDDEF(DH%NCID)
1552 call netcdf_err(stat,Status)
1553 if(Status /= WRF_NO_ERR) then
1554 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
1555 call wrf_debug ( WARN , TRIM(msg))
1558 DH%FileStatus = WRF_FILE_OPENED_FOR_WRITE
1560 end subroutine ext_pnc_enddef
1562 subroutine ext_pnc_ioinit(SysDepInfo, Status)
1565 include 'wrf_status_codes.h'
1566 CHARACTER*(*), INTENT(IN) :: SysDepInfo
1567 INTEGER ,INTENT(INOUT) :: Status
1569 WrfIOnotInitialized = .false.
1570 WrfDataHandles(1:WrfDataHandleMax)%Free = .true.
1571 WrfDataHandles(1:WrfDataHandleMax)%TimesName = 'Times'
1572 WrfDataHandles(1:WrfDataHandleMax)%DimUnlimName = 'Time'
1573 WrfDataHandles(1:WrfDataHandleMax)%FileStatus = WRF_FILE_NOT_OPENED
1576 end subroutine ext_pnc_ioinit
1579 subroutine ext_pnc_inquiry (Inquiry, Result, Status)
1582 include 'wrf_status_codes.h'
1583 character *(*), INTENT(IN) :: Inquiry
1584 character *(*), INTENT(OUT) :: Result
1585 integer ,INTENT(INOUT) :: Status
1586 SELECT CASE (Inquiry)
1587 CASE ("RANDOM_WRITE","RANDOM_READ","SEQUENTIAL_WRITE","SEQUENTIAL_READ")
1589 CASE ("OPEN_READ","OPEN_COMMIT_WRITE")
1591 CASE ("OPEN_WRITE","OPEN_COMMIT_READ","PARALLEL_IO")
1593 CASE ("SELF_DESCRIBING","SUPPORT_METADATA","SUPPORT_3D_FIELDS")
1598 Result = 'No Result for that inquiry!'
1602 end subroutine ext_pnc_inquiry
1607 subroutine ext_pnc_ioexit(Status)
1609 use ext_pnc_support_routines
1611 include 'wrf_status_codes.h'
1612 # include "pnetcdf.inc"
1613 integer , INTENT(INOUT) ::Status
1615 type(wrf_data_handle),pointer :: DH
1618 if(WrfIOnotInitialized) then
1619 Status = WRF_IO_NOT_INITIALIZED
1620 write(msg,*) 'ext_pnc_ioinit was not called ',__FILE__,', line', __LINE__
1621 call wrf_debug ( FATAL , msg)
1624 do i=1,WrfDataHandleMax
1625 CALL deallocHandle( i , stat )
1628 end subroutine ext_pnc_ioexit
1630 subroutine ext_pnc_get_dom_ti_real(DataHandle,Element,Data,Count,OutCount,Status)
1631 #define ROUTINE_TYPE 'REAL'
1632 #define TYPE_DATA real,intent(out) :: Data(*)
1633 #define TYPE_COUNT integer,intent(in) :: Count
1634 #define TYPE_OUTCOUNT integer,intent(out) :: OutCOunt
1635 #define TYPE_BUFFER real,allocatable :: Buffer(:)
1636 #define NF_TYPE NF_FLOAT
1637 #define NF_ROUTINE NFMPI_GET_ATT_REAL
1638 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1639 #include "ext_pnc_get_dom_ti.code"
1640 end subroutine ext_pnc_get_dom_ti_real
1642 subroutine ext_pnc_get_dom_ti_integer(DataHandle,Element,Data,Count,OutCount,Status)
1649 #define ROUTINE_TYPE 'INTEGER'
1650 #define TYPE_DATA integer,intent(out) :: Data(*)
1651 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
1652 #define NF_TYPE NF_INT
1653 #define NF_ROUTINE NFMPI_GET_ATT_INT
1654 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1655 #include "ext_pnc_get_dom_ti.code"
1656 end subroutine ext_pnc_get_dom_ti_integer
1658 subroutine ext_pnc_get_dom_ti_double(DataHandle,Element,Data,Count,OutCount,Status)
1665 #define ROUTINE_TYPE 'DOUBLE'
1666 #define TYPE_DATA real*8,intent(out) :: Data(*)
1667 #define TYPE_BUFFER real*8,allocatable :: Buffer(:)
1668 #define NF_TYPE NF_DOUBLE
1669 #define NF_ROUTINE NFMPI_GET_ATT_DOUBLE
1670 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))
1671 #include "ext_pnc_get_dom_ti.code"
1672 end subroutine ext_pnc_get_dom_ti_double
1674 subroutine ext_pnc_get_dom_ti_logical(DataHandle,Element,Data,Count,OutCount,Status)
1681 #define ROUTINE_TYPE 'LOGICAL'
1682 #define TYPE_DATA logical,intent(out) :: Data(*)
1683 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
1684 #define NF_TYPE NF_INT
1685 #define NF_ROUTINE NFMPI_GET_ATT_INT
1686 #define COPY Data(1:min(Len,Count)) = Buffer(1:min(Len,Count))==1
1687 #include "ext_pnc_get_dom_ti.code"
1688 end subroutine ext_pnc_get_dom_ti_logical
1690 subroutine ext_pnc_get_dom_ti_char(DataHandle,Element,Data,Status)
1694 #undef TYPE_OUTCOUNT
1697 #define ROUTINE_TYPE 'CHAR'
1698 #define TYPE_DATA character*(*),intent(out) :: Data
1700 #define TYPE_OUTCOUNT
1702 #define NF_TYPE NF_CHAR
1704 #include "ext_pnc_get_dom_ti.code"
1706 end subroutine ext_pnc_get_dom_ti_char
1708 subroutine ext_pnc_put_dom_ti_real(DataHandle,Element,Data,Count,Status)
1715 #define ROUTINE_TYPE 'REAL'
1716 #define TYPE_DATA real ,intent(in) :: Data(*)
1717 #define TYPE_COUNT integer,intent(in) :: Count
1718 #define NF_ROUTINE NFMPI_PUT_ATT_REAL
1719 #define ARGS NF_FLOAT,i2offset(Count),Data
1720 #include "ext_pnc_put_dom_ti.code"
1721 end subroutine ext_pnc_put_dom_ti_real
1723 subroutine ext_pnc_put_dom_ti_integer(DataHandle,Element,Data,Count,Status)
1730 #define ROUTINE_TYPE 'INTEGER'
1731 #define TYPE_DATA integer,intent(in) :: Data(*)
1732 #define TYPE_COUNT integer,intent(in) :: Count
1733 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1734 #define ARGS NF_INT,i2offset(Count),Data
1735 #include "ext_pnc_put_dom_ti.code"
1736 end subroutine ext_pnc_put_dom_ti_integer
1738 subroutine ext_pnc_put_dom_ti_double(DataHandle,Element,Data,Count,Status)
1745 #define ROUTINE_TYPE 'DOUBLE'
1746 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1747 #define TYPE_COUNT integer,intent(in) :: Count
1748 #define NF_ROUTINE NFMPI_PUT_ATT_DOUBLE
1749 #define ARGS NF_DOUBLE,i2offset(Count),Data
1750 #include "ext_pnc_put_dom_ti.code"
1751 end subroutine ext_pnc_put_dom_ti_double
1753 subroutine ext_pnc_put_dom_ti_logical(DataHandle,Element,Data,Count,Status)
1759 #define ROUTINE_TYPE 'LOGICAL'
1760 #define TYPE_DATA logical,intent(in) :: Data(*)
1761 #define TYPE_COUNT integer,intent(in) :: Count
1762 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1763 #define ARGS NF_INT,i2offset(Count),Buffer
1765 #include "ext_pnc_put_dom_ti.code"
1766 end subroutine ext_pnc_put_dom_ti_logical
1768 subroutine ext_pnc_put_dom_ti_char(DataHandle,Element,Data,Status)
1775 #define ROUTINE_TYPE 'CHAR'
1776 #define TYPE_DATA character*(*),intent(in) :: Data
1777 #define TYPE_COUNT integer,parameter :: Count=1
1778 #define NF_ROUTINE NFMPI_PUT_ATT_TEXT
1779 #define ARGS i2offset(len_trim(Data)),Data
1780 #include "ext_pnc_put_dom_ti.code"
1781 end subroutine ext_pnc_put_dom_ti_char
1783 subroutine ext_pnc_put_var_ti_real(DataHandle,Element,Var,Data,Count,Status)
1790 #define ROUTINE_TYPE 'REAL'
1791 #define TYPE_DATA real ,intent(in) :: Data(*)
1792 #define TYPE_COUNT integer ,intent(in) :: Count
1793 #define NF_ROUTINE NFMPI_PUT_ATT_REAL
1794 #define ARGS NF_FLOAT,i2offset(Count),Data
1795 #include "ext_pnc_put_var_ti.code"
1796 end subroutine ext_pnc_put_var_ti_real
1798 subroutine ext_pnc_put_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,Status)
1807 #define ROUTINE_TYPE 'REAL'
1808 #define TYPE_DATA real ,intent(in) :: Data(*)
1809 #define TYPE_COUNT integer ,intent(in) :: Count
1810 #define NF_ROUTINE NFMPI_PUT_VARA_REAL_ALL
1811 #define NF_TYPE NF_FLOAT
1812 #define LENGTH Count
1814 #include "ext_pnc_put_var_td.code"
1815 end subroutine ext_pnc_put_var_td_real
1817 subroutine ext_pnc_put_var_ti_double(DataHandle,Element,Var,Data,Count,Status)
1824 #define ROUTINE_TYPE 'DOUBLE'
1825 #define TYPE_DATA real*8 ,intent(in) :: Data(*)
1826 #define TYPE_COUNT integer ,intent(in) :: Count
1827 #define NF_ROUTINE NFMPI_PUT_ATT_DOUBLE
1828 #define ARGS NF_DOUBLE,i2offset(Count),Data
1829 #include "ext_pnc_put_var_ti.code"
1830 end subroutine ext_pnc_put_var_ti_double
1832 subroutine ext_pnc_put_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,Status)
1841 #define ROUTINE_TYPE 'DOUBLE'
1842 #define TYPE_DATA real*8,intent(in) :: Data(*)
1843 #define TYPE_COUNT integer ,intent(in) :: Count
1844 #define NF_ROUTINE NFMPI_PUT_VARA_DOUBLE_ALL
1845 #define NF_TYPE NF_DOUBLE
1846 #define LENGTH Count
1848 #include "ext_pnc_put_var_td.code"
1849 end subroutine ext_pnc_put_var_td_double
1851 subroutine ext_pnc_put_var_ti_integer(DataHandle,Element,Var,Data,Count,Status)
1858 #define ROUTINE_TYPE 'INTEGER'
1859 #define TYPE_DATA integer ,intent(in) :: Data(*)
1860 #define TYPE_COUNT integer ,intent(in) :: Count
1861 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1862 #define ARGS NF_INT,i2offset(Count),Data
1863 #include "ext_pnc_put_var_ti.code"
1864 end subroutine ext_pnc_put_var_ti_integer
1866 subroutine ext_pnc_put_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,Status)
1875 #define ROUTINE_TYPE 'INTEGER'
1876 #define TYPE_DATA integer ,intent(in) :: Data(*)
1877 #define TYPE_COUNT integer ,intent(in) :: Count
1878 #define NF_ROUTINE NFMPI_PUT_VARA_INT_ALL
1879 #define NF_TYPE NF_INT
1880 #define LENGTH Count
1882 #include "ext_pnc_put_var_td.code"
1883 end subroutine ext_pnc_put_var_td_integer
1885 subroutine ext_pnc_put_var_ti_logical(DataHandle,Element,Var,Data,Count,Status)
1891 #define ROUTINE_TYPE 'LOGICAL'
1892 #define TYPE_DATA logical ,intent(in) :: Data(*)
1893 #define TYPE_COUNT integer ,intent(in) :: Count
1894 #define NF_ROUTINE NFMPI_PUT_ATT_INT
1896 #define ARGS NF_INT,i2offset(Count),Buffer
1897 #include "ext_pnc_put_var_ti.code"
1898 end subroutine ext_pnc_put_var_ti_logical
1900 subroutine ext_pnc_put_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,Status)
1908 #define ROUTINE_TYPE 'LOGICAL'
1909 #define TYPE_DATA logical ,intent(in) :: Data(*)
1910 #define TYPE_COUNT integer ,intent(in) :: Count
1911 #define NF_ROUTINE NFMPI_PUT_VARA_INT_ALL
1912 #define NF_TYPE NF_INT
1914 #define LENGTH Count
1916 #include "ext_pnc_put_var_td.code"
1917 end subroutine ext_pnc_put_var_td_logical
1919 subroutine ext_pnc_put_var_ti_char(DataHandle,Element,Var,Data,Status)
1926 #define ROUTINE_TYPE 'CHAR'
1927 #define TYPE_DATA character*(*) ,intent(in) :: Data
1929 #define NF_ROUTINE NFMPI_PUT_ATT_TEXT
1930 #define ARGS i2offset(len_trim(Data)),trim(Data)
1932 #include "ext_pnc_put_var_ti.code"
1934 end subroutine ext_pnc_put_var_ti_char
1936 subroutine ext_pnc_put_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
1945 #define ROUTINE_TYPE 'CHAR'
1946 #define TYPE_DATA character*(*) ,intent(in) :: Data
1948 #define NF_ROUTINE NFMPI_PUT_VARA_TEXT_ALL
1949 #define NF_TYPE NF_CHAR
1950 #define LENGTH len(Data)
1951 #include "ext_pnc_put_var_td.code"
1952 end subroutine ext_pnc_put_var_td_char
1954 subroutine ext_pnc_get_var_ti_real(DataHandle,Element,Var,Data,Count,OutCount,Status)
1959 #undef TYPE_OUTCOUNT
1963 #define ROUTINE_TYPE 'REAL'
1964 #define TYPE_DATA real ,intent(out) :: Data(*)
1965 #define TYPE_BUFFER real ,allocatable :: Buffer(:)
1966 #define TYPE_COUNT integer,intent(in) :: Count
1967 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
1968 #define NF_TYPE NF_FLOAT
1969 #define NF_ROUTINE NFMPI_GET_ATT_REAL
1970 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
1971 #include "ext_pnc_get_var_ti.code"
1972 end subroutine ext_pnc_get_var_ti_real
1974 subroutine ext_pnc_get_var_td_real(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
1979 #undef TYPE_OUTCOUNT
1984 #define ROUTINE_TYPE 'REAL'
1985 #define TYPE_DATA real ,intent(out) :: Data(*)
1986 #define TYPE_BUFFER real
1987 #define TYPE_COUNT integer,intent(in) :: Count
1988 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
1989 #define NF_TYPE NF_FLOAT
1990 #define NF_ROUTINE NFMPI_GET_VARA_REAL_ALL
1991 #define LENGTH min(Count,Len1)
1992 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
1993 #include "ext_pnc_get_var_td.code"
1994 end subroutine ext_pnc_get_var_td_real
1996 subroutine ext_pnc_get_var_ti_double(DataHandle,Element,Var,Data,Count,OutCount,Status)
2001 #undef TYPE_OUTCOUNT
2005 #define ROUTINE_TYPE 'DOUBLE'
2006 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2007 #define TYPE_BUFFER real*8 ,allocatable :: Buffer(:)
2008 #define TYPE_COUNT integer,intent(in) :: Count
2009 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2010 #define NF_TYPE NF_DOUBLE
2011 #define NF_ROUTINE NFMPI_GET_ATT_DOUBLE
2012 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2013 #include "ext_pnc_get_var_ti.code"
2014 end subroutine ext_pnc_get_var_ti_double
2016 subroutine ext_pnc_get_var_td_double(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2021 #undef TYPE_OUTCOUNT
2026 #define ROUTINE_TYPE 'DOUBLE'
2027 #define TYPE_DATA real*8 ,intent(out) :: Data(*)
2028 #define TYPE_BUFFER real*8
2029 #define TYPE_COUNT integer,intent(in) :: Count
2030 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2031 #define NF_TYPE NF_DOUBLE
2032 #define NF_ROUTINE NFMPI_GET_VARA_DOUBLE_ALL
2033 #define LENGTH min(Count,Len1)
2034 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2035 #include "ext_pnc_get_var_td.code"
2036 end subroutine ext_pnc_get_var_td_double
2038 subroutine ext_pnc_get_var_ti_integer(DataHandle,Element,Var,Data,Count,OutCount,Status)
2043 #undef TYPE_OUTCOUNT
2047 #define ROUTINE_TYPE 'INTEGER'
2048 #define TYPE_DATA integer,intent(out) :: Data(*)
2049 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2050 #define TYPE_COUNT integer,intent(in) :: Count
2051 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2052 #define NF_TYPE NF_INT
2053 #define NF_ROUTINE NFMPI_GET_ATT_INT
2054 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))
2055 #include "ext_pnc_get_var_ti.code"
2056 end subroutine ext_pnc_get_var_ti_integer
2058 subroutine ext_pnc_get_var_td_integer(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2063 #undef TYPE_OUTCOUNT
2068 #define ROUTINE_TYPE 'INTEGER'
2069 #define TYPE_DATA integer,intent(out) :: Data(*)
2070 #define TYPE_BUFFER integer
2071 #define TYPE_COUNT integer,intent(in) :: Count
2072 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2073 #define NF_TYPE NF_INT
2074 #define NF_ROUTINE NFMPI_GET_VARA_INT_ALL
2075 #define LENGTH min(Count,Len1)
2076 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))
2077 #include "ext_pnc_get_var_td.code"
2078 end subroutine ext_pnc_get_var_td_integer
2080 subroutine ext_pnc_get_var_ti_logical(DataHandle,Element,Var,Data,Count,OutCount,Status)
2085 #undef TYPE_OUTCOUNT
2089 #define ROUTINE_TYPE 'LOGICAL'
2090 #define TYPE_DATA logical,intent(out) :: Data(*)
2091 #define TYPE_BUFFER integer,allocatable :: Buffer(:)
2092 #define TYPE_COUNT integer,intent(in) :: Count
2093 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2094 #define NF_TYPE NF_INT
2095 #define NF_ROUTINE NFMPI_GET_ATT_INT
2096 #define COPY Data(1:min(XLen,Count)) = Buffer(1:min(XLen,Count))==1
2097 #include "ext_pnc_get_var_ti.code"
2098 end subroutine ext_pnc_get_var_ti_logical
2100 subroutine ext_pnc_get_var_td_logical(DataHandle,Element,DateStr,Var,Data,Count,OutCount,Status)
2105 #undef TYPE_OUTCOUNT
2110 #define ROUTINE_TYPE 'LOGICAL'
2111 #define TYPE_DATA logical,intent(out) :: Data(*)
2112 #define TYPE_BUFFER integer
2113 #define TYPE_COUNT integer,intent(in) :: Count
2114 #define TYPE_OUTCOUNT integer,intent(out) :: OutCount
2115 #define NF_TYPE NF_INT
2116 #define NF_ROUTINE NFMPI_GET_VARA_INT_ALL
2117 #define LENGTH min(Count,Len1)
2118 #define COPY Data(1:min(Len1,Count)) = Buffer(1:min(Len1,Count))==1
2119 #include "ext_pnc_get_var_td.code"
2120 end subroutine ext_pnc_get_var_td_logical
2122 subroutine ext_pnc_get_var_ti_char(DataHandle,Element,Var,Data,Status)
2127 #undef TYPE_OUTCOUNT
2131 #define ROUTINE_TYPE 'CHAR'
2132 #define TYPE_DATA character*(*) ,intent(out) :: Data
2134 #define TYPE_COUNT integer :: Count = 1
2135 #define TYPE_OUTCOUNT
2136 #define NF_TYPE NF_CHAR
2137 #define NF_ROUTINE NFMPI_GET_ATT_TEXT
2140 #include "ext_pnc_get_var_ti.code"
2142 end subroutine ext_pnc_get_var_ti_char
2144 subroutine ext_pnc_get_var_td_char(DataHandle,Element,DateStr,Var,Data,Status)
2149 #undef TYPE_OUTCOUNT
2153 #define ROUTINE_TYPE 'CHAR'
2154 #define TYPE_DATA character*(*) ,intent(out) :: Data
2155 #define TYPE_BUFFER character (80)
2156 #define TYPE_COUNT integer :: Count = 1
2157 #define TYPE_OUTCOUNT
2158 #define NF_TYPE NF_CHAR
2159 #define NF_ROUTINE NFMPI_GET_VARA_TEXT_ALL
2162 #include "ext_pnc_get_var_td.code"
2164 end subroutine ext_pnc_get_var_td_char
2166 subroutine ext_pnc_put_dom_td_real(DataHandle,Element,DateStr,Data,Count,Status)
2167 integer ,intent(in) :: DataHandle
2168 character*(*) ,intent(in) :: Element
2169 character*(*) ,intent(in) :: DateStr
2170 real ,intent(in) :: Data(*)
2171 integer ,intent(in) :: Count
2172 integer ,intent(out) :: Status
2174 call ext_pnc_put_var_td_real(DataHandle,Element,DateStr, &
2175 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2177 end subroutine ext_pnc_put_dom_td_real
2179 subroutine ext_pnc_put_dom_td_integer(DataHandle,Element,DateStr,Data,Count,Status)
2180 integer ,intent(in) :: DataHandle
2181 character*(*) ,intent(in) :: Element
2182 character*(*) ,intent(in) :: DateStr
2183 integer ,intent(in) :: Data(*)
2184 integer ,intent(in) :: Count
2185 integer ,intent(out) :: Status
2187 call ext_pnc_put_var_td_integer(DataHandle,Element,DateStr, &
2188 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2190 end subroutine ext_pnc_put_dom_td_integer
2192 subroutine ext_pnc_put_dom_td_double(DataHandle,Element,DateStr,Data,Count,Status)
2193 integer ,intent(in) :: DataHandle
2194 character*(*) ,intent(in) :: Element
2195 character*(*) ,intent(in) :: DateStr
2196 real*8 ,intent(in) :: Data(*)
2197 integer ,intent(in) :: Count
2198 integer ,intent(out) :: Status
2200 call ext_pnc_put_var_td_double(DataHandle,Element,DateStr, &
2201 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2203 end subroutine ext_pnc_put_dom_td_double
2205 subroutine ext_pnc_put_dom_td_logical(DataHandle,Element,DateStr,Data,Count,Status)
2206 integer ,intent(in) :: DataHandle
2207 character*(*) ,intent(in) :: Element
2208 character*(*) ,intent(in) :: DateStr
2209 logical ,intent(in) :: Data(*)
2210 integer ,intent(in) :: Count
2211 integer ,intent(out) :: Status
2213 call ext_pnc_put_var_td_logical(DataHandle,Element,DateStr, &
2214 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,Status)
2216 end subroutine ext_pnc_put_dom_td_logical
2218 subroutine ext_pnc_put_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2219 integer ,intent(in) :: DataHandle
2220 character*(*) ,intent(in) :: Element
2221 character*(*) ,intent(in) :: DateStr
2222 character*(*) ,intent(in) :: Data
2223 integer ,intent(out) :: Status
2225 call ext_pnc_put_var_td_char(DataHandle,Element,DateStr, &
2226 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2228 end subroutine ext_pnc_put_dom_td_char
2230 subroutine ext_pnc_get_dom_td_real(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2231 integer ,intent(in) :: DataHandle
2232 character*(*) ,intent(in) :: Element
2233 character*(*) ,intent(in) :: DateStr
2234 real ,intent(out) :: Data(*)
2235 integer ,intent(in) :: Count
2236 integer ,intent(out) :: OutCount
2237 integer ,intent(out) :: Status
2238 call ext_pnc_get_var_td_real(DataHandle,Element,DateStr, &
2239 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2241 end subroutine ext_pnc_get_dom_td_real
2243 subroutine ext_pnc_get_dom_td_integer(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2244 integer ,intent(in) :: DataHandle
2245 character*(*) ,intent(in) :: Element
2246 character*(*) ,intent(in) :: DateStr
2247 integer ,intent(out) :: Data(*)
2248 integer ,intent(in) :: Count
2249 integer ,intent(out) :: OutCount
2250 integer ,intent(out) :: Status
2251 call ext_pnc_get_var_td_integer(DataHandle,Element,DateStr, &
2252 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2254 end subroutine ext_pnc_get_dom_td_integer
2256 subroutine ext_pnc_get_dom_td_double(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2257 integer ,intent(in) :: DataHandle
2258 character*(*) ,intent(in) :: Element
2259 character*(*) ,intent(in) :: DateStr
2260 real*8 ,intent(out) :: Data(*)
2261 integer ,intent(in) :: Count
2262 integer ,intent(out) :: OutCount
2263 integer ,intent(out) :: Status
2264 call ext_pnc_get_var_td_double(DataHandle,Element,DateStr, &
2265 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2267 end subroutine ext_pnc_get_dom_td_double
2269 subroutine ext_pnc_get_dom_td_logical(DataHandle,Element,DateStr,Data,Count,OutCount,Status)
2270 integer ,intent(in) :: DataHandle
2271 character*(*) ,intent(in) :: Element
2272 character*(*) ,intent(in) :: DateStr
2273 logical ,intent(out) :: Data(*)
2274 integer ,intent(in) :: Count
2275 integer ,intent(out) :: OutCount
2276 integer ,intent(out) :: Status
2277 call ext_pnc_get_var_td_logical(DataHandle,Element,DateStr, &
2278 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Count,OutCount,Status)
2280 end subroutine ext_pnc_get_dom_td_logical
2282 subroutine ext_pnc_get_dom_td_char(DataHandle,Element,DateStr,Data,Status)
2283 integer ,intent(in) :: DataHandle
2284 character*(*) ,intent(in) :: Element
2285 character*(*) ,intent(in) :: DateStr
2286 character*(*) ,intent(out) :: Data
2287 integer ,intent(out) :: Status
2288 call ext_pnc_get_var_td_char(DataHandle,Element,DateStr, &
2289 'E_X_T_D_O_M_A_I_N_M_E_T_A_DATA_' ,Data,Status)
2291 end subroutine ext_pnc_get_dom_td_char
2294 subroutine ext_pnc_write_field(DataHandle,DateStr,Var,Field,FieldType,Comm, &
2295 IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2296 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2298 use ext_pnc_support_routines
2300 include 'wrf_status_codes.h'
2301 # include "pnetcdf.inc"
2302 integer ,intent(in) :: DataHandle
2303 character*(*) ,intent(in) :: DateStr
2304 character*(*) ,intent(in) :: Var
2305 integer ,intent(inout) :: Field(*)
2306 integer ,intent(in) :: FieldType
2307 integer ,intent(inout) :: Comm
2308 integer ,intent(inout) :: IOComm
2309 integer ,intent(in) :: DomainDesc
2310 character*(*) ,intent(in) :: MemoryOrdIn
2311 character*(*) ,intent(in) :: Stagger ! Dummy for now
2312 character*(*) ,dimension(*) ,intent(in) :: DimNames
2313 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2314 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2315 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2316 integer ,intent(out) :: Status
2317 character (3) :: MemoryOrder
2318 type(wrf_data_handle) ,pointer :: DH
2321 character (VarNameLen) :: VarName
2322 character (3) :: MemO
2323 character (3) :: UCMemO
2325 integer ,dimension(NVarDims) :: Length_global, Length_native
2326 integer ,dimension(NVarDims) :: Length
2327 integer ,dimension(NVarDims) :: VDimIDs
2328 character(80),dimension(NVarDims) :: RODimNames
2329 integer ,dimension(NVarDims) :: StoredStart
2330 integer ,dimension(:,:,:,:),allocatable :: XField
2334 integer :: i1,i2,j1,j2,k1,k2
2335 integer :: x1,x2,y1,y2,z1,z2
2336 integer :: p1,p2,q1,q2,r1,r2
2337 integer :: l1,l2,m1,m2,n1,n2
2340 character (80) :: NullName
2343 ! Local, possibly adjusted, copies of MemoryStart and MemoryEnd
2344 integer ,dimension(NVarDims) :: lMemoryStart, lMemoryEnd
2345 MemoryOrder = trim(adjustl(MemoryOrdIn))
2347 call GetDim(MemoryOrder,NDim,Status)
2348 if(Status /= WRF_NO_ERR) then
2349 write(msg,*) 'Warning BAD MEMORY ORDER |',MemoryOrder,'| in ',__FILE__,', line', __LINE__
2350 call wrf_debug ( WARN , TRIM(msg))
2353 call DateCheck(DateStr,Status)
2354 if(Status /= WRF_NO_ERR) then
2355 write(msg,*) 'Warning DATE STRING ERROR |',DateStr,'| in ',__FILE__,', line', __LINE__
2356 call wrf_debug ( WARN , TRIM(msg))
2360 call GetDH(DataHandle,DH,Status)
2361 if(Status /= WRF_NO_ERR) then
2362 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2363 call wrf_debug ( WARN , TRIM(msg))
2368 write(msg,*)'ext_pnc_write_field: called for ',TRIM(Var)
2369 CALL wrf_debug( 100, msg )
2372 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2373 Length_native(1:NDim) = Length(1:NDim)
2374 Length_global(1:NDim) = DomainEnd(1:NDim)-DomainStart(1:NDim)+1
2376 call ExtOrder(MemoryOrder,Length,Status)
2377 call ExtOrder(MemoryOrder,Length_global,Status)
2379 call ExtOrderStr(MemoryOrder,DimNames,RODimNames,Status)
2381 ! Magic number to identify call from IO server when doing quilting
2382 quilting = (MemoryStart(1) == -998899 .AND. MemoryEnd(1) == -998899)
2384 lMemoryStart(1:NDim) = 1
2385 lMemoryEnd(1:NDim) = Length(1:NDim)
2387 lMemoryStart(1:NDim) = MemoryStart(1:NDim)
2388 lMemoryEnd(1:NDim) = MemoryEnd(1:NDim)
2391 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2392 Status = WRF_WARN_FILE_NOT_OPENED
2393 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2394 call wrf_debug ( WARN , TRIM(msg))
2395 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
2396 Status = WRF_WARN_WRITE_RONLY_FILE
2397 write(msg,*) 'Warning WRITE READ ONLY FILE in ',__FILE__,', line', __LINE__
2398 call wrf_debug ( WARN , TRIM(msg))
2399 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2401 if(DH%VarNames(NVar) == VarName ) then
2402 Status = WRF_WARN_2DRYRUNS_1VARIABLE
2403 write(msg,*) 'Warning 2 DRYRUNS 1 VARIABLE (',TRIM(VarName),') in ',__FILE__,', line', __LINE__
2404 call wrf_debug ( WARN , TRIM(msg))
2406 elseif(DH%VarNames(NVar) == NO_NAME) then
2407 DH%VarNames(NVar) = VarName
2410 elseif(NVar == MaxVars) then
2411 Status = WRF_WARN_TOO_MANY_VARIABLES
2412 write(msg,*) 'Warning TOO MANY VARIABLES in ',__FILE__,', line', __LINE__
2413 call wrf_debug ( WARN , TRIM(msg))
2418 if(RODimNames(j) == NullName .or. RODimNames(j) == '') then
2420 if(DH%DimLengths(i) == Length_global(j)) then
2422 elseif(DH%DimLengths(i) == NO_DIM) then
2423 stat = NFMPI_DEF_DIM(NCID,DH%DimNames(i),i2offset(Length_global(j)),DH%DimIDs(i))
2424 call netcdf_err(stat,Status)
2425 if(Status /= WRF_NO_ERR) then
2426 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2427 call wrf_debug ( WARN , TRIM(msg))
2430 DH%DimLengths(i) = Length_global(j)
2432 elseif(i == MaxDims) then
2433 Status = WRF_WARN_TOO_MANY_DIMS
2434 write(msg,*) 'Warning TOO MANY DIMENSIONS (',i,') in (',TRIM(VarName),') ',__FILE__,', line', __LINE__
2435 call wrf_debug ( WARN , TRIM(msg))
2439 else !look for input name and check if already defined
2442 if (DH%DimNames(i) == RODimNames(j)) then
2443 if (DH%DimLengths(i) == Length_global(j)) then
2447 Status = WRF_WARN_DIMNAME_REDEFINED
2448 write(msg,*) 'Warning DIM ',i,', NAME ',TRIM(DH%DimNames(i)),' REDEFINED by var ', &
2449 TRIM(Var),' ',DH%DimLengths(i),Length_global(j) ,' in ', __FILE__ ,' line', __LINE__
2450 call wrf_debug ( WARN , TRIM(msg))
2457 if (DH%DimLengths(i) == NO_DIM) then
2458 DH%DimNames(i) = RODimNames(j)
2459 stat = NFMPI_DEF_DIM(NCID,DH%DimNames(i),i2offset(Length_global(j)),DH%DimIDs(i))
2460 call netcdf_err(stat,Status)
2461 if(Status /= WRF_NO_ERR) then
2462 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2463 call wrf_debug ( WARN , TRIM(msg))
2466 DH%DimLengths(i) = Length_global(j)
2468 elseif(i == MaxDims) then
2469 Status = WRF_WARN_TOO_MANY_DIMS
2470 write(msg,*) 'Warning TOO MANY DIMENSIONS in ',__FILE__,', line', __LINE__
2471 call wrf_debug ( WARN , TRIM(msg))
2477 VDimIDs(j) = DH%DimIDs(i)
2478 DH%VarDimLens(j,NVar) = Length_global(j)
2480 VDimIDs(NDim+1) = DH%DimUnlimID
2481 select case (FieldType)
2491 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
2492 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
2493 call wrf_debug ( WARN , TRIM(msg))
2498 stat = NFMPI_DEF_VAR(NCID,VarName,XType,NDim+1,VDimIDs,VarID)
2499 call netcdf_err(stat,Status)
2500 if(Status /= WRF_NO_ERR) then
2501 write(msg,*) 'ext_pnc_write_field: NetCDF error for ',TRIM(VarName),' in ',__FILE__,', line', __LINE__
2502 call wrf_debug ( WARN , TRIM(msg))
2505 DH%VarIDs(NVar) = VarID
2506 stat = NFMPI_PUT_ATT_INT(NCID,VarID,'FieldType',NF_INT,i2offset(1),FieldType)
2507 call netcdf_err(stat,Status)
2508 if(Status /= WRF_NO_ERR) then
2509 write(msg,*) 'ext_pnc_write_field: NetCDF error in ',__FILE__,', line', __LINE__
2510 call wrf_debug ( WARN , TRIM(msg))
2513 call reorder(MemoryOrder,MemO)
2514 call uppercase(MemO,UCMemO)
2515 stat = NFMPI_PUT_ATT_TEXT(NCID,VarID,'MemoryOrder',i2offset(3),UCMemO)
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 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
2523 do NVar=1,DH%NumVars
2524 if(DH%VarNames(NVar) == VarName) then
2526 elseif(NVar == DH%NumVars) then
2527 Status = WRF_WARN_VAR_NF
2528 write(msg,*) 'Warning VARIABLE NOT FOUND in ',__FILE__,', line', __LINE__
2529 call wrf_debug ( WARN , TRIM(msg))
2533 VarID = DH%VarIDs(NVar)
2535 if(Length_global(j) /= DH%VarDimLens(j,NVar) .AND. DH%FileStatus /= WRF_FILE_OPENED_FOR_UPDATE ) then
2536 Status = WRF_WARN_WRTLEN_NE_DRRUNLEN
2537 write(msg,*) 'Warning LENGTH != DRY RUN LENGTH for |', &
2538 VarName,'| dim ',j,' in ',__FILE__,', line', __LINE__
2539 call wrf_debug ( WARN , TRIM(msg))
2540 write(msg,*) ' LENGTH ',Length_global(j),' DRY RUN LENGTH ',DH%VarDimLens(j,NVar)
2541 call wrf_debug ( WARN , TRIM(msg))
2543 !jm 061024 elseif(PatchStart(j) < MemoryStart(j)) then
2544 !jm elseif(DomainStart(j) < MemoryStart(j)) then
2545 elseif(PatchStart(j) < lMemoryStart(j)) then
2546 Status = WRF_WARN_DIMENSION_ERROR
2547 write(msg,*) 'Warning DIMENSION ERROR for |',VarName, &
2548 '| in ',__FILE__,', line', __LINE__
2549 call wrf_debug ( WARN , TRIM(msg))
2554 call GetIndices(NDim,lMemoryStart,lMemoryEnd,l1,l2,m1,m2,n1,n2)
2555 call GetIndices(NDim,StoredStart,Length ,x1,x2,y1,y2,z1,z2)
2556 call GetIndices(NDim,StoredStart,Length_native ,p1,p2,q1,q2,r1,r2)
2557 call GetIndices(NDim,PatchStart, PatchEnd ,i1,i2,j1,j2,k1,k2)
2559 if(FieldType == WRF_DOUBLE) di=2
2560 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
2562 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
2563 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
2564 call wrf_debug ( FATAL , TRIM(msg))
2569 WRITE(msg,*) 'ARPDBG: MemoryStart = ',lMemoryStart(1:NDim)
2570 CALL wrf_message(msg)
2571 WRITE(msg,*) 'ARPDBG: lMemoryEnd = ',lMemoryEnd(1:NDim)
2572 CALL wrf_message(msg)
2573 WRITE(msg,*) 'ARPDBG: Length = ',Length(1:NDim)
2574 CALL wrf_message(msg)
2578 ! Don't pass in PatchStart and PatchEnd here since we want to
2579 ! take transpose of whole patch of data which has been sent to
2580 ! the IO server and passed down to us.
2581 ! JM: the field and patch dimensions must be reordered or xpose is a noop
2582 call Transpose('write',MemoryOrder,di, Field,p1,p2,q1,q2,r1,r2 &
2583 ,XField,x1,x2,y1,y2,z1,z2 &
2584 ,p1,p2,q1,q2,r1,r2 )
2586 call Transpose('write',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2587 ,XField,x1,x2,y1,y2,z1,z2 &
2588 ,i1,i2,j1,j2,k1,k2 )
2591 StoredStart(1:NDim) = PatchStart(1:NDim)
2592 call ExtOrder(MemoryOrder,StoredStart,Status)
2593 call FieldIO('write',DataHandle,DateStr,StoredStart,Length,MemoryOrder, &
2594 FieldType,NCID,VarID,XField,Status)
2595 if(Status /= WRF_NO_ERR) then
2596 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2597 call wrf_debug ( WARN , TRIM(msg))
2600 deallocate(XField, STAT=stat)
2602 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
2603 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
2604 call wrf_debug ( FATAL , TRIM(msg))
2608 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2609 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2610 call wrf_debug ( FATAL , TRIM(msg))
2612 DH%first_operation = .FALSE.
2614 end subroutine ext_pnc_write_field
2616 subroutine ext_pnc_read_field(DataHandle,DateStr,Var,Field,FieldType,Comm, &
2617 IOComm, DomainDesc, MemoryOrdIn, Stagger, DimNames, &
2618 DomainStart,DomainEnd,MemoryStart,MemoryEnd,PatchStart,PatchEnd,Status)
2620 use ext_pnc_support_routines
2622 include 'wrf_status_codes.h'
2623 # include "pnetcdf.inc"
2624 integer ,intent(in) :: DataHandle
2625 character*(*) ,intent(in) :: DateStr
2626 character*(*) ,intent(in) :: Var
2627 integer ,intent(out) :: Field(*)
2628 integer ,intent(in) :: FieldType
2629 integer ,intent(inout) :: Comm
2630 integer ,intent(inout) :: IOComm
2631 integer ,intent(in) :: DomainDesc
2632 character*(*) ,intent(in) :: MemoryOrdIn
2633 character*(*) ,intent(in) :: Stagger ! Dummy for now
2634 character*(*) , dimension (*) ,intent(in) :: DimNames
2635 integer ,dimension(*) ,intent(in) :: DomainStart, DomainEnd
2636 integer ,dimension(*) ,intent(in) :: MemoryStart, MemoryEnd
2637 integer ,dimension(*) ,intent(in) :: PatchStart, PatchEnd
2638 integer ,intent(out) :: Status
2639 character (3) :: MemoryOrder
2640 character (NF_MAX_NAME) :: dimname
2641 type(wrf_data_handle) ,pointer :: DH
2644 character (VarNameLen) :: VarName
2646 integer ,dimension(NVarDims) :: VCount
2647 integer ,dimension(NVarDims) :: VStart
2648 integer ,dimension(NVarDims) :: Length
2649 integer ,dimension(NVarDims) :: VDimIDs
2650 integer ,dimension(NVarDims) :: MemS
2651 integer ,dimension(NVarDims) :: MemE
2652 integer ,dimension(NVarDims) :: StoredStart
2653 integer ,dimension(NVarDims) :: StoredLen
2654 integer(KIND=MPI_OFFSET_KIND) ,dimension(NVarDims) :: StoredLen_okind
2655 integer ,dimension(:,:,:,:) ,allocatable :: XField
2658 integer :: i1,i2,j1,j2,k1,k2
2659 integer :: x1,x2,y1,y2,z1,z2
2660 integer :: l1,l2,m1,m2,n1,n2
2661 character (VarNameLen) :: Name
2663 integer :: StoredDim
2665 integer(KIND=MPI_OFFSET_KIND) :: Len
2670 MemoryOrder = trim(adjustl(MemoryOrdIn))
2671 call GetDim(MemoryOrder,NDim,Status)
2672 if(Status /= WRF_NO_ERR) then
2673 write(msg,*) 'Warning BAD MEMORY ORDER |',TRIM(MemoryOrder),'| for |', &
2674 TRIM(Var),'| in ext_pnc_read_field ',__FILE__,', line', __LINE__
2675 call wrf_debug ( WARN , TRIM(msg))
2678 call DateCheck(DateStr,Status)
2679 if(Status /= WRF_NO_ERR) then
2680 write(msg,*) 'Warning DATE STRING ERROR |',TRIM(DateStr),'| for |',TRIM(Var), &
2681 '| in ext_pnc_read_field ',__FILE__,', line', __LINE__
2682 call wrf_debug ( WARN , TRIM(msg))
2686 call GetDH(DataHandle,DH,Status)
2687 if(Status /= WRF_NO_ERR) then
2688 write(msg,*) 'Warning Status = ',Status,' in ext_pnc_read_field ',__FILE__,', line', __LINE__
2689 call wrf_debug ( WARN , TRIM(msg))
2692 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2693 Status = WRF_WARN_FILE_NOT_OPENED
2694 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2695 call wrf_debug ( WARN , TRIM(msg))
2696 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2697 ! jm it is okay to have a dry run read. means read is called between ofrb and ofrc. Just return.
2698 ! Status = WRF_WARN_DRYRUN_READ
2699 ! write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
2700 ! call wrf_debug ( WARN , TRIM(msg))
2703 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
2704 Status = WRF_WARN_READ_WONLY_FILE
2705 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
2706 call wrf_debug ( WARN , TRIM(msg))
2707 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
2710 Length(1:NDim) = PatchEnd(1:NDim)-PatchStart(1:NDim)+1
2711 StoredStart(1:NDim) = PatchStart(1:NDim)
2713 call ExtOrder(MemoryOrder,Length,Status)
2715 stat = NFMPI_INQ_VARID(NCID,VarName,VarID)
2716 call netcdf_err(stat,Status)
2717 if(Status /= WRF_NO_ERR) then
2718 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__,' Varname ',Varname
2719 call wrf_debug ( WARN , TRIM(msg))
2722 stat = NFMPI_INQ_VAR(NCID,VarID,Name,XType,StoredDim,VDimIDs,NAtts)
2723 call netcdf_err(stat,Status)
2724 if(Status /= WRF_NO_ERR) then
2725 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2726 call wrf_debug ( WARN , TRIM(msg))
2729 stat = NFMPI_GET_ATT_INT(NCID,VarID,'FieldType',FType)
2730 call netcdf_err(stat,Status)
2731 if(Status /= WRF_NO_ERR) then
2732 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2733 call wrf_debug ( WARN , TRIM(msg))
2736 ! allow coercion between double and single prec real
2737 !jm if(FieldType /= Ftype) then
2738 if( (FieldType == WRF_REAL .OR. FieldType == WRF_DOUBLE) ) then
2739 if ( .NOT. (Ftype == WRF_REAL .OR. Ftype == WRF_DOUBLE )) then
2740 Status = WRF_WARN_TYPE_MISMATCH
2741 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
2742 call wrf_debug ( WARN , TRIM(msg))
2745 else if(FieldType /= Ftype) then
2746 Status = WRF_WARN_TYPE_MISMATCH
2747 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
2748 call wrf_debug ( WARN , TRIM(msg))
2751 select case (FieldType)
2753 ! allow coercion between double and single prec real
2754 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
2755 Status = WRF_WARN_TYPE_MISMATCH
2756 write(msg,*) 'Warning REAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
2759 ! allow coercion between double and single prec real
2760 if(.NOT. (XType == NF_FLOAT .OR. XType == NF_DOUBLE) ) then
2761 Status = WRF_WARN_TYPE_MISMATCH
2762 write(msg,*) 'Warning DOUBLE TYPE MISMATCH in ',__FILE__,', line', __LINE__
2765 if(XType /= NF_INT) then
2766 Status = WRF_WARN_TYPE_MISMATCH
2767 write(msg,*) 'Warning INTEGER TYPE MISMATCH in ',__FILE__,', line', __LINE__
2770 if(XType /= NF_INT) then
2771 Status = WRF_WARN_TYPE_MISMATCH
2772 write(msg,*) 'Warning LOGICAL TYPE MISMATCH in ',__FILE__,', line', __LINE__
2775 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
2776 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
2778 if(Status /= WRF_NO_ERR) then
2779 call wrf_debug ( WARN , TRIM(msg))
2782 ! NDim=0 for scalars. Handle read of old NDim=1 files. TBH: 20060502
2783 IF ( ( NDim == 0 ) .AND. ( StoredDim == 2 ) ) THEN
2784 stat = NFMPI_INQ_DIMNAME(NCID,VDimIDs(1),dimname)
2785 call netcdf_err(stat,Status)
2786 if(Status /= WRF_NO_ERR) then
2787 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2788 call wrf_debug ( WARN , TRIM(msg))
2791 IF ( dimname(1:10) == 'ext_scalar' ) THEN
2796 if(StoredDim /= NDim+1) then
2797 Status = WRF_ERR_FATAL_BAD_VARIABLE_DIM
2798 write(msg,*) 'Fatal error BAD VARIABLE DIMENSION in ext_pnc_read_field ',TRIM(Var),TRIM(DateStr)
2799 call wrf_debug ( FATAL , msg)
2800 write(msg,*) ' StoredDim ', StoredDim, ' .NE. NDim+1 ', NDim+1
2801 call wrf_debug ( FATAL , msg)
2805 stat = NFMPI_INQ_DIMLEN(NCID,VDimIDs(j),StoredLen_okind(j))
2806 call netcdf_err(stat,Status)
2807 if(Status /= WRF_NO_ERR) then
2808 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
2809 call wrf_debug ( WARN , TRIM(msg))
2812 StoredLen(j) = StoredLen_okind(j)
2813 if(Length(j) > StoredLen(j)) then
2814 Status = WRF_WARN_READ_PAST_EOF
2815 write(msg,*) 'Warning READ PAST EOF in ext_pnc_read_field of ',TRIM(Var),Length(j),'>',StoredLen(j)
2816 call wrf_debug ( WARN , TRIM(msg))
2818 elseif(Length(j) <= 0) then
2819 Status = WRF_WARN_ZERO_LENGTH_READ
2820 write(msg,*) 'Warning ZERO LENGTH READ in ',__FILE__,', line', __LINE__
2821 call wrf_debug ( WARN , TRIM(msg))
2827 call GetIndices(NDim,MemoryStart,MemoryEnd,l1,l2,m1,m2,n1,n2)
2828 call GetIndices(NDim,StoredStart,Length,x1,x2,y1,y2,z1,z2)
2829 !jm call GetIndices(NDim,DomainStart,DomainEnd,i1,i2,j1,j2,k1,k2)
2830 call GetIndices(NDim,PatchStart,PatchEnd,i1,i2,j1,j2,k1,k2)
2832 StoredStart(1:NDim) = PatchStart(1:NDim)
2833 call ExtOrder(MemoryOrder,StoredStart,Status)
2836 if(FieldType == WRF_DOUBLE) di=2
2837 allocate(XField(di,x1:x2,y1:y2,z1:z2), STAT=stat)
2839 Status = WRF_ERR_FATAL_ALLOCATION_ERROR
2840 write(msg,*) 'Fatal ALLOCATION ERROR in ',__FILE__,', line', __LINE__
2841 call wrf_debug ( FATAL , msg)
2844 call FieldIO('read',DataHandle,DateStr,StoredStart,Length,MemoryOrder, &
2845 FieldType,NCID,VarID,XField,Status)
2846 if(Status /= WRF_NO_ERR) then
2847 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2848 call wrf_debug ( WARN , TRIM(msg))
2851 call Transpose('read',MemoryOrder,di, Field,l1,l2,m1,m2,n1,n2 &
2852 ,XField,x1,x2,y1,y2,z1,z2 &
2853 ,i1,i2,j1,j2,k1,k2 )
2854 deallocate(XField, STAT=stat)
2856 Status = WRF_ERR_FATAL_DEALLOCATION_ERR
2857 write(msg,*) 'Fatal DEALLOCATION ERROR in ',__FILE__,', line', __LINE__
2858 call wrf_debug ( FATAL , msg)
2862 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2863 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2864 call wrf_debug ( FATAL , msg)
2866 DH%first_operation = .FALSE.
2868 end subroutine ext_pnc_read_field
2870 subroutine ext_pnc_inquire_opened( DataHandle, FileName , FileStatus, Status )
2872 use ext_pnc_support_routines
2874 include 'wrf_status_codes.h'
2875 integer ,intent(in) :: DataHandle
2876 character*(*) ,intent(in) :: FileName
2877 integer ,intent(out) :: FileStatus
2878 integer ,intent(out) :: Status
2879 type(wrf_data_handle) ,pointer :: DH
2881 call GetDH(DataHandle,DH,Status)
2882 if(Status /= WRF_NO_ERR) then
2883 FileStatus = WRF_FILE_NOT_OPENED
2886 if(FileName /= DH%FileName) then
2887 FileStatus = WRF_FILE_NOT_OPENED
2889 FileStatus = DH%FileStatus
2893 end subroutine ext_pnc_inquire_opened
2895 subroutine ext_pnc_inquire_filename( Datahandle, FileName, FileStatus, Status )
2897 use ext_pnc_support_routines
2899 include 'wrf_status_codes.h'
2900 integer ,intent(in) :: DataHandle
2901 character*(*) ,intent(out) :: FileName
2902 integer ,intent(out) :: FileStatus
2903 integer ,intent(out) :: Status
2904 type(wrf_data_handle) ,pointer :: DH
2905 FileStatus = WRF_FILE_NOT_OPENED
2906 call GetDH(DataHandle,DH,Status)
2907 if(Status /= WRF_NO_ERR) then
2908 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2909 call wrf_debug ( WARN , TRIM(msg))
2912 FileName = DH%FileName
2913 FileStatus = DH%FileStatus
2916 end subroutine ext_pnc_inquire_filename
2918 subroutine ext_pnc_set_time(DataHandle, DateStr, Status)
2920 use ext_pnc_support_routines
2922 include 'wrf_status_codes.h'
2923 integer ,intent(in) :: DataHandle
2924 character*(*) ,intent(in) :: DateStr
2925 integer ,intent(out) :: Status
2926 type(wrf_data_handle) ,pointer :: DH
2929 call DateCheck(DateStr,Status)
2930 if(Status /= WRF_NO_ERR) then
2931 write(msg,*) 'Warning DATE STRING ERROR in ',__FILE__,', line', __LINE__
2932 call wrf_debug ( WARN , TRIM(msg))
2935 call GetDH(DataHandle,DH,Status)
2936 if(Status /= WRF_NO_ERR) then
2937 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2938 call wrf_debug ( WARN , TRIM(msg))
2941 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2942 Status = WRF_WARN_FILE_NOT_OPENED
2943 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2944 call wrf_debug ( WARN , TRIM(msg))
2945 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2946 Status = WRF_WARN_FILE_NOT_COMMITTED
2947 write(msg,*) 'Warning FILE NOT COMMITTED in ',__FILE__,', line', __LINE__
2948 call wrf_debug ( WARN , TRIM(msg))
2949 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
2950 Status = WRF_WARN_READ_WONLY_FILE
2951 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
2952 call wrf_debug ( WARN , TRIM(msg))
2953 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
2955 if(DH%Times(i)==DateStr) then
2959 if(i==MaxTimes) then
2960 Status = WRF_WARN_TIME_NF
2964 DH%CurrentVariable = 0
2967 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
2968 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
2969 call wrf_debug ( FATAL , msg)
2972 end subroutine ext_pnc_set_time
2974 subroutine ext_pnc_get_next_time(DataHandle, DateStr, Status)
2976 use ext_pnc_support_routines
2978 include 'wrf_status_codes.h'
2979 integer ,intent(in) :: DataHandle
2980 character*(*) ,intent(out) :: DateStr
2981 integer ,intent(out) :: Status
2982 type(wrf_data_handle) ,pointer :: DH
2984 call GetDH(DataHandle,DH,Status)
2985 if(Status /= WRF_NO_ERR) then
2986 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
2987 call wrf_debug ( WARN , TRIM(msg))
2990 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
2991 Status = WRF_WARN_FILE_NOT_OPENED
2992 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
2993 call wrf_debug ( WARN , TRIM(msg))
2994 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
2995 Status = WRF_WARN_DRYRUN_READ
2996 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
2997 call wrf_debug ( WARN , TRIM(msg))
2998 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
2999 Status = WRF_WARN_READ_WONLY_FILE
3000 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3001 call wrf_debug ( WARN , TRIM(msg))
3002 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE ) then
3003 if(DH%CurrentTime >= DH%NumberTimes) then
3004 write(msg,*) 'Warning ext_pnc_get_next_time: DH%CurrentTime >= DH%NumberTimes ',DH%CurrentTime,DH%NumberTimes
3005 call wrf_debug ( WARN , TRIM(msg))
3006 Status = WRF_WARN_TIME_EOF
3009 DH%CurrentTime = DH%CurrentTime +1
3010 DateStr = DH%Times(DH%CurrentTime)
3011 DH%CurrentVariable = 0
3014 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3015 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3016 call wrf_debug ( FATAL , msg)
3019 end subroutine ext_pnc_get_next_time
3021 subroutine ext_pnc_get_previous_time(DataHandle, DateStr, Status)
3023 use ext_pnc_support_routines
3025 include 'wrf_status_codes.h'
3026 integer ,intent(in) :: DataHandle
3027 character*(*) ,intent(out) :: DateStr
3028 integer ,intent(out) :: Status
3029 type(wrf_data_handle) ,pointer :: DH
3031 call GetDH(DataHandle,DH,Status)
3032 if(Status /= WRF_NO_ERR) then
3033 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3034 call wrf_debug ( WARN , TRIM(msg))
3037 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3038 Status = WRF_WARN_FILE_NOT_OPENED
3039 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3040 call wrf_debug ( WARN , TRIM(msg))
3041 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3042 Status = WRF_WARN_DRYRUN_READ
3043 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3044 call wrf_debug ( WARN , TRIM(msg))
3045 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3046 Status = WRF_WARN_READ_WONLY_FILE
3047 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3048 call wrf_debug ( WARN , TRIM(msg))
3049 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
3050 if(DH%CurrentTime.GT.0) then
3051 DH%CurrentTime = DH%CurrentTime -1
3053 DateStr = DH%Times(DH%CurrentTime)
3054 DH%CurrentVariable = 0
3057 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3058 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3059 call wrf_debug ( FATAL , msg)
3062 end subroutine ext_pnc_get_previous_time
3064 subroutine ext_pnc_get_next_var(DataHandle, VarName, Status)
3066 use ext_pnc_support_routines
3068 include 'wrf_status_codes.h'
3069 # include "pnetcdf.inc"
3070 integer ,intent(in) :: DataHandle
3071 character*(*) ,intent(out) :: VarName
3072 integer ,intent(out) :: Status
3073 type(wrf_data_handle) ,pointer :: DH
3075 character (80) :: Name
3077 call GetDH(DataHandle,DH,Status)
3078 if(Status /= WRF_NO_ERR) then
3079 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3080 call wrf_debug ( WARN , TRIM(msg))
3083 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3084 Status = WRF_WARN_FILE_NOT_OPENED
3085 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3086 call wrf_debug ( WARN , TRIM(msg))
3087 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3088 Status = WRF_WARN_DRYRUN_READ
3089 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3090 call wrf_debug ( WARN , TRIM(msg))
3091 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3092 Status = WRF_WARN_READ_WONLY_FILE
3093 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3094 call wrf_debug ( WARN , TRIM(msg))
3095 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3097 DH%CurrentVariable = DH%CurrentVariable +1
3098 if(DH%CurrentVariable > DH%NumVars) then
3099 Status = WRF_WARN_VAR_EOF
3102 VarName = DH%VarNames(DH%CurrentVariable)
3105 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3106 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3107 call wrf_debug ( FATAL , msg)
3110 end subroutine ext_pnc_get_next_var
3112 subroutine ext_pnc_end_of_frame(DataHandle, Status)
3114 use ext_pnc_support_routines
3116 # include "pnetcdf.inc"
3117 include 'wrf_status_codes.h'
3118 integer ,intent(in) :: DataHandle
3119 integer ,intent(out) :: Status
3120 type(wrf_data_handle) ,pointer :: DH
3122 call GetDH(DataHandle,DH,Status)
3124 end subroutine ext_pnc_end_of_frame
3126 ! NOTE: For scalar variables NDim is set to zero and DomainStart and
3127 ! NOTE: DomainEnd are left unmodified.
3128 subroutine ext_pnc_get_var_info(DataHandle,Name,NDim,MemoryOrder,Stagger,DomainStart,DomainEnd,WrfType,Status)
3130 use ext_pnc_support_routines
3132 # include "pnetcdf.inc"
3133 include 'wrf_status_codes.h'
3134 integer ,intent(in) :: DataHandle
3135 character*(*) ,intent(in) :: Name
3136 integer ,intent(out) :: NDim
3137 character*(*) ,intent(out) :: MemoryOrder
3138 character*(*) :: Stagger ! Dummy for now
3139 integer ,dimension(*) ,intent(out) :: DomainStart, DomainEnd
3140 integer ,intent(out) :: WrfType
3141 integer ,intent(out) :: Status
3142 type(wrf_data_handle) ,pointer :: DH
3144 integer ,dimension(NVarDims) :: VDimIDs
3149 call GetDH(DataHandle,DH,Status)
3150 if(Status /= WRF_NO_ERR) then
3151 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
3152 call wrf_debug ( WARN , TRIM(msg))
3155 if(DH%FileStatus == WRF_FILE_NOT_OPENED) then
3156 Status = WRF_WARN_FILE_NOT_OPENED
3157 write(msg,*) 'Warning FILE NOT OPENED in ',__FILE__,', line', __LINE__
3158 call wrf_debug ( WARN , TRIM(msg))
3160 elseif(DH%FileStatus == WRF_FILE_OPENED_NOT_COMMITTED) then
3161 Status = WRF_WARN_DRYRUN_READ
3162 write(msg,*) 'Warning DRYRUN READ in ',__FILE__,', line', __LINE__
3163 call wrf_debug ( WARN , TRIM(msg))
3165 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_WRITE) then
3166 Status = WRF_WARN_READ_WONLY_FILE
3167 write(msg,*) 'Warning READ WRITE ONLY FILE in ',__FILE__,', line', __LINE__
3168 call wrf_debug ( WARN , TRIM(msg))
3170 elseif(DH%FileStatus == WRF_FILE_OPENED_FOR_READ .OR. DH%FileStatus == WRF_FILE_OPENED_FOR_UPDATE) then
3171 stat = NFMPI_INQ_VARID(DH%NCID,Name,VarID)
3172 call netcdf_err(stat,Status)
3173 if(Status /= WRF_NO_ERR) then
3174 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3175 call wrf_debug ( WARN , TRIM(msg))
3178 stat = NFMPI_INQ_VARTYPE(DH%NCID,VarID,XType)
3179 call netcdf_err(stat,Status)
3180 if(Status /= WRF_NO_ERR) then
3181 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3182 call wrf_debug ( WARN , TRIM(msg))
3185 stat = NFMPI_GET_ATT_INT(DH%NCID,VarID,'FieldType',WrfType)
3186 call netcdf_err(stat,Status)
3187 if(Status /= WRF_NO_ERR) then
3188 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3189 call wrf_debug ( WARN , TRIM(msg))
3194 Status = WRF_WARN_BAD_DATA_TYPE
3195 write(msg,*) 'Warning BYTE IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3196 call wrf_debug ( WARN , TRIM(msg))
3199 Status = WRF_WARN_BAD_DATA_TYPE
3200 write(msg,*) 'Warning CHAR IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3201 call wrf_debug ( WARN , TRIM(msg))
3204 Status = WRF_WARN_BAD_DATA_TYPE
3205 write(msg,*) 'Warning SHORT IS BAD DATA TYPE in ',__FILE__,', line', __LINE__
3206 call wrf_debug ( WARN , TRIM(msg))
3209 if(WrfType /= WRF_INTEGER .and. WrfType /= WRF_LOGICAL) then
3210 Status = WRF_WARN_BAD_DATA_TYPE
3211 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3212 call wrf_debug ( WARN , TRIM(msg))
3216 if(WrfType /= WRF_REAL) then
3217 Status = WRF_WARN_BAD_DATA_TYPE
3218 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3219 call wrf_debug ( WARN , TRIM(msg))
3223 if(WrfType /= WRF_DOUBLE) then
3224 Status = WRF_WARN_BAD_DATA_TYPE
3225 write(msg,*) 'Warning BAD DATA TYPE in ',__FILE__,', line', __LINE__
3226 call wrf_debug ( WARN , TRIM(msg))
3230 Status = WRF_WARN_DATA_TYPE_NOT_FOUND
3231 write(msg,*) 'Warning DATA TYPE NOT FOUND in ',__FILE__,', line', __LINE__
3232 call wrf_debug ( WARN , TRIM(msg))
3236 stat = NFMPI_GET_ATT_TEXT(DH%NCID,VarID,'MemoryOrder',MemoryOrder)
3237 call netcdf_err(stat,Status)
3238 if(Status /= WRF_NO_ERR) then
3239 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3240 call wrf_debug ( WARN , TRIM(msg))
3243 call GetDim(MemoryOrder,NDim,Status)
3244 if(Status /= WRF_NO_ERR) then
3245 write(msg,*) 'Warning BAD MEMORY ORDER ',TRIM(MemoryOrder),' in ',__FILE__,', line', __LINE__
3246 call wrf_debug ( WARN , TRIM(msg))
3249 stat = NFMPI_INQ_VARDIMID(DH%NCID,VarID,VDimIDs)
3250 call netcdf_err(stat,Status)
3251 if(Status /= WRF_NO_ERR) then
3252 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3253 call wrf_debug ( WARN , TRIM(msg))
3258 stat = NFMPI_INQ_DIMLEN(DH%NCID,VDimIDs(j),DomainEnd(j))
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 Status = WRF_ERR_FATAL_BAD_FILE_STATUS
3268 write(msg,*) 'Fatal error BAD FILE STATUS in ',__FILE__,', line', __LINE__
3269 call wrf_debug ( FATAL , msg)
3272 end subroutine ext_pnc_get_var_info
3274 subroutine ext_pnc_warning_str( Code, ReturnString, Status)
3276 use ext_pnc_support_routines
3278 # include "pnetcdf.inc"
3279 include 'wrf_status_codes.h'
3281 integer , intent(in) ::Code
3282 character *(*), intent(out) :: ReturnString
3283 integer, intent(out) ::Status
3287 ReturnString='No error'
3291 ReturnString= 'File not found (or file is incomplete)'
3295 ReturnString='Metadata not found'
3299 ReturnString= 'Timestamp not found'
3303 ReturnString= 'No more timestamps'
3307 ReturnString= 'Variable not found'
3311 ReturnString= 'No more variables for the current time'
3315 ReturnString= 'Too many open files'
3319 ReturnString= 'Data type mismatch'
3323 ReturnString= 'Attempt to write read-only file'
3327 ReturnString= 'Attempt to read write-only file'
3331 ReturnString= 'Attempt to access unopened file'
3335 ReturnString= 'Attempt to do 2 trainings for 1 variable'
3339 ReturnString= 'Attempt to read past EOF'
3343 ReturnString= 'Bad data handle'
3347 ReturnString= 'Write length not equal to training length'
3351 ReturnString= 'More dimensions requested than training'
3355 ReturnString= 'Attempt to read more data than exists'
3359 ReturnString= 'Input dimensions inconsistent'
3363 ReturnString= 'Input MemoryOrder not recognized'
3367 ReturnString= 'A dimension name with 2 different lengths'
3371 ReturnString= 'String longer than provided storage'
3375 ReturnString= 'Function not supportable'
3379 ReturnString= 'Package implements this routine as NOOP'
3383 !netcdf-specific warning messages
3385 ReturnString= 'Bad data type'
3389 ReturnString= 'File not committed'
3393 ReturnString= 'File is opened for reading'
3397 ReturnString= 'Attempt to write metadata after open commit'
3401 ReturnString= 'I/O not initialized'
3405 ReturnString= 'Too many variables requested'
3409 ReturnString= 'Attempt to close file during a dry run'
3413 ReturnString= 'Date string not 19 characters in length'
3417 ReturnString= 'Attempt to read zero length words'
3421 ReturnString= 'Data type not found'
3425 ReturnString= 'Badly formatted date string'
3429 ReturnString= 'Attempt at read during a dry run'
3433 ReturnString= 'Attempt to get zero words'
3437 ReturnString= 'Attempt to put zero length words'
3441 ReturnString= 'NetCDF error'
3445 ReturnString= 'Requested length <= 1'
3449 ReturnString= 'More data available than requested'
3453 ReturnString= 'New date less than previous date'
3458 ReturnString= 'This warning code is not supported or handled directly by WRF and NetCDF. &
3459 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
3460 & to be calling a package-specific routine to return a message for this warning code.'
3465 end subroutine ext_pnc_warning_str
3468 !returns message string for all WRF and netCDF warning/error status codes
3469 !Other i/o packages must provide their own routines to return their own status messages
3470 subroutine ext_pnc_error_str( Code, ReturnString, Status)
3472 use ext_pnc_support_routines
3474 # include "pnetcdf.inc"
3475 include 'wrf_status_codes.h'
3477 integer , intent(in) ::Code
3478 character *(*), intent(out) :: ReturnString
3479 integer, intent(out) ::Status
3483 ReturnString= 'Allocation Error'
3487 ReturnString= 'Deallocation Error'
3491 ReturnString= 'Bad File Status'
3495 ReturnString= 'Variable on disk is not 3D'
3499 ReturnString= 'Metadata on disk is not 1D'
3503 ReturnString= 'Time dimension too small'
3507 ReturnString= 'This error code is not supported or handled directly by WRF and NetCDF. &
3508 & Might be an erroneous number, or specific to an i/o package other than NetCDF; you may need &
3509 & to be calling a package-specific routine to return a message for this error code.'
3514 end subroutine ext_pnc_error_str
3518 subroutine ext_pnc_end_independent_mode(DataHandle, Status)
3520 use ext_pnc_support_routines
3521 include 'wrf_status_codes.h'
3522 # include "pnetcdf.inc"
3523 integer ,intent(in) :: DataHandle
3524 integer ,intent(out) :: Status
3525 type(wrf_data_handle) ,pointer :: DH
3528 DH => WrfDataHandles(DataHandle)
3529 DH%Collective = .TRUE.
3530 stat = NFMPI_END_INDEP_DATA(DH%NCID)
3532 call netcdf_err(stat,Status)
3533 if(Status /= WRF_NO_ERR) then
3534 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3535 call wrf_debug ( WARN , TRIM(msg))
3539 end subroutine ext_pnc_end_independent_mode
3541 subroutine ext_pnc_start_independent_mode(DataHandle, Status)
3543 use ext_pnc_support_routines
3544 include 'wrf_status_codes.h'
3545 # include "pnetcdf.inc"
3546 integer ,intent(in) :: DataHandle
3547 integer ,intent(out) :: Status
3548 type(wrf_data_handle) ,pointer :: DH
3551 DH => WrfDataHandles(DataHandle)
3552 DH%Collective = .FALSE.
3553 stat = NFMPI_BEGIN_INDEP_DATA(DH%NCID)
3554 call netcdf_err(stat,Status)
3555 if(Status /= WRF_NO_ERR) then
3556 write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
3557 call wrf_debug ( WARN , TRIM(msg))
3562 end subroutine ext_pnc_start_independent_mode