/[MITgcm]/MITgcm/pkg/zonal_filt/zonal_filt_init.F
ViewVC logotype

Contents of /MITgcm/pkg/zonal_filt/zonal_filt_init.F

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


Revision 1.7 - (show annotations) (download)
Fri Feb 18 00:33:05 2011 UTC (13 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint63g, checkpoint63, checkpoint63p, checkpoint63q, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.6: +33 -21 lines
stop if trying to use several tiles in X direction

1 C $Header: /u/gcmpack/MITgcm/pkg/zonal_filt/zonal_filt_init.F,v 1.6 2010/01/31 20:24:16 jmc Exp $
2 C $Name: $
3
4 #include "ZONAL_FILT_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: ZONAL_FILT_INIT
8 C !INTERFACE:
9 SUBROUTINE ZONAL_FILT_INIT( myThid )
10
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | S/R ZONAL_FILT_INIT
14 C | o Initialise FFT filter for latitude circle.
15 C *==========================================================*
16 C | The details of particular FFT libraries may differ.
17 C | Changing to a different library may entail modifying the
18 C | code here. However, the broad process is usually the
19 C | same.
20 C | Note - Fourier modes for sNx and sNx+1 are damped in the
21 C | same way. This is because we have not implemented
22 C | a scheme that sets the damping factor for the
23 C | highest wave number for odd sNx. Instead the
24 C | highest wave number for odd sNx. Instead only
25 C | wave numbers 1:INT(sNx/2) are partially damped.
26 C | Wave number sNx/2 (if it exists) is removed
27 C | altogether.
28 C *==========================================================*
29 C \ev
30
31 C !USES:
32 IMPLICIT NONE
33
34 C === Global variables ===
35 #include "SIZE.h"
36 #include "EEPARAMS.h"
37 #include "PARAMS.h"
38 #include "GRID.h"
39 #include "ZONAL_FILT.h"
40 #include "FFTPACK.h"
41
42 C !INPUT/OUTPUT PARAMETERS:
43 C == Routine arguments ==
44 C myThid :: my Thread Id number
45 INTEGER myThid
46 CEOP
47
48 #ifdef ALLOW_ZONAL_FILT
49 C !LOCAL VARIABLES:
50 C == Local variables ==
51 C alpha :: Used to evaluate frequency and latitude dependent
52 C amplitude damping factor.
53 C wvNum :: Wave number
54 C lat :: Temporary holding latitude
55 C nWv :: No. of waves that fit on grid.
56 C msgBuf :: Informational/error message buffer
57 c _RL alpha, wvNum
58 c INTEGER nWv
59 INTEGER i, j, bi, bj
60 CHARACTER*(MAX_LEN_MBUF) msgBuf
61
62 C !FUNCTIONS:
63 _RL oneRL
64 _RL ampfact
65 PARAMETER( oneRL = 1. _d 0 )
66 _RS lat
67 ampfact(lat,i) = MIN( oneRL,
68 & ( COS( ABS(lat)*deg2rad )
69 & /COS( zonal_filt_lat*deg2rad ) )**zonal_filt_cospow
70 & /(SIN( PI*FLOAT(i)/FLOAT(Nx) ) )**zonal_filt_sinpow
71 & )
72
73 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
74
75 IF ( sNx.NE.Nx ) THEN
76 WRITE(msgBuf,'(A,I3,A)')
77 & 'S/R ZONAL_FILT_INIT: Multi-tiles ( nSx*nPx=', nSx*nPx, ' )'
78 CALL PRINT_ERROR( msgBuf, myThid )
79 WRITE(msgBuf,'(A)')
80 & ' in Zonal (X) dir. not implemented in Zonal-Filter code'
81 CALL PRINT_ERROR( msgBuf, myThid )
82 STOP 'ABNORMAL END: S/R ZONAL_FILT_INIT'
83 ENDIF
84
85 _BEGIN_MASTER(myThid)
86 C o Initialise specific library FFT package
87 DO bj=1,nSy
88 c CALL R8FFTI( Nx, FFTPACKWS(1,bj) )
89 CALL R8FFTI1( Nx, FFTPACKWS2(1,bj), FFTPACKWS3(1,bj) )
90 ENDDO
91
92 C o Set amplitude scale factor as function of latitude and mode number
93 DO bj=1,nSy
94 DO bi=1,nSx
95 DO j=1-oLy,sNy+Oly
96 ampFactor(1,j,bi,bj) = oneRL
97 ampFactorV(1,j,bi,bj) = oneRL
98 DO i=1,Nx/2-1
99 ampFactor(2*i,j,bi,bj) = ampfact( yC(1,j,bi,bj) , I )
100 c IF (ampFactor(2*i,j,bi,bj).LE..9) ampFactor(2*i,j,bi,bj)=0.
101 ampFactor(2*I+1,j,bi,bj) = ampFactor(2*i,j,bi,bj)
102 ampFactorV(2*i,j,bi,bj) = ampfact( yG(1,j,bi,bj) , I )
103 c IF (ampFactorV(2*i,j,bi,bj).LE..9) ampFactorV(2*i,j,bi,bj)=0.
104 ampFactorV(2*I+1,j,bi,bj) = ampFactorV(2*i,j,bi,bj)
105 ENDDO
106
107 i=Nx/2
108 IF ( zonal_filt_mode2dx.EQ.0 ) THEN
109 ampFactor(Nx,j,bi,bj) = ampfact( yC(1,j,bi,bj) , i )
110 ampFactorV(Nx,j,bi,bj) = ampfact( yG(1,j,bi,bj) , i )
111 ELSE
112 ampFactor(Nx,j,bi,bj) = 0.
113 ampFactorV(Nx,j,bi,bj) = 0.
114 ENDIF
115
116 ENDDO
117 ENDDO
118 ENDDO
119 _END_MASTER(myThid)
120 CALL BAR2(myThid)
121
122 CALL WRITE_REC_XY_RL( 'ampFactor', ampFactor, 1, 0, myThid )
123
124 #endif /* ALLOW_ZONAL_FILT */
125
126 RETURN
127 END

  ViewVC Help
Powered by ViewVC 1.1.22