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

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

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


Revision 1.3 - (show annotations) (download)
Tue Jan 20 20:47:42 2015 UTC (9 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, HEAD
Changes since 1.2: +14 -12 lines
- move ALLOW_EDDYPSI block out of DYNVARS.h and merge it into FFIELDS.h
- rename uMean,vMean --> uEulerMean,vEulerMean

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_eddy_stress.F,v 1.2 2014/01/01 23:25:59 m_bates Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_GMREDI
8 # include "GMREDI_OPTIONS.h"
9 #endif
10
11 SUBROUTINE CALC_EDDY_STRESS( bi, bj, myThid )
12
13 C !DESCRIPTION: \bv
14 C *==========================================================*
15 C | S/R CALC_EDDY_STRESS
16 C | o Calculates the eddy stress when running a residual
17 C | ocean model
18 C *==========================================================*
19 C | Calculates the eddy stress. Later this will be added to
20 C | gU the same as external sources (e.g. wind stress, bottom
21 C | friction, etc.
22 C *==========================================================*
23 C \ev
24
25 IMPLICIT NONE
26
27 #include "SIZE.h"
28 #include "EEPARAMS.h"
29 #include "PARAMS.h"
30 #include "GRID.h"
31 #include "FFIELDS.h"
32 #ifdef ALLOW_GMREDI
33 # include "GMREDI.h"
34 #endif
35
36 C !INPUT/OUTPUT PARAMETERS:
37 C == Routine arguments ==
38 INTEGER bi,bj
39 INTEGER myThid
40 #ifdef ALLOW_EDDYPSI
41
42 C !LOCAL VARIABLES:
43 C == Local variables ==
44 C Loop counters
45 INTEGER i,j,k
46 C Interpolated stream function and coriolis
47 _RL psix, psiy, coriU, coriV
48
49 C Calculate the eddy stress from the eddy induced streamfunction
50 #ifdef ALLOW_GMREDI
51 IF ( GM_InMomAsStress ) THEN
52 #endif
53 DO k=1,Nr
54
55 DO j=1-OLy,sNy+OLy-1
56 DO i=1-OLx+1,sNx+OLx
57 #ifdef ALLOW_GMREDI
58 psiy = op25*(GM_PsiY(i, j ,k,bi,bj)
59 & +GM_PsiY(i, j+1,k,bi,bj)
60 & +GM_PsiY(i-1,j ,k,bi,bj)
61 & +GM_PsiY(i-1,j+1,k,bi,bj))
62 #else
63 psiy = op25*(eddyPsiY(i, j ,k,bi,bj)
64 & +eddyPsiY(i, j+1,k,bi,bj)
65 & +eddyPsiY(i-1,j ,k,bi,bj)
66 & +eddyPsiY(i-1,j+1,k,bi,bj))
67 #endif
68 coriU = op5*(fcori(i-1,j,bi,bj)
69 & +fCori(i ,j,bi,bj))
70 tauxEddy(i,j,k,bi,bj) = rhoConst*coriU*psiy
71 ENDDO
72 ENDDO
73
74 DO j=1-OLy+1,sNy+OLy
75 DO i=1-OLx,sNx+OLx-1
76 #ifdef ALLOW_GMREDI
77 psix = op25*(GM_PsiX(i, j ,k,bi,bj)
78 & +GM_PsiX(i+1,j ,k,bi,bj)
79 & +GM_PsiX(i ,j-1,k,bi,bj)
80 & +GM_PsiX(i+1,j-1,k,bi,bj))
81 #else
82 psix = op25*(eddyPsiX(i, j ,k,bi,bj)
83 & +eddyPsiX(i+1,j ,k,bi,bj)
84 & +eddyPsiX(i ,j-1,k,bi,bj)
85 & +eddyPsiX(i+1,j-1,k,bi,bj))
86 #endif
87 coriV = op5*(fcori(i,j-1,bi,bj)
88 & +fCori(i,j ,bi,bj))
89 tauyEddy(i,j,k,bi,bj) = -rhoConst*coriV*psix
90 ENDDO
91 ENDDO
92
93 C- end k loop
94 ENDDO
95
96 #ifdef ALLOW_DIAGNOSTICS
97 IF ( useDiagnostics ) THEN
98 CALL DIAGNOSTICS_FILL( tauxEddy, 'TAUXEDDY',
99 & 0, Nr, 1, bi, bj, myThid )
100 CALL DIAGNOSTICS_FILL( tauyEddy, 'TAUYEDDY',
101 & 0, Nr, 1, bi, bj, myThid )
102 ENDIF
103 #endif
104 #ifdef ALLOW_GMREDI
105 ENDIF
106 #endif
107 #endif /* ALLOW_EDDYPSI */
108 RETURN
109 END

  ViewVC Help
Powered by ViewVC 1.1.22