/[MITgcm]/MITgcm/pkg/seaice/seaice_model.F
ViewVC logotype

Contents of /MITgcm/pkg/seaice/seaice_model.F

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


Revision 1.8 - (show annotations) (download)
Tue Jun 24 16:09:37 2003 UTC (20 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51, checkpoint51f_post, checkpoint51d_post, checkpoint51j_post, checkpoint51b_pre, checkpoint51h_pre, branchpoint-genmake2, checkpoint51b_post, checkpoint51c_post, checkpoint51i_pre, checkpoint51e_post, checkpoint51f_pre, checkpoint51g_post, checkpoint51a_post
Branch point for: branch-genmake2
Changes since 1.7: +21 -4 lines
Merging for c51 vs. e34

1 C $Header:
2
3 #include "SEAICE_OPTIONS.h"
4
5 CStartOfInterface
6 SUBROUTINE seaice_model( myTime, myIter, myThid )
7 C /===========================================================\
8 C | SUBROUTINE SEAICE_MODEL |
9 C | o Time stepping of a dynamic/thermodynamic sea ice model. |
10 C | Dynamics solver: Zhang/Hibler, JGR, 102, 8691-8702, 1997 |
11 C | and Zhang/Rothrock, JGR, 105, 3325-3338, 2000 |
12 C | Thermodynamics: Hibler, MWR, 108, 1943-1973, 1980 |
13 C | Rheology: Hibler, JPO, 9, 815- 846, 1979 |
14 C | Snow: Zhang et al. , JPO, 28, 191- 217, 1998 |
15 C | Parallel forward ice model written by Jinlun Zhang PSC/UW|
16 C | & coupled into MITgcm by Dimitris Menemenlis (JPL) 2/2001|
17 C | zhang@apl.washington.edu / menemenlis@jpl.nasa.gov |
18 C |===========================================================|
19 C \===========================================================/
20 IMPLICIT NONE
21
22 C === Global variables ===
23 #include "SIZE.h"
24 #include "EEPARAMS.h"
25 #include "DYNVARS.h"
26 #include "PARAMS.h"
27 #include "FFIELDS.h"
28 #include "SEAICE.h"
29 #include "SEAICE_PARAMS.h"
30 #include "SEAICE_EXTERNAL.h"
31
32 #ifdef SEAICE_EXTERNAL_FORCING
33 # include "SEAICE_FFIELDS.h"
34 #endif
35
36 #ifdef ALLOW_AUTODIFF_TAMC
37 # include "tamc.h"
38 #endif
39
40 C === Routine arguments ===
41 C myTime - Simulation time
42 C myIter - Simulation timestep number
43 C myThid - Thread no. that called this routine.
44 _RL myTime
45 INTEGER myIter
46 INTEGER myThid
47 CEndOfInterface
48
49 #ifdef ALLOW_SEAICE
50
51 C === Local variables ===
52 C i,j,bi,bj - Loop counters
53
54 INTEGER i, j, bi, bj
55
56 #ifdef SEAICE_EXTERNAL_FORCING
57 C-- Atmospheric state and runoff are from
58 C pkg/exf, which does not update edges.
59 _EXCH_XY_R8( uwind, myThid )
60 _EXCH_XY_R8( vwind, myThid )
61 _EXCH_XY_R8( atemp, myThid )
62 _EXCH_XY_R8( aqh, myThid )
63 _EXCH_XY_R8( lwdown, myThid )
64 _EXCH_XY_R8( swdown, mythid )
65 _EXCH_XY_R8( precip, myThid )
66 _EXCH_XY_R8( evap, myThid )
67 _EXCH_XY_R8( runoff, myThid )
68 #else /* SEAICE_EXTERNAL_FORCING */
69 C-- Load atmospheric state and runoff.
70 CALL SEAICE_GET_FORCING ( myTime, myIter, myThid )
71 #endif /* SEAICE_EXTERNAL_FORCING */
72
73 C-- Third level model velocity is used as proxy for geostrophic velocity
74 DO bj=myByLo(myThid),myByHi(myThid)
75 DO bi=myBxLo(myThid),myBxHi(myThid)
76 DO j=0,sNy+1
77 DO i=0,sNx+1
78 GWATX(I,J,bi,bj)=HALF*(uVel(i,j,3,bi,bj)
79 & +uVel(i,j-1,3,bi,bj))
80 GWATY(I,J,bi,bj)=HALF*(vVel(i,j,3,bi,bj)
81 & +vVel(i-1,j,3,bi,bj))
82 #ifdef SEAICE_DEBUG
83 c write(*,'(2i4,2i2,f7.1,7f12.3)')
84 c & ,i,j,bi,bj,UVM(I,J,bi,bj)
85 c & ,GWATX(I,J,bi,bj),GWATY(I,J,bi,bj)
86 c & ,uVel(i+1,j,3,bi,bj),uVel(i+1,j+1,3,bi,bj)
87 c & ,vVel(i,j+1,3,bi,bj),vVel(i+1,j+1,3,bi,bj)
88 #endif
89 ENDDO
90 ENDDO
91 ENDDO
92 ENDDO
93
94 #ifdef ALLOW_AUTODIFF_TAMC
95 CADJ STORE uwind = comlev1, key = ikey_dynamics
96 CADJ STORE vwind = comlev1, key = ikey_dynamics
97 # ifdef SEAICE_ALLOW_DYNAMICS
98 CADJ STORE heff = comlev1, key = ikey_dynamics
99 CADJ STORE area = comlev1, key = ikey_dynamics
100 # endif
101 #endif /* ALLOW_AUTODIFF_TAMC */
102
103 C solve ice momentum equations and calculate ocean surface stress
104 CALL DYNSOLVER ( myTime, myIter, myThid )
105
106 #ifdef ALLOW_AUTODIFF_TAMC
107 # ifdef SEAICE_ALLOW_DYNAMICS
108 CADJ STORE heff = comlev1, key = ikey_dynamics
109 CADJ STORE area = comlev1, key = ikey_dynamics
110 CADJ STORE uice = comlev1, key = ikey_dynamics
111 CADJ STORE vice = comlev1, key = ikey_dynamics
112 # endif
113 #endif /* ALLOW_AUTODIFF_TAMC */
114 C NOW DO ADVECTION
115 CALL ADVECT( UICE, VICE, HEFF, HEFFM, myThid )
116 CALL ADVECT( UICE, VICE, AREA, HEFFM, myThid )
117
118 C NOW DO GROWTH
119 C MUST CALL GROWTH ONLY AFTER CALLING ADVECTION
120 CALL GROWTH( myTime, myIter, myThid)
121
122 C-- Update overlap regions for a bunch of stuff
123 _BARRIER
124 CALL SEAICE_EXCH( HEFF, myThid )
125 CALL SEAICE_EXCH( AREA, myThid )
126 _EXCH_XY_R4(fu , myThid )
127 _EXCH_XY_R4(fv , myThid )
128 _EXCH_XY_R4(EmPmR, myThid )
129 _EXCH_XY_R4(Qnet , myThid )
130 _EXCH_XY_R4(surfaceTendencyTice, myThid )
131 #ifdef SHORTWAVE_HEATING
132 _EXCH_XY_R4(Qsw , myThid )
133 #endif
134 _EXCH_XYZ_R8(theta , myThid )
135
136 C-- Sea ice diagnostics.
137 CALL SEAICE_DO_DIAGS( myTime, myIter, myThid )
138
139 C-- Write sea ice restart files
140 CALL SEAICE_WRITE_PICKUP ( .FALSE.,
141 & myTime+deltaTClock, myIter+1, myThid )
142
143 C---------------------------------------------------
144 C OOH NOOOO we need to move the whole stuff
145 C---------------------------------------------------
146 #ifdef ALLOW_AUTODIFF_TAMC
147 CRG CADJ store UICE,VICE,AREA,HEFF,fu,fv,EmPmR,Qnet,Qsw = comlev1_bibj
148 #endif
149
150 C-- Call sea-ice cost function routine
151 CRG CALL SEAICE_COST( myTime, myIter, myThid )
152
153 #endif /* ALLOW_SEAICE */
154
155 RETURN
156 END

  ViewVC Help
Powered by ViewVC 1.1.22