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

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

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


Revision 1.6 - (hide annotations) (download)
Mon Jun 22 15:26:25 1998 UTC (25 years, 11 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint11, checkpoint10, checkpoint13, checkpoint9, checkpoint8, checkpoint12, branch-point-rdot
Branch point for: branch-rdot
Changes since 1.5: +3 -2 lines
Various changes including time-dependant forcing:
 o logic for controlling external forcing fields now allows
   for time-dependant forcing: load_external_fields.F
 o genmake.dec needed a special line for the above file.
 o theta* and salt* time-stepping algorithm were re-implemented.
The 4x4 global configuration has been "double-checked" against
CNH's version. However, we do not assume any responsibility for
the correctness of this code ...  8-)

1 adcroft 1.6 C $Header: /u/gcmpack/models/MITgcmUV/model/src/correction_step.F,v 1.5 1998/06/16 17:07:11 adcroft Exp $
2 adcroft 1.1
3     #include "CPP_EEOPTIONS.h"
4    
5     C /==========================================================\
6     C | S/R CORRECTION_STEP |
7     C | o Corrects the horizontal flow fields with the surface |
8     C | pressure gradient. |
9     C \==========================================================/
10     SUBROUTINE CORRECTION_STEP( bi, bj, iMin, iMax, jMin, jMax,
11     I K, pSurfX, pSurfY,
12 adcroft 1.6 I myCurrentTime, myThid )
13 adcroft 1.1 implicit none
14     ! Common
15     #include "SIZE.h"
16     #include "DYNVARS.h"
17 cnh 1.3 #include "EEPARAMS.h"
18 adcroft 1.1 #include "PARAMS.h"
19     #include "GRID.h"
20     #include "CG2D.h"
21     C == Routine Arguments ==
22     INTEGER bi,bj,iMin,iMax,jMin,jMax
23     INTEGER K
24     _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
25     _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
26     INTEGER myThid
27 adcroft 1.6 _RL myCurrentTime
28 adcroft 1.1 C == Local variables ==
29     INTEGER i,j
30 adcroft 1.5 _RL hxFac,hyFac,rRhoNil
31 adcroft 1.1
32     C On/off scaling paramters
33     hxFac = pfFacMom
34     hyFac = pfFacMom
35    
36     rRhoNil=1. / rhonil
37    
38     C Step forward zonal velocity
39     DO j=jMin,jMax
40     DO i=iMin,iMax
41     uVel(i,j,k,bi,bj)=( gUNm1(i,j,k,bi,bj)
42     & -deltaTmom*hxFac*rRhonil *pSurfX(i,j)
43     & )*_maskW(i,j,k,bi,bj)
44     gUNm1(i,j,k,bi,bj)=gU(i,j,k,bi,bj)
45     ENDDO
46     ENDDO
47    
48     C Step forward meridional velocity
49     DO j=jMin,jMax
50     DO i=iMin,iMax
51     vVel(i,j,k,bi,bj)=( gVNm1(i,j,k,bi,bj)
52     & -deltaTmom*hyFac*rRhonil *pSurfY(i,j)
53     & )*_maskS(i,j,k,bi,bj)
54     gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj)
55 adcroft 1.4 ENDDO
56     ENDDO
57    
58     C Rotate theta/gT/gTnm1
59     DO j=jMin,jMax
60     DO i=iMin,iMax
61     theta(i,j,k,bi,bj)=gTNm1(i,j,k,bi,bj)
62     gTNm1(i,j,k,bi,bj)=gT(i,j,k,bi,bj)
63     ENDDO
64     ENDDO
65    
66     C Rotate salt/gS/gSnm1
67     DO j=jMin,jMax
68     DO i=iMin,iMax
69     salt(i,j,k,bi,bj)=gSNm1(i,j,k,bi,bj)
70     gSNm1(i,j,k,bi,bj)=gS(i,j,k,bi,bj)
71 adcroft 1.1 ENDDO
72     ENDDO
73    
74     RETURN
75     END

  ViewVC Help
Powered by ViewVC 1.1.22