/[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.2 - (show annotations) (download)
Mon Nov 13 16:32:58 2000 UTC (23 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: branch-atmos-merge-freeze, branch-atmos-merge-start, branch-atmos-merge-shapiro, checkpoint33, checkpoint32, checkpoint34, branch-atmos-merge-zonalfilt, branch-atmos-merge-phase5, branch-atmos-merge-phase4, branch-atmos-merge-phase7, branch-atmos-merge-phase6, branch-atmos-merge-phase1, branch-atmos-merge-phase3, branch-atmos-merge-phase2
Branch point for: branch-atmos-merge
Changes since 1.1: +7 -3 lines
Rescaling of forcing fields done immediately after reading fields.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/find_alpha.F,v 1.1 1999/05/18 18:01:13 adcroft 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 cph( write statement confuses TAMC
80 cph write(0,*) 'FIND_ALPHA: eqn = ',eqn
81 cph)
82 stop 'FIND_ALPHA: argument "eqn" has illegal value'
83 endif
84
85 return
86 end
87
88 subroutine FIND_BETA (
89 I bi, bj, iMin, iMax, jMin, jMax, k, kRef, eqn,
90 O betaloc )
91 C /==========================================================\
92 C | o SUBROUTINE FIND_BETA |
93 C | Calculates [drho(S,T,z) / dS] of a horizontal slice |
94 C |==========================================================|
95 C | |
96 C | k - is the Theta/Salt level |
97 C | kRef - determines pressure reference level |
98 C | (not used in 'LINEAR' mode) |
99 C | eqn - determines the eqn. of state: 'LINEAR' or 'POLY3' |
100 C | |
101 C | betaloc - drho / dS (kg/m^3/PSU) |
102 C | |
103 C \==========================================================/
104 implicit none
105
106 c Common
107 #include "SIZE.h"
108 #include "DYNVARS.h"
109 #include "EEPARAMS.h"
110 #include "PARAMS.h"
111
112 c Arguments
113 integer bi,bj,iMin,iMax,jMin,jMax
114 integer k ! Level of Theta/Salt slice
115 integer kRef ! Pressure reference level
116 character*(*) eqn
117 _RL betaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
118
119 c Local
120 integer i,j
121 _RL refTemp,refSalt,tP,sP
122
123 if (eqn.eq.'LINEAR') then
124
125 do j=jMin,jMax
126 do i=iMin,iMax
127 betaloc(i,j) = rhonil * sBeta
128 enddo
129 enddo
130
131 elseif (eqn.eq.'POLY3') then
132
133 refTemp=eosRefT(kRef)
134 refSalt=eosRefS(kRef)
135
136 do j=jMin,jMax
137 do i=iMin,iMax
138 tP=theta(i,j,k,bi,bj)-refTemp
139 sP=salt(i,j,k,bi,bj)-refSalt
140 #ifdef USE_FACTORIZED_POLY
141 betaloc(i,j) =
142 & ( eosC(9,kRef)*sP*3. + eosC(5,kRef)*2. )*sP + eosC(2,kRef)
143 & + ( eosC(7,kRef)*tP
144 & +eosC(8,kRef)*sP*2. + eosC(4,kRef)
145 & )*tP
146 #else
147 betaloc(i,j) =
148 & eosC(2,kRef) +
149 & eosC(4,kRef)*tP +
150 & eosC(5,kRef) *sP*2. +
151 & eosC(7,kRef)*tP*tP +
152 & eosC(8,kRef)*tP *sP*2. +
153 & eosC(9,kRef) *sP*sP*3.
154 #endif
155 enddo
156 enddo
157
158 else
159 cph( write statement confuses TAMC
160 cph write(0,*) 'FIND_BETA: eqn = ',eqn
161 cph)
162 stop 'FIND_BETA: argument "eqn" has illegal value'
163 endif
164
165 return
166 end

  ViewVC Help
Powered by ViewVC 1.1.22