/[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.2 - (show annotations) (download)
Mon Jun 1 22:27:14 1998 UTC (26 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint5, checkpoint6
Changes since 1.1: +1 -9 lines
Implemented implicit vertical diffusion (tracers only).
Involved introducing a "total" diffusivity array (local 3D)
calculated by calc_diffusivity().
Made some small changes to time-stepping algorithm.
Switched on by setting implicitZdiffusion.
(note: *Not* fully tested with topography. But when switched off
this does produce identical results)

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/correction_step.F,v 1.1 1998/06/01 20:36:13 adcroft Exp $
2
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 I myThid )
13 implicit none
14 ! Common
15 #include "SIZE.h"
16 #include "DYNVARS.h"
17 #include "PARAMS.h"
18 #include "GRID.h"
19 #include "EEPARAMS.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 C == Local variables ==
28 INTEGER i,j
29 _RL ab15,ab05,hxFac,hyFac,rRhoNil
30
31 C Adams-Bashforth timestepping weights
32 ab15=1.5+abeps
33 ab05=-0.5-abeps
34
35 C On/off scaling paramters
36 hxFac = pfFacMom
37 hyFac = pfFacMom
38
39 rRhoNil=1. / rhonil
40
41 C Step forward zonal velocity
42 DO j=jMin,jMax
43 DO i=iMin,iMax
44 uVel(i,j,k,bi,bj)=( gUNm1(i,j,k,bi,bj)
45 & -deltaTmom*hxFac*rRhonil *pSurfX(i,j)
46 & )*_maskW(i,j,k,bi,bj)
47 gUNm1(i,j,k,bi,bj)=gU(i,j,k,bi,bj)
48 ENDDO
49 ENDDO
50
51 C Step forward meridional velocity
52 DO j=jMin,jMax
53 DO i=iMin,iMax
54 vVel(i,j,k,bi,bj)=( gVNm1(i,j,k,bi,bj)
55 & -deltaTmom*hyFac*rRhonil *pSurfY(i,j)
56 & )*_maskS(i,j,k,bi,bj)
57 gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj)
58 ENDDO
59 ENDDO
60
61 RETURN
62 END

  ViewVC Help
Powered by ViewVC 1.1.22