1 | ! |
---|
2 | ! $Id: read_inp.F,v 1.21 2005/10/10 13:40:18 pmarches Exp $ |
---|
3 | ! |
---|
4 | #include "cppdefs.h" |
---|
5 | ! Read/report model input |
---|
6 | subroutine read_inp (ierr) ! parameters from startup |
---|
7 | ! script file using keywords |
---|
8 | ! to recognize variables. |
---|
9 | ! implicit none |
---|
10 | #include "param.h" |
---|
11 | #include "scalars.h" |
---|
12 | #include "ncscrum.h" |
---|
13 | #include "sediment.h" |
---|
14 | #ifdef FLOATS |
---|
15 | # include "ncscrum_floats.h" |
---|
16 | #endif |
---|
17 | #ifdef STATIONS |
---|
18 | # include "nc_sta.h" |
---|
19 | #endif |
---|
20 | #ifdef PSOURCE |
---|
21 | # include "sources.h" |
---|
22 | #endif |
---|
23 | #ifdef MPI |
---|
24 | include 'mpif.h' |
---|
25 | #endif |
---|
26 | integer kwsize, testunit, input |
---|
27 | parameter (kwsize=32, testunit=40, input=15) |
---|
28 | character end_signal*3, keyword*32, fname*64 |
---|
29 | parameter (end_signal='end') |
---|
30 | integer ierr, iargc, is,ie, kwlen, lstr, lenstr |
---|
31 | #ifdef SOLVE3D |
---|
32 | & , itrc |
---|
33 | #endif |
---|
34 | #ifdef AGRIF |
---|
35 | integer Agrif_lev_sedim |
---|
36 | #endif |
---|
37 | ! |
---|
38 | ! Use pre-set default startup filename for known applications, |
---|
39 | ! or get it as an argument from command line via iargc-getarg |
---|
40 | ! (override default). NOTE: The usage of the executable should |
---|
41 | ! be either |
---|
42 | ! roms |
---|
43 | ! or |
---|
44 | ! roms startup_file_name |
---|
45 | ! |
---|
46 | ! WITHOUT the UNIX redirection (<): roms<startup_file like it |
---|
47 | ! used to be. |
---|
48 | ! |
---|
49 | #if defined SOLITON |
---|
50 | fname='roms.in.Soliton' |
---|
51 | #elif defined OVERFLOW |
---|
52 | fname='roms.in.Overflow' |
---|
53 | #elif defined SEAMOUNT |
---|
54 | fname='roms.in.Seamount' |
---|
55 | #elif defined UPWELLING |
---|
56 | fname='roms.in.Upwelling' |
---|
57 | #elif defined REGIONAL |
---|
58 | fname='roms.in' |
---|
59 | #elif defined BASIN |
---|
60 | fname='roms.in.Basin' |
---|
61 | #elif defined INNERSHELF |
---|
62 | fname='roms.in.Innershelf' |
---|
63 | #elif defined SHELFRONT |
---|
64 | fname='roms.in.Shelfront' |
---|
65 | #elif defined CANYON_A || defined CANYON_B |
---|
66 | fname='roms.in.Canyon' |
---|
67 | #elif defined GRAV_ADJ |
---|
68 | fname='roms.in.Grav_adj' |
---|
69 | #else |
---|
70 | fname='roms.in' |
---|
71 | #endif |
---|
72 | #ifdef MPI |
---|
73 | if (mynode.eq.0 .and. iargc().GT.0) call getarg(1,fname) |
---|
74 | call MPI_Bcast(fname,64,MPI_BYTE, 0, MPI_COMM_WORLD,ierr) |
---|
75 | #else |
---|
76 | if (iargc().eq.1) call getarg(1,fname) |
---|
77 | #endif |
---|
78 | ! |
---|
79 | ! if child grid, use an input name: fname.1, .2, .3, ... |
---|
80 | ! |
---|
81 | #ifdef AGRIF |
---|
82 | if (.Not.Agrif_Root()) then |
---|
83 | lstr=lenstr(fname) |
---|
84 | fname=fname(1:lstr) / / '.' / / Agrif_Cfixed() |
---|
85 | lstr=lenstr(fname) |
---|
86 | endif |
---|
87 | #endif |
---|
88 | ! |
---|
89 | wrthis(indxTime)=.false. |
---|
90 | #ifdef AVERAGES |
---|
91 | wrtavg(indxTime)=.false. |
---|
92 | #endif |
---|
93 | ! |
---|
94 | ! Read in keyword: keep trying, until keyword is found. |
---|
95 | ! ==== == ======== ==== ======= ===== ======= == ====== |
---|
96 | ! |
---|
97 | ierr=0 ! <-- reset error counter |
---|
98 | call setup_kwds (ierr) |
---|
99 | open (input,file=fname,status='old',form='formatted',err=97) |
---|
100 | 1 keyword=' ' |
---|
101 | read(input,'(A)',err=1,end=99) keyword |
---|
102 | if (ichar(keyword(1:1)).eq.33) goto 1 |
---|
103 | is=1 |
---|
104 | 2 if (is.eq.kwsize) then |
---|
105 | goto 1 |
---|
106 | elseif (keyword(is:is).eq.' ') then |
---|
107 | is=is+1 |
---|
108 | goto 2 |
---|
109 | endif |
---|
110 | ie=is |
---|
111 | 3 if (keyword(ie:ie).eq.':') then |
---|
112 | keyword(ie:ie)=' ' |
---|
113 | goto 4 !--> recognized keyword. |
---|
114 | elseif (keyword(ie:ie).ne.' ' .and. ie.lt.kwsize) then |
---|
115 | ie=ie+1 |
---|
116 | goto 3 |
---|
117 | endif |
---|
118 | goto 1 |
---|
119 | 4 kwlen=ie-is |
---|
120 | if (is.gt.1) keyword(1:kwlen)=keyword(is:is+kwlen-1) |
---|
121 | ! |
---|
122 | ! Read input parameters according to the keyword: |
---|
123 | ! ==== ===== ========== ========= == === ======== |
---|
124 | ! |
---|
125 | ! Title |
---|
126 | ! |
---|
127 | if (keyword(1:kwlen).eq.'title') then |
---|
128 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
129 | read(input,'(A)',err=95) title |
---|
130 | lstr=lenstr(title) |
---|
131 | MPI_master_only write(stdout,'(/1x,A)') title(1:lstr) |
---|
132 | ! |
---|
133 | ! Time-stepping parameters |
---|
134 | ! |
---|
135 | elseif (keyword(1:kwlen).eq.'time_stepping') then |
---|
136 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
137 | read(input,*,err=95) ntimes,dt,ndtfast, ninfo |
---|
138 | MPI_master_only write(stdout, |
---|
139 | & '(I10,2x,A,1x,A /F10.2,2x,A,2(/I10,2x,A,1x,A)/F10.4,2x,A)' |
---|
140 | & ) ntimes, 'ntimes Total number of timesteps for', |
---|
141 | & '3D equations.', |
---|
142 | & dt, 'dt Timestep [sec] for 3D equations', |
---|
143 | & ndtfast, 'ndtfast Number of 2D timesteps within each', |
---|
144 | & '3D step.', |
---|
145 | & ninfo, 'ninfo Number of timesteps between', |
---|
146 | & 'runtime diagnostics.' |
---|
147 | dtfast=dt/float(ndtfast) ! set barotropic time step. |
---|
148 | |
---|
149 | if (NWEIGHT.lt.(2*ndtfast-1)) then |
---|
150 | write(stdout,'(a,i3)') |
---|
151 | & ' Error - Number of 2D timesteps (2*ndtfast-1): ', |
---|
152 | & 2*ndtfast-1 |
---|
153 | write(stdout,'(a,i3)') |
---|
154 | & ' exceeds barotopic weight dimension: ',NWEIGHT |
---|
155 | goto 95 |
---|
156 | endif |
---|
157 | |
---|
158 | #ifdef SOLVE3D |
---|
159 | ! |
---|
160 | ! Vertical S-coordinates parameters. |
---|
161 | ! |
---|
162 | elseif (keyword(1:kwlen).eq.'S-coord') then |
---|
163 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
164 | read(input,*,err=95) theta_s, theta_b, Tcline |
---|
165 | MPI_master_only write(stdout, |
---|
166 | & '(3(1pe10.3,2x,A,1x,A/),32x,A)') |
---|
167 | & theta_s, 'theta_s S-coordinate surface control', |
---|
168 | & 'parameter.', |
---|
169 | & theta_b, 'theta_b S-coordinate bottom control', |
---|
170 | & 'parameter.', |
---|
171 | & Tcline, 'Tcline S-coordinate surface/bottom layer', |
---|
172 | & 'width used in', 'vertical coordinate stretching, meters.' |
---|
173 | #endif |
---|
174 | ! |
---|
175 | ! Initial conditions file name. Check its availability (in the case |
---|
176 | ! of analytical initial conditions and nrrec=0 initial conditions are |
---|
177 | ! created internally and no file is needed). |
---|
178 | ! |
---|
179 | elseif (keyword(1:kwlen).eq.'initial') then |
---|
180 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
181 | read(input,*,err=95) nrrec |
---|
182 | #ifdef ANA_INITIAL |
---|
183 | if (nrrec.gt.0) then |
---|
184 | #endif |
---|
185 | read(input,'(A)',err=95) fname |
---|
186 | lstr=lenstr(fname) |
---|
187 | #if defined MPI && defined PARALLEL_FILES |
---|
188 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
189 | #endif |
---|
190 | open (testunit, file=fname(1:lstr), status='old', err=97) |
---|
191 | close(testunit) |
---|
192 | ininame=fname(1:lstr) |
---|
193 | MPI_master_only write(stdout,'(1x,A,2x,A,4x,A,I3)') |
---|
194 | & 'Initial State File:', ininame(1:lstr), 'Record:',nrrec |
---|
195 | #ifdef ANA_INITIAL |
---|
196 | endif |
---|
197 | #endif |
---|
198 | #ifndef ANA_GRID |
---|
199 | ! |
---|
200 | ! Grid file name. Check its availability. |
---|
201 | ! |
---|
202 | elseif (keyword(1:kwlen).eq.'grid') then |
---|
203 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
204 | read(input,'(A)',err=95) fname |
---|
205 | lstr=lenstr(fname) |
---|
206 | # if defined MPI && defined PARALLEL_FILES |
---|
207 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
208 | # endif |
---|
209 | open(testunit,file=fname(1:lstr), status='old', err=97) |
---|
210 | close(testunit) |
---|
211 | grdname=fname(1:lstr) |
---|
212 | MPI_master_only write(stdout,'(10x,A,2x,A)') |
---|
213 | & 'Grid File:', grdname(1:lstr) |
---|
214 | #endif |
---|
215 | #if !defined ANA_SMFLUX || defined SOLVE3D && (\ |
---|
216 | !defined ANA_STFLUX || !defined ANA_BTFLUX \ |
---|
217 | ||(defined BBL && !defined ANA_BSEDIM && !defined SEDIMENT) \ |
---|
218 | ||(defined BBL && !defined ANA_WWAVE) \ |
---|
219 | ||(defined QCORRECTION && !defined ANA_SST) \ |
---|
220 | ||(defined SALINITY && !defined ANA_SSFLUX) \ |
---|
221 | ||(defined LMD_SKPP && !defined ANA_SRFLUX) \ |
---|
222 | ||(defined SALINITY && defined QCORRECTION && \ |
---|
223 | defined SFLX_CORR && !defined ANA_SSS) \ |
---|
224 | || defined BULK_FLUX) |
---|
225 | ! |
---|
226 | ! Forcing file name. Check its availability. |
---|
227 | ! |
---|
228 | elseif (keyword(1:kwlen).eq.'forcing') then |
---|
229 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
230 | read(input,'(A)',err=95) fname |
---|
231 | lstr=lenstr(fname) |
---|
232 | # if defined MPI && defined PARALLEL_FILES |
---|
233 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
234 | # endif |
---|
235 | open (testunit, file=fname(1:lstr), status='old', err=97) |
---|
236 | close(testunit) |
---|
237 | frcname=fname(1:lstr) |
---|
238 | MPI_master_only write(stdout,'(2x,A,2x,A)') |
---|
239 | & 'Forcing Data File:', frcname(1:lstr) |
---|
240 | ! |
---|
241 | ! Biology forcing: iron dust deposition. |
---|
242 | ! |
---|
243 | # ifdef BIOLOGY |
---|
244 | elseif (keyword(1:kwlen).eq.'biology') then |
---|
245 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
246 | read(input,'(A)',err=95) fname |
---|
247 | lstr=lenstr(fname) |
---|
248 | # if defined MPI && defined PARALLEL_FILES |
---|
249 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
250 | # endif |
---|
251 | c open (testunit, file=fname(1:lstr), status='old', err=97) |
---|
252 | c close(testunit) |
---|
253 | bioname=fname(1:lstr) |
---|
254 | MPI_master_only write(stdout,'(2x,A,2x,A)') |
---|
255 | & 'Biology Forcing Data File:', bioname(1:lstr) |
---|
256 | # endif |
---|
257 | ! |
---|
258 | ! Bulk file name. Check its availability. |
---|
259 | ! |
---|
260 | # if defined BULK_FLUX |
---|
261 | elseif (keyword(1:kwlen).eq.'bulk_forcing') then |
---|
262 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
263 | read(input,'(A)',err=95) fname |
---|
264 | lstr=lenstr(fname) |
---|
265 | # if defined MPI && defined PARALLEL_FILES |
---|
266 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
267 | # endif |
---|
268 | open (testunit, file=fname(1:lstr), status='old', err=97) |
---|
269 | close(testunit) |
---|
270 | bulkname=fname(1:lstr) |
---|
271 | MPI_master_only write(stdout,'(2x,A,2x,A)') |
---|
272 | & ' Bulk Data File:', bulkname(1:lstr) |
---|
273 | # endif |
---|
274 | #endif |
---|
275 | #if (defined TCLIMATOLOGY && !defined ANA_TCLIMA) || \ |
---|
276 | (defined ZCLIMATOLOGY && !defined ANA_SSH) || \ |
---|
277 | (defined M2CLIMATOLOGY && !defined ANA_M2CLIMA) || \ |
---|
278 | (defined M3CLIMATOLOGY && !defined ANA_M3CLIMA) |
---|
279 | ! |
---|
280 | ! Climatology file name. Check availability. |
---|
281 | ! |
---|
282 | elseif (keyword(1:kwlen).eq.'climatology') then |
---|
283 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
284 | # ifdef AGRIF |
---|
285 | if (Agrif_Root()) then |
---|
286 | # endif |
---|
287 | read(input,'(A)',err=95) fname |
---|
288 | lstr=lenstr(fname) |
---|
289 | # if defined MPI && defined PARALLEL_FILES |
---|
290 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
291 | # endif |
---|
292 | open (testunit, file=fname(1:lstr), status='old', err=97) |
---|
293 | close(testunit) |
---|
294 | clmname=fname(1:lstr) |
---|
295 | MPI_master_only write(stdout,'(3x,A,2x,A)') |
---|
296 | & 'Climatology File:', clmname(1:lstr) |
---|
297 | # ifdef AGRIF |
---|
298 | endif |
---|
299 | # endif |
---|
300 | #endif |
---|
301 | #if defined T_FRC_BRY || defined M2_FRC_BRY || \ |
---|
302 | defined M3_FRC_BRY || defined Z_FRC_BRY |
---|
303 | # ifndef ANA_BRY |
---|
304 | ! |
---|
305 | ! Boundary file name. Check availability. |
---|
306 | ! |
---|
307 | elseif (keyword(1:kwlen).eq.'boundary') then |
---|
308 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
309 | # ifdef AGRIF |
---|
310 | if (Agrif_Root()) then |
---|
311 | # endif |
---|
312 | read(input,'(A)',err=95) fname |
---|
313 | lstr=lenstr(fname) |
---|
314 | # if defined MPI && defined PARALLEL_FILES |
---|
315 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
316 | # endif |
---|
317 | open (testunit, file=fname(1:lstr), status='old', err=97) |
---|
318 | close(testunit) |
---|
319 | bry_file=fname(1:lstr) |
---|
320 | MPI_master_only write(stdout,'(6x,A,2x,A)') |
---|
321 | & 'Boundary File:', bry_file(1:lstr) |
---|
322 | # ifdef AGRIF |
---|
323 | endif |
---|
324 | # endif |
---|
325 | # endif |
---|
326 | #endif |
---|
327 | ! |
---|
328 | ! Restart file name. |
---|
329 | ! |
---|
330 | elseif (keyword(1:kwlen).eq.'restart') then |
---|
331 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
332 | read(input,*,err=95) nrst, nrpfrst |
---|
333 | read(input,'(A)',err=95) fname |
---|
334 | lstr=lenstr(fname) |
---|
335 | # if defined MPI && defined PARALLEL_FILES |
---|
336 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
337 | # endif |
---|
338 | rstname=fname(1:lstr) |
---|
339 | MPI_master_only write(stdout, |
---|
340 | & '(7x,A,2x,A,4x,A,I6,4x,A,I4)') |
---|
341 | & 'Restart File:', rstname(1:lstr), |
---|
342 | & 'nrst =', nrst, 'rec/file: ', nrpfrst |
---|
343 | ! |
---|
344 | ! History file name. |
---|
345 | ! |
---|
346 | elseif (keyword(1:kwlen).eq.'history') then |
---|
347 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
348 | read(input,*,err=95) ldefhis, nwrt, nrpfhis |
---|
349 | read(input,'(A)',err=95) fname |
---|
350 | lstr=lenstr(fname) |
---|
351 | # if defined MPI && defined PARALLEL_FILES |
---|
352 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
353 | # endif |
---|
354 | hisname=fname(1:lstr) |
---|
355 | MPI_master_only write(stdout, |
---|
356 | & '(7x,A,2x,A,2x,A,1x,L1,2x,A,I4,2x,A,I3)') |
---|
357 | & 'History File:', hisname(1:lstr), 'Create new:', |
---|
358 | & ldefhis, 'nwrt =', nwrt, 'rec/file =', nrpfhis |
---|
359 | #ifdef AVERAGES |
---|
360 | ! |
---|
361 | ! Averages file name. |
---|
362 | ! |
---|
363 | elseif (keyword(1:kwlen).eq.'averages') then |
---|
364 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
365 | read(input,*,err=95) ntsavg, navg, nrpfavg |
---|
366 | read(input,'(A)',err=95) fname |
---|
367 | lstr=lenstr(fname) |
---|
368 | # if defined MPI && defined PARALLEL_FILES |
---|
369 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
370 | # endif |
---|
371 | avgname=fname(1:lstr) |
---|
372 | MPI_master_only write(stdout, |
---|
373 | & '(2(I10,2x,A,1x,A/32x,A/),6x,A,2x,A,1x,A,I3)') |
---|
374 | & ntsavg, 'ntsavg Starting timestep for the', |
---|
375 | & 'accumulation of output', 'time-averaged data.', |
---|
376 | & navg, 'navg Number of timesteps between', |
---|
377 | & 'writing of time-averaged','data into averages file.', |
---|
378 | & 'Averages File:', avgname(1:lstr), |
---|
379 | & 'rec/file =', nrpfavg |
---|
380 | #endif |
---|
381 | |
---|
382 | #if defined DIAGNOSTICS_TS |
---|
383 | ! |
---|
384 | ! Diagnostics file name. |
---|
385 | ! |
---|
386 | elseif (keyword(1:kwlen).eq.'diagnostics') then |
---|
387 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
388 | read(input,*,err=95) ldefdia, nwrtdia, nrpfdia |
---|
389 | read(input,'(A)',err=95) fname |
---|
390 | lstr=lenstr(fname) |
---|
391 | # if defined MPI && defined PARALLEL_FILES |
---|
392 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
393 | # endif |
---|
394 | dianame=fname(1:lstr) |
---|
395 | MPI_master_only write(stdout, |
---|
396 | & '(9x,A,2x,A/,32x,A,1x,L1,2x,A,I4,2x,A,I3)') |
---|
397 | & 'Tracer Diag File:',dianame(1:lstr),'Create new:', |
---|
398 | & ldefdia,'nwrt =',nwrtdia,'rec/file =',nrpfdia |
---|
399 | ! |
---|
400 | # ifdef AVERAGES |
---|
401 | elseif (keyword(1:kwlen).eq.'diag_avg') then |
---|
402 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
403 | read(input,*,err=95) ldefdia_avg, ntsdia_avg, nwrtdia_avg, |
---|
404 | & nrpfdia_avg |
---|
405 | read(input,'(A)',err=95) fname |
---|
406 | lstr=lenstr(fname) |
---|
407 | # if defined MPI && defined PARALLEL_FILES |
---|
408 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
409 | # endif |
---|
410 | dianame_avg=fname(1:lstr) |
---|
411 | MPI_master_only write(stdout, |
---|
412 | & '(5x,A,2x,A/,32x,A,1x,L1,2x,A,I4,2x,A,I3,/32x,A,I10)') |
---|
413 | & 'Tracer AVG Diag File:',dianame_avg(1:lstr),'Create new:', |
---|
414 | & ldefdia_avg,'nwrt =',nwrtdia_avg,'rec/file =',nrpfdia_avg, |
---|
415 | & 'Starting timestep = ',ntsdia_avg |
---|
416 | # endif |
---|
417 | #endif /* DIAGOSTICS_TS */ |
---|
418 | #if defined DIAGNOSTICS_UV |
---|
419 | ! |
---|
420 | ! Diagnostics Momentum file name. |
---|
421 | ! |
---|
422 | elseif (keyword(1:kwlen).eq.'diagnosticsM') then |
---|
423 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
424 | read(input,*,err=95) ldefdiaM, nwrtdiaM, nrpfdiaM |
---|
425 | read(input,'(A)',err=95) fname |
---|
426 | lstr=lenstr(fname) |
---|
427 | # if defined MPI && defined PARALLEL_FILES |
---|
428 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
429 | # endif |
---|
430 | dianameM=fname(1:lstr) |
---|
431 | MPI_master_only write(stdout, |
---|
432 | & '(7x,A,2x,A/,32x,A,1x,L1,2x,A,I4,2x,A,I3)') |
---|
433 | & 'Momemtum Diag File:', dianameM(1:lstr), 'Create new:', |
---|
434 | & ldefdiaM, 'nwrt =', nwrtdiaM, 'rec/file =', nrpfdiaM |
---|
435 | ! |
---|
436 | # ifdef AVERAGES |
---|
437 | elseif (keyword(1:kwlen).eq.'diagM_avg') then |
---|
438 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
439 | read(input,*,err=95) ldefdiaM_avg, ntsdiaM_avg, nwrtdiaM_avg, |
---|
440 | & nrpfdiaM_avg |
---|
441 | read(input,'(A)',err=95) fname |
---|
442 | lstr=lenstr(fname) |
---|
443 | # if defined MPI && defined PARALLEL_FILES |
---|
444 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
445 | # endif |
---|
446 | dianameM_avg=fname(1:lstr) |
---|
447 | MPI_master_only write(stdout, |
---|
448 | & '(3x,A,2x,A/,32x,A,1x,L1,2x,A,I4,2x,A,I3,/32x,A,I10)') |
---|
449 | & 'Momemtum AVG Diag File:',dianameM_avg(1:lstr),'Create new:', |
---|
450 | & ldefdiaM_avg,'nwrt =',nwrtdiaM_avg,'rec/file =',nrpfdiaM_avg, |
---|
451 | & 'Starting timestep = ',ntsdiaM_avg |
---|
452 | |
---|
453 | # endif |
---|
454 | #endif /*DIAGNOSTICS_UV */ |
---|
455 | #ifdef DIAGNOSTICS_BIO |
---|
456 | ! |
---|
457 | ! Diagnostics Biology file name. |
---|
458 | ! |
---|
459 | elseif (keyword(1:kwlen).eq.'diagnostics_bio') then |
---|
460 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
461 | read(input,*,err=95) ldefdiabio, nwrtdiabio, nrpfdiabio |
---|
462 | read(input,'(A)',err=95) fname |
---|
463 | lstr=lenstr(fname) |
---|
464 | # if defined MPI && defined PARALLEL_FILES |
---|
465 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
466 | # endif |
---|
467 | dianamebio=fname(1:lstr) |
---|
468 | MPI_master_only write(stdout, |
---|
469 | & '(8x,A,2x,A/,32x,A,1x,L1,2x,A,I4,2x,A,I3)') |
---|
470 | & 'Biology Diag File:', dianamebio(1:lstr), 'Create new:', |
---|
471 | & ldefdiabio, 'nwrt =', nwrtdiabio, 'rec/file =', nrpfdiabio |
---|
472 | ! |
---|
473 | # ifdef AVERAGES |
---|
474 | elseif (keyword(1:kwlen).eq.'diagbio_avg') then |
---|
475 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
476 | read(input,*,err=95) ldefdiabio_avg, ntsdiabio_avg, |
---|
477 | & nwrtdiabio_avg, nrpfdiabio_avg |
---|
478 | read(input,'(A)',err=95) fname |
---|
479 | lstr=lenstr(fname) |
---|
480 | # if defined MPI && defined PARALLEL_FILES |
---|
481 | call insert_node (fname, lstr, mynode, NNODES, ierr) |
---|
482 | # endif |
---|
483 | dianamebio_avg=fname(1:lstr) |
---|
484 | MPI_master_only write(stdout, |
---|
485 | & '(4x,A,2x,A/,32x,A,1x,L1,2x,A,I4,2x,A,I3,/32x,A,I10)') |
---|
486 | & 'Biology AVG Diag File:',dianamebio_avg(1:lstr),'Create new:', |
---|
487 | & ldefdiabio_avg,'nwrt =',nwrtdiabio_avg,'rec/file =', |
---|
488 | & nrpfdiabio_avg,'Starting timestep = ',ntsdiabio_avg |
---|
489 | # endif |
---|
490 | #endif /* DIAGOSTICS_BIO */ |
---|
491 | |
---|
492 | #ifdef FLOATS |
---|
493 | ! |
---|
494 | ! Floats file name. |
---|
495 | ! |
---|
496 | elseif (keyword(1:kwlen).eq.'floats') then |
---|
497 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
498 | #ifdef AGRIF |
---|
499 | if (Agrif_Root()) then |
---|
500 | #endif |
---|
501 | read(input,*,err=95) ldefflt, nflt, nrpfflt |
---|
502 | read(input,'(A)',err=95) fposnam |
---|
503 | read(input,'(A)',err=95) fname |
---|
504 | lstr=lenstr(fname) |
---|
505 | # if defined MPI && defined PARALLEL_FILES |
---|
506 | call insert_node (fposnam, lstr, mynode, NNODES, ierr) |
---|
507 | # endif |
---|
508 | fltname=fname(1:lstr) |
---|
509 | MPI_master_only write(stdout, |
---|
510 | & '(9x,A,2x,A,2x,A,1x,L1,2x,A,I4,2x,A,I3)') |
---|
511 | & 'Float File:',fltname(1:lstr), 'Create new:', |
---|
512 | & ldefflt, 'nflt =', nflt, 'rec/file =', nrpfflt |
---|
513 | #ifdef AGRIF |
---|
514 | else |
---|
515 | ldefflt=Agrif_Parent("ldefflt",ldefflt) |
---|
516 | nflt=Agrif_Parent("nflt",nflt) |
---|
517 | nrpfflt=Agrif_Parent("nrpfflt",nrpfflt) |
---|
518 | endif |
---|
519 | #endif |
---|
520 | #endif /* FLOATS */ |
---|
521 | |
---|
522 | #ifdef STATIONS |
---|
523 | ! |
---|
524 | ! Stations file name. |
---|
525 | ! |
---|
526 | elseif (keyword(1:kwlen).eq.'stations') then |
---|
527 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
528 | # ifdef AGRIF |
---|
529 | if (Agrif_Root()) then |
---|
530 | # endif |
---|
531 | read(input,*,err=95) ldefsta, nsta, nrpfsta |
---|
532 | read(input,'(A)',err=95) staposname |
---|
533 | read(input,'(A)',err=95) fname |
---|
534 | lstr=lenstr(fname) |
---|
535 | # if defined MPI && defined PARALLEL_FILES |
---|
536 | call insert_node (staposname, lstr, mynode, NNODES, ierr) |
---|
537 | # endif |
---|
538 | staname=fname(1:lstr) |
---|
539 | MPI_master_only write(stdout, |
---|
540 | & '(9x,A,2x,A,2x,A,1x,L1,2x,A,I4,2x,A,I3)') |
---|
541 | & 'Station File:',staname(1:lstr), 'Create new:', |
---|
542 | & ldefsta, 'nsta =', nsta, 'rec/file =', nrpfsta |
---|
543 | # ifdef AGRIF |
---|
544 | else |
---|
545 | ldefsta=Agrif_Parent("ldefsta",ldefsta) |
---|
546 | nsta=Agrif_Parent("nsta",nsta) |
---|
547 | nrpfsta=Agrif_Parent("nrpfsta",nrpfsta) |
---|
548 | endif |
---|
549 | # endif |
---|
550 | #endif /* STATIONS */ |
---|
551 | |
---|
552 | #ifdef ASSIMILATION |
---|
553 | ! |
---|
554 | ! Assimilation input/output file names. |
---|
555 | ! |
---|
556 | elseif (keyword(1:kwlen).eq.'assimilation') then |
---|
557 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
558 | read(input,'(A)',err=95) aparnam |
---|
559 | read(input,'(A)',err=95) assname |
---|
560 | fname=aparnam |
---|
561 | lstr=lenstr(aparnam) |
---|
562 | open (testunit,file=aparnam(1:lstr),status='old',err=97) |
---|
563 | close(testunit) |
---|
564 | MPI_master_only write(stdout,'(1x,A,2x,A)') |
---|
565 | & 'Assimilation Parameters File:', aparnam(1:lstr) |
---|
566 | fname=assname |
---|
567 | lstr=lenstr(assname) |
---|
568 | open (testunit,file=assname(1:lstr),status='old',err=97) |
---|
569 | close(testunit) |
---|
570 | MPI_master_only write(stdout,'(12x,A,2x,A)') |
---|
571 | & 'Assimilation File:', assname(1:lstr) |
---|
572 | #endif |
---|
573 | ! |
---|
574 | ! Switches for fields to be saved into history file. |
---|
575 | ! |
---|
576 | elseif (keyword(1:kwlen).eq.'primary_history_fields') then |
---|
577 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
578 | read(input,*,err=95) wrthis(indxZ), wrthis(indxUb) |
---|
579 | & , wrthis(indxVb) |
---|
580 | #ifdef SOLVE3D |
---|
581 | & , wrthis(indxU), wrthis(indxV) |
---|
582 | & , (wrthis(itrc), itrc=indxT,indxT+NT-1) |
---|
583 | #endif |
---|
584 | if ( wrthis(indxZ) .or. wrthis(indxUb) .or. wrthis(indxVb) |
---|
585 | #ifdef SOLVE3D |
---|
586 | & .or. wrthis(indxU) .or. wrthis(indxV) |
---|
587 | #endif |
---|
588 | & ) wrthis(indxTime)=.true. |
---|
589 | |
---|
590 | MPI_master_only write(stdout,'(/1x,A,5(/6x,l1,2x,A,1x,A))') |
---|
591 | & 'Fields to be saved in history file: (T/F)' |
---|
592 | & , wrthis(indxZ), 'write zeta ', 'free-surface.' |
---|
593 | & , wrthis(indxUb), 'write UBAR ', '2D U-momentum component.' |
---|
594 | & , wrthis(indxVb), 'write VBAR ', '2D V-momentum component.' |
---|
595 | #ifdef SOLVE3D |
---|
596 | & , wrthis(indxU), 'write U ', '3D U-momentum component.' |
---|
597 | & , wrthis(indxV), 'write V ', '3D V-momentum component.' |
---|
598 | do itrc=1,NT |
---|
599 | if (wrthis(indxT+itrc-1)) wrthis(indxTime)=.true. |
---|
600 | MPI_master_only write(stdout, '(6x,L1,2x,A,I1,A,I2,A)') |
---|
601 | & wrthis(indxT+itrc-1), 'write T(', itrc, |
---|
602 | & ') Tracer of index', itrc,'.' |
---|
603 | enddo |
---|
604 | |
---|
605 | elseif (keyword(1:kwlen).eq.'auxiliary_history_fields') then |
---|
606 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
607 | read(input,*,err=95) |
---|
608 | & wrthis(indxR), wrthis(indxO) |
---|
609 | & , wrthis(indxW), wrthis(indxAkv), wrthis(indxAkt) |
---|
610 | # ifdef SALINITY |
---|
611 | & , wrthis(indxAks) |
---|
612 | # endif |
---|
613 | # ifdef LMD_SKPP |
---|
614 | & , wrthis(indxHbl) |
---|
615 | # endif |
---|
616 | # ifdef LMD_BKPP |
---|
617 | & , wrthis(indxHbbl) |
---|
618 | # endif |
---|
619 | # ifdef BIOLOGY |
---|
620 | & , wrthis(indxHel) |
---|
621 | # ifdef BIO_NPZD |
---|
622 | & , wrthis(indxChC) |
---|
623 | # ifdef OXYGEN |
---|
624 | & , wrthis(indxU10) |
---|
625 | & , wrthis(indxKvO2) |
---|
626 | & , wrthis(indxO2sat) |
---|
627 | # endif |
---|
628 | # elif defined BIO_N2P2Z2D2 && defined VAR_CHL_C |
---|
629 | & , wrthis(indxChC1) |
---|
630 | & , wrthis(indxChC2) |
---|
631 | # endif |
---|
632 | # endif |
---|
633 | & , wrthis(indxBostr) |
---|
634 | |
---|
635 | if ( wrthis(indxR) .or. wrthis(indxO) .or. wrthis(indxW) |
---|
636 | & .or. wrthis(indxAkv) .or. wrthis(indxAkt) |
---|
637 | # ifdef SALINITY |
---|
638 | & .or. wrthis(indxAks) |
---|
639 | # endif |
---|
640 | # ifdef LMD_SKPP |
---|
641 | & .or. wrthis(indxHbl) |
---|
642 | # endif |
---|
643 | # ifdef LMD_BKPP |
---|
644 | & .or. wrthis(indxHbbl) |
---|
645 | # endif |
---|
646 | # ifdef BIOLOGY |
---|
647 | & .or. wrthis(indxHel) |
---|
648 | # ifdef BIO_NPZD |
---|
649 | & .or. wrthis(indxChC) |
---|
650 | # ifdef OXYGEN |
---|
651 | & .or. wrthis(indxU10) |
---|
652 | & .or. wrthis(indxKvO2) |
---|
653 | & .or. wrthis(indxO2sat) |
---|
654 | # endif |
---|
655 | # elif defined BIO_N2P2Z2D2 && defined VAR_CHL_C |
---|
656 | & .or. wrthis(indxChC1) |
---|
657 | & .or. wrthis(indxChC2) |
---|
658 | # endif |
---|
659 | # endif |
---|
660 | & .or. wrthis(indxBostr) |
---|
661 | & ) wrthis(indxTime)=.true. |
---|
662 | |
---|
663 | # if (defined UV_VIS2 || defined TS_DIF2) && defined SMAGORINSKY |
---|
664 | wrthis(indxVisc)=.true. |
---|
665 | # endif |
---|
666 | |
---|
667 | MPI_master_only write(stdout,'(8(/6x,l1,2x,A,1x,A))') |
---|
668 | & wrthis(indxR), 'write RHO ', 'Density anomaly.' |
---|
669 | & , wrthis(indxO), 'write Omega', 'Omega vertical velocity.' |
---|
670 | & , wrthis(indxW), 'write W ', 'True vertical velocity.' |
---|
671 | & , wrthis(indxAkv), 'write Akv ', 'Vertical viscosity.' |
---|
672 | & , wrthis(indxAkt), 'write Akt ', |
---|
673 | & 'Vertical diffusivity for temperature.' |
---|
674 | # ifdef SALINITY |
---|
675 | & , wrthis(indxAks), 'write Aks ', |
---|
676 | & 'Vertical diffusivity for salinity.' |
---|
677 | # endif |
---|
678 | # ifdef LMD_SKPP |
---|
679 | & , wrthis(indxHbl), 'write Hbl ', |
---|
680 | & 'Depth of KPP-model boundary layer.' |
---|
681 | # endif |
---|
682 | # ifdef LMD_BKPP |
---|
683 | & , wrthis(indxHbbl), 'write Hbbl ', |
---|
684 | & 'Depth of bottom planetary boundary layer.' |
---|
685 | # endif |
---|
686 | # ifdef BIOLOGY |
---|
687 | & , wrthis(indxHel), 'write Hel ', |
---|
688 | & 'Depth of the euphotic layer' |
---|
689 | # ifdef BIO_NPZD |
---|
690 | & , wrthis(indxChC), 'write ChC ', |
---|
691 | & 'Chlorophyll to Carbon ratio' |
---|
692 | # ifdef OXYGEN |
---|
693 | & , wrthis(indxU10), 'write u10 ', |
---|
694 | & 'Wind speed at 10 m height' |
---|
695 | & , wrthis(indxKvO2), 'write Kv_O2 ', |
---|
696 | & 'Gas transfer coefficient for O2' |
---|
697 | & , wrthis(indxO2sat), 'write O2sat ', |
---|
698 | & 'Saturation concentration of O2' |
---|
699 | # endif |
---|
700 | # elif defined BIO_N2P2Z2D2 && defined VAR_CHL_C |
---|
701 | & , wrthis(indxChC1), 'write ChC1 ', |
---|
702 | & 'Chlorophyll to Carbon ratio for Phy1' |
---|
703 | & , wrthis(indxChC2), 'write ChC2 ', |
---|
704 | & 'Chlorophyll to Carbon ratio for Phy2' |
---|
705 | # endif |
---|
706 | # endif |
---|
707 | & , wrthis(indxBostr), 'write Bostr', 'Bottom Stress.' |
---|
708 | # if (defined UV_VIS2 || defined TS_DIF2) && defined SMAGORINSKY |
---|
709 | & , wrthis(indxVisc), 'write Visc3d', 'Horizontal viscosity.' |
---|
710 | # endif |
---|
711 | #endif /* SOLVE3D */ |
---|
712 | #ifdef AVERAGES |
---|
713 | ! |
---|
714 | ! Switches for fields to be saved into averages file. |
---|
715 | ! |
---|
716 | elseif (keyword(1:kwlen).eq.'primary_averages') then |
---|
717 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
718 | read(input,*,err=95) wrtavg(indxZ), wrtavg(indxUb) |
---|
719 | & , wrtavg(indxVb) |
---|
720 | # ifdef SOLVE3D |
---|
721 | & , wrtavg(indxU), wrtavg(indxV) |
---|
722 | & , (wrtavg(itrc), itrc=indxT,indxT+NT-1) |
---|
723 | # endif |
---|
724 | if ( wrtavg(indxZ) .or. wrtavg(indxUb) .or. wrtavg(indxVb) |
---|
725 | # ifdef SOLVE3D |
---|
726 | & .or. wrtavg(indxU) .or. wrtavg(indxV) |
---|
727 | # endif |
---|
728 | & ) wrtavg(indxTime)=.true. |
---|
729 | |
---|
730 | MPI_master_only write(stdout,'(/1x,A,5(/6x,l1,2x,A,1x,A))') |
---|
731 | & 'Fields to be saved in averages file: (T/F)' |
---|
732 | & , wrtavg(indxZ), 'write zeta ', 'free-surface.' |
---|
733 | & , wrtavg(indxUb), 'write UBAR ', '2D U-momentum component.' |
---|
734 | & , wrtavg(indxVb), 'write VBAR ', '2D V-momentum component.' |
---|
735 | # ifdef SOLVE3D |
---|
736 | & , wrtavg(indxU), 'write U ', '3D U-momentum component.' |
---|
737 | & , wrtavg(indxV), 'write V ', '3D V-momentum component.' |
---|
738 | do itrc=1,NT |
---|
739 | if (wrtavg(indxT+itrc-1)) wrtavg(indxTime)=.true. |
---|
740 | MPI_master_only write(stdout, |
---|
741 | & '(6x,L1,2x,A,I1,A,2x,A,I2,A)') |
---|
742 | & wrtavg(indxT+itrc-1), 'write T(', |
---|
743 | & itrc,')', 'Tracer of index', itrc,'.' |
---|
744 | enddo |
---|
745 | |
---|
746 | elseif (keyword(1:kwlen).eq.'auxiliary_averages') then |
---|
747 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
748 | read(input,*,err=95) wrtavg(indxR), wrtavg(indxO) |
---|
749 | & , wrtavg(indxW), wrtavg(indxAkv), wrtavg(indxAkt) |
---|
750 | # ifdef SALINITY |
---|
751 | & , wrtavg(indxAks) |
---|
752 | # endif |
---|
753 | # ifdef LMD_SKPP |
---|
754 | & , wrtavg(indxHbl) |
---|
755 | # endif |
---|
756 | # ifdef LMD_BKPP |
---|
757 | & , wrtavg(indxHbbl) |
---|
758 | # endif |
---|
759 | # ifdef BIOLOGY |
---|
760 | & , wrtavg(indxHel) |
---|
761 | # ifdef BIO_NPZD |
---|
762 | & , wrtavg(indxChC) |
---|
763 | # ifdef OXYGEN |
---|
764 | & , wrtavg(indxU10) |
---|
765 | & , wrtavg(indxKvO2) |
---|
766 | & , wrtavg(indxO2sat) |
---|
767 | # endif |
---|
768 | # elif defined BIO_N2P2Z2D2 && defined VAR_CHL_C |
---|
769 | & , wrtavg(indxChC1) |
---|
770 | & , wrtavg(indxChC2) |
---|
771 | # endif |
---|
772 | # endif |
---|
773 | & , wrtavg(indxBostr) |
---|
774 | |
---|
775 | if ( wrtavg(indxR) .or. wrtavg(indxO) .or. wrtavg(indxW) |
---|
776 | & .or. wrtavg(indxAkv) .or. wrtavg(indxAkt) |
---|
777 | # ifdef SALINITY |
---|
778 | & .or. wrtavg(indxAks) |
---|
779 | # endif |
---|
780 | # ifdef LMD_SKPP |
---|
781 | & .or. wrtavg(indxHbl) |
---|
782 | # endif |
---|
783 | # ifdef LMD_BKPP |
---|
784 | & .or. wrtavg(indxHbbl) |
---|
785 | # endif |
---|
786 | # ifdef BIOLOGY |
---|
787 | & .or. wrtavg(indxHel) |
---|
788 | # ifdef BIO_NPZD |
---|
789 | & .or. wrtavg(indxChC) |
---|
790 | # ifdef OXYGEN |
---|
791 | & .or. wrtavg(indxU10) |
---|
792 | & .or. wrtavg(indxKvO2) |
---|
793 | & .or. wrtavg(indxO2sat) |
---|
794 | # endif |
---|
795 | # elif defined BIO_N2P2Z2D2 && defined VAR_CHL_C |
---|
796 | & .or. wrtavg(indxChC1) |
---|
797 | & .or. wrtavg(indxChC2) |
---|
798 | # endif |
---|
799 | # endif |
---|
800 | & .or. wrtavg(indxBostr) |
---|
801 | & ) wrtavg(indxTime)=.true. |
---|
802 | |
---|
803 | # if (defined UV_VIS2 || defined TS_DIF2) && defined SMAGORINSKY |
---|
804 | wrtavg(indxVisc)=.true. |
---|
805 | # endif |
---|
806 | |
---|
807 | MPI_master_only write(stdout,'(8(/6x,l1,2x,A,1x,A))') |
---|
808 | & wrtavg(indxR), 'write RHO ', 'Density anomaly' |
---|
809 | & , wrtavg(indxO), 'write Omega', 'Omega vertical velocity.' |
---|
810 | & , wrtavg(indxW), 'write W ', 'True vertical velocity.' |
---|
811 | & , wrtavg(indxAkv), 'write Akv ', 'Vertical viscosity' |
---|
812 | & , wrtavg(indxAkt), 'write Akt ', |
---|
813 | & 'Vertical diffusivity for temperature.' |
---|
814 | # ifdef SALINITY |
---|
815 | & , wrtavg(indxAks), 'write Aks ', |
---|
816 | & 'Vertical diffusivity for salinity.' |
---|
817 | # endif |
---|
818 | # ifdef LMD_SKPP |
---|
819 | & , wrtavg(indxHbl), 'write Hbl ', |
---|
820 | & 'Depth of KPP-model boundary layer' |
---|
821 | # endif |
---|
822 | # ifdef LMD_BKPP |
---|
823 | & , wrtavg(indxHbbl), 'write Hbbl ', |
---|
824 | & 'Depth of the bottom planetary boundary layer' |
---|
825 | # endif |
---|
826 | # ifdef BIOLOGY |
---|
827 | & , wrtavg(indxHel),'write Hel ', |
---|
828 | & 'Depth of the euphotic layer' |
---|
829 | # ifdef BIO_NPZD |
---|
830 | & , wrtavg(indxChC),'write ChC ', |
---|
831 | & 'Chlorophyll to Carbon ratio' |
---|
832 | # ifdef OXYGEN |
---|
833 | & , wrtavg(indxU10),'write u10 ', |
---|
834 | & 'Wind speed at 10 m height' |
---|
835 | & , wrtavg(indxKvO2),'write Kv_O2 ', |
---|
836 | & 'Gas transfer coefficient for O2' |
---|
837 | & , wrtavg(indxO2sat),'write O2sat ', |
---|
838 | & 'Saturation concentration of O2' |
---|
839 | # endif |
---|
840 | # elif defined BIO_N2P2Z2D2 && defined VAR_CHL_C |
---|
841 | & , wrtavg(indxChC1),'write ChC1 ', |
---|
842 | & 'Chlorophyll to Carbon ratio for Phy1' |
---|
843 | & , wrtavg(indxChC2),'write ChC2 ', |
---|
844 | & 'Chlorophyll to Carbon ratio for Phy2' |
---|
845 | # endif |
---|
846 | # endif |
---|
847 | & , wrtavg(indxBostr),'write Bostr', 'Bottom Stress.' |
---|
848 | # if (defined UV_VIS2 || defined TS_DIF2) && defined SMAGORINSKY |
---|
849 | & , wrtavg(indxVisc),'write Visc3d', 'Horizontal viscosity' |
---|
850 | # endif |
---|
851 | # endif /* SOLVE3D */ |
---|
852 | #endif /* AVERAGES */ |
---|
853 | |
---|
854 | |
---|
855 | #ifdef FLOATS |
---|
856 | ! |
---|
857 | ! Switches for fields to be saved into floats output file. |
---|
858 | ! |
---|
859 | elseif (keyword(1:kwlen).eq.'float_fields') then |
---|
860 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
861 | #ifdef AGRIF |
---|
862 | if (Agrif_Root()) then |
---|
863 | #endif |
---|
864 | read(input,*,err=95) wrtflt(indxfltGrd), |
---|
865 | & wrtflt(indxfltTemp), wrtflt(indxfltSalt), |
---|
866 | & wrtflt(indxfltRho), wrtflt(indxfltVel) |
---|
867 | MPI_master_only write(stdout,'(/1x,A,5(/6x,l1,2x,A))') |
---|
868 | & 'Fields to be saved in floats output file: (T/F)' |
---|
869 | & , wrtflt(indxfltGrd), 'write Grid location variables' |
---|
870 | & , wrtflt(indxfltTemp), 'write temperature.' |
---|
871 | & , wrtflt(indxfltSalt), 'write salinity.' |
---|
872 | & , wrtflt(indxfltRho), 'write density.' |
---|
873 | & , wrtflt(indxfltVel), 'write mean float velocity' |
---|
874 | #ifdef AGRIF |
---|
875 | endif |
---|
876 | #endif |
---|
877 | #endif /* FLOATS */ |
---|
878 | |
---|
879 | #if defined DIAGNOSTICS_TS |
---|
880 | ! |
---|
881 | ! Switches for fields to be saved into tracer diagnostics file. |
---|
882 | ! |
---|
883 | wrtdia=.true. |
---|
884 | # ifdef AVERAGES |
---|
885 | wrtdia_avg=.true. |
---|
886 | # endif |
---|
887 | #endif /*DIAGNOSTICS_TS */ |
---|
888 | #if defined DIAGNOSTICS_UV |
---|
889 | ! |
---|
890 | ! Switches for fields to be saved into momentum diagnostics file. |
---|
891 | ! |
---|
892 | wrtdiaM=.true. |
---|
893 | # ifdef AVERAGES |
---|
894 | wrtdiaM_avg=.true. |
---|
895 | # endif |
---|
896 | #endif /*DIAGNOSTICS_UV */ |
---|
897 | #ifdef DIAGNOSTICS_BIO |
---|
898 | ! |
---|
899 | ! Switches for fields to be saved into biology diagnostics file. |
---|
900 | ! |
---|
901 | wrtdiabio=.true. |
---|
902 | # ifdef AVERAGES |
---|
903 | wrtdiabio_avg=.true. |
---|
904 | # endif /* AVERAGES */ |
---|
905 | #endif /* DIAGNOSTICS_BIO */ |
---|
906 | |
---|
907 | #ifdef STATIONS |
---|
908 | ! |
---|
909 | ! Switches for fields to be saved into stations output file. |
---|
910 | ! |
---|
911 | elseif (keyword(1:kwlen).eq.'station_fields') then |
---|
912 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
913 | # ifdef AGRIF |
---|
914 | if (Agrif_Root()) then |
---|
915 | # endif |
---|
916 | read(input,*,err=95) wrtsta(indxstaGrd), |
---|
917 | & wrtsta(indxstaTemp), wrtsta(indxstaSalt), |
---|
918 | & wrtsta(indxstaRho), wrtsta(indxstaVel) |
---|
919 | MPI_master_only write(stdout,'(/1x,A,5(/6x,l1,2x,A))') |
---|
920 | & 'Fields to be saved in stations output (T/F)' |
---|
921 | & , wrtsta(indxstaGrd), 'write Grid location variables' |
---|
922 | & , wrtsta(indxstaTemp), 'write temperature.' |
---|
923 | & , wrtsta(indxstaSalt), 'write salinity.' |
---|
924 | & , wrtsta(indxstaRho), 'write density.' |
---|
925 | & , wrtsta(indxstaVel), 'write mean station velocity' |
---|
926 | # ifdef AGRIF |
---|
927 | endif |
---|
928 | # endif |
---|
929 | #endif /* STATIONS */ |
---|
930 | |
---|
931 | ! |
---|
932 | ! Boussinesq Approximation mean density. |
---|
933 | ! |
---|
934 | elseif (keyword(1:kwlen).eq.'rho0') then |
---|
935 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
936 | read(input,*,err=95) rho0 |
---|
937 | MPI_master_only write(stdout,'(F10.4,2x,A,1x,A)') |
---|
938 | & rho0, 'rho0 Boussinesq approximation', |
---|
939 | & 'mean density, kg/m3.' |
---|
940 | #if defined UV_VIS2 || defined UV_VIS4 |
---|
941 | ! |
---|
942 | ! Horizontal viscosity coefficients. |
---|
943 | ! |
---|
944 | elseif (keyword(1:kwlen).eq.'lateral_visc') then |
---|
945 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
946 | read(input,*,err=95) visc2, visc4 |
---|
947 | #endif |
---|
948 | #ifdef UV_VIS2 |
---|
949 | MPI_master_only write(stdout,9) visc2 |
---|
950 | 9 format(1pe10.3,2x,'visc2 Horizontal Laplacian ', |
---|
951 | & 'mixing coefficient [m2/s]',/,32x,'for momentum.') |
---|
952 | #endif |
---|
953 | #ifdef UV_VIS4 |
---|
954 | MPI_master_only write(stdout,10) visc4 |
---|
955 | 10 format(1pe10.3,2x,'visc4 Horizontal biharmonic ', |
---|
956 | & 'mixing coefficient [m4/s]',/,32x,'for momentum.') |
---|
957 | #endif |
---|
958 | ! |
---|
959 | ! Bottom drag coefficients. |
---|
960 | ! |
---|
961 | elseif (keyword(1:kwlen).eq.'bottom_drag') then |
---|
962 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
963 | read(input,*,err=95) rdrg, rdrg2, Zob, Cdb_min, Cdb_max |
---|
964 | MPI_master_only write(stdout,'(5(1pe10.3,2x,A/))') |
---|
965 | & rdrg, 'rdrg Linear bottom drag coefficient (m/si).', |
---|
966 | & rdrg2, 'rdrg2 Quadratic bottom drag coefficient.', |
---|
967 | & Zob, 'Zob Bottom roughness for logarithmic law (m).', |
---|
968 | & Cdb_min, 'Cdb_min Minimum bottom drag coefficient.', |
---|
969 | & Cdb_max, 'Cdb_max Maximum bottom drag coefficient.' |
---|
970 | ! |
---|
971 | ! Lateral boundary slipperness. |
---|
972 | ! |
---|
973 | elseif (keyword(1:kwlen).eq.'gamma2') then |
---|
974 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
975 | read(input,*,err=95) gamma2 |
---|
976 | MPI_master_only write(stdout,'(f10.2,2x,A,1x,A)') |
---|
977 | & gamma2, 'gamma2 Slipperiness parameter:', |
---|
978 | & 'free-slip +1, or no-slip -1.' |
---|
979 | #ifdef SOLVE3D |
---|
980 | # ifdef TS_DIF2 |
---|
981 | ! |
---|
982 | ! Horizontal Laplacian mixing coefficients for tracers. |
---|
983 | ! |
---|
984 | elseif (keyword(1:kwlen).eq.'tracer_diff2') then |
---|
985 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
986 | read(input,*,err=95) (tnu2(itrc),itrc=1,NT) |
---|
987 | do itrc=1,NT |
---|
988 | MPI_master_only write(stdout,7) tnu2(itrc), itrc, itrc |
---|
989 | 7 format(1pe10.3,' tnu2(',i1,') Horizontal Laplacian ' |
---|
990 | & ,'mixing coefficient (m2/s)',/,32x,'for tracer ',i1,'.') |
---|
991 | enddo |
---|
992 | # endif |
---|
993 | # ifdef TS_DIF4 |
---|
994 | ! |
---|
995 | ! Horizontal biharmonic mixing coefficients for tracer. |
---|
996 | ! |
---|
997 | elseif (keyword(1:kwlen).eq.'tracer_diff4') then |
---|
998 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
999 | read(input,*,err=95) (tnu4(itrc),itrc=1,NT) |
---|
1000 | do itrc=1,NT |
---|
1001 | MPI_master_only write(stdout,8) tnu4(itrc), itrc, itrc |
---|
1002 | 8 format(1pe10.3,' tnu4(',i1,') Horizontal biharmonic' |
---|
1003 | & ,' mixing coefficient [m4/s]',/,32x,'for tracer ',i1,'.') |
---|
1004 | enddo |
---|
1005 | |
---|
1006 | # endif |
---|
1007 | # if !defined LMD_MIXING && !defined BVF_MIXING\ |
---|
1008 | && !defined MY2_MIXING && !defined MY25_MIXING |
---|
1009 | ! |
---|
1010 | ! Background vertical viscosity and mixing coefficients for tracers. |
---|
1011 | ! |
---|
1012 | elseif (keyword(1:kwlen).eq.'vertical_mixing') then |
---|
1013 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1014 | read(input,*,err=95) Akv_bak,(Akt_bak(itrc),itrc=1,NT) |
---|
1015 | MPI_master_only write(stdout,'(1pe10.3,2x,A,1x,A)') |
---|
1016 | & Akv_bak, 'Akv_bak Background vertical viscosity', |
---|
1017 | & 'coefficient, m2/s.' |
---|
1018 | do itrc=1,NT |
---|
1019 | MPI_master_only write(stdout, |
---|
1020 | & '(1pe10.3,2x,A,I1,A,1x,A/32x,A,I3,A)') |
---|
1021 | & Akt_bak(itrc), 'Akt_bak(', itrc, ')', |
---|
1022 | & 'Background vertical mixing coefficient, m2/s,', |
---|
1023 | & 'for tracer ', itrc, '.' |
---|
1024 | enddo |
---|
1025 | # endif |
---|
1026 | # ifdef MY25_MIXING |
---|
1027 | ! |
---|
1028 | ! Mellor-Yamada Level 2.5 turbulent closure parameters. |
---|
1029 | ! |
---|
1030 | elseif (keyword(1:kwlen).eq.'MY_bak_mixing') then |
---|
1031 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1032 | read(input,*,err=95) Akq_bak, q2nu2, q2nu4 |
---|
1033 | MPI_master_only write(stdout,13) Akq_bak |
---|
1034 | 13 format(1pe10.3,2x,'Akq_bak Background vertical mixing', |
---|
1035 | & ' coefficient [m2/s]',/,32x,'for turbulent energy.') |
---|
1036 | # ifdef Q_DIF2 |
---|
1037 | MPI_master_only write(stdout,14) q2nu2 |
---|
1038 | 14 format(1pe10.3,2x,'q2nu2 Horizontal Laplacian ', |
---|
1039 | & 'mixing coefficient [m2/s]',/,32x,'for turbulent energy.') |
---|
1040 | # endif |
---|
1041 | # ifdef Q_DIF4 |
---|
1042 | MPI_master_only write(stdout,15) q2nu4 |
---|
1043 | 15 format(1pe10.3,2x,'q2nu4 Horizontal, biharmonic ', |
---|
1044 | & 'mixing coefficient, m2/s,',/,32x,'for turbulent energy.') |
---|
1045 | # endif |
---|
1046 | # endif |
---|
1047 | # ifdef BODYFORCE |
---|
1048 | elseif (keyword(1:kwlen).eq.'bodyforce') then |
---|
1049 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1050 | read(input,*,err=95) levsfrc,levbfrc |
---|
1051 | if (levsfrc.lt.1 .or. levsfrc.gt.N) then |
---|
1052 | MPI_master_only write(stdout,19) 'LEVSFRC = ',levsfrc |
---|
1053 | 19 format(' READ_INP - Illegal bodyforce level, ',A,i4) |
---|
1054 | ierr=ierr+1 |
---|
1055 | endif |
---|
1056 | if (levbfrc.lt.1 .or. levbfrc.gt.N) then |
---|
1057 | MPI_master_only write(stdout,19) 'LEVBFRC = ',levbfrc |
---|
1058 | ierr=ierr+1 |
---|
1059 | endif |
---|
1060 | MPI_master_only write(stdout,20) levsfrc, levbfrc |
---|
1061 | 20 format(4x,i6,2x,'levsfrc ', |
---|
1062 | & 'Deepest level to apply surface stress as a ', |
---|
1063 | & 'bodyforce.',/, |
---|
1064 | & 4x,i6,2x,'levbfrc ', |
---|
1065 | & 'Shallowest level to apply bottom stress as a ', |
---|
1066 | & 'bodyforce.') |
---|
1067 | # endif |
---|
1068 | #endif |
---|
1069 | #if defined SPONGE || \ |
---|
1070 | defined TNUDGING || defined M2NUDGING || \ |
---|
1071 | defined M3NUDGING || defined ZNUDGING |
---|
1072 | ! |
---|
1073 | ! Parameters for sponge layers |
---|
1074 | ! |
---|
1075 | elseif (keyword(1:kwlen).eq.'sponge') then |
---|
1076 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1077 | read(input,*,err=95) x_sponge, v_sponge |
---|
1078 | MPI_master_only write(stdout,'(1pe10.2,2x,A,1x,A)') |
---|
1079 | & x_sponge,'x_sponge Thickness of sponge', |
---|
1080 | & 'and/or nudging layer (m)' |
---|
1081 | MPI_master_only write(stdout,'(f10.2,2x,A)') |
---|
1082 | & v_sponge,'v_sponge Viscosity in sponge layer (m2/s)' |
---|
1083 | #endif |
---|
1084 | #if defined T_FRC_BRY || defined M2_FRC_BRY || \ |
---|
1085 | defined M3_FRC_BRY || defined Z_FRC_BRY || \ |
---|
1086 | defined TNUDGING || defined M2NUDGING || \ |
---|
1087 | defined M3NUDGING || defined ZNUDGING |
---|
1088 | ! |
---|
1089 | ! Nudging parameters for OBC and nudging layers |
---|
1090 | ! (converted from [days] to [sec^-1] |
---|
1091 | ! |
---|
1092 | elseif (keyword(1:kwlen).eq.'nudg_cof') then |
---|
1093 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1094 | # if defined AGRIF && !defined AGRIF_OBC_M2ORLANSKI && \ |
---|
1095 | !defined AGRIF_OBC_M3ORLANSKI && !defined AGRIF_OBC_TORLANSKI |
---|
1096 | if (Agrif_Root()) then |
---|
1097 | # endif |
---|
1098 | read(input,*,err=95) tauT_in,tauT_out,tauM_in,tauM_out |
---|
1099 | tauT_in =1./(tauT_in *86400.) |
---|
1100 | tauT_out=1./(tauT_out*86400.) |
---|
1101 | tauM_in =1./(tauM_in *86400.) |
---|
1102 | tauM_out=1./(tauM_out*86400.) |
---|
1103 | MPI_master_only write(stdout,'(1pe10.3,2x,A)') |
---|
1104 | & tauT_in,'tauT_in Nudging coefficients [sec^-1]' |
---|
1105 | MPI_master_only write(stdout,'(1pe10.3,2x,A)') |
---|
1106 | & tauT_out,'tauT_out Nudging coefficients [sec^-1]' |
---|
1107 | MPI_master_only write(stdout,'(1pe10.3,2x,A)') |
---|
1108 | & tauM_in,'tauM_in Nudging coefficients [sec^-1]' |
---|
1109 | MPI_master_only write(stdout,'(1pe10.3,2x,A/)') |
---|
1110 | & tauM_out,'tauM_out Nudging coefficients [sec^-1]' |
---|
1111 | # if defined AGRIF && !defined AGRIF_OBC_M2ORLANSKI && \ |
---|
1112 | !defined AGRIF_OBC_M3ORLANSKI && !defined AGRIF_OBC_TORLANSKI |
---|
1113 | endif |
---|
1114 | # endif |
---|
1115 | #endif |
---|
1116 | #ifdef SOLVE3D |
---|
1117 | # ifndef NONLIN_EOS |
---|
1118 | ! |
---|
1119 | ! Parameters for linear equations of state. |
---|
1120 | ! |
---|
1121 | elseif (keyword(1:kwlen).eq.'lin_EOS_cff') then |
---|
1122 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1123 | read(input,*,err=95) R0, T0, S0, Tcoef, Scoef |
---|
1124 | MPI_master_only write(stdout,'(5(f10.4,2x,A,1x,A/))') |
---|
1125 | & T0, 'T0 Background value for potential', |
---|
1126 | & 'temperature (Celsius).', |
---|
1127 | & S0, 'S0 Background salinity (PSU),', 'constant.', |
---|
1128 | & R0, 'R0 Background density (kg/m3) used in', |
---|
1129 | & 'linear EOS.', |
---|
1130 | & Tcoef, 'Tcoef Thermal expansion coefficient', |
---|
1131 | & '(kg/m3/Celsius).', |
---|
1132 | & Scoef, 'Scoef Saline contraction coefficient', |
---|
1133 | & '(kg/m3/PSU).' |
---|
1134 | # endif |
---|
1135 | # ifdef SEDIMENT |
---|
1136 | ! |
---|
1137 | ! Sediments input file name. |
---|
1138 | ! |
---|
1139 | elseif (keyword(1:kwlen).eq.'sediments') then |
---|
1140 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1141 | read(input,'(A)',err=95) sedname |
---|
1142 | lstr=lenstr(sedname) |
---|
1143 | |
---|
1144 | MPI_master_only write(stdout, |
---|
1145 | & '(/9x,A,2x,A)') |
---|
1146 | & 'Sediment input file:',sedname |
---|
1147 | |
---|
1148 | elseif (keyword(1:kwlen).eq.'sediment_history_fields') then |
---|
1149 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1150 | read(input,*,err=95) |
---|
1151 | & (wrthis(itrc), itrc=indxSed,indxSed+NST+1) |
---|
1152 | |
---|
1153 | MPI_master_only write(stdout,'(2(/6x,L1,2x,A,1x,A))') |
---|
1154 | & wrthis(indxBTHK), 'write bed_thick ', |
---|
1155 | & 'thickness of sediment bed layer.' |
---|
1156 | & , wrthis(indxBPOR), 'write bed_poros ', |
---|
1157 | & 'porosity of sediment bed layer.' |
---|
1158 | do itrc=1,NST |
---|
1159 | MPI_master_only write(stdout, '(6x,L1,2x,A,I1,A,I2,A)') |
---|
1160 | & wrthis(indxBFRA+itrc-1), 'write bed_frac(',itrc, |
---|
1161 | & ') Sediment fraction of index ',itrc,'.' |
---|
1162 | enddo |
---|
1163 | # endif /* SEDIMENT */ |
---|
1164 | #endif |
---|
1165 | #ifdef BBL |
---|
1166 | elseif (keyword(1:kwlen).eq.'bbl_history_fields') then |
---|
1167 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1168 | read(input,*,err=95) |
---|
1169 | & (wrthis(itrc), itrc=indxBBL,indxBBL+5) |
---|
1170 | |
---|
1171 | MPI_master_only write(stdout,'(6(/6x,L1,2x,A,1x,A))') |
---|
1172 | & wrthis(indxAbed), 'write Abed ', |
---|
1173 | & 'bed wave excursion amplitude.' |
---|
1174 | & , wrthis(indxHrip), 'write Hripple ', |
---|
1175 | & 'Bed ripple height.' |
---|
1176 | & , wrthis(indxLrip), 'write Lripple ', |
---|
1177 | & 'Bed ripple length.' |
---|
1178 | & , wrthis(indxZbnot), 'write Zbnot ', |
---|
1179 | & 'Physical bottom roughness.' |
---|
1180 | & , wrthis(indxZbapp), 'write Zbapp ', |
---|
1181 | & 'Apparent bottom roughness.' |
---|
1182 | & , wrthis(indxBostrw), 'write Bostrw ', |
---|
1183 | & 'Wave-induced bottom stress.' |
---|
1184 | |
---|
1185 | #endif /* BBL */ |
---|
1186 | #ifdef ANA_PSOURCE |
---|
1187 | ! |
---|
1188 | ! Set-up point Sources/Sink number (Nsrc), direction (Dsrc), I- and |
---|
1189 | ! J-grid locations (Isrc,Jsrc), and logical switch for type of tracer |
---|
1190 | ! to apply (Lsrc). Currently, the direction can be along XI-direction |
---|
1191 | ! (Dsrc = 0) or along ETA-direction (Dsrc > 0). The mass sources are |
---|
1192 | ! located at U- or V-points so the grid locations should range from |
---|
1193 | ! 1 =< Isrc =< L and 1 =< Jsrc =< M. |
---|
1194 | ! |
---|
1195 | elseif (keyword(1:kwlen).eq.'psource') then |
---|
1196 | call cancel_kwd (keyword(1:kwlen), ierr) |
---|
1197 | read(input,*,err=95) Nsrc |
---|
1198 | MPI_master_only write(stdout,'(/6x,i6,2x,A,1x,A)') |
---|
1199 | & Nsrc, 'Number of point sources' |
---|
1200 | do is=1,Nsrc |
---|
1201 | read(input,*,err=95) Isrc(is), Jsrc(is), Dsrc(is), Qbar(is), |
---|
1202 | & (Lsrc(is,itrc), itrc=1,NT), |
---|
1203 | & (Tsrc0(is,itrc), itrc=1,NT) |
---|
1204 | MPI_master_only write(stdout,'(3(/6x,i6,2x,A))') |
---|
1205 | & Isrc(is), 'I point source indice' |
---|
1206 | & , Jsrc(is), 'J point source indice' |
---|
1207 | & , Dsrc(is), 'Direction of point source flow' |
---|
1208 | MPI_master_only write(stdout,'(1pe10.3,2x,A)') |
---|
1209 | & Qbar(is), 'Total transport at point source' |
---|
1210 | do itrc=1,NT |
---|
1211 | MPI_master_only write(stdout, |
---|
1212 | & '(6x,L1,2x,A,I1,A,2x,A,I2,A)') |
---|
1213 | & Lsrc(is,itrc), 'write Lsrc(', |
---|
1214 | & itrc,')', 'Tracer of index', itrc,'.' |
---|
1215 | enddo |
---|
1216 | do itrc=1,NT |
---|
1217 | MPI_master_only write(stdout, |
---|
1218 | & '(6x,1pe10.3,2x,A,I1,A,2x,A,I2,A)') |
---|
1219 | & Tsrc0(is,itrc), 'write Tsrc(', |
---|
1220 | & itrc,')', 'Tracer of index', itrc,'.' |
---|
1221 | enddo |
---|
1222 | enddo |
---|
1223 | #endif /* ANA_SOURCES */ |
---|
1224 | |
---|
1225 | else |
---|
1226 | MPI_master_only write(stdout,'(/3(1x,A)/)') |
---|
1227 | & 'WARNING: Unrecognized keyword:', |
---|
1228 | & keyword(1:kwlen),' --> DISREGARDED.' |
---|
1229 | endif |
---|
1230 | if (keyword(1:kwlen) .eq. end_signal) goto 99 |
---|
1231 | goto 1 |
---|
1232 | ! |
---|
1233 | ! Error while reading input parameters. |
---|
1234 | ! |
---|
1235 | 95 write(stdout,'(/1x,4A/)') 'READ_INP ERROR while reading block', |
---|
1236 | & ' with keyword ''', keyword(1:kwlen), '''.' |
---|
1237 | ierr=ierr+1 |
---|
1238 | goto 99 |
---|
1239 | 97 write(stdout,'(/1x,4A/)') 'READ_INP ERROR: Cannot find input ', |
---|
1240 | & 'file ''', fname(1:lstr), '''.' |
---|
1241 | ierr=ierr+1 |
---|
1242 | 99 close (input) |
---|
1243 | ! |
---|
1244 | ! Check that all keywords were canceled |
---|
1245 | ! Complain if some of them are left |
---|
1246 | ! |
---|
1247 | if (ierr.eq.0) then |
---|
1248 | call check_kwds (ierr) |
---|
1249 | ! |
---|
1250 | ! Check CPP-switches for consistency. This operation is split into |
---|
1251 | ! two stages because the first subroutine, "check_switches1", is |
---|
1252 | ! generated by special program cppcheck (file cppcheck.F) by |
---|
1253 | ! examining and documention all available switches in cppdefs.h. |
---|
1254 | ! This subroutine creates log of all switches defined in "cppdefs.h", |
---|
1255 | ! as well as traps multiply defined global configurations (project |
---|
1256 | ! switches, such as REGIONAL, etc). |
---|
1257 | ! The second routine, "check_switches2" is hand written and it |
---|
1258 | ! contains traps for mutually exclussive definition of all other |
---|
1259 | ! CPP-switches (i.e. those which are NOT project selection switches, |
---|
1260 | ! for example, it traps multiply defined vertical mixing schemes or |
---|
1261 | ! lateral boundary conditions). |
---|
1262 | ! |
---|
1263 | ! Both codes are written in transparent mode: they assumed that error |
---|
1264 | ! variable (ierr) is initialized at entry and they add 1 for each |
---|
1265 | ! error discovered. |
---|
1266 | ! |
---|
1267 | call check_srcs |
---|
1268 | call check_switches1 (ierr) |
---|
1269 | call check_switches2 (ierr) |
---|
1270 | endif |
---|
1271 | if (ierr.ne.0) then |
---|
1272 | write(stdout,'(/1x,2A,I3,1x,A/)') 'READ_INP ERROR: ', |
---|
1273 | & 'A total of', ierr, 'configuration errors discovered.' |
---|
1274 | return |
---|
1275 | endif |
---|
1276 | #ifdef MPI |
---|
1277 | ! call MPI_Barrier (MPI_COMM_WORLD, ierr) |
---|
1278 | #endif |
---|
1279 | return |
---|
1280 | end |
---|
1281 | |
---|
1282 | c |
---|