File indexing completed on 2018-03-09 06:09:17 UTC
view on githubraw file Latest commit 84f265d4 on 2018-03-08 18:49:29 UTC
108a00eab9 Ryan*0001 #include "LAYERS_OPTIONS.h"
0002
0003
0004
0005 SUBROUTINE LAYERS_READPARMS( myThid )
0006
0007
0008
0009 IMPLICIT NONE
0010 #include "SIZE.h"
0011 #include "EEPARAMS.h"
0012 #include "PARAMS.h"
0013 #include "LAYERS_SIZE.h"
0014 #include "LAYERS.h"
0015
0016
0017 INTEGER myThid
0018
0019 #ifdef ALLOW_LAYERS
ae4c29e0db Jean*0020
0021
0022
0023
0024 CHARACTER*(MAX_LEN_MBUF) msgBuf
0025 INTEGER iUnit, k, iLa
2c2df907b2 Jean*0026 INTEGER errCount
0027
0028
0029
0030 INTEGER LAYER_nb, layers_kref
0031 LOGICAL useBOLUS
0032 _RL layers_G(Nlayers+1)
0033
0034
108a00eab9 Ryan*0035
0036 NAMELIST /LAYERS_PARM01/
a60f4dd950 Davi*0037 & layers_G, layers_taveFreq, layers_diagFreq,
406891c1a3 Gael*0038 & LAYER_nb, layers_kref, useBOLUS, layers_bolus,
0039 & layers_name, layers_bounds, layers_krho
108a00eab9 Ryan*0040
ae4c29e0db Jean*0041 IF ( .NOT.useLayers ) THEN
0042
0043 _BEGIN_MASTER(myThid)
0044
0045
0046 CALL PACKAGES_UNUSED_MSG( 'useLayers', ' ', ' ' )
0047 _END_MASTER(myThid)
0048 RETURN
0049 ENDIF
108a00eab9 Ryan*0050
0051 _BEGIN_MASTER(myThid)
2c2df907b2 Jean*0052 errCount = 0
108a00eab9 Ryan*0053
0054
0055
2c2df907b2 Jean*0056 layers_taveFreq = taveFreq
0057 layers_diagFreq = dumpFreq
108a00eab9 Ryan*0058
0059 layers_MNC = .FALSE.
0060 layers_MDSIO = .TRUE.
0061
406891c1a3 Gael*0062 DO iLa=1,layers_maxNum
0063 layers_name(iLa) = ' '
0064 layers_krho(iLa)= 1
20c24d88bc Jean*0065 layers_bolus(iLa) = useGMRedi
406891c1a3 Gael*0066 DO k=1,Nlayers+1
0067 layers_bounds(k,iLa) = UNSET_RL
0068 ENDDO
0069 ENDDO
20c24d88bc Jean*0070
2c2df907b2 Jean*0071
406891c1a3 Gael*0072 LAYER_nb = 0
17ce8d85dd Davi*0073 layers_kref = 1
20c24d88bc Jean*0074 useBOLUS = useGMRedi
2c2df907b2 Jean*0075 DO k=1,Nlayers+1
0076 layers_G(k) = UNSET_RL
0077 ENDDO
108a00eab9 Ryan*0078
0079 WRITE(msgBuf,'(A)') 'LAYERS_READPARMS: opening data.layers'
0080 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
0081 & SQUEEZE_RIGHT , 1)
0082 CALL OPEN_COPY_DATA_FILE(
0083 I 'data.layers', 'LAYERS_READPARMS',
0084 O iUnit,
0085 I myThid )
0086
0087
0088 READ(UNIT=iUnit,NML=LAYERS_PARM01)
0089 WRITE(msgBuf,'(A)')
0090 & 'LAYERS_READPARMS: finished reading data.layers'
0091 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
0092 & SQUEEZE_RIGHT , 1)
0093
7a77863887 Mart*0094 #ifdef SINGLE_DISK_IO
108a00eab9 Ryan*0095 CLOSE(iUnit)
7a77863887 Mart*0096 #else
0097 CLOSE(iUnit,STATUS='DELETE')
0098 #endif /* SINGLE_DISK_IO */
108a00eab9 Ryan*0099
2c2df907b2 Jean*0100
0101 IF ( LAYER_nb.LT.0 .OR. LAYER_nb.GT.3 ) THEN
0102 WRITE(msgBuf,'(2A,I2,A,I9)') 'LAYERS_READPARMS: ',
0103 & 'Invalid LAYER_nb=', LAYER_nb
0104 CALL PRINT_ERROR( msgBuf, myThid )
0105 errCount = errCount + 1
0106 ENDIF
0107 IF ( LAYER_nb.EQ.0 ) THEN
0108 IF ( layers_kref.NE.1 ) errCount = errCount + 1
0109 DO k=1,Nlayers+1
0110 IF ( layers_G(k).NE.UNSET_RL ) errCount = errCount + 1
0111 ENDDO
0112 ELSE
0113 DO iLa=1,layers_maxNum
0114 IF ( layers_name(iLa).NE.' ' ) errCount = errCount + 1
0115 IF ( layers_krho(iLa).NE.1 ) errCount = errCount + 1
0116 DO k=1,Nlayers+1
0117 IF ( layers_bounds(k,iLa).NE.UNSET_RL ) errCount = errCount+1
0118 ENDDO
0119 ENDDO
0120
0121 IF ( LAYER_nb.EQ.1 ) layers_name(1) = 'TH '
0122 IF ( LAYER_nb.EQ.2 ) layers_name(1) = 'SLT'
0123 IF ( LAYER_nb.EQ.3 ) layers_name(1) = 'RHO'
406891c1a3 Gael*0124 layers_krho(1) = layers_kref
d30f939fae Davi*0125 layers_bolus(1) = useBOLUS
406891c1a3 Gael*0126 DO k=1,Nlayers+1
0127 layers_bounds(k,1) = layers_G(k)
0128 ENDDO
2c2df907b2 Jean*0129 ENDIF
0130 IF ( errCount.GE.1 ) THEN
0131 WRITE(msgBuf,'(2A)') 'LAYERS_READPARMS: ',
0132 & 'Cannot mix old params setting (LAYER_nb > 0)'
0133 CALL PRINT_ERROR( msgBuf, myThid )
0134 WRITE(msgBuf,'(2A)') 'LAYERS_READPARMS: ',
0135 & ' with new params setting (layer_name(#)= ...)'
0136 CALL PRINT_ERROR( msgBuf, myThid )
0137 WRITE(msgBuf,'(2A,I4,A)') 'LAYERS_READPARMS: ',
0138 & 'Detected', errCount,' fatal error/conflict(s)'
0139 CALL PRINT_ERROR( msgBuf, myThid )
0140 CALL ALL_PROC_DIE( 0 )
0141 STOP 'ABNORMAL END: S/R LAYERS_READPARMS'
406891c1a3 Gael*0142 ENDIF
0143
2c2df907b2 Jean*0144
406891c1a3 Gael*0145 DO iLa=1,layers_maxNum
2c2df907b2 Jean*0146 layers_num(iLa) = 0
0147 IF ( layers_name(iLa).EQ.'TH ' ) layers_num(iLa) = 1
0148 IF ( layers_name(iLa).EQ.'SLT' ) layers_num(iLa) = 2
0149 IF ( layers_name(iLa).EQ.'RHO' ) layers_num(iLa) = 3
84f265d4b3 dfer*0150 IF ( layers_name(iLa).EQ.'MSE' ) layers_num(iLa) = 4
2c2df907b2 Jean*0151 IF ( layers_name(iLa).NE.' ' .AND.
0152 & layers_num(iLa).EQ.0 ) THEN
0153 WRITE(msgBuf,'(2A,I2,3A)') 'LAYERS_READPARMS: ',
0154 & 'invalid layers_name(',iLa,')= "',layers_name(iLa),'"'
0155 CALL PRINT_ERROR( msgBuf, myThid )
0156 errCount = errCount + 1
20c24d88bc Jean*0157 ENDIF
0158
0159 layers_bolus(iLa) = layers_bolus(iLa) .AND. useGMRedi
406891c1a3 Gael*0160 ENDDO
20c24d88bc Jean*0161
406891c1a3 Gael*0162
0163 DO iLa=1,layers_maxNum
20c24d88bc Jean*0164 IF ( layers_num(iLa).NE.0 ) THEN
0165 DO k=1,Nlayers+1
2c2df907b2 Jean*0166 IF ( layers_bounds(k,iLa).EQ.UNSET_RL ) THEN
0167 WRITE(msgBuf,'(2A,I4,A,I3,A)') 'LAYERS_READPARMS: ',
0168 & 'No value for layers_bounds(k=',k,', iLa=', iLa, ')'
108a00eab9 Ryan*0169 CALL PRINT_ERROR( msgBuf, myThid )
2c2df907b2 Jean*0170 errCount = errCount + 1
20c24d88bc Jean*0171 ENDIF
0172 ENDDO
0173 ENDIF
84f265d4b3 dfer*0174 #ifndef LAYERS_PRHO_REF
0175 IF ( layers_num(iLa).EQ.3 ) THEN
0176 WRITE(msgBuf,'(2A,I4)')
0177 & 'S/R LAYERS_READPARMS: ',
0178 & 'Layering in density requires to define LAYERS_PRHO_REF'
0179 CALL PRINT_ERROR( msgBuf, myThid )
0180 errCount = errCount + 1
0181 ENDIF
0182 #endif
0183 #ifndef LAYERS_MSE
0184 IF ( layers_num(iLa).EQ.4 ) THEN
0185 WRITE(msgBuf,'(2A,I4)')
0186 & 'S/R LAYERS_READPARMS: ',
0187 & 'Layering in MSE requires to define LAYERS_MSE'
0188 CALL PRINT_ERROR( msgBuf, myThid )
0189 errCount = errCount + 1
0190 ENDIF
0191 #endif
0192
0193 IF ( layers_num(iLa).EQ.4 ) THEN
0194 IF ( useAtm_Phys .OR. useFizhi ) THEN
0195 WRITE(msgBuf,'(2A,I4)')
0196 & 'S/R LAYERS_READPARMS: ',
0197 & 'Layering in MSE only work with Aim_v23'
0198 CALL PRINT_ERROR( msgBuf, myThid )
0199 errCount = errCount + 1
0200 ENDIF
0201 ENDIF
406891c1a3 Gael*0202 ENDDO
108a00eab9 Ryan*0203
0204
0205 layers_MNC = layers_MNC .AND. useMNC
0206 #ifndef ALLOW_MNC
0207
0208 layers_MNC = .FALSE.
0209 #endif
0210 layers_MDSIO = (.NOT. layers_MNC) .OR. outputTypesInclusive
0211
2c2df907b2 Jean*0212 IF ( errCount.GE.1 ) THEN
0213 WRITE(msgBuf,'(A,I3,A)')
0214 & 'LAYERS_READPARMS: detected', errCount,' fatal error(s)'
0215 CALL PRINT_ERROR( msgBuf, myThid )
0216 CALL ALL_PROC_DIE( 0 )
0217 STOP 'ABNORMAL END: S/R LAYERS_READPARMS'
0218 ENDIF
0219
108a00eab9 Ryan*0220 _END_MASTER(myThid)
0221
0222
0223 _BARRIER
0224
0225 #endif /* ALLOW_MYPACKAGE */
0226
0227 RETURN
0228 END