1 |
C $Header: /u/gcmpack/MITgcm/pkg/obcs/OBCS_PARAMS.h,v 1.3 2012/03/09 20:13:03 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#ifdef ALLOW_OBCS |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: OBCS_PARAMS.h |
8 |
C !INTERFACE: |
9 |
C #include "OBCS_PARAMS.h" |
10 |
|
11 |
C !DESCRIPTION: |
12 |
C *==========================================================* |
13 |
C | OBCS_PARAMS.h |
14 |
C | o Header file containing OBCS parameters |
15 |
C *==========================================================* |
16 |
C | o Note: does not (and should not) contain any conditional |
17 |
C | statement that depends on OBCS options ; therefore |
18 |
C | can be safely included without OBCS_OPTIONS.h |
19 |
C *==========================================================* |
20 |
CEOP |
21 |
|
22 |
C tidalComponents :: number of tidal components to be applied |
23 |
INTEGER tidalComponents |
24 |
PARAMETER ( tidalComponents = 10 ) |
25 |
|
26 |
C-- COMMON /OBC_PARM_I/ OBCS integer-type parameter |
27 |
C OBCS_u1_adv_T :: >0: use 1rst O. upwind adv-scheme @ OB (=1: only if outflow) |
28 |
C OBCS_u1_adv_S :: >0: use 1rst O. upwind adv-scheme @ OB (=1: only if outflow) |
29 |
C OBCS_monSelect :: select group of variables to monitor |
30 |
C spongeThickness :: number grid points that make up the sponge layer (def=0) |
31 |
COMMON /OBC_PARM_I/ |
32 |
& OBCS_u1_adv_T, OBCS_u1_adv_S, |
33 |
& OBCS_monSelect, |
34 |
& spongeThickness |
35 |
INTEGER OBCS_u1_adv_T, OBCS_u1_adv_S |
36 |
INTEGER OBCS_monSelect |
37 |
INTEGER spongeThickness |
38 |
|
39 |
C-- COMMON /OBC_PARM_L/ OBCS logical-type parameter |
40 |
C useOrlanskiNorth/South/East/West |
41 |
C :: specify Orlanski boundary conditions for northern/ |
42 |
C southern/eastern/Western |
43 |
C useStevensNorth/South/East/West |
44 |
C :: use open boundary computations following Stevens (1990) |
45 |
C useStevensPhaseVel |
46 |
C :: use phase velocity contribution for open boundary |
47 |
C computations following Stevens (1990), default = true |
48 |
C useStevensAdvection |
49 |
C :: use advective contribution for open boundary |
50 |
C computations following Stevens (1990), default = true |
51 |
C useOBCSsponge :: turns on sponge layer along boundary (def=false) |
52 |
C useOBCSbalance :: balance the volume flux through boundary |
53 |
C at every time step |
54 |
C useOBCStides :: modify OB normal flow to add tidal forcing |
55 |
C useOBCSprescribe :: read boundary conditions from a file |
56 |
C (overrides Orlanski and other boundary values) |
57 |
C OBCSprintDiags :: print boundary values to STDOUT (def=true) |
58 |
C OBCSfixTopo :: check and adjust topography for problematic gradients |
59 |
C across boundaries (def=true) |
60 |
|
61 |
COMMON /OBC_PARM_L/ |
62 |
& useOrlanskiNorth,useOrlanskiSouth, |
63 |
& useOrlanskiEast,useOrlanskiWest, |
64 |
& useStevensNorth,useStevensSouth, |
65 |
& useStevensEast,useStevensWest, |
66 |
& useStevensPhaseVel, useStevensAdvection, |
67 |
& useOBCSsponge, useOBCSbalance, useOBCStides, useOBCSprescribe, |
68 |
& OBCSprintDiags, |
69 |
& OBCSfixTopo |
70 |
LOGICAL useOrlanskiNorth |
71 |
LOGICAL useOrlanskiSouth |
72 |
LOGICAL useOrlanskiEast |
73 |
LOGICAL useOrlanskiWest |
74 |
LOGICAL useStevensNorth |
75 |
LOGICAL useStevensSouth |
76 |
LOGICAL useStevensEast |
77 |
LOGICAL useStevensWest |
78 |
LOGICAL useStevensPhaseVel |
79 |
LOGICAL useStevensAdvection |
80 |
LOGICAL useOBCSsponge |
81 |
LOGICAL useOBCSbalance |
82 |
LOGICAL useOBCStides |
83 |
LOGICAL useOBCSprescribe |
84 |
LOGICAL OBCSprintDiags |
85 |
LOGICAL OBCSfixTopo |
86 |
|
87 |
C-- COMMON /OBC_PARM_R/ OBCS real-type parameter |
88 |
C OBCS_balanceFacN/S/E/W :: weighting factor for balancing OB normal flow |
89 |
C OBCS_uvApplyFac :: multiplying factor to U,V normal comp. when applying |
90 |
C OBC to 2nd column/row (for backward compatibility). |
91 |
C OBCS_monitorFreq :: monitor output frequency (s) for OB statistics |
92 |
C U/Vrelaxobcsinner/bound :: relaxation time scale (in seconds) on the boundary |
93 |
C (bound) and at the innermost grid point of the sponge |
94 |
C layer (inner); relaxation time scales in-between |
95 |
C are linearly interpolated from these values |
96 |
C T/SrelaxStevens :: relaxation time scale (in seconds) for T/S-points |
97 |
C for Stevens boundary conditions |
98 |
C tidalPeriod :: tidal period (s) |
99 |
COMMON /OBC_PARM_R/ |
100 |
& OBCS_balanceFacN, OBCS_balanceFacS, |
101 |
& OBCS_balanceFacE, OBCS_balanceFacW, |
102 |
& OBCS_uvApplyFac, |
103 |
& OBCS_monitorFreq, |
104 |
& tidalPeriod, |
105 |
& Urelaxobcsinner,Urelaxobcsbound, |
106 |
& Vrelaxobcsinner,Vrelaxobcsbound, |
107 |
& TrelaxStevens, SrelaxStevens |
108 |
_RL OBCS_balanceFacN, OBCS_balanceFacS |
109 |
_RL OBCS_balanceFacE, OBCS_balanceFacW |
110 |
_RL OBCS_uvApplyFac |
111 |
_RL OBCS_monitorFreq |
112 |
_RL tidalPeriod(tidalComponents) |
113 |
_RS Urelaxobcsinner |
114 |
_RS Urelaxobcsbound |
115 |
_RS Vrelaxobcsinner |
116 |
_RS Vrelaxobcsbound |
117 |
_RS TrelaxStevens |
118 |
_RS SrelaxStevens |
119 |
|
120 |
C-- COMMON /OBC_FILES/ OBCS character-type parameter |
121 |
C insideOBmaskFile :: File to specify Inside OB region mask (zero beyond OB). |
122 |
C OB[N,S,E,W][u,v,w,t,s,eta,am,ph]File :: Files with boundary conditions, |
123 |
C the letter combinations mean: |
124 |
C N/S/E/W :: northern/southern/eastern/western boundary |
125 |
C u/v/w/t/s :: ocean u/v/w velocities, temperature/salinity |
126 |
C eta :: sea surface height |
127 |
C am/ph :: tidal amplitude (m/s) / phase (s) |
128 |
COMMON /OBC_FILES/ |
129 |
& OBNuFile, OBSuFile, OBEuFile, OBWuFile, |
130 |
& OBNvFile, OBSvFile, OBEvFile, OBWvFile, |
131 |
& OBNwFile, OBSwFile, OBEwFile, OBWwFile, |
132 |
& OBNtFile, OBStFile, OBEtFile, OBWtFile, |
133 |
& OBNsFile, OBSsFile, OBEsFile, OBWsFile, |
134 |
& OBNetaFile,OBSetaFile,OBEetaFile,OBWetaFile, |
135 |
& OBNamFile, OBSamFile, OBEamFile, OBWamFile, |
136 |
& OBNphFile, OBSphFile, OBEphFile, OBWphFile, |
137 |
& insideOBmaskFile |
138 |
CHARACTER*(MAX_LEN_FNAM) |
139 |
& OBNuFile, OBSuFile, OBEuFile, OBWuFile, |
140 |
& OBNvFile, OBSvFile, OBEvFile, OBWvFile, |
141 |
& OBNwFile, OBSwFile, OBEwFile, OBWwFile, |
142 |
& OBNtFile, OBStFile, OBEtFile, OBWtFile, |
143 |
& OBNsFile, OBSsFile, OBEsFile, OBWsFile, |
144 |
& OBNetaFile,OBSetaFile,OBEetaFile,OBWetaFile, |
145 |
& OBNamFile, OBSamFile, OBEamFile, OBWamFile, |
146 |
& OBNphFile, OBSphFile, OBEphFile, OBWphFile, |
147 |
& insideOBmaskFile |
148 |
|
149 |
#endif /* ALLOW_OBCS */ |