1 | module eofdump_mod |
---|
2 | ! this module controls the production of data for EOFs |
---|
3 | ! it won't work if run in parallel (but it's OK, we don't use it anymore...) |
---|
4 | ! Mainly kept for reference. |
---|
5 | implicit none |
---|
6 | ! Dump profiles for EOFs every ieofs physics timesteps, |
---|
7 | ! starting at first call; |
---|
8 | integer :: ieofs |
---|
9 | ! Dump profiles every eofskip points in each direction |
---|
10 | ! on the model grid. |
---|
11 | integer, parameter :: eofskip = 4 |
---|
12 | ! Units for writing EOF header and data: |
---|
13 | integer, parameter :: uehead = 82, uedata = 83 |
---|
14 | |
---|
15 | contains |
---|
16 | |
---|
17 | subroutine eofdump(ngrid,nlayer,u,v,t,rho,ps) |
---|
18 | |
---|
19 | use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat |
---|
20 | implicit none |
---|
21 | ! |
---|
22 | ! Dumps profiles for calculation of variability EOFs |
---|
23 | ! Modified to include rho, FF 09/2004 |
---|
24 | ! Corrected small bug in sampling rate/count, EM 11/2007 |
---|
25 | ! |
---|
26 | ! |
---|
27 | |
---|
28 | integer,intent(in) :: ngrid ! total number of physics grid points |
---|
29 | integer,intent(in) :: nlayer ! number of atmospheric layers |
---|
30 | real,intent(in) :: u(ngrid,nlayer) |
---|
31 | real,intent(in) :: v(ngrid,nlayer) |
---|
32 | real,intent(in) :: t(ngrid,nlayer) |
---|
33 | real,intent(in) :: rho(ngrid,nlayer) |
---|
34 | real,intent(in) :: ps(ngrid) |
---|
35 | integer,save :: count=0 |
---|
36 | integer i,j,l, ig |
---|
37 | |
---|
38 | LOGICAL,SAVE :: firstcall=.true. |
---|
39 | |
---|
40 | !$OMP THREADPRIVATE(count,firstcall) |
---|
41 | |
---|
42 | !------------------------------------------------------- |
---|
43 | ! Initialization at first call: |
---|
44 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
45 | IF (firstcall) THEN |
---|
46 | write(*,*) 'CALL ineofdump' |
---|
47 | CALL ineofdump(ngrid,nlayer) |
---|
48 | firstcall=.false. |
---|
49 | END IF |
---|
50 | |
---|
51 | !------------------------------------------------------- |
---|
52 | ! Dumps every ieofs physics timesteps |
---|
53 | ! |
---|
54 | ! write(*,*)'eofdump:count=',count,' ps(1)=',ps(1) |
---|
55 | ! if ((ieofs.gt.0).and.(mod(count,ieofs).eq.0)) then |
---|
56 | if (mod(count+1,ieofs).eq.0) then |
---|
57 | ! write(*,*)'eofdump: dump --> ps(1)=',ps(1) |
---|
58 | do i=1,nbp_lon,eofskip |
---|
59 | do j=1+eofskip/2,nbp_lat,eofskip |
---|
60 | ig = 1+ (j-2)*nbp_lon +i |
---|
61 | #ifdef NC_DOUBLE |
---|
62 | write(uedata) (real(u(ig,l)),l=1,nlayer) |
---|
63 | write(uedata) (real(v(ig,l)),l=1,nlayer) |
---|
64 | write(uedata) (real(t(ig,l)),l=1,nlayer) |
---|
65 | write(uedata) (real(rho(ig,l)),l=1,nlayer) |
---|
66 | write(uedata) real(ps(ig)) |
---|
67 | #else |
---|
68 | write(uedata) (u(ig,l),l=1,nlayer) |
---|
69 | write(uedata) (v(ig,l),l=1,nlayer) |
---|
70 | write(uedata) (t(ig,l),l=1,nlayer) |
---|
71 | write(uedata) (rho(ig,l),l=1,nlayer) |
---|
72 | write(uedata) ps(ig) |
---|
73 | #endif |
---|
74 | enddo |
---|
75 | enddo |
---|
76 | endif |
---|
77 | count=count+1 |
---|
78 | |
---|
79 | end subroutine eofdump |
---|
80 | |
---|
81 | |
---|
82 | subroutine ineofdump(ngrid,nlayer) |
---|
83 | |
---|
84 | use geometry_mod, only: longitude, latitude |
---|
85 | use nrtype, only: pi |
---|
86 | use time_phylmdz_mod, only: daysec, dtphys |
---|
87 | USE vertical_layers_mod, ONLY: aps,bps |
---|
88 | use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat |
---|
89 | implicit none |
---|
90 | ! |
---|
91 | ! Initialise dumping of profiles for EOF calculations |
---|
92 | ! |
---|
93 | |
---|
94 | integer,intent(in) :: ngrid ! total number of physics grid points |
---|
95 | integer,intent(in) :: nlayer ! number of atmospheric layers |
---|
96 | integer ig,i,j,l |
---|
97 | logical,save :: firstcall=.true. |
---|
98 | integer,save :: npgrid |
---|
99 | |
---|
100 | !$OMP THREADPRIVATE(firstcall,npgrid) |
---|
101 | |
---|
102 | |
---|
103 | if (firstcall) then |
---|
104 | npgrid=ngrid+2*(nbp_lon-1) |
---|
105 | firstcall=.false. |
---|
106 | endif |
---|
107 | |
---|
108 | ! |
---|
109 | ! Set frequency for dumping at once per day |
---|
110 | ! |
---|
111 | ieofs=nint(daysec/dtphys) |
---|
112 | if (abs(float(ieofs)-daysec/dtphys).gt.1.e-8*daysec) then |
---|
113 | call abort_physic("ineofdump",' 1 day .ne. n physics timesteps',1) |
---|
114 | endif |
---|
115 | ! |
---|
116 | ! Header |
---|
117 | ! |
---|
118 | open(uehead,file='profiles.hdr',form='formatted') |
---|
119 | write(uehead,*) 0.E+0,0,0,ieofs,1,0 |
---|
120 | write(uehead,*) nbp_lon,npgrid/nbp_lon,npgrid,nlayer |
---|
121 | |
---|
122 | do i=1,nbp_lon,eofskip |
---|
123 | do j=1+eofskip/2,nbp_lat,eofskip |
---|
124 | ig = 1+ (j-2)*nbp_lon +i |
---|
125 | if(j.eq.1) then |
---|
126 | call abort_physic("ineofdump",'Error: j==1',1) |
---|
127 | endif |
---|
128 | if(j.eq.nbp_lat) then |
---|
129 | call abort_physic("ineofdump",'Error: j==nbp_lat',1) |
---|
130 | endif |
---|
131 | #ifdef NC_DOUBLE |
---|
132 | write(uehead,*) real(longitude(ig)*180./pi),real(latitude(ig)*180./pi) |
---|
133 | #else |
---|
134 | write(uehead,*) longitude(ig)*180./pi, latitude(ig)*180./pi |
---|
135 | #endif |
---|
136 | ! write(*,*) 'eof grid j=',j,' lat= ', latitude(ig)*180./pi |
---|
137 | enddo |
---|
138 | enddo |
---|
139 | |
---|
140 | #ifdef NC_DOUBLE |
---|
141 | write(uehead,*) real(aps) |
---|
142 | write(uehead,*) real(bps) |
---|
143 | #else |
---|
144 | write(uehead,*) aps |
---|
145 | write(uehead,*) bps |
---|
146 | #endif |
---|
147 | close(uehead) |
---|
148 | ! |
---|
149 | ! Main profile file |
---|
150 | ! |
---|
151 | open(uedata,file='profiles.dat',form='unformatted') |
---|
152 | end subroutine ineofdump |
---|
153 | |
---|
154 | end module eofdump_mod |
---|