1 |
jmc |
1.8 |
C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_check.F,v 1.7 2012/09/19 18:48:18 gforget 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_CHECK( myThid ) |
9 |
|
|
|
10 |
|
|
C Check dependances with other packages |
11 |
|
|
|
12 |
|
|
IMPLICIT NONE |
13 |
|
|
#include "SIZE.h" |
14 |
|
|
#include "EEPARAMS.h" |
15 |
|
|
#include "PARAMS.h" |
16 |
|
|
#include "EOS.h" |
17 |
|
|
#include "LAYERS_SIZE.h" |
18 |
|
|
#include "LAYERS.h" |
19 |
|
|
|
20 |
|
|
C myThid :: my Thread Id number |
21 |
|
|
INTEGER myThid |
22 |
|
|
|
23 |
|
|
C LOCAL VARIABLES: |
24 |
jmc |
1.3 |
C msgBuf :: Informational/error message buffer |
25 |
rpa |
1.1 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
26 |
jmc |
1.8 |
INTEGER iLa, k |
27 |
|
|
_RL tmpVar |
28 |
rpa |
1.1 |
|
29 |
|
|
#ifdef ALLOW_LAYERS |
30 |
|
|
_BEGIN_MASTER(myThid) |
31 |
|
|
|
32 |
|
|
WRITE(msgBuf,'(A)') 'LAYERS_CHECK: #define LAYERS' |
33 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
34 |
|
|
& SQUEEZE_RIGHT , 1) |
35 |
|
|
|
36 |
jmc |
1.8 |
C-- Print out some key parameters : |
37 |
rpa |
1.1 |
CALL WRITE_0D_I( NZZ, INDEX_NONE, 'NZZ =', |
38 |
|
|
& ' /* number of levels in the fine vertical grid */') |
39 |
rpa |
1.2 |
CALL WRITE_1D_RL( dZZf, NZZ, INDEX_K, 'dZZf =', |
40 |
|
|
& ' /* fine vertical grid spacing for isopycnal interp */') |
41 |
|
|
|
42 |
gforget |
1.7 |
DO iLa=1,layers_maxNum |
43 |
jmc |
1.8 |
IF ( layers_num(iLa).NE.0 ) THEN |
44 |
rpa |
1.1 |
|
45 |
jmc |
1.8 |
CALL WRITE_0D_I( layers_num(iLa), INDEX_NONE, 'layers_num =', |
46 |
gforget |
1.7 |
& '/* (1) theta; (2) salt; (3) prho; as averaging field */' ) |
47 |
jmc |
1.8 |
CALL WRITE_0D_C( layers_name(iLa),-1,INDEX_NONE,'layers_name =', |
48 |
|
|
& '/* (TH) theta; (SLT) salt; (RHO) prho; as averaging field */') |
49 |
|
|
CALL WRITE_0D_L ( layers_bolus(iLa), INDEX_NONE, |
50 |
|
|
& 'layers_bolus =',' /* include potential GM bolus velocity */') |
51 |
|
|
IF ( layers_num(iLa).EQ.3 ) |
52 |
gforget |
1.7 |
& CALL WRITE_0D_I( layers_krho(iLa), INDEX_NONE, 'layers_krho =', |
53 |
dfer |
1.6 |
& ' /* model level to reference potential density to */' ) |
54 |
jmc |
1.8 |
CALL WRITE_1D_RL( layers_bounds(1,iLa), Nlayers+1, INDEX_K, |
55 |
|
|
& 'layers_bounds =', ' /* boundaries of tracer-averaging bins */') |
56 |
dfer |
1.6 |
|
57 |
jmc |
1.8 |
ENDIF !IF ( layers_num(iLa).NE.0 ) THEN |
58 |
gforget |
1.7 |
ENDDO !DO iLa=1,layers_maxNum |
59 |
dfer |
1.6 |
|
60 |
jmc |
1.8 |
C-- Check parameters: |
61 |
|
|
DO iLa=1,layers_maxNum |
62 |
rpa |
1.1 |
|
63 |
jmc |
1.8 |
IF ( layers_num(iLa).NE.0 ) THEN |
64 |
|
|
C- check for inconsistent density layers_bounds specification |
65 |
|
|
C make sure layers_bounds is increasing: |
66 |
|
|
DO k=1,Nlayers |
67 |
|
|
IF ( layers_bounds(k,iLa).GE.layers_bounds(k+1,iLa) ) THEN |
68 |
|
|
WRITE(msgBuf,'(A,I2,A,I4)') 'LAYERS_CHECK(iLa=', iLa, |
69 |
|
|
& '): layers_bounds k -> k+1 not increasing at k=', k |
70 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
71 |
|
|
STOP 'ABNORMAL END: S/R LAYERS_CHECK' |
72 |
|
|
ENDIF |
73 |
|
|
ENDDO |
74 |
|
|
ENDIF |
75 |
|
|
|
76 |
|
|
IF ( layers_num(iLa).EQ.3 ) THEN |
77 |
|
|
C Pot.Density is now expressed as rho-1000 (previously just rho): |
78 |
|
|
C check for realistic layers_bounds values: |
79 |
|
|
tmpVar = layers_bounds(Nlayers+1,iLa) - layers_bounds(1,iLa) |
80 |
|
|
IF ( tmpVar.LE.50. .AND. layers_bounds(1,iLa).GE.950. ) THEN |
81 |
|
|
WRITE(msgBuf,'(A,I2,A)') 'LAYERS_CHECK(iLa=', iLa, |
82 |
|
|
& '): layers_bounds seems to be expressed as "rho"' |
83 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
84 |
|
|
WRITE(msgBuf,'(A,I2,A)') 'LAYERS_CHECK(iLa=', iLa, |
85 |
|
|
& '): while it should be expressed as "rho - 1000"' |
86 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
87 |
|
|
STOP 'ABNORMAL END: S/R LAYERS_CHECK' |
88 |
|
|
ENDIF |
89 |
|
|
ENDIF |
90 |
rpa |
1.1 |
|
91 |
jmc |
1.8 |
ENDDO |
92 |
rpa |
1.1 |
|
93 |
|
|
_END_MASTER(myThid) |
94 |
jmc |
1.3 |
#endif /* ALLOW_LAYERS */ |
95 |
rpa |
1.1 |
|
96 |
|
|
RETURN |
97 |
|
|
END |