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

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

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


Revision 1.33 - (hide annotations) (download)
Thu Nov 13 06:35:14 2003 UTC (20 years, 6 months ago) by dimitri
Branch: MAIN
CVS Tags: hrcube4, checkpoint52d_pre, checkpoint52j_pre, checkpoint52k_post, checkpoint52f_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint52e_pre, checkpoint52e_post, checkpoint52b_pre, checkpoint52b_post, checkpoint52c_post, checkpoint52f_pre, checkpoint52d_post, checkpoint52i_post, checkpoint52h_pre, checkpoint52j_post, branch-netcdf, checkpoint52a_post
Branch point for: netcdf-sm0
Changes since 1.32: +3 -1 lines
o modifications to make FREEZE flux visible to pkg/kpp
  - moved surfaceTendencyTice from pkg/seaice to main code
  - FREEZE moved to FORWARD_STEP
  - subroutine FREEZE now limits only surface temperature
    this means new output.txt for global_ocean.90x40x15,
    global_ocean.cs32x15, and global_with_exf, but note
    that results for these three experiments remain
    bit-identical to before if allowFreezing=.FALSE.
o added surface flux output variables to TIMEAVE_STATVARS
o time-averaged output for pkg/ptracers

1 dimitri 1.33 C $Header: /usr/local/gcmpack/MITgcm/model/src/initialise_varia.F,v 1.32 2003/10/09 04:19:18 edhill Exp $
2 jmc 1.22 C $Name: $
3 adcroft 1.1
4 edhill 1.32 #include "PACKAGES_CONFIG.h"
5 adcroft 1.1 #include "CPP_OPTIONS.h"
6    
7 cnh 1.20 CBOP
8     C !ROUTINE: INITIALISE_VARIA
9     C !INTERFACE:
10 adcroft 1.1 SUBROUTINE INITIALISE_VARIA(myThid)
11 cnh 1.20 C !DESCRIPTION: \bv
12     C *==========================================================*
13     C | SUBROUTINE INITIALISE_VARIA
14     C | o Set the initial conditions for dynamics variables
15     C | and time dependent arrays
16     C *==========================================================*
17     C | This routine reads/writes data from an input file and
18     C | from various binary files.
19     C | Each thread invokes an instance of this routine as does
20     C | each process in a multi-process parallel environment like
21     C | MPI.
22     C *==========================================================*
23     C \ev
24    
25     C !CALLING SEQUENCE:
26     C INITIALISE_VARIA
27     C |
28     C |-- INI_LINEAR_PHISURF
29     C |
30     C |-- INI_CORI
31     C |
32     C |-- INI_CG2D
33     C |
34     C |-- INI_CG3D
35     C |
36     C |-- INI_MIXING
37     C |
38     C |-- INI_DYNVARS
39     C |
40     C |-- INI_FIELDS
41     C |
42     C |-- INI_TR1
43     C |
44 dimitri 1.28 C |-- INI_AUTODIFF
45     C |
46 adcroft 1.23 C |-- PACKAGES_INIT_VARIABLES
47     C |
48 jmc 1.30 C |-- THE_CORRECTION_STEP (if restart from old pickup files)
49 cnh 1.20 C |
50     C |-- CALL CONVECTIVE_ADJUSTMENT_INI
51     C |
52     C |-- CALC_SURF_DR
53     C |
54     C |-- UPDATE_SURF_DR
55     C |
56     C |-- UPDATE_CG2D
57     C |
58 jmc 1.26 C |-- INTEGR_CONTINUITY
59 cnh 1.20 C |
60     C |-- TIMEAVE_STATVARS
61 dimitri 1.33 C |
62     C |-- PTRACERS_STATVARS
63 jmc 1.26 C |
64     C |-- STATE_SUMMARY
65 cnh 1.20
66     C !USES:
67 adcroft 1.1 IMPLICIT NONE
68     C == Global variables ==
69     #include "SIZE.h"
70     #include "EEPARAMS.h"
71     #include "PARAMS.h"
72 adcroft 1.5 #include "DYNVARS.h"
73 adcroft 1.1
74 cnh 1.20 C !INPUT/OUTPUT PARAMETERS:
75 adcroft 1.1 C == Routine arguments ==
76     INTEGER myThid
77    
78 cnh 1.20 C !LOCAL VARIABLES:
79 adcroft 1.1 C == Local variables ==
80 adcroft 1.5 INTEGER bi,bj,K,iMin,iMax,jMin,jMax
81 cnh 1.20 CEOP
82 heimbach 1.14
83     #ifdef ALLOW_TAMC_CHECKPOINTING
84 heimbach 1.15
85     nIter0 = INT( startTime/deltaTClock )
86    
87     C-- Set Bo_surf => define the Linear Relation: Phi_surf(eta)
88     CALL INI_LINEAR_PHISURF( myThid )
89    
90     C-- Set coriolis operators
91     CALL INI_CORI( myThid )
92    
93     C-- Set laplace operators for use in 2D conjugate gradient solver.
94     CALL INI_CG2D( myThid )
95    
96     #ifdef ALLOW_NONHYDROSTATIC
97     C-- Set laplace operators for use in 3D conjugate gradient solver.
98 edhill 1.32 ceh3 should add an IF ( useNONHYDROSTATIC ) THEN
99 heimbach 1.15 CALL INI_CG3D( myThid )
100 heimbach 1.14 #endif
101 adcroft 1.1
102 heimbach 1.15 #endif /* ALLOW_TAMC_CHECKPOINTING */
103 heimbach 1.16 _BARRIER
104 heimbach 1.15
105 heimbach 1.16 C-- Initialise 3-dim. diffusivities
106     CALL INI_MIXING( myThid )
107 adcroft 1.5 _BARRIER
108 heimbach 1.16
109 adcroft 1.12 C-- Initialize DYNVARS arrays (state fields + G terms: Gu,Gv,...) to zero [always]
110     CALL INI_DYNVARS( myThid )
111    
112 adcroft 1.1 C-- Initialise model fields.
113     C Starting values of U, V, W, temp., salt. and tendency terms
114     C are set here. Fields are either set to default or read from
115     C stored files.
116     CALL INI_FIELDS( myThid )
117     _BARRIER
118 jmc 1.8
119 heimbach 1.15 #ifdef ALLOW_PASSIVE_TRACER
120 heimbach 1.13 C-- Initialise passive tracer(s)
121     CALL INI_TR1( myThid )
122 heimbach 1.21 _BARRIER
123     #endif
124    
125     #ifdef ALLOW_AUTODIFF_TAMC
126     C-- Initialise active fields to help TAMC
127     CALL INI_AUTODIFF( myThid )
128 heimbach 1.13 _BARRIER
129 heimbach 1.15 #endif
130 heimbach 1.13
131 heimbach 1.25 C-- Initialize variable data for packages
132     CALL PACKAGES_INIT_VARIABLES( myThid )
133    
134 heimbach 1.18 #ifndef ALLOW_AUTODIFF_TAMC
135 heimbach 1.11 IF ( usePickupBeforeC35 ) THEN
136 jmc 1.8 C-- IMPORTANT : Need to activate the following call to restart from
137     C a pickup file written by MITgcmUV_checkpoint34 or earlier.
138 heimbach 1.11 IF ( startTime .NE. 0. ) THEN
139     CALL THE_CORRECTION_STEP(startTime, nIter0, myThid)
140     ENDIF
141     ENDIF
142 heimbach 1.18 #endif
143 adcroft 1.1
144 adcroft 1.5 C-- Initial conditions are convectively adjusted (for historical reasons)
145 edhill 1.32 IF ( startTime .EQ. 0. .AND. cAdjFreq .NE. 0. ) THEN
146 heimbach 1.13 CADJ loop = parallel
147 adcroft 1.1 DO bj = myByLo(myThid), myByHi(myThid)
148 heimbach 1.13 CADJ loop = parallel
149 adcroft 1.1 DO bi = myBxLo(myThid), myBxHi(myThid)
150 adcroft 1.5 iMin=1-Olx
151     iMax=sNx+Olx
152     jMin=1-Oly
153     jMax=sNy+Oly
154 heimbach 1.11 CALL CONVECTIVE_ADJUSTMENT_INI(
155 adcroft 1.5 I bi, bj, iMin, iMax, jMin, jMax,
156 jmc 1.7 I startTime, nIter0, myThid )
157 adcroft 1.1 ENDDO
158     ENDDO
159 adcroft 1.5 _BARRIER
160     END IF
161 jmc 1.17
162     #ifdef NONLIN_FRSURF
163 jmc 1.19 C-- Compute the surface level thickness <-- function of etaH(n)
164 jmc 1.17 C and modify hFac(C,W,S) accordingly :
165 jmc 1.29 IF ( select_rStar.NE.0 )
166     & CALL CALC_R_STAR(etaH, startTime, -1 , myThid )
167 jmc 1.17 IF (nonlinFreeSurf.GT.0) THEN
168 jmc 1.29 IF ( select_rStar.GT.0 ) THEN
169     CALL UPDATE_R_STAR( startTime, nIter0, myThid )
170     ELSE
171 jmc 1.22 CALL CALC_SURF_DR(etaH, startTime, -1 , myThid )
172 jmc 1.17 CALL UPDATE_SURF_DR( startTime, nIter0, myThid )
173 jmc 1.29 ENDIF
174 jmc 1.17 ENDIF
175     C- update also CG2D matrix (and preconditioner)
176     IF ( nonlinFreeSurf.GT.2) THEN
177     CALL UPDATE_CG2D( startTime, nIter0, myThid )
178     ENDIF
179 jmc 1.29 #endif /* NONLIN_FRSURF */
180 adcroft 1.1
181 jmc 1.26 DO bj=myByLo(myThid),myByHi(myThid)
182     DO bi=myBxLo(myThid),myBxHi(myThid)
183 adcroft 1.1
184 jmc 1.26 C-- Integrate continuity vertically for vertical velocity
185     CALL INTEGR_CONTINUITY( bi, bj, uVel, vVel,
186     I startTime, nIter0, myThid )
187    
188 jmc 1.9 #ifdef ALLOW_TIMEAVE
189 jmc 1.7 C-- initialise time-average arrays with initial state values
190 jmc 1.26 IF (taveFreq.GT.0.) THEN
191 jmc 1.9 CALL TIMEAVE_STATVARS(startTime, nIter0, bi, bj, myThid)
192 jmc 1.26 ENDIF
193 stephd 1.31 cswdptr -- add ---
194     #ifdef ALLOW_PTRACERS
195 edhill 1.32 che3 needs an IF ( usePTRACERS ) THEN
196 stephd 1.31 CALL PTRACERS_STATVARS(startTime, nIter0, bi, bj, myThid)
197     #endif
198     cswdptr -- end add ---
199 jmc 1.9 #endif /* ALLOW_TIMEAVE */
200 jmc 1.26
201     C-- end bi,bj loop.
202     ENDDO
203     ENDDO
204    
205     C-- Fill in overlap regions for wVel :
206     _EXCH_XYZ_R8(wVel,myThid)
207    
208     C-- Finally summarise the model state
209     CALL STATE_SUMMARY( myThid )
210 jmc 1.7
211     RETURN
212 adcroft 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22