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

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

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


Revision 1.3 - (hide annotations) (download)
Wed Jun 7 01:55:12 2006 UTC (18 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58l_post, mitgcm_mapl_00, checkpoint58u_post, checkpoint58w_post, checkpoint60, checkpoint61, checkpoint58r_post, checkpoint58n_post, checkpoint58x_post, checkpoint58t_post, checkpoint58h_post, checkpoint58q_post, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint58j_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint59, checkpoint58i_post, checkpoint58o_post, checkpoint58y_post, checkpoint58k_post, checkpoint58v_post, checkpoint58s_post, checkpoint61b, checkpoint61c, checkpoint58p_post, checkpoint61a, checkpoint58m_post
Changes since 1.2: +3 -3 lines
Modifications for bottom topography control
o replace hFacC by _hFacC at various places
o replace ALLOW_HFACC_CONTROL by ALLOW_DEPTH_CONTROL
o add non-self-adjoint cg2d_nsa
o update autodiff support routines
o re-initialise hfac after ctrl_depth_ini
o works for 5x5 box, doesnt work for global_ocean.90x40x15

1 heimbach 1.3 C $Header: /u/gcmpack/MITgcm/model/src/diags_rho.F,v 1.2 2005/11/05 01:00:57 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: DIAGS_RHO
9     C !INTERFACE:
10     SUBROUTINE DIAGS_RHO(
11     I k, bi, bj,
12     I rhoK, rhoKm1,
13     I myTime, myIter, myThid)
14     C !DESCRIPTION: \bv
15     C *==========================================================*
16     C | S/R DIAGS_RHO
17     C | o Buoyancy Flux diagnostics
18     C *==========================================================*
19     C *==========================================================*
20     C \ev
21    
22     C !USES:
23     IMPLICIT NONE
24     C == Global variables ==
25     #include "SIZE.h"
26     #include "EEPARAMS.h"
27     #include "PARAMS.h"
28     #include "GRID.h"
29     c #include "SURFACE.h"
30     #include "DYNVARS.h"
31    
32     C !INPUT/OUTPUT PARAMETERS:
33     C == Routine Arguments ==
34     C k, bi,bj :: level & tile indices
35     C phiHydC :: hydrostatic potential anomaly at cell center
36     C (atmos: =Geopotential ; ocean-z: =Pressure/rho)
37     C myTime :: Current time
38     C myIter :: Current iteration number
39     C myThid :: Instance number for this call of the routine.
40     INTEGER k, bi,bj
41     _RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
42     _RL rhoKm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
43     _RL myTime
44     INTEGER myIter, myThid
45     CEOP
46    
47     #ifdef ALLOW_DIAGNOSTICS
48    
49     C !LOCAL VARIABLES:
50     C == Local variables ==
51     C i,j :: Loop counters
52     INTEGER i,j
53 jmc 1.2 _RL tmpFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
54 jmc 1.1 LOGICAL DIAGNOSTICS_IS_ON
55     EXTERNAL DIAGNOSTICS_IS_ON
56    
57     IF ( DIAGNOSTICS_IS_ON('RHOANOSQ',myThid) ) THEN
58     DO j=1,sNy
59     DO i=1,sNx
60     tmpFld(i,j) = rhoK(i,j)*rhoK(i,j)
61     ENDDO
62     ENDDO
63     CALL DIAGNOSTICS_FILL(tmpFld,'RHOANOSQ',k,1,2,bi,bj,myThid)
64     ENDIF
65    
66     IF ( DIAGNOSTICS_IS_ON('URHOMASS',myThid) ) THEN
67     DO j=1,sNy
68     DO i=1,sNx+1
69 heimbach 1.3 tmpFld(i,j) = uVel(i,j,k,bi,bj)*_hFacW(i,j,k,bi,bj)
70 jmc 1.1 & *(rhoK(i-1,j)+rhoK(i,j))*0.5 _d 0
71     ENDDO
72     ENDDO
73     CALL DIAGNOSTICS_FILL(tmpFld,'URHOMASS',k,1,2,bi,bj,myThid)
74     ENDIF
75    
76     IF ( DIAGNOSTICS_IS_ON('VRHOMASS',myThid) ) THEN
77     DO j=1,sNy+1
78     DO i=1,sNx
79 heimbach 1.3 tmpFld(i,j) = vVel(i,j,k,bi,bj)*_hFacS(i,j,k,bi,bj)
80 jmc 1.1 & *(rhoK(i,j-1)+rhoK(i,j))*0.5 _d 0
81     ENDDO
82     ENDDO
83     CALL DIAGNOSTICS_FILL(tmpFld,'VRHOMASS',k,1,2,bi,bj,myThid)
84     ENDIF
85    
86     IF ( DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) THEN
87     IF ( k.EQ.1 ) THEN
88     DO j=1,sNy
89     DO i=1,sNx
90     tmpFld(i,j) = wVel(i,j,k,bi,bj)*rhoK(i,j)
91     ENDDO
92     ENDDO
93     ELSE
94     DO j=1,sNy
95     DO i=1,sNx
96     tmpFld(i,j) = wVel(i,j,k,bi,bj)
97     & *(rhoKm1(i,j)+rhoK(i,j))*0.5 _d 0
98     ENDDO
99     ENDDO
100     ENDIF
101     CALL DIAGNOSTICS_FILL(tmpFld,'WRHOMASS',k,1,2,bi,bj,myThid)
102     ENDIF
103    
104     #endif /* ALLOW_DIAGNOSTICS */
105    
106     RETURN
107     END

  ViewVC Help
Powered by ViewVC 1.1.22