/[MITgcm]/MITgcm/model/src/find_alpha.F
ViewVC logotype

Contents of /MITgcm/model/src/find_alpha.F

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


Revision 1.1 - (show annotations) (download)
Tue May 18 18:01:13 1999 UTC (25 years, 1 month ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint28, checkpoint29, checkpoint22, checkpoint23, checkpoint24, checkpoint25, checkpoint27, checkpoint26, checkpoint31, checkpoint30
Modifications/additions for KPP mixing scheme. Instigated by Dimitri.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/find_rho.F,v 1.9 1998/11/06 22:44:46 cnh Exp $
2
3 #include "CPP_OPTIONS.h"
4 #define USE_FACTORIZED_POLY
5
6 subroutine FIND_ALPHA (
7 I bi, bj, iMin, iMax, jMin, jMax, k, kRef, eqn,
8 O alphaloc )
9 C /==========================================================\
10 C | o SUBROUTINE FIND_ALPHA |
11 C | Calculates [drho(S,T,z) / dT] of a horizontal slice |
12 C |==========================================================|
13 C | |
14 C | k - is the Theta/Salt level |
15 C | kRef - determines pressure reference level |
16 C | (not used in 'LINEAR' mode) |
17 C | eqn - determines the eqn. of state: 'LINEAR' or 'POLY3' |
18 C | |
19 C | alphaloc - drho / dT (kg/m^3/C) |
20 C | |
21 C \==========================================================/
22 implicit none
23
24 c Common
25 #include "SIZE.h"
26 #include "DYNVARS.h"
27 #include "EEPARAMS.h"
28 #include "PARAMS.h"
29
30 c Arguments
31 integer bi,bj,iMin,iMax,jMin,jMax
32 integer k ! Level of Theta/Salt slice
33 integer kRef ! Pressure reference level
34 character*(*) eqn
35 _RL alphaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
36
37 c Local
38 integer i,j
39 _RL refTemp,refSalt,tP,sP
40
41 if (eqn.eq.'LINEAR') then
42
43 do j=jMin,jMax
44 do i=iMin,iMax
45 alphaloc(i,j) = -rhonil * tAlpha
46 enddo
47 enddo
48
49 elseif (eqn.eq.'POLY3') then
50
51 refTemp=eosRefT(kRef)
52 refSalt=eosRefS(kRef)
53
54 do j=jMin,jMax
55 do i=iMin,iMax
56 tP=theta(i,j,k,bi,bj)-refTemp
57 sP=salt(i,j,k,bi,bj)-refSalt
58 #ifdef USE_FACTORIZED_POLY
59 alphaloc(i,j) =
60 & ( eosC(6,kRef)
61 & *tP*3.
62 & +(eosC(7,kRef)*sP + eosC(3,kRef))*2.
63 & )*tP
64 & +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef)
65 &
66 #else
67 alphaloc(i,j) =
68 & eosC(1,kRef) +
69 & eosC(3,kRef)*tP*2. +
70 & eosC(4,kRef) *sP +
71 & eosC(6,kRef)*tP*tP*3. +
72 & eosC(7,kRef)*tP*2. *sP +
73 & eosC(8,kRef) *sP*sP
74 #endif
75 enddo
76 enddo
77
78 else
79 write(0,*) 'FIND_ALPHA: eqn = ',eqn
80 stop 'FIND_ALPHA: argument "eqn" has illegal value'
81 endif
82
83 return
84 end
85
86 subroutine FIND_BETA (
87 I bi, bj, iMin, iMax, jMin, jMax, k, kRef, eqn,
88 O betaloc )
89 C /==========================================================\
90 C | o SUBROUTINE FIND_BETA |
91 C | Calculates [drho(S,T,z) / dS] of a horizontal slice |
92 C |==========================================================|
93 C | |
94 C | k - is the Theta/Salt level |
95 C | kRef - determines pressure reference level |
96 C | (not used in 'LINEAR' mode) |
97 C | eqn - determines the eqn. of state: 'LINEAR' or 'POLY3' |
98 C | |
99 C | betaloc - drho / dS (kg/m^3/PSU) |
100 C | |
101 C \==========================================================/
102 implicit none
103
104 c Common
105 #include "SIZE.h"
106 #include "DYNVARS.h"
107 #include "EEPARAMS.h"
108 #include "PARAMS.h"
109
110 c Arguments
111 integer bi,bj,iMin,iMax,jMin,jMax
112 integer k ! Level of Theta/Salt slice
113 integer kRef ! Pressure reference level
114 character*(*) eqn
115 _RL betaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
116
117 c Local
118 integer i,j
119 _RL refTemp,refSalt,tP,sP
120
121 if (eqn.eq.'LINEAR') then
122
123 do j=jMin,jMax
124 do i=iMin,iMax
125 betaloc(i,j) = rhonil * sBeta
126 enddo
127 enddo
128
129 elseif (eqn.eq.'POLY3') then
130
131 refTemp=eosRefT(kRef)
132 refSalt=eosRefS(kRef)
133
134 do j=jMin,jMax
135 do i=iMin,iMax
136 tP=theta(i,j,k,bi,bj)-refTemp
137 sP=salt(i,j,k,bi,bj)-refSalt
138 #ifdef USE_FACTORIZED_POLY
139 betaloc(i,j) =
140 & ( eosC(9,kRef)*sP*3. + eosC(5,kRef)*2. )*sP + eosC(2,kRef)
141 & + ( eosC(7,kRef)*tP
142 & +eosC(8,kRef)*sP*2. + eosC(4,kRef)
143 & )*tP
144 #else
145 betaloc(i,j) =
146 & eosC(2,kRef) +
147 & eosC(4,kRef)*tP +
148 & eosC(5,kRef) *sP*2. +
149 & eosC(7,kRef)*tP*tP +
150 & eosC(8,kRef)*tP *sP*2. +
151 & eosC(9,kRef) *sP*sP*3.
152 #endif
153 enddo
154 enddo
155
156 else
157 write(0,*) 'FIND_BETA: eqn = ',eqn
158 stop 'FIND_BETA: argument "eqn" has illegal value'
159 endif
160
161 return
162 end

  ViewVC Help
Powered by ViewVC 1.1.22