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

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

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


Revision 1.3 - (hide annotations) (download)
Sun Feb 4 14:38:47 2001 UTC (23 years, 5 months ago) by cnh
Branch: MAIN
CVS Tags: pre38tag1, c37_adj, pre38-close, checkpoint37, checkpoint36, checkpoint35
Branch point for: pre38
Changes since 1.2: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

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

  ViewVC Help
Powered by ViewVC 1.1.22