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

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

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


Revision 1.3 - (hide annotations) (download)
Fri Jul 7 20:10:36 2006 UTC (17 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58l_post, mitgcm_mapl_00, checkpoint58u_post, checkpoint58r_post, checkpoint58n_post, checkpoint58t_post, checkpoint58q_post, checkpoint58o_post, checkpoint58v_post, checkpoint58s_post, checkpoint58p_post, checkpoint58m_post
Changes since 1.2: +7 -15 lines
take bi,bj loops outside of calc_gw & timestep_wvel.F
 and move the call inside DYNAMICS bi,bj loops

1 jmc 1.3 C $Header: /u/gcmpack/MITgcm/model/src/timestep_wvel.F,v 1.2 2006/02/23 20:55:49 jmc Exp $
2 jmc 1.1 C $Name: $
3    
4     #include "PACKAGES_CONFIG.h"
5     #include "CPP_OPTIONS.h"
6    
7     CBOP
8     C !ROUTINE: TIMESTEP_WVEL
9     C !INTERFACE:
10     SUBROUTINE TIMESTEP_WVEL(
11 jmc 1.3 I bi,bj, myTime, myIter, myThid )
12 jmc 1.1 C !DESCRIPTION: \bv
13     C *==========================================================*
14     C | S/R TIMESTEP_WVEL
15     C | o Step model vertical velocity forward in time
16     C *==========================================================*
17     C \ev
18    
19     C !USES:
20     IMPLICIT NONE
21     C == Global variables ==
22     #include "SIZE.h"
23     #include "EEPARAMS.h"
24     #include "PARAMS.h"
25     #include "GRID.h"
26     #include "DYNVARS.h"
27     #include "NH_VARS.h"
28     c #include "SURFACE.h"
29    
30     C !INPUT/OUTPUT PARAMETERS:
31     C == Routine Arguments ==
32 jmc 1.3 INTEGER bi,bj
33 jmc 1.1 _RL myTime
34     INTEGER myIter, myThid
35    
36     #ifdef ALLOW_NONHYDROSTATIC
37     C !LOCAL VARIABLES:
38     C == Local variables ==
39     INTEGER iMin,iMax,jMin,jMax
40 jmc 1.3 INTEGER i,j,k
41 jmc 1.1 _RL gWtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
42 jmc 1.2 _RL tmpFac, nh_Fac, igwFac
43 jmc 1.1 CEOP
44    
45     iMin = 1
46     iMax = sNx
47     jMin = 1
48     jMax = sNy
49    
50 jmc 1.2 igwFac = 0.
51     IF ( implicitIntGravWave ) igwFac = horiVertRatio*horiVertRatio
52    
53 jmc 1.1 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
54    
55 jmc 1.3 IF ( nonHydrostatic ) THEN
56 jmc 1.2 nh_Fac = 0.
57     IF ( nh_Am2.NE.0. ) nh_Fac = 1. _d 0 / nh_Am2
58    
59     k = 1
60     DO j=1-Oly,sNy+Oly
61     DO i=1-Olx,sNx+Olx
62     gW(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
63     ENDDO
64     ENDDO
65    
66 jmc 1.1 DO k=2,Nr
67    
68 jmc 1.2 C apply mask to gW and keep a copy of wVel in gW:
69     DO j=1-Oly,sNy+Oly
70     DO i=1-Olx,sNx+Olx
71     gWtmp(i,j) = gW(i,j,k,bi,bj)
72 jmc 1.1 & *maskC(i,j,k,bi,bj)*maskC(i,j,k-1,bi,bj)
73 jmc 1.2 gW(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
74 jmc 1.1 ENDDO
75     ENDDO
76 jmc 1.2 C Step forward vertical velocity
77     tmpFac = nh_Fac + igwFac*dBdrRef(k)*deltaTMom*dTtracerLev(k)
78     IF (tmpFac.GT.0. ) tmpFac = 1. _d 0 / tmpFac
79 jmc 1.1 DO j=jMin,jMax
80     DO i=iMin,iMax
81     wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
82 jmc 1.2 & + deltaTmom*tmpFac*gWtmp(i,j)
83 jmc 1.1 ENDDO
84     ENDDO
85    
86     C- End of k loop
87     ENDDO
88    
89     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
90    
91     #ifdef ALLOW_OBCS
92     C-- This call is aesthetic: it makes the W field
93     C consistent with the OBs but this has no algorithmic
94     C impact. This is purely for diagnostic purposes.
95     IF (useOBCS) THEN
96     DO k=1,Nr
97     CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
98     ENDDO
99     ENDIF
100     #endif /* ALLOW_OBCS */
101    
102 jmc 1.3 ELSEIF ( implicitIntGravWave ) THEN
103 jmc 1.2 C keep a copy of wVel in gW:
104     DO k=1,Nr
105     DO j=1-Oly,sNy+Oly
106     DO i=1-Olx,sNx+Olx
107     gW(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
108     ENDDO
109     ENDDO
110     ENDDO
111    
112     C- End if nonHydrostatic / elseif implicitIntGravWave
113 jmc 1.3 ENDIF
114 jmc 1.1
115     #endif /* ALLOW_NONHYDROSTATIC */
116    
117     RETURN
118     END

  ViewVC Help
Powered by ViewVC 1.1.22