/[MITgcm]/MITgcm/pkg/mdsio/mdsio_wr_rec_rl.F
ViewVC logotype

Contents of /MITgcm/pkg/mdsio/mdsio_wr_rec_rl.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (show annotations) (download)
Tue Dec 30 02:00:51 2008 UTC (15 years, 4 months ago) by jmc
Branch: MAIN
Read/Write one reccord from already opened io-unit (S/R:MDS_RD/WR_REC_RL/RS)
- replace set of 4 S/R MDS_READ/WRITE_RL/RS_VEC
- honor filePrec option ; _BYTESWAPIO implemented

1 C $Header: /u/gcmpack/MITgcm/pkg/mdsio/mdsio_read_rs_vec.F,v 1.3 2004/11/30 16:11:10 heimbach Exp $
2 C $Name: $
3
4 #include "MDSIO_OPTIONS.h"
5 #define REPRODUCE_BOTTOM_CTRL_OLDRESULTS
6
7 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
8 CBOP 0
9 C !ROUTINE: MDS_WR_REC_RL
10
11 C !INTERFACE:
12 SUBROUTINE MDS_WR_REC_RL(
13 I arr,
14 O r4Buf, r8Buf,
15 I fPrec, dUnit, iRec, nArr, myThid )
16
17 C !DESCRIPTION:
18 C Write one reccord to already opened io-unit "dUnit", from RL array "arr"
19
20 C !USES:
21 IMPLICIT NONE
22 #include "EEPARAMS.h"
23 #include "SIZE.h"
24 #include "PARAMS.h"
25
26 C !INPUT PARAMETERS:
27 C arr RL :: vector array to write
28 C fPrec integer :: file precision
29 C dUnit integer :: 'Opened' I/O channel
30 C iRec integer :: record number to WRITE
31 C nArr integer :: dimension off array "arr"
32 C myThid integer :: my Thread Id number
33 C !OUTPUT PARAMETERS:
34 C r4Buf real*4 :: buffer array
35 C r8Buf real*8 :: buffer array
36 INTEGER fPrec
37 INTEGER dUnit
38 INTEGER iRec
39 INTEGER nArr
40 INTEGER myThid
41 _RL arr(nArr)
42 Real*4 r4Buf(nArr)
43 Real*8 r8Buf(nArr)
44 CEOP
45
46 C !LOCAL VARIABLES:
47 CHARACTER*(MAX_LEN_MBUF) msgBuf
48 INTEGER k
49
50 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
51 IF ( debugLevel.GT.debLevB ) THEN
52 WRITE(msgBuf,'(A,I9,2x,I9)')
53 & ' MDS_WR_REC_RL: iRec,Dim = ', iRec, nArr
54 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
55 & SQUEEZE_RIGHT , myThid )
56 ENDIF
57
58 IF ( fPrec.EQ.precFloat32 ) THEN
59 DO k=1,nArr
60 r4Buf(k) = arr(k)
61 ENDDO
62 #ifdef _BYTESWAPIO
63 CALL MDS_BYTESWAPR4( nArr, r4Buf )
64 #endif
65 WRITE( dUnit, rec=iRec ) r4Buf
66 ELSEIF ( fPrec.EQ.precFloat64 ) THEN
67 #ifdef REPRODUCE_BOTTOM_CTRL_OLDRESULTS
68 WRITE( dUnit, rec=iRec ) arr
69 #else /* REPRODUCE_BOTTOM_CTRL_OLDRESULTS */
70 DO k=1,nArr
71 r8Buf(k) = arr(k)
72 ENDDO
73 #ifdef _BYTESWAPIO
74 CALL MDS_BYTESWAPR8( nArr, r8Buf )
75 #endif
76 WRITE( dUnit, rec=iRec ) r8Buf
77 #endif /* REPRODUCE_BOTTOM_CTRL_OLDRESULTS */
78 ELSE
79 WRITE(msgBuf,'(A,I9)')
80 & ' MDS_WR_REC_RL: illegal value for fPrec=',fPrec
81 CALL PRINT_ERROR( msgBuf, myThid )
82 STOP 'ABNORMAL END: S/R MDS_WR_REC_RL'
83 ENDIF
84
85 RETURN
86 END

  ViewVC Help
Powered by ViewVC 1.1.22