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

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

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


Revision 1.6 - (hide annotations) (download)
Wed Oct 19 01:28:45 2011 UTC (12 years, 7 months ago) by dfer
Branch: MAIN
CVS Tags: checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g
Changes since 1.5: +7 -1 lines
Include potential density as new coordinate (Thanks to David Munday)

1 dfer 1.6 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_check.F,v 1.5 2010/12/04 23:50:32 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_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    
27     #ifdef ALLOW_LAYERS
28    
29     _BEGIN_MASTER(myThid)
30    
31     WRITE(msgBuf,'(A)') 'LAYERS_CHECK: #define LAYERS'
32     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
33     & SQUEEZE_RIGHT , 1)
34    
35     C-- Print out some key parameters :
36     CALL WRITE_0D_I( NZZ, INDEX_NONE, 'NZZ =',
37     & ' /* number of levels in the fine vertical grid */')
38    
39 rpa 1.2 CALL WRITE_1D_RL( dZZf, NZZ, INDEX_K, 'dZZf =',
40     & ' /* fine vertical grid spacing for isopycnal interp */')
41    
42 rpa 1.1 CALL WRITE_1D_RL(layers_G,Nlayers+1, INDEX_K,'layers_G =',
43     & ' /* boundaries of isopycnal-averaging bins */')
44    
45 dfer 1.6 CALL WRITE_0D_I( layers_kref, INDEX_NONE, 'layers_kref =',
46     & ' /* model level to reference potential density to */' )
47    
48     CALL WRITE_0D_I( LAYER_nb, INDEX_NONE, 'LAYER_nb =',
49     & '/* (1) theta; (2) salt; (3) prho; as averaging field */' )
50    
51 rpa 1.1 C-- Check parameters:
52    
53    
54     C For now the package will only work if density ~ temperature
55 dfer 1.5 c IF ( (eosType .EQ. 'LINEAR')
56     c & .AND. (sBeta .EQ. 0.0 _d 0) ) THEN
57 jmc 1.4 C we are good
58 dfer 1.5 c ELSE
59     c WRITE(msgBuf,'(2A)') 'eosType must be eosType=''LINEAR''',
60     c & ' and sBeta must = 0.0'
61     c CALL PRINT_ERROR( msgBuf , 1)
62     c STOP 'ABNORMAL END: S/R LAYERS_CHECK'
63     c ENDIF
64 rpa 1.1
65     _END_MASTER(myThid)
66    
67 jmc 1.3 #endif /* ALLOW_LAYERS */
68 rpa 1.1
69     RETURN
70     END

  ViewVC Help
Powered by ViewVC 1.1.22