/[MITgcm]/MITgcm/pkg/layers/layers_readparms.F
ViewVC logotype

Annotation of /MITgcm/pkg/layers/layers_readparms.F

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


Revision 1.8 - (hide annotations) (download)
Wed Jan 9 00:03:08 2013 UTC (11 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64x, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64c, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f
Changes since 1.7: +19 -17 lines
ensure that "layers_bolus" switch is False when useGMRedi=F

1 jmc 1.8 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_readparms.F,v 1.7 2013/01/07 02:07:55 dfer Exp $
2 rpa 1.1 C $Name: $
3    
4     #include "LAYERS_OPTIONS.h"
5    
6     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
7    
8     SUBROUTINE LAYERS_READPARMS( myThid )
9    
10     C Read LAYERS parameters from data file.
11    
12     IMPLICIT NONE
13     #include "SIZE.h"
14     #include "EEPARAMS.h"
15     #include "PARAMS.h"
16     #include "LAYERS_SIZE.h"
17     #include "LAYERS.h"
18    
19     C INPUT PARAMETERS:
20     INTEGER myThid
21    
22     #ifdef ALLOW_LAYERS
23    
24     NAMELIST /LAYERS_PARM01/
25 dfer 1.3 & layers_G, layers_taveFreq, layers_diagFreq,
26 gforget 1.6 & LAYER_nb, layers_kref, useBOLUS, layers_bolus,
27     & layers_name, layers_bounds, layers_krho
28 rpa 1.1
29     C === Local variables ===
30 jmc 1.8 C msgBuf - Informational/error message buffer
31 rpa 1.1 C iUnit - Work variable for IO unit number
32     C k - index
33     CHARACTER*(MAX_LEN_MBUF) msgBuf
34 gforget 1.6 INTEGER iUnit, k, iLa
35 rpa 1.1
36     _BEGIN_MASTER(myThid)
37    
38     C-- Default values for LAYERS
39    
40     C The MNC stuff is not working yet
41     layers_MNC = .FALSE.
42     layers_MDSIO = .TRUE.
43    
44 gforget 1.6 DO iLa=1,layers_maxNum
45     layers_name(iLa) = ' '
46     layers_num(iLa) = 0
47     layers_krho(iLa)= 1
48 jmc 1.8 layers_bolus(iLa) = useGMRedi
49 gforget 1.6 DO k=1,Nlayers+1
50     layers_bounds(k,iLa) = UNSET_RL
51     ENDDO
52     ENDDO
53 jmc 1.8
54 rpa 1.1 DO k=1,Nlayers+1
55     layers_G(k) = UNSET_RL
56     ENDDO
57 dfer 1.2 layers_taveFreq = taveFreq
58     layers_diagFreq = dumpFreq
59 gforget 1.6 LAYER_nb = 0
60 dfer 1.5 layers_kref = 1
61 jmc 1.8 useBOLUS = useGMRedi
62 rpa 1.1
63     WRITE(msgBuf,'(A)') 'LAYERS_READPARMS: opening data.layers'
64     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
65     & SQUEEZE_RIGHT , 1)
66     CALL OPEN_COPY_DATA_FILE(
67     I 'data.layers', 'LAYERS_READPARMS',
68     O iUnit,
69     I myThid )
70    
71 gforget 1.6
72 rpa 1.1 C Read parameters from open data file
73     READ(UNIT=iUnit,NML=LAYERS_PARM01)
74     WRITE(msgBuf,'(A)')
75     & 'LAYERS_READPARMS: finished reading data.layers'
76     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
77     & SQUEEZE_RIGHT , 1)
78     C Close the open data file
79     CLOSE(iUnit)
80    
81 gforget 1.6 c revert to old approach
82 jmc 1.8 c (unless I simply retore the old ways, I need a
83 gforget 1.6 c print statement or two as I override the new way)
84     IF ( LAYER_nb.NE.0 ) THEN
85     layers_num(1) = LAYER_nb
86     layers_krho(1) = layers_kref
87 dfer 1.7 layers_bolus(1) = useBOLUS
88 gforget 1.6 DO k=1,Nlayers+1
89     layers_bounds(k,1) = layers_G(k)
90     ENDDO
91     DO iLa=2,layers_maxNum
92     layers_num(iLa) = 0
93     layers_name(iLa) = ' '
94     layers_krho(iLa) = 1
95 jmc 1.8 layers_bolus(iLa) = useGMRedi
96 gforget 1.6 DO k=1,Nlayers+1
97     layers_bounds(k,iLa) = UNSET_RL
98     ENDDO
99     ENDDO
100     ENDIF
101    
102     C-- ensure layers_name/layers_num setup consistency
103     DO iLa=1,layers_maxNum
104     IF ( ( layers_name(iLa).EQ.'TH' ).OR.
105     & ( layers_num(iLa).EQ.1 ) ) THEN
106     layers_name(iLa)='TH'
107     layers_num(iLa)=1
108     ELSEIF ( ( layers_name(iLa).EQ.'SLT' ).OR.
109     & ( layers_num(iLa).EQ.2 ) ) THEN
110     layers_name(iLa)='SLT'
111     layers_num(iLa)=2
112     ELSEIF ( ( layers_name(iLa).EQ.'RHO' ).OR.
113     & ( layers_num(iLa).EQ.3 ) ) THEN
114     layers_name(iLa)='RHO'
115     layers_num(iLa)=3
116     ELSE
117     layers_name(iLa)=' '
118     layers_num(iLa)=0
119 jmc 1.8 ENDIF
120     C-- bolus contribution only available if using GMRedi
121     layers_bolus(iLa) = layers_bolus(iLa) .AND. useGMRedi
122 gforget 1.6 ENDDO
123 jmc 1.8
124 gforget 1.6 C-- Make sure the layers_bounds we just read is big enough
125     DO iLa=1,layers_maxNum
126 jmc 1.8 IF ( layers_num(iLa).NE.0 ) THEN
127     DO k=1,Nlayers+1
128     IF ( layers_bounds(k,iLa) .EQ. UNSET_RL ) THEN
129 gforget 1.6 WRITE(msgBuf,'(2A,I4)')
130     & 'S/R LAYERS_READPARMS: ',
131     & 'No value for layers_bounds at k =', k
132 rpa 1.1 CALL PRINT_ERROR( msgBuf, myThid )
133     STOP 'ABNORMAL END: S/R LAYERS_READPARMS'
134 jmc 1.8 ELSEIF ( k .EQ. 1 ) THEN
135 rpa 1.1 C Do nothing
136 jmc 1.8 ELSEIF ( layers_bounds(k,iLa) .LE.
137 gforget 1.6 & layers_bounds(k-1,iLa) ) THEN
138     C Check to make sure layers_bounds is increasing
139     WRITE(msgBuf,'(2A,I4)')
140     & 'S/R LAYERS_READPARMS: ',
141     & 'layers_bounds is not increasing at k =', k
142 rpa 1.1 CALL PRINT_ERROR( msgBuf, myThid )
143     STOP 'ABNORMAL END: S/R LAYERS_READPARMS'
144 jmc 1.8 ENDIF
145     ENDDO
146     ENDIF
147 gforget 1.6 ENDDO
148 rpa 1.1
149     C-- Make sure that we locally honor the global MNC on/off flag
150     layers_MNC = layers_MNC .AND. useMNC
151     #ifndef ALLOW_MNC
152     C Fix to avoid running without getting any output:
153     layers_MNC = .FALSE.
154     #endif
155     layers_MDSIO = (.NOT. layers_MNC) .OR. outputTypesInclusive
156    
157     _END_MASTER(myThid)
158    
159     C-- Everyone else must wait for the parameters to be loaded
160     _BARRIER
161    
162     #endif /* ALLOW_MYPACKAGE */
163    
164     RETURN
165     END

  ViewVC Help
Powered by ViewVC 1.1.22