1 |
jmc |
1.7 |
C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_fixed.F,v 1.6 2009/07/08 21:59:01 jmc Exp $ |
2 |
heimbach |
1.1 |
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "SEAICE_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
CStartOfInterface |
7 |
|
|
SUBROUTINE SEAICE_INIT_FIXED( myThid ) |
8 |
jmc |
1.6 |
C *==========================================================* |
9 |
|
|
C | SUBROUTINE SEAICE_INIT_FIXED |
10 |
|
|
C | o Initialization of sea ice model. |
11 |
|
|
C *==========================================================* |
12 |
|
|
C *==========================================================* |
13 |
heimbach |
1.1 |
IMPLICIT NONE |
14 |
jmc |
1.6 |
|
15 |
heimbach |
1.1 |
C === Global variables === |
16 |
|
|
#include "SIZE.h" |
17 |
|
|
#include "EEPARAMS.h" |
18 |
|
|
#include "PARAMS.h" |
19 |
|
|
#include "GRID.h" |
20 |
jmc |
1.7 |
#include "FFIELDS.h" |
21 |
|
|
#include "SEAICE_PARAMS.h" |
22 |
heimbach |
1.1 |
#include "SEAICE.h" |
23 |
|
|
CML#include "SEAICE_GRID.h" |
24 |
|
|
|
25 |
|
|
C === Routine arguments === |
26 |
|
|
C myThid - Thread no. that called this routine. |
27 |
|
|
INTEGER myThid |
28 |
|
|
CEndOfInterface |
29 |
jmc |
1.6 |
|
30 |
heimbach |
1.1 |
C === Local variables === |
31 |
|
|
C i,j,k,bi,bj - Loop counters |
32 |
|
|
|
33 |
jmc |
1.7 |
INTEGER i, j, bi, bj |
34 |
mlosch |
1.4 |
INTEGER kSurface |
35 |
jmc |
1.6 |
#ifndef SEAICE_CGRID |
36 |
mlosch |
1.4 |
_RS mask_uice |
37 |
jmc |
1.6 |
#endif |
38 |
heimbach |
1.3 |
cif( |
39 |
|
|
cif Helper variable for determining the fraction of sw radiation |
40 |
|
|
cif penetrating the model's shallowest layer |
41 |
|
|
INTEGER dummyIter |
42 |
|
|
_RL dummyTime |
43 |
|
|
_RL swfracba(2) |
44 |
|
|
_RL FACTORM |
45 |
|
|
INTEGER IMAX |
46 |
|
|
cif) |
47 |
heimbach |
1.1 |
|
48 |
mlosch |
1.4 |
IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN |
49 |
|
|
kSurface = Nr |
50 |
|
|
ELSE |
51 |
|
|
kSurface = 1 |
52 |
heimbach |
1.1 |
ENDIF |
53 |
|
|
|
54 |
jmc |
1.6 |
C Initialize MNC variable information for SEAICE |
55 |
|
|
IF ( useMNC .AND. |
56 |
|
|
& (seaice_tave_mnc.OR.seaice_dump_mnc.OR.SEAICE_mon_mnc) |
57 |
|
|
& ) THEN |
58 |
|
|
CALL SEAICE_MNC_INIT( myThid ) |
59 |
|
|
ENDIF |
60 |
|
|
|
61 |
heimbach |
1.3 |
cif( |
62 |
|
|
#ifdef SHORTWAVE_HEATING |
63 |
|
|
IMAX = 2 |
64 |
|
|
FACTORM = -1.0 |
65 |
|
|
dummyTime = 1.0 |
66 |
|
|
dummyIter = 0 |
67 |
jmc |
1.6 |
swfracba(1) = abs(rF(1)) |
68 |
heimbach |
1.3 |
swfracba(2) = abs(rF(2)) |
69 |
|
|
CALL SWFRAC( |
70 |
|
|
I IMAX,FACTORM, |
71 |
|
|
U swfracba, |
72 |
|
|
I dummyTime,dummyIter,myThid) |
73 |
|
|
SWFRACB = swfracba(2) |
74 |
|
|
#endif |
75 |
|
|
cif) |
76 |
|
|
|
77 |
mlosch |
1.4 |
C-- Initialize grid info |
78 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
79 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
80 |
|
|
DO j=1-OLy,sNy+OLy |
81 |
|
|
DO i=1-OLx,sNx+OLx |
82 |
|
|
HEFFM(i,j,bi,bj) = 0. _d 0 |
83 |
|
|
ENDDO |
84 |
|
|
ENDDO |
85 |
|
|
DO j=1-OLy,sNy+OLy |
86 |
|
|
DO i=1-OLx,sNx+OLx |
87 |
|
|
HEFFM(i,j,bi,bj)= 1. _d 0 |
88 |
|
|
IF (_hFacC(i,j,kSurface,bi,bj).eq.0.) |
89 |
|
|
& HEFFM(i,j,bi,bj)= 0. _d 0 |
90 |
|
|
ENDDO |
91 |
|
|
ENDDO |
92 |
|
|
DO j=1-OLy+1,sNy+OLy |
93 |
|
|
DO i=1-OLx+1,sNx+OLx |
94 |
heimbach |
1.5 |
#ifndef SEAICE_CGRID |
95 |
mlosch |
1.4 |
UVM(i,j,bi,bj)=0. _d 0 |
96 |
|
|
mask_uice=HEFFM(i,j, bi,bj)+HEFFM(i-1,j-1,bi,bj) |
97 |
|
|
& +HEFFM(i,j-1,bi,bj)+HEFFM(i-1,j, bi,bj) |
98 |
|
|
IF(mask_uice.GT.3.5 _d 0) UVM(i,j,bi,bj)=1. _d 0 |
99 |
|
|
#endif /* SEAICE_CGRID */ |
100 |
|
|
ENDDO |
101 |
|
|
ENDDO |
102 |
|
|
ENDDO |
103 |
|
|
ENDDO |
104 |
|
|
|
105 |
|
|
#ifdef SEAICE_CGRID |
106 |
jmc |
1.6 |
C coefficients for metric terms |
107 |
mlosch |
1.4 |
DO bj=myByLo(myThid),myByHi(myThid) |
108 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
109 |
|
|
DO j=1-OLy,sNy+OLy |
110 |
|
|
DO i=1-OLx,sNx+OLx |
111 |
|
|
k1AtC(I,J,bi,bj) = 0.0 _d 0 |
112 |
|
|
k1AtZ(I,J,bi,bj) = 0.0 _d 0 |
113 |
|
|
k2AtC(I,J,bi,bj) = 0.0 _d 0 |
114 |
|
|
k2AtZ(I,J,bi,bj) = 0.0 _d 0 |
115 |
|
|
ENDDO |
116 |
|
|
ENDDO |
117 |
|
|
IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN |
118 |
jmc |
1.6 |
C This is the only case where tan(phi) is not zero. In this case |
119 |
|
|
C C and U points, and Z and V points have the same phi, so that we |
120 |
mlosch |
1.4 |
C only need a copy here. Do not use tan(YC) and tan(YG), because these |
121 |
jmc |
1.6 |
C can be the geographical coordinates and not the correct grid |
122 |
mlosch |
1.4 |
C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0) |
123 |
|
|
DO j=1-OLy,sNy+OLy |
124 |
|
|
DO i=1-OLx,sNx+OLx |
125 |
|
|
k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
126 |
|
|
k2AtZ(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere |
127 |
|
|
ENDDO |
128 |
|
|
ENDDO |
129 |
|
|
ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN |
130 |
|
|
C compute metric term coefficients from finite difference approximation |
131 |
|
|
DO j=1-OLy,sNy+OLy |
132 |
|
|
DO i=1-OLx,sNx+OLx-1 |
133 |
|
|
k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj) |
134 |
|
|
& * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) ) |
135 |
|
|
& * _recip_dxF(I,J,bi,bj) |
136 |
|
|
ENDDO |
137 |
|
|
ENDDO |
138 |
|
|
DO j=1-OLy,sNy+OLy |
139 |
|
|
DO i=1-OLx+1,sNx+OLx |
140 |
jmc |
1.6 |
k1AtZ(I,J,bi,bj) = _recip_dyU(I,J,bi,bj) |
141 |
mlosch |
1.4 |
& * ( _dyC(I,J,bi,bj) - _dyC(I-1,J,bi,bj) ) |
142 |
|
|
& * _recip_dxV(I,J,bi,bj) |
143 |
|
|
ENDDO |
144 |
|
|
ENDDO |
145 |
|
|
DO j=1-OLy,sNy+OLy-1 |
146 |
|
|
DO i=1-OLx,sNx+OLx |
147 |
jmc |
1.6 |
k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj) |
148 |
mlosch |
1.4 |
& * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) ) |
149 |
|
|
& * _recip_dyF(I,J,bi,bj) |
150 |
|
|
ENDDO |
151 |
|
|
ENDDO |
152 |
|
|
DO j=1-OLy+1,sNy+OLy |
153 |
|
|
DO i=1-OLx,sNx+OLx |
154 |
jmc |
1.6 |
k2AtC(I,J,bi,bj) = _recip_dxV(I,J,bi,bj) |
155 |
mlosch |
1.4 |
& * ( _dxC(I,J,bi,bj) - _dxC(I,J-1,bi,bj) ) |
156 |
|
|
& * _recip_dyU(I,J,bi,bj) |
157 |
|
|
ENDDO |
158 |
|
|
ENDDO |
159 |
|
|
ENDIF |
160 |
|
|
ENDDO |
161 |
|
|
ENDDO |
162 |
|
|
#endif /* SEAICE_CGRID */ |
163 |
|
|
|
164 |
|
|
#ifndef SEAICE_CGRID |
165 |
|
|
C-- Choose a proxy level for geostrophic velocity, |
166 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
167 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
168 |
|
|
DO j=1-OLy,sNy+OLy |
169 |
|
|
DO i=1-OLx,sNx+OLx |
170 |
|
|
KGEO(i,j,bi,bj) = 0 |
171 |
|
|
ENDDO |
172 |
|
|
ENDDO |
173 |
|
|
DO j=1-OLy,sNy+OLy |
174 |
|
|
DO i=1-OLx,sNx+OLx |
175 |
|
|
#ifdef SEAICE_BICE_STRESS |
176 |
|
|
KGEO(i,j,bi,bj) = 1 |
177 |
|
|
#else /* SEAICE_BICE_STRESS */ |
178 |
|
|
IF (klowc(i,j,bi,bj) .LT. 2) THEN |
179 |
|
|
KGEO(i,j,bi,bj) = 1 |
180 |
|
|
ELSE |
181 |
|
|
KGEO(i,j,bi,bj) = 2 |
182 |
|
|
DO WHILE ( abs(rC(KGEO(i,j,bi,bj))) .LT. 50.0 _d 0 .AND. |
183 |
|
|
& KGEO(i,j,bi,bj) .LT. (klowc(i,j,bi,bj)-1) ) |
184 |
|
|
KGEO(i,j,bi,bj) = KGEO(i,j,bi,bj) + 1 |
185 |
|
|
ENDDO |
186 |
|
|
ENDIF |
187 |
|
|
#endif /* SEAICE_BICE_STRESS */ |
188 |
|
|
ENDDO |
189 |
|
|
ENDDO |
190 |
|
|
ENDDO |
191 |
|
|
ENDDO |
192 |
|
|
#endif /* SEAICE_CGRID */ |
193 |
|
|
|
194 |
|
|
#ifdef ALLOW_DIAGNOSTICS |
195 |
|
|
IF ( useDiagnostics ) THEN |
196 |
|
|
CALL SEAICE_DIAGNOSTICS_INIT( myThid ) |
197 |
|
|
ENDIF |
198 |
|
|
#endif |
199 |
|
|
|
200 |
heimbach |
1.1 |
RETURN |
201 |
|
|
END |