1 | !$Id: sulfint.F90 104 2008-12-23 10:28:51Z acosce $ |
---|
2 | !! ========================================================================= |
---|
3 | !! INCA - INteraction with Chemistry and Aerosols |
---|
4 | !! |
---|
5 | !! Copyright Laboratoire des Sciences du Climat et de l'Environnement (LSCE) |
---|
6 | !! Unite mixte CEA-CNRS-UVSQ |
---|
7 | !! |
---|
8 | !! Contributors to this INCA subroutine: |
---|
9 | !! |
---|
10 | !! Didier Hauglustaine, LSCE, hauglustaine@cea.fr |
---|
11 | !! |
---|
12 | !! Anne Cozic, LSCE, anne.cozic@cea.fr |
---|
13 | !! Yann Meurdesoif, LSCE, yann.meurdesoif@cea.fr |
---|
14 | !! |
---|
15 | !! This software is a computer program whose purpose is to simulate the |
---|
16 | !! atmospheric gas phase and aerosol composition. The model is designed to be |
---|
17 | !! used within a transport model or a general circulation model. This version |
---|
18 | !! of INCA was designed to be coupled to the LMDz GCM. LMDz-INCA accounts |
---|
19 | !! for emissions, transport (resolved and sub-grid scale), photochemical |
---|
20 | !! transformations, and scavenging (dry deposition and washout) of chemical |
---|
21 | !! species and aerosols interactively in the GCM. Several versions of the INCA |
---|
22 | !! model are currently used depending on the envisaged applications with the |
---|
23 | !! chemistry-climate model. |
---|
24 | !! |
---|
25 | !! This software is governed by the CeCILL license under French law and |
---|
26 | !! abiding by the rules of distribution of free software. You can use, |
---|
27 | !! modify and/ or redistribute the software under the terms of the CeCILL |
---|
28 | !! license as circulated by CEA, CNRS and INRIA at the following URL |
---|
29 | !! "http://www.cecill.info". |
---|
30 | !! |
---|
31 | !! As a counterpart to the access to the source code and rights to copy, |
---|
32 | !! modify and redistribute granted by the license, users are provided only |
---|
33 | !! with a limited warranty and the software's author, the holder of the |
---|
34 | !! economic rights, and the successive licensors have only limited |
---|
35 | !! liability. |
---|
36 | !! |
---|
37 | !! In this respect, the user's attention is drawn to the risks associated |
---|
38 | !! with loading, using, modifying and/or developing or reproducing the |
---|
39 | !! software by the user in light of its specific status of free software, |
---|
40 | !! that may mean that it is complicated to manipulate, and that also |
---|
41 | !! therefore means that it is reserved for developers and experienced |
---|
42 | !! professionals having in-depth computer knowledge. Users are therefore |
---|
43 | !! encouraged to load and test the software's suitability as regards their |
---|
44 | !! requirements in conditions enabling the security of their systems and/or |
---|
45 | !! data to be ensured and, more generally, to use and operate it in the |
---|
46 | !! same conditions as regards security. |
---|
47 | !! |
---|
48 | !! The fact that you are presently reading this means that you have had |
---|
49 | !! knowledge of the CeCILL license and that you accept its terms. |
---|
50 | !! ========================================================================= |
---|
51 | |
---|
52 | #include <inca_define.h> |
---|
53 | SUBROUTINE SULF_INTI (filename) |
---|
54 | !---------------------------------------------------------------------- |
---|
55 | ! ... SULFATE climatologies initialize |
---|
56 | ! Didier Hauglustaine, IPSL, 2000. |
---|
57 | !---------------------------------------------------------------------- |
---|
58 | |
---|
59 | USE AC_SULF |
---|
60 | USE INCA_DIM |
---|
61 | USE MOD_INCA_PARA |
---|
62 | IMPLICIT NONE |
---|
63 | |
---|
64 | INCLUDE 'netcdf.inc' |
---|
65 | |
---|
66 | !---------------------------------------------------------------------- |
---|
67 | ! ... Dummy args |
---|
68 | !---------------------------------------------------------------------- |
---|
69 | CHARACTER(len=*), INTENT(in) :: filename |
---|
70 | |
---|
71 | !---------------------------------------------------------------------- |
---|
72 | ! ... Local variables |
---|
73 | !---------------------------------------------------------------------- |
---|
74 | INTEGER :: iret |
---|
75 | INTEGER :: varid |
---|
76 | INTEGER :: error |
---|
77 | |
---|
78 | ! ... Open the file |
---|
79 | !$OMP MASTER |
---|
80 | IF (is_mpi_root) THEN |
---|
81 | iret = nf_open(filename, 0, ncid) |
---|
82 | CALL check_err(iret, 'SULF_INTI', 'problem to open file') |
---|
83 | |
---|
84 | ! ... Get lengths |
---|
85 | iret = nf_inq_dimid(ncid, 'time', varid) |
---|
86 | CALL check_err(iret, 'SULF_INTI', 'problem to check dimid time') |
---|
87 | iret = nf_inq_dimlen(ncid, varid, ntimes) |
---|
88 | CALL check_err(iret, 'SULF_INTI', 'problem to check dimlen time') |
---|
89 | iret = nf_inq_dimid(ncid, 'lev', varid) |
---|
90 | CALL check_err(iret, 'SULF_INTI', 'problem to check dimid lev') |
---|
91 | iret = nf_inq_dimlen(ncid, varid, nlevs) |
---|
92 | CALL check_err(iret, 'SULF_INTI', 'problem to check dimlen lev') |
---|
93 | iret = nf_inq_dimid(ncid, 'vector', varid) |
---|
94 | CALL check_err(iret, 'SULF_INTI', 'problem to check dimid vector') |
---|
95 | iret = nf_inq_dimlen(ncid, varid, klons) |
---|
96 | CALL check_err(iret, 'SULF_INTI', 'problem to check dimlen vector') |
---|
97 | ENDIF |
---|
98 | !$OMP END MASTER |
---|
99 | |
---|
100 | CALL bcast(ntimes) |
---|
101 | CALL bcast(nlevs) |
---|
102 | CALL bcast(klons) |
---|
103 | |
---|
104 | IF (klons/=nbp_glo) THEN |
---|
105 | write(lunout,*) 'sulf_inti : [klons] [nbp_glo]', klons, nbp_glo |
---|
106 | call print_err(3, 'SULF_INTI','there a problem of vector size', 'check klons and nbp_glo', '') |
---|
107 | ELSE |
---|
108 | klons=PLON |
---|
109 | ENDIF |
---|
110 | |
---|
111 | ! ... Allocate variables |
---|
112 | ALLOCATE(times(ntimes), STAT=error) |
---|
113 | IF (error /= 0) call print_err(3, 'SULF_INTI', 'Space requested not possible for variable times ', '', '') |
---|
114 | ALLOCATE(levs(nlevs), STAT=error) |
---|
115 | IF (error /= 0) CALL print_err(3, 'SULF_INTI', 'Space requested not possible for variable levs', '', '') |
---|
116 | ALLOCATE(so4bd(klons,nlevs,2), STAT=error) |
---|
117 | IF (error /= 0) CALL print_err(3, 'SULF_INTI', 'Space requested not possible for variable so4bd ', '', '') |
---|
118 | |
---|
119 | ! ... Get the coordinates |
---|
120 | !$OMP MASTER |
---|
121 | IF (is_mpi_root) THEN |
---|
122 | iret = nf_inq_varid(ncid, 'time', varid) |
---|
123 | CALL check_err(iret, 'SULF_INTI', 'problem to check varid time') |
---|
124 | iret = nf_get_var_double(ncid, varid, times) |
---|
125 | CALL check_err(iret, 'SULF_INTI', 'problem to read time') |
---|
126 | |
---|
127 | |
---|
128 | iret = nf_inq_varid(ncid, 'lev', varid) |
---|
129 | CALL check_err(iret, 'SULF_INTI', 'problem to check varid lev') |
---|
130 | iret = nf_get_var_double(ncid, varid, levs) |
---|
131 | CALL check_err(iret, 'SULF_INTI', 'problem to read lev') |
---|
132 | |
---|
133 | |
---|
134 | |
---|
135 | ! ... Variable ids |
---|
136 | iret = nf_inq_varid(ncid, 'so4', so4_id) |
---|
137 | CALL check_err(iret, 'SULF_INTI', 'problem to check varid so4') |
---|
138 | ENDIF |
---|
139 | !$OMP END MASTER |
---|
140 | CALL bcast(times) |
---|
141 | CALL bcast(so4_id) |
---|
142 | call bcast(levs) |
---|
143 | |
---|
144 | END SUBROUTINE SULF_INTI |
---|
145 | |
---|
146 | SUBROUTINE SULF_READ(calday) |
---|
147 | !---------------------------------------------------------------------- |
---|
148 | ! ... Read SO4 distributions |
---|
149 | ! Didier Hauglustaine, IPSL, 2000. |
---|
150 | !---------------------------------------------------------------------- |
---|
151 | |
---|
152 | USE AC_SULF |
---|
153 | USE INCA_DIM |
---|
154 | USE MOD_INCA_PARA |
---|
155 | USE PRINT_INCA |
---|
156 | IMPLICIT NONE |
---|
157 | |
---|
158 | INCLUDE 'netcdf.inc' |
---|
159 | |
---|
160 | !---------------------------------------------------------------------- |
---|
161 | ! ... Dummy args |
---|
162 | !---------------------------------------------------------------------- |
---|
163 | REAL, INTENT(IN) :: calday |
---|
164 | |
---|
165 | !---------------------------------------------------------------------- |
---|
166 | ! ... Local variables |
---|
167 | !---------------------------------------------------------------------- |
---|
168 | INTEGER :: iret |
---|
169 | INTEGER :: varid |
---|
170 | INTEGER :: oldlotime, oldhitime |
---|
171 | INTEGER, DIMENSION(3) :: start3, count3 |
---|
172 | INTEGER, DIMENSION(2) :: start2, count2 |
---|
173 | REAL :: so4bd_glo(nbp_glo,nlevs,2) |
---|
174 | |
---|
175 | ! ... Check to see if model time still bounded by data set |
---|
176 | oldlotime = lotime |
---|
177 | oldhitime = hitime |
---|
178 | CALL FINDPLB (times, ntimes, calday, lotime) |
---|
179 | IF ( cyclical ) THEN |
---|
180 | hitime = MOD(lotime,ntimes) + 1 |
---|
181 | ELSE |
---|
182 | hitime = lotime + 1 |
---|
183 | END IF |
---|
184 | |
---|
185 | ! ... Read new data if needed |
---|
186 | IF ( hitime /= oldhitime ) THEN |
---|
187 | !$OMP MASTER |
---|
188 | IF (is_mpi_root) THEN |
---|
189 | loind = hiind |
---|
190 | hiind = MOD( loind, 2) + 1 |
---|
191 | |
---|
192 | count2(1) = nbp_glo |
---|
193 | count2(2) = 1 |
---|
194 | start2(1) = 1 |
---|
195 | count3(1) = nbp_glo |
---|
196 | count3(2) = nlevs |
---|
197 | count3(3) = 1 |
---|
198 | start3(1) = 1 |
---|
199 | start3(2) = 1 |
---|
200 | |
---|
201 | WRITE(lunout, *) 'SULF_READ : read new data for time ', times(hitime) |
---|
202 | start3(3) = hitime |
---|
203 | start2(2) = hitime |
---|
204 | iret = nf_get_vara_double(ncid, so4_id, start3, count3, so4bd_glo(1,1,hiind)) |
---|
205 | CALL check_err(iret, 'SULF_READ', 'problem to read so4bd_glo hiind') |
---|
206 | |
---|
207 | ENDIF |
---|
208 | !$OMP END MASTER |
---|
209 | CALL scatter(so4bd_glo(:,:,hiind),so4bd(:,:,hiind)) |
---|
210 | |
---|
211 | IF ( lotime /= oldlotime ) THEN |
---|
212 | !$OMP MASTER |
---|
213 | IF (is_mpi_root) THEN |
---|
214 | WRITE(lunout, *) 'SULF_READ : read new data for time ', times(lotime) |
---|
215 | start3(3) = lotime |
---|
216 | start2(2) = lotime |
---|
217 | iret = nf_get_vara_double(ncid, so4_id, start3, count3, so4bd_glo(1,1,loind)) |
---|
218 | CALL check_err(iret, 'SULF_READ', 'problem to read so4bd_glo loind') |
---|
219 | |
---|
220 | ENDIF |
---|
221 | !$OMP END MASTER |
---|
222 | CALL scatter(so4bd_glo(:,:,loind),so4bd(:,:,loind)) |
---|
223 | |
---|
224 | END IF |
---|
225 | |
---|
226 | END IF |
---|
227 | |
---|
228 | END SUBROUTINE SULF_READ |
---|
229 | |
---|
230 | SUBROUTINE SULF_INTERP(calday, pmid, zmid) |
---|
231 | !---------------------------------------------------------------------- |
---|
232 | ! ... Interpolate SO4 on model time |
---|
233 | ! Didier Hauglustaine, IPSL, 2000. |
---|
234 | !---------------------------------------------------------------------- |
---|
235 | |
---|
236 | USE AC_SULF |
---|
237 | USE INCA_DIM |
---|
238 | IMPLICIT NONE |
---|
239 | |
---|
240 | !---------------------------------------------------------------------- |
---|
241 | ! ... Dummy args |
---|
242 | !---------------------------------------------------------------------- |
---|
243 | REAL, INTENT(IN) :: calday |
---|
244 | REAL, INTENT(IN) :: pmid(PLON,PLEV) |
---|
245 | REAL, INTENT(IN) :: zmid(PLON,PLEV) |
---|
246 | |
---|
247 | !---------------------------------------------------------------------- |
---|
248 | ! ... Local variables |
---|
249 | !---------------------------------------------------------------------- |
---|
250 | REAL, PARAMETER :: dzs = 1. !specific to Considine's input file |
---|
251 | INTEGER :: k, i, m, j, l |
---|
252 | INTEGER :: lonlev |
---|
253 | REAL :: dt, dt1, tint |
---|
254 | REAL :: wght |
---|
255 | REAL, DIMENSION(PLON,PLEV) :: so4wrk |
---|
256 | |
---|
257 | lonlev = PLON * PLEV |
---|
258 | |
---|
259 | !---------------------------------------------------------------------- |
---|
260 | ! Note: 365 day year inconsistent with LMDz (360 days) !!! |
---|
261 | !---------------------------------------------------------------------- |
---|
262 | ! ... interpolate linearly on current model time |
---|
263 | |
---|
264 | IF ( times(hitime) < times(lotime) ) THEN |
---|
265 | dt = 365. - times(lotime) + times(hitime) |
---|
266 | IF (calday <= times(hitime) ) THEN |
---|
267 | dt1 = 365. - times(lotime) + calday |
---|
268 | ELSE |
---|
269 | dt1 = calday - times(lotime) |
---|
270 | END IF |
---|
271 | ELSE |
---|
272 | dt = times(hitime) - times(lotime) |
---|
273 | dt1 = calday - times(lotime) |
---|
274 | END IF |
---|
275 | |
---|
276 | tint = dt1/dt |
---|
277 | |
---|
278 | ! ... Interpolate to local time |
---|
279 | |
---|
280 | CALL linintp (& |
---|
281 | lonlev, & |
---|
282 | 0., & |
---|
283 | 1., & |
---|
284 | tint, & |
---|
285 | so4bd(1,1,loind), & |
---|
286 | so4bd(1,1,hiind), & |
---|
287 | so4wrk) |
---|
288 | |
---|
289 | !-------------------------------------------------------- |
---|
290 | ! ... Change units (from kg/m3 to g/cm3) |
---|
291 | !-------------------------------------------------------- |
---|
292 | |
---|
293 | DO i = 1, PLON |
---|
294 | DO l = 1, PLEV |
---|
295 | |
---|
296 | so4(i,l) = so4wrk(i,l) *1.e-3 |
---|
297 | |
---|
298 | END DO |
---|
299 | END DO |
---|
300 | |
---|
301 | |
---|
302 | END SUBROUTINE SULF_INTERP |
---|
303 | |
---|