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

Contents of /MITgcm/model/src/correction_step.F

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


Revision 1.10 - (show annotations) (download)
Fri Nov 6 22:44:45 1998 UTC (25 years, 7 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint19, checkpoint18
Changes since 1.9: +16 -12 lines
Changes to allow for atmospheric integration builds of the code

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/correction_step.F,v 1.9 1998/08/20 20:05:00 cnh Exp $
2
3 #include "CPP_OPTIONS.h"
4
5 C /==========================================================\
6 C | S/R CORRECTION_STEP |
7 C | o Corrects the horizontal flow fields with the surface |
8 C | slope. |
9 C \==========================================================/
10 SUBROUTINE CORRECTION_STEP( bi, bj, iMin, iMax, jMin, jMax,
11 I K, etaSurfX, etaSurfY,
12 I myCurrentTime, myThid )
13 IMPLICIT NONE
14 ! Common
15 #include "SIZE.h"
16 #include "DYNVARS.h"
17 #include "EEPARAMS.h"
18 #include "PARAMS.h"
19 #include "GRID.h"
20 #include "CG2D.h"
21 C == Routine Arguments ==
22 C etaSurfX, etaSurfY - Surface slope
23 C bi,bj,iMin,iMax,jMin,jMax, K - Loop counters
24 C myThid - Instance number for
25 C this call to S/R CORRECTION_STEP
26 C myCurrentTime - Current simulation time for this instance.
27 _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
28 _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
29 INTEGER bi,bj,iMin,iMax,jMin,jMax
30 INTEGER K
31 INTEGER myThid
32 _RL myCurrentTime
33
34 C == Local variables ==
35 INTEGER i,j
36 _RL hxFac,hyFac,rRhoNil
37
38 C On/off scaling paramters
39 hxFac = pfFacMom
40 hyFac = pfFacMom
41
42 C Step forward zonal velocity
43 DO j=jMin,jMax
44 DO i=iMin,iMax
45 uVel(i,j,k,bi,bj)=( gUNm1(i,j,k,bi,bj)
46 & -deltaTmom*hxFac*gBaro*etaSurfX(i,j)
47 & )*_maskW(i,j,k,bi,bj)
48 gUNm1(i,j,k,bi,bj)=gU(i,j,k,bi,bj)
49 ENDDO
50 ENDDO
51
52 C Step forward meridional velocity
53 DO j=jMin,jMax
54 DO i=iMin,iMax
55 vVel(i,j,k,bi,bj)=( gVNm1(i,j,k,bi,bj)
56 & -deltaTmom*hyFac*gBaro*etaSurfY(i,j)
57 & )*_maskS(i,j,k,bi,bj)
58 gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj)
59 ENDDO
60 ENDDO
61
62 C Rotate theta/gT/gTnm1
63 IF ( tempStepping ) THEN
64 DO j=jMin,jMax
65 DO i=iMin,iMax
66 theta(i,j,k,bi,bj)=gTNm1(i,j,k,bi,bj)
67 gTNm1(i,j,k,bi,bj)=gT(i,j,k,bi,bj)
68 ENDDO
69 ENDDO
70 ENDIF
71
72 C Rotate salt/gS/gSnm1
73 IF ( saltStepping ) THEN
74 DO j=jMin,jMax
75 DO i=iMin,iMax
76 salt(i,j,k,bi,bj)=gSNm1(i,j,k,bi,bj)
77 gSNm1(i,j,k,bi,bj)=gS(i,j,k,bi,bj)
78 ENDDO
79 ENDDO
80 ENDIF
81
82 RETURN
83 END

  ViewVC Help
Powered by ViewVC 1.1.22