source: palm/trunk/SOURCE/chemistry_model_mod.f90 @ 3173

Last change on this file since 3173 was 3173, checked in by suehring, 6 years ago

Further revision of 2D surface output for radiation and chemistry quantities; bugfix for commit 3170

  • Property svn:keywords set to Id
File size: 76.2 KB
Line 
1!> @file chemistry_model_mod.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2018 Leibniz Universitaet Hannover
18! Copyright 2017-2018 Karlsruhe Institute of Technology
19! Copyright 2017-2018 Freie Universitaet Berlin
20!------------------------------------------------------------------------------!
21!
22! Current revisions:
23! -----------------
24!
25!
26! Former revisions:
27! -----------------
28! $Id: chemistry_model_mod.f90 3173 2018-07-26 12:55:23Z suehring $
29! Revise output of surface quantities in case of overhanging structures
30!
31! 3045 2018-05-28 07:55:41Z Giersch
32! error messages revised
33!
34! 3014 2018-05-09 08:42:38Z maronga
35! Bugfix: nzb_do and nzt_do were not used for 3d data output
36!
37! 3004 2018-04-27 12:33:25Z Giersch
38! Comment concerning averaged data output added
39!
40! 2932 2018-03-26 09:39:22Z maronga
41! renamed chemistry_par to chemistry_parameters
42!
43! 2894 2018-03-15 09:17:58Z Giersch
44! Calculations of the index range of the subdomain on file which overlaps with
45! the current subdomain are already done in read_restart_data_mod,
46! chem_last_actions was renamed to chem_wrd_local, chem_read_restart_data was
47! renamed to chem_rrd_local, chem_write_var_list was renamed to
48! chem_wrd_global, chem_read_var_list was renamed to chem_rrd_global,
49! chem_skip_var_list has been removed, variable named found has been
50! introduced for checking if restart data was found, reading of restart strings
51! has been moved completely to read_restart_data_mod, chem_rrd_local is already
52! inside the overlap loop programmed in read_restart_data_mod, todo list has
53! bees extended, redundant characters in chem_wrd_local have been removed,
54! the marker *** end chemistry *** is not necessary anymore, strings and their
55! respective lengths are written out and read now in case of restart runs to
56! get rid of prescribed character lengths
57!
58! 2815 2018-02-19 11:29:57Z suehring
59! Bugfix in restart mechanism,
60! rename chem_tendency to chem_prognostic_equations,
61! implement vector-optimized version of chem_prognostic_equations,
62! some clean up (incl. todo list)
63!
64! 2773 2018-01-30 14:12:54Z suehring
65! Declare variables required for nesting as public
66!
67! 2772 2018-01-29 13:10:35Z suehring
68! Bugfix in string handling
69!
70! 2768 2018-01-24 15:38:29Z kanani
71! Shorten lines to maximum length of 132 characters
72!
73! 2766 2018-01-22 17:17:47Z kanani
74! Removed preprocessor directive __chem
75!
76! 2756 2018-01-16 18:11:14Z suehring
77! Fill values in 3D output introduced.
78!
79! 2718 2018-01-02 08:49:38Z maronga
80! Initial revision
81!
82!
83!
84!
85! Authors:
86! --------
87! @author Renate Forkel
88! @author Farah Kanani-Suehring
89! @author Klaus Ketelsen
90! @author Basit Khan
91!
92!
93!------------------------------------------------------------------------------!
94! Description:
95! ------------
96!> Chemistry model for PALM-4U
97!> @todo Adjust chem_rrd_local to CASE structure of others modules. It is not
98!>       allowed to use the chemistry model in a precursor run and additionally
99!>       not using it in a main run
100!> @todo Update/clean-up todo list! (FK)
101!> @todo Set proper fill values (/= 0) for chem output arrays! (FK)
102!> @todo Add routine chem_check_parameters, add checks for inconsistent or
103!>       unallowed parameter settings.
104!>       CALL of chem_check_parameters from check_parameters. (FK)
105!> @todo Make routine chem_header available, CALL from header.f90
106!>       (see e.g. how it is done in routine lsm_header in
107!>        land_surface_model_mod.f90). chem_header should include all setup
108!>        info about chemistry parameter settings. (FK)
109!> @todo Implement turbulent inflow of chem spcs in inflow_turbulence.
110!> @todo Separate boundary conditions for each chem spcs to be implemented
111!> @todo Currently only total concentration are calculated. Resolved, parameterized
112!>       and chemistry fluxes although partially and some completely coded but
113!>       are not operational/activated in this version. bK.
114!> @todo slight differences in passive scalar and chem spcs when chem reactions
115!>       turned off. Need to be fixed. bK
116!> @todo test nesting for chem spcs, was implemented by suehring (kanani)
117!> @todo subroutine set_const_initial_values to be taken out from chemistry_model_mod !bK.
118!> @todo chemistry error messages
119!> @todo Format this module according to PALM coding standard (see e.g. module
120!>       template under http://palm.muk.uni-hannover.de/mosaik/downloads/8 or
121!>       D3_coding_standard.pdf under https://palm.muk.uni-hannover.de/trac/downloads/16)
122!
123!------------------------------------------------------------------------------!
124
125MODULE chemistry_model_mod
126
127   USE kinds,              ONLY: wp, iwp
128   USE indices,            ONLY: nz, nzb,nzt,nysg,nyng,nxlg,nxrg,nys,nyn
129   USE pegrid,             ONLY: myid, threads_per_task
130   USE control_parameters, ONLY: dt_3d, ws_scheme_sca, initializing_actions, message_string, &
131                                 omega, tsc, intermediate_timestep_count,      &
132                                 intermediate_timestep_count_max,              &
133                                 timestep_scheme, use_prescribed_profile_data 
134   USE arrays_3d,          ONLY: hyp, pt, rdf_sc, tend, zu                     
135   USE chem_gasphase_mod,  ONLY: nspec, spc_names, nkppctrl, nmaxfixsteps,     &
136                                 t_steps, fill_temp, chem_gasphase_integrate,  &
137                                 nvar, atol, rtol, nphot, phot_names
138   USE cpulog,             ONLY: cpu_log, log_point
139
140   USE chem_modules
141 
142
143   IMPLICIT   NONE
144   PRIVATE
145   SAVE
146
147!- Define chemical variables
148
149   TYPE   species_def
150      CHARACTER(LEN=8)                                   :: name
151      CHARACTER(LEN=16)                                  :: unit
152      REAL(kind=wp),POINTER,DIMENSION(:,:,:)             :: conc
153      REAL(kind=wp),POINTER,DIMENSION(:,:,:)             :: conc_av
154      REAL(kind=wp),POINTER,DIMENSION(:,:,:)             :: conc_p
155      REAL(kind=wp),POINTER,DIMENSION(:,:,:)             :: tconc_m
156      REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:)           :: cssws_av           
157      REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:)           :: flux_s_cs
158      REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:)           :: diss_s_cs
159      REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:,:)         :: flux_l_cs
160      REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:,:)         :: diss_l_cs
161      REAL(kind=wp),ALLOCATABLE,DIMENSION(:)             :: conc_pr_init
162   END TYPE species_def
163
164
165   TYPE   photols_def                                                           
166      CHARACTER(LEN=8)                                   :: name
167      CHARACTER(LEN=16)                                  :: unit
168      REAL(kind=wp),POINTER,DIMENSION(:,:,:)             :: freq
169   END TYPE photols_def
170
171
172   PUBLIC  species_def
173   PUBLIC  photols_def
174
175
176   TYPE(species_def),ALLOCATABLE,DIMENSION(:),TARGET, PUBLIC     :: chem_species
177   TYPE(photols_def),ALLOCATABLE,DIMENSION(:),TARGET, PUBLIC     :: phot_frequen 
178
179   REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:,:,:),TARGET   :: spec_conc_1
180   REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:,:,:),TARGET   :: spec_conc_2
181   REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:,:,:),TARGET   :: spec_conc_3
182   REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:,:,:),TARGET   :: spec_conc_av       
183   REAL(kind=wp),ALLOCATABLE,DIMENSION(:,:,:,:),TARGET   :: freq_1
184
185   INTEGER,DIMENSION(nkppctrl)                           :: icntrl                            ! Fine tuning kpp
186   REAL(kind=wp),DIMENSION(nkppctrl)                     :: rcntrl                            ! Fine tuning kpp
187
188   CHARACTER(10), PUBLIC :: photolysis_scheme
189                                         ! 'constant',
190                                         ! 'simple' (Simple parameterisation from MCM, Saunders et al., 2003, Atmos. Chem. Phys., 3, 161-180
191                                         ! 'fastj'  (Wild et al., 2000, J. Atmos. Chem., 37, 245-282) STILL NOT IMPLEMENTED
192
193   PUBLIC nspec
194   PUBLIC nvar       
195   PUBLIC spc_names
196   PUBLIC spec_conc_2 
197
198!- Interface section
199  INTERFACE chem_boundary_conds
200      MODULE PROCEDURE chem_boundary_conds
201  END INTERFACE chem_boundary_conds
202
203   INTERFACE chem_init
204      MODULE PROCEDURE chem_init
205   END INTERFACE chem_init
206
207   INTERFACE chem_init_profiles
208      MODULE PROCEDURE chem_init_profiles
209   END INTERFACE chem_init_profiles
210
211   INTERFACE chem_parin
212      MODULE PROCEDURE chem_parin
213   END INTERFACE chem_parin
214
215   INTERFACE chem_integrate
216      MODULE PROCEDURE chem_integrate_ij
217   END INTERFACE chem_integrate
218
219   INTERFACE chem_swap_timelevel
220      MODULE PROCEDURE chem_swap_timelevel
221   END INTERFACE chem_swap_timelevel
222
223   INTERFACE chem_define_netcdf_grid
224      MODULE PROCEDURE chem_define_netcdf_grid
225   END INTERFACE chem_define_netcdf_grid
226
227   INTERFACE chem_data_output_3d
228      MODULE PROCEDURE chem_data_output_3d
229   END INTERFACE chem_data_output_3d
230
231   INTERFACE chem_check_data_output
232      MODULE PROCEDURE chem_check_data_output
233   END INTERFACE chem_check_data_output
234
235   INTERFACE chem_check_data_output_pr
236      MODULE PROCEDURE chem_check_data_output_pr
237   END INTERFACE chem_check_data_output_pr
238
239   INTERFACE chem_3d_data_averaging
240      MODULE PROCEDURE chem_3d_data_averaging 
241   END INTERFACE chem_3d_data_averaging
242
243   INTERFACE chem_wrd_local
244      MODULE PROCEDURE chem_wrd_local 
245   END INTERFACE chem_wrd_local
246
247   INTERFACE chem_rrd_local
248      MODULE PROCEDURE chem_rrd_local
249   END INTERFACE chem_rrd_local
250
251   INTERFACE chem_prognostic_equations
252      MODULE PROCEDURE chem_prognostic_equations
253      MODULE PROCEDURE chem_prognostic_equations_ij
254   END INTERFACE chem_prognostic_equations
255
256   INTERFACE chem_header
257      MODULE PROCEDURE chem_header
258   END INTERFACE chem_header
259
260   INTERFACE chem_emissions
261      MODULE PROCEDURE chem_emissions
262   END INTERFACE chem_emissions
263
264!    INTERFACE chem_wrd_global
265!       MODULE PROCEDURE chem_wrd_global
266!    END INTERFACE chem_wrd_global
267!
268!    INTERFACE chem_rrd_global
269!       MODULE PROCEDURE chem_rrd_global
270!    END INTERFACE chem_rrd_global
271
272
273   PUBLIC chem_3d_data_averaging, chem_boundary_conds, chem_check_data_output, &
274          chem_check_data_output_pr, chem_data_output_3d,                      &
275          chem_define_netcdf_grid, chem_emissions, chem_header, chem_init,     &
276          chem_init_profiles, chem_integrate, chem_wrd_local,                  &
277          chem_parin, chem_prognostic_equations,                               &
278          chem_rrd_local, chem_swap_timelevel
279         
280
281 CONTAINS
282
283!------------------------------------------------------------------------------!
284!
285! Description:
286! ------------
287!> Subroutine to initialize and set all boundary conditions for chemical species
288!------------------------------------------------------------------------------!
289
290 SUBROUTINE chem_boundary_conds( mode )                                           
291                                                                                 
292    USE control_parameters,                                                    & 
293        ONLY: air_chemistry, outflow_l, outflow_n, outflow_r, outflow_s                 
294    USE indices,                                                               & 
295        ONLY: nxl, nxr,  nxlg, nxrg, nyng, nysg, nzt                             
296                                                                                 
297!    USE prognostic_equations_mod,                                             &
298
299    USE arrays_3d,                                                             &     
300        ONLY: dzu                                               
301    USE surface_mod,                                                           &
302        ONLY: bc_h                                                           
303
304    CHARACTER (len=*), INTENT(IN) :: mode
305    INTEGER(iwp) ::  i                                                            !< grid index x direction.
306    INTEGER(iwp) ::  j                                                            !< grid index y direction.
307    INTEGER(iwp) ::  k                                                            !< grid index z direction.
308    INTEGER(iwp) ::  kb                                                           !< variable to set respective boundary value, depends on facing.
309    INTEGER(iwp) ::  l                                                            !< running index boundary type, for up- and downward-facing walls.
310    INTEGER(iwp) ::  m                                                            !< running index surface elements.
311    INTEGER(iwp) ::  lsp                                                          !< running index for chem spcs.
312    INTEGER(iwp) ::  lph                                                          !< running index for photolysis frequencies
313
314
315    SELECT CASE ( TRIM( mode ) )       
316       CASE ( 'init' )     
317          DO lsp = 1, nspec           
318            IF ( surface_csflux(lsp) == 9999999.9_wp )  THEN
319                 constant_csflux(lsp) = .FALSE.           
320            ENDIF
321          ENDDO
322
323          IF ( bc_cs_b == 'dirichlet' )    THEN
324             ibc_cs_b = 0 
325          ELSEIF ( bc_cs_b == 'neumann' )  THEN
326             ibc_cs_b = 1 
327          ELSE
328!             message_string = 'unknown boundary condition: bc_cs_b ="' // TRIM( bc_cs_b ) // '"'  ! bK commented
329             CALL message( 'chem_boundary_conds', 'CM0010', 1, 2, 0, 6, 0 )     !< chemistry_model_mod should have special error numbers --> "CHEM###",
330          ENDIF                                                                 
331!
332!--       Set Integer flags and check for possible erroneous settings for top
333!--       boundary condition. bK added *_cs_* here.
334          IF ( bc_cs_t == 'dirichlet' )             THEN
335             ibc_cs_t = 0 
336          ELSEIF ( bc_cs_t == 'neumann' )           THEN
337             ibc_cs_t = 1
338          ELSEIF ( bc_cs_t == 'initial_gradient' )  THEN
339             ibc_cs_t = 2
340          ELSEIF ( bc_cs_t == 'nested' )            THEN         
341             ibc_cs_t = 3
342          ELSE
343!            message_string = 'unknown boundary condition: bc_c_t ="' // TRIM( bc_cs_t ) // '"' 
344             CALL message( 'check_parameters', 'CM0011', 1, 2, 0, 6, 0 )
345          ENDIF
346
347     
348       CASE ( 'set_bc_bottomtop' )                   
349!--      Bottom boundary condtions for chemical species     
350          DO lsp = 1, nspec                                                     
351             IF ( ibc_cs_b == 0 ) THEN                   
352                DO l = 0, 1 
353!--             Set index kb: For upward-facing surfaces (l=0), kb=-1, i.e.
354!--             the chem_species(nsp)%conc_p value at the topography top (k-1)
355!--             is set; for downward-facing surfaces (l=1), kb=1, i.e. the
356!--             value at the topography bottom (k+1) is set.
357
358                   kb = MERGE( -1, 1, l == 0 )
359                   !$OMP PARALLEL DO PRIVATE( i, j, k )
360                   DO m = 1, bc_h(l)%ns
361                      i = bc_h(l)%i(m)           
362                      j = bc_h(l)%j(m)
363                      k = bc_h(l)%k(m)
364                      chem_species(lsp)%conc_p(k+kb,j,i) = chem_species(lsp)%conc(k+kb,j,i) 
365                   ENDDO                                       
366                ENDDO                                       
367         
368             ELSEIF ( ibc_cs_b == 1 ) THEN
369         ! in boundary_conds there is som extra loop over m here for passive tracer
370                DO l = 0, 1
371                   kb = MERGE( -1, 1, l == 0 )
372                   !$OMP PARALLEL DO PRIVATE( i, j, k )                                           
373                   DO m = 1, bc_h(l)%ns
374                      i = bc_h(l)%i(m)           
375                      j = bc_h(l)%j(m)
376                      k = bc_h(l)%k(m)
377                      chem_species(lsp)%conc_p(k+kb,j,i) = chem_species(lsp)%conc_p(k,j,i)
378
379!< @todo: chem_species(nsp)%conc_p(k+kb,j,i) = chem_species(nsp)%conc(k+kb,j,i),    &
380!<  pls loop over nsp=1, NSPEC.
381!< @todo: We should also think about the possibility to have &
382!< individual boundary conditions for each species? This means, bc_cs_b,       &
383!< bc_cs_t, ibc_cs_b, ibc_cs_t would need to be added to TYPE chem_species(nsp)%,   &
384!< and the loop over nsp would have to be put outside of this IF-clause.
385!< i think its better we keep the same bonundary cond i.e. dirichlet or neumann
386!< for all chem spcs. ... !bK
387
388                   ENDDO
389                ENDDO
390             ENDIF
391          ENDDO    ! end lsp loop 
392
393!--    Top boundary conditions for chemical species - Should this not be done for all species?
394          IF ( ibc_cs_t == 0 )  THEN                   
395             DO lsp = 1, nspec
396                chem_species(lsp)%conc_p(nzt+1,:,:) = chem_species(lsp)%conc(nzt+1,:,:)       
397             ENDDO
398          ELSEIF ( ibc_cs_t == 1 )  THEN
399             DO lsp = 1, nspec
400                chem_species(lsp)%conc_p(nzt+1,:,:) = chem_species(lsp)%conc_p(nzt,:,:)
401             ENDDO
402          ELSEIF ( ibc_cs_t == 2 )  THEN
403             DO lsp = 1, nspec
404                chem_species(lsp)%conc_p(nzt+1,:,:) = chem_species(lsp)%conc_p(nzt,:,:) + bc_cs_t_val(lsp) * dzu(nzt+1)
405             ENDDO
406                                    !@todo: bc_cs_t_val needs to be calculated,   
407                                    !see bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1)
408                                    !(in time_integration). pt_init would be the counterpart to
409                                    !chem_species(i)%conc_pr_init (see kchem_driver_FKa1408.f90 of my
410                                    !"Hints: prescribing initial concentration.
411          ENDIF
412!
413       CASE ( 'set_bc_lateral' )                       ! bK commented it
414!--          Lateral boundary conditions for chem species at inflow boundary
415!--          are automatically set when chem_species concentration is
416!--          initialized. The initially set value at the inflow boundary is not
417!--          touched during time integration, hence, this boundary value remains
418!--          at a constant value, which is the concentration that flows into the
419!--          domain.                                                           
420!--          Lateral boundary conditions for chem species at outflow boundary
421
422          IF ( outflow_s )  THEN
423             DO lsp = 1, nspec
424                chem_species(lsp)%conc_p(:,nys-1,:) = chem_species(lsp)%conc_p(:,nys,:)
425             ENDDO
426          ELSEIF ( outflow_n )  THEN
427             DO lsp = 1, nspec
428                chem_species(lsp)%conc_p(:,nyn+1,:) = chem_species(lsp)%conc_p(:,nyn,:)
429             ENDDO
430          ELSEIF ( outflow_l )  THEN
431             DO lsp = 1, nspec
432                chem_species(lsp)%conc_p(:,:,nxl-1) = chem_species(lsp)%conc_p(:,:,nxl)
433             ENDDO
434          ELSEIF ( outflow_r )  THEN
435             DO lsp = 1, nspec
436                chem_species(lsp)%conc_p(:,:,nxr+1) = chem_species(lsp)%conc_p(:,:,nxr)
437             ENDDO
438          ENDIF
439         
440    END SELECT
441
442 END SUBROUTINE chem_boundary_conds
443!
444!------------------------------------------------------------------------------!
445!
446! Description:
447! ------------
448!> Subroutine defining initial vertical profiles of chemical species (given by
449!> namelist parameters chem_profiles and chem_heights)  --> which should work
450!> analogue to parameters u_profile, v_profile and uv_heights)
451!------------------------------------------------------------------------------!
452 SUBROUTINE chem_init_profiles              !< SUBROUTINE is called from chem_init in case of
453                                            !< TRIM( initializing_actions ) /= 'read_restart_data'
454                                            !< We still need to see what has to be done in case of restart run
455    USE chem_modules
456
457    IMPLICIT NONE
458
459!-- Local variables
460    INTEGER ::  lsp        !< running index for number of species in derived data type species_def
461    INTEGER ::  lsp_pr     !< running index for number of species in cs_names, cs_profiles etc
462    INTEGER ::  lpr_lev    !< running index for profile level for each chem spcs.
463    INTEGER ::  npr_lev    !< the next available profile lev
464
465!-----------------
466!-- To prescribe/initialize a vertically constant  'cs_profile', another parameter
467!-- "cs_surface" is introduced. If "cs_profile" and "cs_heights" are prescribed,
468!-- their values will override the constant profile given by "cs_surface".
469
470!     IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
471       lsp_pr = 1
472       DO  WHILE ( TRIM( cs_name( lsp_pr ) ) /= 'novalue' )   !'novalue' is the default
473          DO  lsp = 1, nspec                                !
474!--          Create initial profile (conc_pr_init) for each chemical species
475             IF ( TRIM( chem_species(lsp)%name ) == TRIM( cs_name(lsp_pr) ) )  THEN   !
476!                IF ( cs_profile(1,1) == 9999999.9_wp ) THEN
477!--               Set a vertically constant profile based on the surface conc (cs_surface(lsp_pr)) of each species
478                  DO lpr_lev = 0, nzt+1
479                     chem_species(lsp)%conc_pr_init(lpr_lev) = cs_surface(lsp_pr)
480                  ENDDO
481
482!                ELSE
483!                   IF ( cs_heights(lsp,1) /= 0.0_wp )  THEN
484!                      message_string = 'cs_heights(1,1) must be 0.0'
485!                      CALL message( 'chem_check_parameters', 'CM0012', 1, 2, 0, 6, 0 )
486!                   ENDIF
487
488!                   IF ( omega /= 0.0_wp )  THEN
489!                      message_string = 'Coriolis force must be switched off (by setting omega=0.0)' //  &
490!                                       ' when prescribing the forcing by u_profile and v_profile'
491!                      CALL message( 'check_parameters', 'PA0347', 1, 2, 0, 6, 0 )
492!                   ENDIF
493
494!                   use_prescribed_profile_data = .TRUE.
495!
496!                   npr_lev = 1
497! !                  chem_sddpecies(lsp)%conc_pr_init(0) = 0.0_wp
498!                   DO  lpr_lev = 1, nz+1
499!                      IF ( npr_lev < 100 )  THEN
500!                         DO  WHILE ( cs_heights(lsp, npr_lev+1) <= zu(lpr_lev) )
501!                            npr_lev = npr_lev + 1
502!                            IF ( npr_lev == 100 )  EXIT
503!                         ENDDO
504!                      ENDIF
505
506!                      IF ( npr_lev < 100  .AND.  cs_heights(lsp, npr_lev + 1) /= 9999999.9_wp )  THEN
507!                         chem_species(lsp)%conc_pr_init(lpr_lev) = cs_profile(lsp, npr_lev) +         &
508!                                                 ( zu(lpr_lev) - cs_heights(lsp, npr_lev) ) /         &
509!                                                 ( cs_heights(lsp, npr_lev + 1) - cs_heights(lsp, npr_lev ) ) * &
510!                                                 ( cs_profile(lsp, npr_lev + 1) - cs_profile(lsp, npr_lev ) )
511!                      ELSE
512!                         chem_species(lsp)%conc_pr_init(lpr_lev) = cs_profile(lsp, npr_lev)
513!                      ENDIF
514!                   ENDDO
515!                ENDIF
516!-- If a profile is prescribed explicity using cs_profiles and cs_heights, we then have to fill the
517!-- chem_species(lsp)%conc_pr_init for the specific "lsp" based on the cs_profiles(lsp_pr,:)
518!-- and cs_heights(lsp_pr,:).
519             ENDIF
520          ENDDO
521          lsp_pr = lsp_pr + 1
522       ENDDO
523!     ELSE
524!       IF (chem_debug1 ) print*,'code to be added for initializing_actions == read_restart_data'   !bK
525!     ENDIF
526
527!-- Now, go back to chem_init and use the contents of chem_species(lsp)%conc_pr_init to
528!-- initialize the 3D conc arrays, as it is currently taken care of in chem_set_constant_values.
529!-- After initializing the 3D arrays, these can be used to set the boundary conditions in the
530!-- subroutine kchem_boundary_conds, which should be called from boundary_conds.f90.
531
532
533 END SUBROUTINE chem_init_profiles
534!
535!------------------------------------------------------------------------------!
536!
537! Description:
538! ------------
539!> Subroutine initializating chemistry_model_mod
540!------------------------------------------------------------------------------!
541   SUBROUTINE chem_init
542
543
544      USE control_parameters,                                                  &
545         ONLY: message_string, io_blocks, io_group, turbulent_inflow
546       
547      USE arrays_3d,                                                           &
548          ONLY: mean_inflow_profiles
549
550      USE pegrid
551
552    IMPLICIT   none
553!-- local variables
554    INTEGER                 :: i,j               !< running index for for horiz numerical grid points
555    INTEGER                 :: lsp               !< running index for chem spcs
556    INTEGER                 :: lpr_lev           !< running index for chem spcs profile level
557!
558!-- NOPOINTER version not implemented yet
559! #if defined( __nopointer )
560!     message_string = 'The chemistry module only runs with POINTER version'
561!     CALL message( 'chemistry_model_mod', 'CM0001', 1, 2, 0, 6, 0 )     
562! #endif
563!
564!-- Allocate memory for chemical species
565    ALLOCATE( chem_species(nspec) )
566    ALLOCATE( spec_conc_1 (nzb:nzt+1,nysg:nyng,nxlg:nxrg,nspec) )
567    ALLOCATE( spec_conc_2 (nzb:nzt+1,nysg:nyng,nxlg:nxrg,nspec) )
568    ALLOCATE( spec_conc_3 (nzb:nzt+1,nysg:nyng,nxlg:nxrg,nspec) )
569    ALLOCATE( spec_conc_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nspec) ) 
570    ALLOCATE( phot_frequen(nphot) ) 
571    ALLOCATE( freq_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nphot) )
572    ALLOCATE( bc_cs_t_val(nspec) )
573!
574!-- Initialize arrays
575    spec_conc_1 (:,:,:,:) = 0.0_wp
576    spec_conc_2 (:,:,:,:) = 0.0_wp
577    spec_conc_3 (:,:,:,:) = 0.0_wp
578    spec_conc_av(:,:,:,:) = 0.0_wp
579
580
581    DO lsp = 1, nspec
582       chem_species(lsp)%name    = spc_names(lsp)
583
584       chem_species(lsp)%conc   (nzb:nzt+1,nysg:nyng,nxlg:nxrg)       => spec_conc_1 (:,:,:,lsp)
585       chem_species(lsp)%conc_p (nzb:nzt+1,nysg:nyng,nxlg:nxrg)       => spec_conc_2 (:,:,:,lsp)
586       chem_species(lsp)%tconc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg)       => spec_conc_3 (:,:,:,lsp)
587       chem_species(lsp)%conc_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg)       => spec_conc_av(:,:,:,lsp)     
588
589       ALLOCATE (chem_species(lsp)%cssws_av(nysg:nyng,nxlg:nxrg))                   
590       chem_species(lsp)%cssws_av    = 0.0_wp
591!
592!--    (todo (FK): This needs to be revised. This block must go somewhere else)                                         
593!      IF ( ws_scheme_sca )  THEN                                               
594          ALLOCATE (chem_species(lsp)%flux_s_cs(nzb+1:nzt,0:threads_per_task-1))   
595          ALLOCATE (chem_species(lsp)%diss_s_cs(nzb+1:nzt,0:threads_per_task-1))   
596          ALLOCATE (chem_species(lsp)%flux_l_cs(nzb+1:nzt,nys:nyn,0:threads_per_task-1)) 
597          ALLOCATE (chem_species(lsp)%diss_l_cs(nzb+1:nzt,nys:nyn,0:threads_per_task-1))   
598          chem_species(lsp)%flux_s_cs = 0.0_wp                                     
599          chem_species(lsp)%flux_l_cs = 0.0_wp                                     
600          chem_species(lsp)%diss_s_cs = 0.0_wp                                     
601          chem_species(lsp)%diss_l_cs = 0.0_wp                                     
602!      ENDIF         
603!
604!--   Allocate memory for initial concentration profiles
605!--   (concentration values come from namelist)
606!--   (todo (FK): Because of this, chem_init is called in palm before
607!--               check_parameters, since conc_pr_init is used there.
608!--               We have to find another solution since chem_init should
609!--               eventually be called from init_3d_model!!)
610      ALLOCATE ( chem_species(lsp)%conc_pr_init(0:nz+1) )
611      chem_species(lsp)%conc_pr_init(:) = 0.0_wp
612
613    ENDDO
614
615!
616!-- Set initial concentration of profiles prescribed by parameters cs_profile
617!-- and cs_heights in the namelist &chemistry_parameters
618!-- (todo (FK): chem_init_profiles not ready yet, has some bugs)
619!     CALL chem_init_profiles
620!
621!-- Initialize model variables
622
623
624    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
625         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
626
627
628!--    First model run of a possible job queue.
629!--    Initial profiles of the variables must be computes.
630       IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
631!            CALL location_message( 'initializing with 1D model profiles', .FALSE. )
632!
633!          CALL init_1d_model    ...... decide to call or not later     !bK
634
635!--        Transfer initial profiles to the arrays of the 3D model
636          DO lsp = 1, nspec
637             DO  i = nxlg, nxrg
638                DO  j = nysg, nyng
639                   DO lpr_lev = 1, nz + 1
640                      chem_species(lsp)%conc(lpr_lev,j,i) = chem_species(lsp)%conc_pr_init(lpr_lev)
641                   ENDDO
642                ENDDO
643             ENDDO   
644          ENDDO
645       
646       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )    &
647       THEN
648!           CALL location_message( 'initializing with constant profiles', .FALSE. )
649
650
651
652!--       Set initial horizontal velocities at the lowest computational grid
653!--       levels to zero in order to avoid too small time steps caused by the
654!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
655!--       in the limiting formula!).
656
657          DO lsp = 1, nspec 
658             DO  i = nxlg, nxrg
659                DO  j = nysg, nyng
660                   chem_species(lsp)%conc(:,j,i) = chem_species(lsp)%conc_pr_init    !ITS THERE bK
661                ENDDO
662             ENDDO
663          ENDDO
664
665!        ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )                  &
666!        THEN
667!           CALL location_message( 'initializing by user', .FALSE. )
668!
669!--       Initialization will completely be done by the user
670!--       (FK: This should be called only once, in init_3d_model, i.e. remove it here)
671!           CALL user_init_3d_model
672!           CALL location_message( 'finished', .TRUE. )
673
674       ENDIF
675
676!-- Store initial chem spcs profile
677!         DO lsp = 1, nvar
678!          hom_cs(:,1,115,:) = SPREAD(  chem_species(lsp)%conc(:,nys,nxl), 2, statistic_regions+1 )
679!         ENDDO
680!
681!-- If required, change the surface chem spcs at the start of the 3D run
682       IF ( cs_surface_initial_change(1) /= 0.0_wp ) THEN           
683          DO lsp = 1, nspec 
684             chem_species(lsp)%conc(nzb,:,:) = chem_species(lsp)%conc(nzb,:,:) +  &
685                                               cs_surface_initial_change(lsp)
686          ENDDO
687       ENDIF 
688!
689!-- Initiale old and new time levels.
690       DO lsp = 1, nvar
691          chem_species(lsp)%tconc_m = 0.0_wp                     
692          chem_species(lsp)%conc_p  = chem_species(lsp)%conc     
693       ENDDO
694
695    ENDIF
696
697
698
699!--- new code add above this line
700    DO lsp = 1, nphot
701       phot_frequen(lsp)%name = phot_names(lsp)
702!        IF( myid == 0 )  THEN
703!           WRITE(6,'(a,i4,3x,a)')  'Photolysis: ',lsp,trim(phot_names(lsp))
704!        ENDIF
705         phot_frequen(lsp)%freq(nzb:nzt+1,nysg:nyng,nxlg:nxrg)    => freq_1(:,:,:,lsp)
706    ENDDO
707
708!--   Set initial values
709!    Not required any more, this can now be done with the namelist  by setting cs_surface
710!    and cs_name without specifying cs_profile (Nevertheless, we still want to keep it for a while)
711!   IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN             
712!         CALL set_const_initial_values
713!   ENDIF
714
715    RETURN
716
717    CONTAINS
718!------------------------------------------------------------------------------!
719!
720! Description:
721! ------------
722!> Subroutine setting constant initial values of chemical species
723!------------------------------------------------------------------------------!
724     SUBROUTINE set_const_initial_values
725!    Not required any more, this can now be done with the namelist  by setting cs_surface
726!    and cs_name without specifying cs_profile (Nevertheless, we still want to keep it for a while)
727         IMPLICIT   none
728
729!--      local variables
730         INTEGER                  :: lsp
731
732         IF(myid == 0)  THEN
733            write(6,'(/,a,/)')  ' chemics >>>> Set constant Initial Values: '
734         ENDIF
735
736!        Default values are taken from smog.def from supplied kpp example
737         DO lsp = 1, nspec
738            IF(trim(chem_species(lsp)%name) == 'NO')           THEN
739!              chem_species(lsp)%conc = 8.725*1.0E+08
740!              chem_species(lsp)%conc = 0.05_wp                          !added by bK
741               chem_species(lsp)%conc = 0.01_wp                          !added by RFo
742            ELSEIF (trim(chem_species(lsp)%name) == 'NO2')     THEN
743!              chem_species(lsp)%conc = 2.240*1.0E+08             
744!              chem_species(lsp)%conc = 0.01_wp                          !added by bK
745               chem_species(lsp)%conc = 0.05_wp                          !added by RFo
746            ELSEIF( trim( chem_species(lsp)%name ) == 'O3' )   THEN
747               chem_species(lsp)%conc = 0.05_wp                          !added by bK
748            ELSEIF(trim(chem_species(lsp)%name) == 'H2O')      THEN
749!              chem_species(lsp)%conc = 5.326*1.0E+11
750               chem_species(lsp)%conc = 1.30*1.0E+4_wp                   !added by bK
751            ELSEIF(trim(chem_species(lsp)%name) == 'O2')       THEN
752               chem_species(lsp)%conc = 2.0*1.0E+5_wp                    !added by bK
753            ELSEIF(trim(chem_species(lsp)%name) == 'RH')       THEN
754               chem_species(lsp)%conc = 0.001_wp                         !added by RFo
755            ELSEIF(trim(chem_species(lsp)%name) == 'CO')       THEN
756               chem_species(lsp)%conc = 0.5_wp                           !added by RFo
757            ELSEIF(trim(chem_species(lsp)%name) == 'RCHO')     THEN
758!              chem_species(lsp)%conc = 2.0_wp                           !added by bK
759               chem_species(lsp)%conc = 0.01_wp                          !added by RFo
760!           ELSEIF(trim(chem_species(lsp)%name) == 'OH')       THEN
761!              chem_species(lsp)%conc = 1.0*1.0E-07_wp                   !added by bK
762!           ELSEIF(trim(chem_species(lsp)%name) == 'HO2')      THEN
763!              chem_species(lsp)%conc = 1*1.0E-7_wp                      !added by bK
764!           ELSEIF(trim(chem_species(lsp)%name) == 'RCOO2')    THEN      ! corrected RFo
765!              chem_species(lsp)%conc = 1.0*1.0E-7_wp                    !added by bK
766!           ELSEIF(trim(chem_species(lsp)%name) == 'RCOO2NO2') THEN
767!              chem_species(lsp)%conc = 1.0*1.0E-7_wp                    !added by bK
768            ELSE
769!              H2O = 2.0e+04;
770               chem_species(lsp)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) = 0.0_wp
771            ENDIF
772
773            IF(myid == 0)  THEN
774               WRITE(6,'(a,3x,a,3x,a,e12.4)')  '   Species:     ',chem_species(lsp)%name(1:7), & 
775                                        'Initial Value = ',chem_species(lsp)%conc(nzb,nysg,nxlg)
776            ENDIF
777         ENDDO
778
779!   #if defined( __nopointer )
780!   !kk      Hier mit message abbrechen
781!            if(myid == 0)  then
782!               write(6,*)  '   KPP does only run with POINTER Version'
783!            end if
784!            stop 'error'
785!   #endif
786
787         RETURN
788      END SUBROUTINE set_const_initial_values
789
790
791   END SUBROUTINE chem_init
792!
793!------------------------------------------------------------------------------!
794!
795! Description:
796! ------------
797!> Subroutine defining parin for &chemistry_parameters for chemistry model
798!------------------------------------------------------------------------------!
799   SUBROUTINE chem_parin
800   
801      USE control_parameters,                                                  &
802          ONLY: air_chemistry
803         
804      USE chem_modules
805     
806      USE kinds
807
808      IMPLICIT   none
809
810      CHARACTER (LEN=80) ::  line   !< dummy string that contains the current line of the parameter file
811     
812      REAL(wp), DIMENSION(nmaxfixsteps) ::   my_steps   !< List of fixed timesteps   my_step(1) = 0.0 automatic stepping
813
814      NAMELIST /chemistry_parameters/   bc_cs_b,                          &
815                                        bc_cs_t,                          &
816                                        call_chem_at_all_substeps,        &
817                                        chem_debug0,                      &
818                                        chem_debug1,                      &
819                                        chem_debug2,                      &
820                                        chem_gasphase_on,                 &
821                                        cs_heights,                       &
822                                        cs_name,                          &
823                                        cs_profile,                       &
824                                        cs_profile_name,                  & 
825                                        cs_surface,                       &
826                                        emiss_factor_main,                &
827                                        emiss_factor_side,                &                     
828                                        icntrl,                           &
829                                        main_street_id,                   &
830                                        max_street_id,                    &
831                                        my_steps,                         &
832                                        rcntrl,                           &
833                                        side_street_id,                   &
834                                        photolysis_scheme,                &
835                                        wall_csflux,                      &
836                                        cs_vertical_gradient,             &
837                                        top_csflux,                       & 
838                                        surface_csflux,                   &
839                                        surface_csflux_name,              &
840                                        cs_surface_initial_change,        &
841                                        cs_vertical_gradient_level
842                             
843!-- analogue to chem_names(nspj) we could invent chem_surfaceflux(nspj) and chem_topflux(nspj)
844!-- so this way we could prescribe a specific flux value for each species
845!>  chemistry_parameters for initial profiles
846!>  cs_names = 'O3', 'NO2', 'NO', ...   to set initial profiles)
847!>  cs_heights(1,:) = 0.0, 100.0, 500.0, 2000.0, .... (height levels where concs will be prescribed for O3)
848!>  cs_heights(2,:) = 0.0, 200.0, 400.0, 1000.0, .... (same for NO2 etc.)
849!>  cs_profiles(1,:) = 10.0, 20.0, 20.0, 30.0, .....  (chem spcs conc at height lvls chem_heights(1,:)) etc.
850!>  If the respective concentration profile should be constant with height, then use "cs_surface( number of spcs)"
851!>  then write these cs_surface values to chem_species(lsp)%conc_pr_init(:)
852
853!
854!--   Read chem namelist   
855!--   (todo: initialize these parameters in declaration part, do this for
856!--          all chemistry_parameters namelist parameters)
857      icntrl    = 0
858      rcntrl    = 0.0_wp
859      my_steps  = 0.0_wp
860      icntrl(2) = 1                                 
861      photolysis_scheme = 'simple'
862      atol = 1.0_wp
863      rtol = 0.01_wp
864!
865!--   Try to find chemistry package
866      REWIND ( 11 )
867      line = ' '
868      DO   WHILE ( INDEX( line, '&chemistry_parameters' ) == 0 )
869         READ ( 11, '(A)', END=10 )  line
870      ENDDO
871      BACKSPACE ( 11 )
872!
873!--   Read chemistry namelist
874      READ ( 11, chemistry_parameters )                         
875!
876!--   Enable chemistry model
877      air_chemistry = .TRUE.                   
878
879     
880 10   CONTINUE
881
882      t_steps = my_steps          !(todo: Why not directly make t_steps a
883                                  !       namelist parameter?)
884
885   END SUBROUTINE chem_parin
886
887 !
888!------------------------------------------------------------------------------!
889!
890! Description:
891! ------------
892!> Subroutine to integrate chemical species in the given chemical mechanism
893!------------------------------------------------------------------------------!
894
895   SUBROUTINE chem_integrate_ij (i, j)
896
897      USE cpulog,                                                              &
898           ONLY: cpu_log, log_point
899      USE statistics,                                                          &   ! ## RFo
900           ONLY:  weight_pres
901       USE control_parameters,                                                 &   ! ## RFo
902           ONLY:  dt_3d, intermediate_timestep_count
903
904      IMPLICIT   none
905      INTEGER,INTENT(IN)       :: i,j
906
907!--   local variables
908      INTEGER(iwp) ::  lsp                                                          !< running index for chem spcs.
909      INTEGER(iwp) ::  lph                                                          !< running index for photolysis frequencies
910      INTEGER                  :: k,m,istatf                                       
911      INTEGER,dimension(20)    :: istatus
912      REAL(kind=wp),dimension(nzb+1:nzt,nspec)                :: tmp_conc           
913      REAL(kind=wp),dimension(nzb+1:nzt)                      :: tmp_temp
914      REAL(kind=wp),dimension(nzb+1:nzt,nphot)                :: tmp_phot
915      REAL(kind=wp),dimension(nzb+1:nzt)                      :: tmp_fact   
916     REAL(kind=wp),dimension(nzb+1:nzt)                      :: tmp_fact_i    !< conversion factor between molecules cm^{-3} and ppm
917     REAL(wp)                         ::  conv                                !< conversion factor
918     REAL(wp), PARAMETER              ::  ppm2fr  = 1.0e-6_wp                 !< Conversion factor ppm to fraction
919     REAL(wp), PARAMETER              ::  pref_i  = 1._wp / 100000.0_wp       !< inverse reference pressure (1/Pa)
920     REAL(wp), PARAMETER              ::  t_std   = 273.15_wp                 !< standard pressure (Pa)
921     REAL(wp), PARAMETER              ::  p_std   = 101325.0_wp               !< standard pressure (Pa)
922     REAL(wp), PARAMETER              ::  r_cp    = 0.286_wp                  !< R / cp (exponent for potential temperature)
923     REAL(wp), PARAMETER              ::  vmolcm  = 22.414e3_wp               !< Mole volume (22.414 l) in cm^{-3}
924     REAL(wp), PARAMETER              ::  xna     = 6.022e23_wp               !< Avogadro number (molecules/mol)
925
926
927      REAL(kind=wp)  :: dt_chem                                             
928
929       CALL cpu_log( log_point(80), '[chem_integrate_ij]', 'start' )
930!<     Set chem_gasphase_on to .FALSE. if you want to skip computation of gas phase chemistry
931    IF (chem_gasphase_on) THEN
932
933       tmp_temp(:) = pt(:,j,i) * ( hyp(nzb+1:nzt) / 100000.0_wp )**0.286_wp
934! ppm to molecules/cm**3
935!      tmp_fact = 1.e-6_wp*6.022e23_wp/(22.414_wp*1000._wp) * 273.15_wp * hyp(nzb+1:nzt)/( 101300.0_wp * tmp_temp )
936       conv = ppm2fr * xna / vmolcm
937       tmp_fact(:) = conv * t_std * hyp(nzb+1:nzt) / (tmp_temp(:) * p_std)
938       tmp_fact_i = 1.0_wp/tmp_fact
939
940       CALL fill_temp (istatf, tmp_temp)                                      !< Load constant temperature into kpp context
941!      CALL fill_temp (istatf, pt(nzb+1:nzt,j,i))                             !< Load temperature into kpp context
942
943       DO lsp = 1,nspec
944          tmp_conc(:,lsp) = chem_species(lsp)%conc(nzb+1:nzt,j,i) * tmp_fact(:) 
945       ENDDO
946
947       DO lph = 1,nphot
948          tmp_phot(:,lph) = phot_frequen(lph)%freq(nzb+1:nzt,j,i)               
949       ENDDO
950
951       IF(myid == 0 .AND. chem_debug0 ) THEN
952          IF (i == 10 .and. j == 10) WRITE(0,*) 'begin chemics step ',dt_3d
953       ENDIF
954
955!--    Compute length of time step     # RFo
956       IF ( call_chem_at_all_substeps )  THEN
957          dt_chem = dt_3d * weight_pres(intermediate_timestep_count)
958       ELSE
959          dt_chem = dt_3d
960       ENDIF
961
962       CALL cpu_log( log_point(81), '{chem_gasphase_integrate}', 'start' )
963
964
965       CALL chem_gasphase_integrate (dt_chem, tmp_conc, tmp_temp, tmp_phot, istatus=istatus)
966
967
968       CALL cpu_log( log_point(81), '{chem_gasphase_integrate}', 'stop' )
969
970       DO lsp = 1,nspec
971          chem_species(lsp)%conc (nzb+1:nzt,j,i) = tmp_conc(:,lsp) * tmp_fact_i(:)  ! RFo
972       ENDDO
973
974!      IF (myid == 0 .AND. chem_debug2 ) THEN
975!         IF (i == 10 .and. j == 10)   WRITE(6,'(a,8i7)') ' KPP Status ',istatus(1:8)
976!         write(6,'(a,8i7)') ' KPP Status ',istatus(1:8)
977!      ENDIF
978
979    ENDIF
980       CALL cpu_log( log_point(80), '[chem_integrate_ij]', 'stop' )
981
982       RETURN
983   END SUBROUTINE chem_integrate_ij
984!
985!------------------------------------------------------------------------------!
986!
987! Description:
988! ------------
989!> Subroutine for swapping of timelevels for chemical species
990!> called out from subroutine swap_timelevel
991!------------------------------------------------------------------------------!
992
993   SUBROUTINE chem_swap_timelevel (level)
994      IMPLICIT   none
995
996      INTEGER,INTENT(IN)                  :: level
997
998!--   local variables
999
1000      INTEGER               :: lsp
1001
1002!        print*,' *** entering chem_swap_timelevel ***) '
1003      if(level == 0)  then
1004         do lsp=1, nvar                                       
1005            chem_species(lsp)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)    => spec_conc_1(:,:,:,lsp)
1006            chem_species(lsp)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)  => spec_conc_2(:,:,:,lsp)
1007         end do
1008      else
1009         do lsp=1, nvar                                       
1010            chem_species(lsp)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)    => spec_conc_2(:,:,:,lsp)
1011            chem_species(lsp)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)  => spec_conc_1(:,:,:,lsp)
1012         end do
1013      end if
1014
1015      RETURN
1016   END SUBROUTINE chem_swap_timelevel
1017
1018!------------------------------------------------------------------------------!
1019!
1020! Description:
1021! ------------
1022!> Subroutine defining appropriate grid for netcdf variables.
1023!> It is called out from subroutine netcdf.
1024!------------------------------------------------------------------------------!
1025   SUBROUTINE chem_define_netcdf_grid( var, found, grid_x, grid_y, grid_z )
1026
1027      IMPLICIT NONE
1028
1029      CHARACTER (LEN=*), INTENT(IN)  ::  var         !<
1030      LOGICAL, INTENT(OUT)           ::  found       !<
1031      CHARACTER (LEN=*), INTENT(OUT) ::  grid_x      !<
1032      CHARACTER (LEN=*), INTENT(OUT) ::  grid_y      !<
1033      CHARACTER (LEN=*), INTENT(OUT) ::  grid_z      !<
1034
1035      found  = .TRUE.
1036
1037      if(var(1:3) == 'kc_')   then                    !< always the same grid for chemistry variables
1038            grid_x = 'x'
1039            grid_y = 'y'
1040            grid_z = 'zu'                             !< kk Use same z axis as u variables. Has to be checked if OK
1041      else
1042            found  = .FALSE.
1043            grid_x = 'none'
1044            grid_y = 'none'
1045            grid_z = 'none'
1046      end if
1047
1048!     write(6,*) 'chem_define_netcdf_grid ',TRIM(var),' ',trim(grid_x),' ',found
1049
1050   END SUBROUTINE chem_define_netcdf_grid
1051!
1052!------------------------------------------------------------------------------!
1053!
1054! Description:
1055! ------------
1056!> Subroutine for checking data output for chemical species
1057!------------------------------------------------------------------------------!
1058
1059   SUBROUTINE chem_check_data_output( var, unit, i, ilen, k )
1060
1061
1062      USE control_parameters,                                                 &
1063         ONLY: data_output, message_string
1064
1065      IMPLICIT NONE
1066
1067      CHARACTER (LEN=*) ::  unit     !<
1068      CHARACTER (LEN=*) ::  var      !<
1069
1070      INTEGER(iwp) :: i, lsp
1071      INTEGER(iwp) :: ilen
1072      INTEGER(iwp) :: k
1073
1074      CHARACTER(len=16)    ::  spec_name
1075
1076      unit = 'illegal'
1077
1078      spec_name = TRIM(var(4:))             !< var 1:3 is 'kc_'.
1079
1080       DO lsp=1,nspec
1081         IF (TRIM(spec_name) == TRIM(chem_species(lsp)%name))   THEN
1082            unit = 'ppm'
1083         ENDIF
1084         ! It is possible  to plant PM10 and PM25 into the gasphase chemistry code
1085!        ! as passive species (e.g. 'passive' in GASPHASE_PREPROC/mechanisms):
1086!        ! set unit to micrograms per m**3 for PM10 and PM25 (PM2.5)
1087         IF (spec_name(1:2) == 'PM')   THEN 
1088            unit = 'ug m-3'
1089         ENDIF
1090       ENDDO
1091
1092       DO lsp=1,nphot                                                     
1093         IF (TRIM(spec_name) == TRIM(phot_frequen(lsp)%name))   THEN
1094            unit = 'sec-1'
1095         ENDIF
1096       ENDDO
1097
1098
1099       RETURN
1100   END SUBROUTINE chem_check_data_output
1101!
1102!------------------------------------------------------------------------------!
1103!
1104! Description:
1105! ------------
1106!> Subroutine for checking data output of profiles for chemistry model
1107!------------------------------------------------------------------------------!
1108
1109   SUBROUTINE chem_check_data_output_pr( variable, var_count, unit, dopr_unit )
1110
1111
1112    USE arrays_3d,                                                             &
1113        ONLY: zu
1114
1115    USE control_parameters,                                                    &
1116        ONLY: data_output_pr, message_string, air_chemistry
1117
1118    USE indices
1119
1120    USE profil_parameter
1121
1122    USE statistics
1123
1124
1125    IMPLICIT NONE
1126
1127    CHARACTER (LEN=*) ::  unit     !<
1128    CHARACTER (LEN=*) ::  variable !<
1129    CHARACTER (LEN=*) ::  dopr_unit
1130    CHARACTER(len=16) ::  spec_name
1131 
1132    INTEGER(iwp) ::  var_count, lsp  !<
1133   
1134
1135    spec_name = TRIM(variable(4:))   
1136!             write(9,*) 'fm #32 .. air_chemistry ', air_chemistry
1137
1138          IF (  .NOT.  air_chemistry )  THEN
1139                 message_string = 'data_output_pr = ' //                        &
1140                 TRIM( data_output_pr(var_count) ) // ' is not imp' // &
1141                          'lemented for air_chemistry = .FALSE.'
1142!                CALL message( 'check_parameters', 'PA0XXX', 1, 2, 0, 6, 0 )
1143          ELSE
1144             DO lsp = 1, nspec
1145                IF (TRIM( spec_name ) == TRIM( chem_species(lsp)%name ) ) THEN
1146                    dopr_index(var_count) = 900 
1147                    dopr_unit  = 'ppm'
1148                    hom(:,2,900,:) = SPREAD( zu, 2, statistic_regions+1 )
1149                ENDIF
1150             ENDDO
1151          ENDIF
1152
1153   END SUBROUTINE chem_check_data_output_pr
1154!
1155!------------------------------------------------------------------------------!
1156!
1157! Description:
1158! ------------
1159!> Subroutine defining 3D output variables for chemical species
1160!------------------------------------------------------------------------------!
1161
1162
1163   SUBROUTINE chem_data_output_3d( av, variable, found, local_pf, fill_value, nzb_do, nzt_do )
1164
1165
1166      USE indices
1167
1168      USE kinds
1169
1170
1171      IMPLICIT NONE
1172
1173      CHARACTER (LEN=*) ::  variable !<
1174
1175      INTEGER(iwp) ::  av    !<
1176      INTEGER(iwp) ::  nzb_do !< lower limit of the data output (usually 0)
1177      INTEGER(iwp) ::  nzt_do !< vertical upper limit of the data output (usually nz_do3d)
1178
1179      LOGICAL      ::  found !<
1180
1181      REAL(wp) ::  fill_value !<
1182      REAL(sp), DIMENSION(nxl:nxr,nys:nyn,nzb_do:nzt_do) ::  local_pf 
1183
1184
1185      !-- local variables
1186
1187      INTEGER              ::  i, j, k, lsp
1188      CHARACTER(len=16)    ::  spec_name
1189
1190
1191      found = .FALSE.
1192
1193      spec_name = TRIM(variable(4:))
1194!av == 0
1195
1196       DO lsp=1,nspec
1197          IF (TRIM(spec_name) == TRIM(chem_species(lsp)%name))   THEN
1198             IF(myid == 0 .AND. chem_debug0 )  WRITE(6,*) 'Output of species ' // TRIM(variable)  //  &
1199                                                          TRIM(chem_species(lsp)%name)       
1200             
1201             IF (av == 0) THEN
1202                DO  i = nxl, nxr
1203                   DO  j = nys, nyn
1204                      DO  k = nzb_do, nzt_do
1205                          local_pf(i,j,k) = MERGE(                             &
1206                                              chem_species(lsp)%conc(k,j,i),   &
1207                                              REAL( fill_value, KIND = wp ),   &
1208                                              BTEST( wall_flags_0(k,j,i), 0 ) )
1209                      ENDDO
1210                   ENDDO
1211                ENDDO
1212           
1213             ELSE
1214                DO  i = nxl, nxr
1215                   DO  j = nys, nyn
1216                      DO  k = nzb_do, nzt_do
1217                          local_pf(i,j,k) = MERGE(                             &
1218                                              chem_species(lsp)%conc_av(k,j,i),&
1219                                              REAL( fill_value, KIND = wp ),   &
1220                                              BTEST( wall_flags_0(k,j,i), 0 ) )
1221                      ENDDO
1222                   ENDDO
1223                ENDDO
1224             ENDIF
1225
1226            found = .TRUE.
1227          ENDIF
1228       ENDDO
1229
1230       RETURN
1231   END SUBROUTINE chem_data_output_3d
1232!
1233!------------------------------------------------------------------------------!
1234!
1235! Description:
1236! ------------
1237!> Subroutine for averaging 3D data of chemical species. Due to the fact that
1238!> the averaged chem arrays are allocated in chem_init, no if-query concerning
1239!> the allocation is required (in any mode). Attention: If you just specify an
1240!> averaged output quantity in the _p3dr file during restarts the first output
1241!> includes the time between the beginning of the restart run and the first
1242!> output time (not necessarily the whole averaging_interval you have
1243!> specified in your _p3d/_p3dr file )
1244!------------------------------------------------------------------------------!
1245
1246    SUBROUTINE chem_3d_data_averaging ( mode, variable )
1247
1248    USE control_parameters
1249
1250    USE indices
1251
1252    USE kinds
1253
1254    USE surface_mod,                                                         &
1255        ONLY: surf_def_h, surf_lsm_h, surf_usm_h   
1256 
1257    IMPLICIT NONE
1258 
1259    CHARACTER (LEN=*) ::  mode    !<
1260    CHARACTER (LEN=*) :: variable !<
1261 
1262    LOGICAL      ::  match_def !< flag indicating natural-type surface
1263    LOGICAL      ::  match_lsm !< flag indicating natural-type surface
1264    LOGICAL      ::  match_usm !< flag indicating urban-type surface
1265
1266    INTEGER(iwp) ::  i                  !< grid index x direction
1267    INTEGER(iwp) ::  j                  !< grid index y direction
1268    INTEGER(iwp) ::  k                  !< grid index z direction
1269    INTEGER(iwp) ::  m                  !< running index surface type
1270    INTEGER(iwp) :: lsp                 !< running index for chem spcs
1271    INTEGER(iwp) :: lsp_2               !< it looks like redundent .. will be delted ..bK
1272 
1273    IF ( mode == 'allocate' )  THEN
1274       DO lsp = 1, nspec
1275          IF (TRIM(variable(4:)) == TRIM(chem_species(lsp)%name)) THEN
1276!                   lsp_2 = lsp
1277             chem_species(lsp)%conc_av = 0.0_wp
1278             
1279          ENDIF
1280       ENDDO
1281
1282    ELSEIF ( mode == 'sum' )  THEN
1283   
1284       DO lsp = 1, nspec
1285          IF (TRIM(variable(4:)) == TRIM(chem_species(lsp)%name)) THEN
1286!                   lsp_2 = lsp
1287             DO  i = nxlg, nxrg
1288                DO  j = nysg, nyng
1289                   DO  k = nzb, nzt+1
1290                       chem_species(lsp)%conc_av(k,j,i) = chem_species(lsp)%conc_av(k,j,i) + &
1291                                                          chem_species(lsp)%conc(k,j,i)
1292                   ENDDO
1293                ENDDO
1294             ENDDO
1295          ELSEIF ( TRIM(variable(4:)) == TRIM('cssws*') )        THEN
1296             DO  i = nxl, nxr
1297                DO  j = nys, nyn
1298                   match_def = surf_def_h(0)%start_index(j,i) <=               &
1299                               surf_def_h(0)%end_index(j,i)
1300                   match_lsm = surf_lsm_h%start_index(j,i) <=                  &
1301                               surf_lsm_h%end_index(j,i)
1302                   match_usm = surf_usm_h%start_index(j,i) <=                  &
1303                               surf_usm_h%end_index(j,i)
1304
1305                   IF ( match_def )  THEN
1306                      m = surf_def_h(0)%end_index(j,i)
1307                      chem_species(lsp)%cssws_av(j,i) =                        &
1308                                             chem_species(lsp)%cssws_av(j,i) + &
1309                                             surf_def_h(0)%cssws(lsp,m)
1310                   ELSEIF ( match_lsm  .AND.  .NOT. match_usm )  THEN
1311                      m = surf_lsm_h%end_index(j,i)
1312                      chem_species(lsp)%cssws_av(j,i) =                        &
1313                                             chem_species(lsp)%cssws_av(j,i) + &
1314                                             surf_lsm_h%cssws(lsp,m)
1315                   ELSEIF ( match_usm )  THEN
1316                      m = surf_usm_h%end_index(j,i)
1317                      chem_species(lsp)%cssws_av(j,i) =                        &
1318                                             chem_species(lsp)%cssws_av(j,i) + &
1319                                             surf_usm_h%cssws(lsp,m)
1320                   ENDIF
1321                ENDDO
1322             ENDDO
1323          ENDIF
1324       ENDDO
1325 
1326    ELSEIF ( mode == 'average' )  THEN
1327       DO lsp = 1, nspec
1328          IF (TRIM(variable(4:)) == TRIM(chem_species(lsp)%name)) THEN
1329!                   lsp_2 = lsp
1330             DO  i = nxlg, nxrg
1331                DO  j = nysg, nyng
1332                   DO  k = nzb, nzt+1
1333                       chem_species(lsp)%conc_av(k,j,i) = chem_species(lsp)%conc_av(k,j,i) / REAL( average_count_3d, KIND=wp )
1334                   ENDDO
1335                ENDDO
1336             ENDDO
1337
1338          ELSEIF (TRIM(variable(4:)) == TRIM('cssws*') )            THEN
1339             DO i = nxlg, nxrg
1340                DO  j = nysg, nyng
1341                     chem_species(lsp)%cssws_av(j,i) = chem_species(lsp)%cssws_av(j,i) / REAL( average_count_3d, KIND=wp )
1342                ENDDO
1343             ENDDO
1344                CALL exchange_horiz_2d( chem_species(lsp)%cssws_av, nbgp )                 
1345          ENDIF
1346       ENDDO
1347       
1348    ENDIF     
1349
1350    END SUBROUTINE chem_3d_data_averaging
1351
1352!------------------------------------------------------------------------------!
1353!
1354! Description:
1355! ------------
1356!> Subroutine to write restart data for chemistry model
1357!------------------------------------------------------------------------------!
1358 SUBROUTINE chem_wrd_local
1359
1360
1361    IMPLICIT NONE
1362
1363    INTEGER(iwp) ::  lsp !<
1364
1365!      REAL(kind=wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg)   :: chems_conc
1366
1367
1368    DO  lsp = 1, nspec
1369
1370       CALL wrd_write_string( TRIM( chem_species(lsp)%name ))
1371       WRITE ( 14 )  chem_species(lsp)%conc
1372
1373       CALL wrd_write_string( TRIM( chem_species(lsp)%name )//'_av' )
1374       WRITE ( 14 )  chem_species(lsp)%conc_av
1375
1376    ENDDO
1377
1378
1379 END SUBROUTINE chem_wrd_local
1380
1381!------------------------------------------------------------------------------!
1382!
1383! Description:
1384! ------------
1385!> Subroutine to read restart data of chemical species
1386!------------------------------------------------------------------------------!
1387
1388 SUBROUTINE chem_rrd_local( i, k, nxlf, nxlc, nxl_on_file, nxrf, nxrc,         &
1389                            nxr_on_file, nynf, nync, nyn_on_file, nysf, nysc,  &
1390                            nys_on_file, tmp_3d, found )   
1391                                     
1392    USE control_parameters
1393           
1394    USE indices
1395   
1396    USE pegrid
1397
1398    IMPLICIT NONE
1399
1400    CHARACTER (LEN=20) :: spc_name_av !<   
1401       
1402    INTEGER(iwp) ::  i, lsp          !<
1403    INTEGER(iwp) ::  k               !<
1404    INTEGER(iwp) ::  nxlc            !<
1405    INTEGER(iwp) ::  nxlf            !<
1406    INTEGER(iwp) ::  nxl_on_file     !<   
1407    INTEGER(iwp) ::  nxrc            !<
1408    INTEGER(iwp) ::  nxrf            !<
1409    INTEGER(iwp) ::  nxr_on_file     !<   
1410    INTEGER(iwp) ::  nync            !<
1411    INTEGER(iwp) ::  nynf            !<
1412    INTEGER(iwp) ::  nyn_on_file     !<   
1413    INTEGER(iwp) ::  nysc            !<
1414    INTEGER(iwp) ::  nysf            !<
1415    INTEGER(iwp) ::  nys_on_file     !<   
1416   
1417    LOGICAL, INTENT(OUT)  :: found 
1418
1419    REAL(wp), DIMENSION(nzb:nzt+1,nys_on_file-nbgp:nyn_on_file+nbgp,nxl_on_file-nbgp:nxr_on_file+nbgp) :: tmp_3d   !< 3D array to temp store data
1420
1421
1422    found = .FALSE. 
1423
1424
1425       IF ( ALLOCATED(chem_species) )  THEN
1426
1427          DO  lsp = 1, nspec
1428
1429              !< for time-averaged chemical conc.
1430             spc_name_av  =  TRIM(chem_species(lsp)%name)//'_av'
1431
1432             IF (restart_string(1:length) == TRIM(chem_species(lsp)%name) )    &
1433             THEN
1434                !< read data into tmp_3d
1435                IF ( k == 1 )  READ ( 13 )  tmp_3d 
1436                !< fill ..%conc in the restart run   
1437                chem_species(lsp)%conc(:,nysc-nbgp:nync+nbgp,                  &
1438                                       nxlc-nbgp:nxrc+nbgp) =                  & 
1439                tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
1440                found = .TRUE.
1441             ELSEIF (restart_string(1:length) == spc_name_av )  THEN
1442                IF ( k == 1 )  READ ( 13 )  tmp_3d
1443                chem_species(lsp)%conc_av(:,nysc-nbgp:nync+nbgp,               &
1444                                          nxlc-nbgp:nxrc+nbgp) =               &
1445                tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
1446                found = .TRUE.
1447             ENDIF
1448
1449          ENDDO
1450
1451       ENDIF
1452
1453
1454 END SUBROUTINE chem_rrd_local
1455
1456
1457!------------------------------------------------------------------------------!
1458!
1459! Description:
1460! ------------
1461!> Subroutine calculating prognostic equations for chemical species
1462!> (cache-optimized).
1463!> Routine is called separately for each chemical species over a loop from
1464!> prognostic_equations.
1465!------------------------------------------------------------------------------!
1466 SUBROUTINE chem_prognostic_equations_ij ( cs_scalar_p, cs_scalar, tcs_scalar_m, pr_init_cs,  &
1467                            i, j, i_omp_start, tn, ilsp, flux_s_cs, diss_s_cs, &
1468                            flux_l_cs, diss_l_cs )
1469    USE pegrid         
1470    USE advec_ws,        ONLY:  advec_s_ws 
1471    USE advec_s_pw_mod,  ONLY:  advec_s_pw
1472    USE advec_s_up_mod,  ONLY:  advec_s_up
1473    USE diffusion_s_mod, ONLY:  diffusion_s
1474    USE indices,         ONLY:  wall_flags_0
1475    USE surface_mod,     ONLY:  surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h,    &
1476                                surf_usm_v
1477
1478
1479    IMPLICIT NONE
1480
1481    REAL(wp), DIMENSION(:,:,:), POINTER   :: cs_scalar_p, cs_scalar, tcs_scalar_m
1482
1483    INTEGER(iwp),INTENT(IN) :: i, j, i_omp_start, tn, ilsp
1484    REAL(wp), DIMENSION(nzb+1:nzt,0:threads_per_task-1)         :: flux_s_cs   !<
1485    REAL(wp), DIMENSION(nzb+1:nzt,0:threads_per_task-1)         :: diss_s_cs   !<
1486    REAL(wp), DIMENSION(nzb+1:nzt,nys:nyn,0:threads_per_task-1) :: flux_l_cs   !<
1487    REAL(wp), DIMENSION(nzb+1:nzt,nys:nyn,0:threads_per_task-1) :: diss_l_cs   !<
1488    REAL(wp), DIMENSION(0:nz+1)                                 :: pr_init_cs  !<
1489
1490!-- local variables
1491
1492    INTEGER :: k
1493    !
1494    !--    Tendency-terms for chem spcs.
1495    tend(:,j,i) = 0.0_wp
1496!   
1497!-- Advection terms
1498    IF ( timestep_scheme(1:5) == 'runge' )  THEN
1499       IF ( ws_scheme_sca )  THEN
1500          CALL advec_s_ws( i, j, cs_scalar, 'kc', flux_s_cs, diss_s_cs,                &
1501             flux_l_cs, diss_l_cs, i_omp_start, tn )
1502       ELSE
1503          CALL advec_s_pw( i, j, cs_scalar )
1504       ENDIF
1505    ELSE
1506         CALL advec_s_up( i, j, cs_scalar )
1507    ENDIF
1508
1509!
1510
1511!-- Diffusion terms (the last three arguments are zero)
1512
1513      CALL diffusion_s( i, j, cs_scalar,                                                 &
1514                        surf_def_h(0)%cssws(ilsp,:), surf_def_h(1)%cssws(ilsp,:),        &
1515                        surf_def_h(2)%cssws(ilsp,:),                                     &
1516                        surf_lsm_h%cssws(ilsp,:), surf_usm_h%cssws(ilsp,:),              &
1517                        surf_def_v(0)%cssws(ilsp,:), surf_def_v(1)%cssws(ilsp,:),        &
1518                        surf_def_v(2)%cssws(ilsp,:), surf_def_v(3)%cssws(ilsp,:),        &
1519                        surf_lsm_v(0)%cssws(ilsp,:), surf_lsm_v(1)%cssws(ilsp,:),        &
1520                        surf_lsm_v(2)%cssws(ilsp,:), surf_lsm_v(3)%cssws(ilsp,:),        &
1521                        surf_usm_v(0)%cssws(ilsp,:), surf_usm_v(1)%cssws(ilsp,:),        &
1522                        surf_usm_v(2)%cssws(ilsp,:), surf_usm_v(3)%cssws(ilsp,:) )
1523 
1524!   
1525!-- Prognostic equation for chem spcs
1526    DO k = nzb+1, nzt
1527       cs_scalar_p(k,j,i) = cs_scalar(k,j,i) + ( dt_3d  *                       &
1528                                               ( tsc(2) * tend(k,j,i) +         &
1529                                                 tsc(3) * tcs_scalar_m(k,j,i) ) & 
1530                                               - tsc(5) * rdf_sc(k)             &
1531                                                        * ( cs_scalar(k,j,i) - pr_init_cs(k) )    &   
1532                                               )                                                  &
1533                                                        * MERGE( 1.0_wp, 0.0_wp,                  &     
1534                                                                BTEST( wall_flags_0(k,j,i), 0 )   &             
1535                                                                 )       
1536
1537       IF ( cs_scalar_p(k,j,i) < 0.0_wp )  cs_scalar_p(k,j,i) = 0.1_wp * cs_scalar(k,j,i)    !FKS6
1538    ENDDO
1539
1540!
1541!-- Calculate tendencies for the next Runge-Kutta step
1542    IF ( timestep_scheme(1:5) == 'runge' )  THEN
1543       IF ( intermediate_timestep_count == 1 )  THEN
1544          DO  k = nzb+1, nzt
1545             tcs_scalar_m(k,j,i) = tend(k,j,i)
1546          ENDDO
1547       ELSEIF ( intermediate_timestep_count < &
1548          intermediate_timestep_count_max )  THEN
1549          DO  k = nzb+1, nzt
1550             tcs_scalar_m(k,j,i) = -9.5625_wp * tend(k,j,i) + &
1551                5.3125_wp * tcs_scalar_m(k,j,i)
1552          ENDDO
1553       ENDIF
1554    ENDIF
1555
1556 END SUBROUTINE chem_prognostic_equations_ij
1557
1558
1559!------------------------------------------------------------------------------!
1560!
1561! Description:
1562! ------------
1563!> Subroutine calculating prognostic equations for chemical species
1564!> (vector-optimized).
1565!> Routine is called separately for each chemical species over a loop from
1566!> prognostic_equations.
1567!------------------------------------------------------------------------------!
1568 SUBROUTINE chem_prognostic_equations ( cs_scalar_p, cs_scalar, tcs_scalar_m,  &
1569                                        pr_init_cs, ilsp )
1570
1571    USE advec_s_pw_mod,                                                        &
1572        ONLY:  advec_s_pw
1573
1574    USE advec_s_up_mod,                                                        &
1575        ONLY:  advec_s_up
1576
1577    USE advec_ws,                                                              &
1578        ONLY:  advec_s_ws 
1579
1580    USE diffusion_s_mod,                                                       &
1581        ONLY:  diffusion_s
1582
1583    USE indices,                                                               &
1584        ONLY:  nxl, nxr, nyn, nys, wall_flags_0
1585
1586    USE pegrid
1587
1588    USE surface_mod,                                                           &
1589        ONLY:  surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h,     &
1590               surf_usm_v
1591
1592    IMPLICIT NONE
1593
1594    INTEGER ::  i   !< running index
1595    INTEGER ::  j   !< running index
1596    INTEGER ::  k   !< running index
1597
1598    INTEGER(iwp),INTENT(IN) ::  ilsp          !<
1599
1600    REAL(wp), DIMENSION(0:nz+1) ::  pr_init_cs   !<
1601
1602    REAL(wp), DIMENSION(:,:,:), POINTER ::  cs_scalar      !<
1603    REAL(wp), DIMENSION(:,:,:), POINTER ::  cs_scalar_p    !<
1604    REAL(wp), DIMENSION(:,:,:), POINTER ::  tcs_scalar_m   !<
1605
1606
1607!
1608!-- Tendency terms for chemical species
1609    tend = 0.0_wp
1610!   
1611!-- Advection terms
1612    IF ( timestep_scheme(1:5) == 'runge' )  THEN
1613       IF ( ws_scheme_sca )  THEN
1614          CALL advec_s_ws( cs_scalar, 'kc' )
1615       ELSE
1616          CALL advec_s_pw( cs_scalar )
1617       ENDIF
1618    ELSE
1619         CALL advec_s_up( cs_scalar )
1620    ENDIF
1621!
1622!-- Diffusion terms  (the last three arguments are zero)
1623    CALL diffusion_s( cs_scalar,                                               &
1624                      surf_def_h(0)%cssws(ilsp,:),                             &
1625                      surf_def_h(1)%cssws(ilsp,:),                             &
1626                      surf_def_h(2)%cssws(ilsp,:),                             &
1627                      surf_lsm_h%cssws(ilsp,:),                                &
1628                      surf_usm_h%cssws(ilsp,:),                                &
1629                      surf_def_v(0)%cssws(ilsp,:),                             &
1630                      surf_def_v(1)%cssws(ilsp,:),                             &
1631                      surf_def_v(2)%cssws(ilsp,:),                             &
1632                      surf_def_v(3)%cssws(ilsp,:),                             &
1633                      surf_lsm_v(0)%cssws(ilsp,:),                             &
1634                      surf_lsm_v(1)%cssws(ilsp,:),                             &
1635                      surf_lsm_v(2)%cssws(ilsp,:),                             &
1636                      surf_lsm_v(3)%cssws(ilsp,:),                             &
1637                      surf_usm_v(0)%cssws(ilsp,:),                             &
1638                      surf_usm_v(1)%cssws(ilsp,:),                             &
1639                      surf_usm_v(2)%cssws(ilsp,:),                             &
1640                      surf_usm_v(3)%cssws(ilsp,:) )
1641!   
1642!-- Prognostic equation for chemical species
1643    DO  i = nxl, nxr
1644       DO  j = nys, nyn   
1645          DO  k = nzb+1, nzt
1646             cs_scalar_p(k,j,i) =   cs_scalar(k,j,i)                           &
1647                                  + ( dt_3d  *                                 &
1648                                      (   tsc(2) * tend(k,j,i)                 &
1649                                        + tsc(3) * tcs_scalar_m(k,j,i)         &
1650                                      )                                        & 
1651                                    - tsc(5) * rdf_sc(k)                       &
1652                                             * ( cs_scalar(k,j,i) - pr_init_cs(k) )    &   
1653                                    )                                          &
1654                                    * MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) )       
1655
1656             IF ( cs_scalar_p(k,j,i) < 0.0_wp )  cs_scalar_p(k,j,i) = 0.1_wp * cs_scalar(k,j,i)
1657          ENDDO
1658       ENDDO
1659    ENDDO
1660!
1661!-- Calculate tendencies for the next Runge-Kutta step
1662    IF ( timestep_scheme(1:5) == 'runge' )  THEN
1663       IF ( intermediate_timestep_count == 1 )  THEN
1664          DO  i = nxl, nxr
1665             DO  j = nys, nyn   
1666                DO  k = nzb+1, nzt
1667                   tcs_scalar_m(k,j,i) = tend(k,j,i)
1668                ENDDO
1669             ENDDO
1670          ENDDO
1671       ELSEIF ( intermediate_timestep_count < &
1672          intermediate_timestep_count_max )  THEN
1673          DO  i = nxl, nxr
1674             DO  j = nys, nyn
1675                DO  k = nzb+1, nzt
1676                   tcs_scalar_m(k,j,i) = - 9.5625_wp * tend(k,j,i)             &
1677                                         + 5.3125_wp * tcs_scalar_m(k,j,i)
1678                ENDDO
1679             ENDDO
1680          ENDDO
1681       ENDIF
1682    ENDIF
1683
1684 END SUBROUTINE chem_prognostic_equations
1685
1686
1687!------------------------------------------------------------------------------!
1688!
1689! Description:
1690! ------------
1691!> Subroutine defining header output for chemistry model
1692!------------------------------------------------------------------------------!
1693 SUBROUTINE chem_header ( io )
1694       
1695    IMPLICIT NONE
1696 
1697       INTEGER(iwp), INTENT(IN) ::  io            !< Unit of the output file
1698
1699!      print*,'the header subroutine is still not operational'     
1700!!
1701!!--    Write chemistry model header
1702!       WRITE( io, 3 )
1703!
1704!       IF ( radiation_scheme == "constant" )  THEN
1705!          WRITE( io, 4 ) net_radiation
1706!       ELSEIF ( radiation_scheme == "clear-sky" )  THEN
1707!          WRITE( io, 5 )
1708!       ELSEIF ( radiation_scheme == "rrtmg" )  THEN
1709!          WRITE( io, 6 )
1710!          IF ( .NOT. lw_radiation )  WRITE( io, 10 )
1711!          IF ( .NOT. sw_radiation )  WRITE( io, 11 )
1712!       ENDIF
1713!
1714!       IF ( albedo_type == 0 )  THEN
1715!          WRITE( io, 7 ) albedo
1716!       ELSE
1717!          WRITE( io, 8 ) TRIM( albedo_type_name(albedo_type) )
1718!       ENDIF
1719!       IF ( constant_albedo )  THEN
1720!          WRITE( io, 9 )
1721!       ENDIF
1722!     
1723!       IF ( radiation .AND. radiation_scheme /= 'constant' )  THEN
1724!          WRITE ( io, 1 )  lambda
1725!          WRITE ( io, 2 )  day_init, time_utc_init
1726!       ENDIF
1727!
1728!       WRITE( io, 12 ) dt_radiation
1729!
1730! 1 FORMAT ('    Geograph. longitude            :   lambda = ',F4.1,' degr')
1731! 2 FORMAT ('    Day of the year at model start :   day_init = ',I3   &
1732!            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
1733! 3 FORMAT (//' Radiation model information:'/                                  &
1734!              ' ----------------------------'/)
1735! 4 FORMAT ('    --> Using constant net radiation: net_radiation = ', F6.2,        &
1736!           // 'W/m**2')
1737! 5 FORMAT ('    --> Simple radiation scheme for clear sky is used (no clouds,',   &
1738!                   ' default)')
1739! 6 FORMAT ('    --> RRTMG scheme is used')
1740! 7 FORMAT (/'    User-specific surface albedo: albedo =', F6.3)
1741! 8 FORMAT (/'    Albedo is set for land surface type: ', A)
1742! 9 FORMAT (/'    --> Albedo is fixed during the run')
1743!10 FORMAT (/'    --> Longwave radiation is disabled')
1744!11 FORMAT (/'    --> Shortwave radiation is disabled.')
1745!12 FORMAT  ('    Timestep: dt_radiation = ', F6.2, '  s')
1746!
1747!
1748 END SUBROUTINE chem_header
1749
1750!------------------------------------------------------------------------------
1751! Description:
1752! ------------
1753!> Subroutine reading restart data for chemistry model input parameters
1754!  (FK: To make restarts work, I had to comment this routine. We actually
1755!       don't need it, since the namelist parameters are always read in,
1756!       also in case of a restart run)
1757!------------------------------------------------------------------------------
1758!  SUBROUTINE chem_rrd_global
1759!
1760!    USE chem_modules
1761!
1762!    USE control_parameters,                                                   &
1763!        ONLY: length, message_string, restart_string
1764!
1765!     
1766!     IMPLICIT NONE
1767!     
1768!     
1769!     
1770!     DO
1771!
1772!        SELECT CASE ( restart_string(1:length) )
1773!       
1774!           CASE ( 'bc_cs_b' )
1775!              READ ( 13 )  bc_cs_b
1776!
1777!           CASE DEFAULT
1778!
1779!              EXIT
1780!           
1781!        END SELECT
1782!       
1783!!
1784!!--     Read next string and its length
1785!        READ ( 13 )  length
1786!        READ ( 13 )  restart_string(1:length)
1787!       
1788!     ENDDO
1789!     
1790!  END SUBROUTINE chem_rrd_global
1791
1792
1793!------------------------------------------------------------------------------!
1794!
1795! Description:
1796! ------------
1797!> Subroutine writing restart data for chemistry model input parameters
1798!  (FK: To make restarts work, I had to comment this routine. We actually
1799!       don't need it, since the namelist parameters are always read in,
1800!       also in case of a restart run)
1801!------------------------------------------------------------------------------!
1802!  SUBROUTINE chem_wrd_global
1803!
1804!     USE chem_modules
1805!     
1806!     USE kinds
1807!
1808!
1809!     IMPLICIT NONE
1810!
1811!     INTEGER(iwp) ::  lsp  !< running index for chem spcs
1812!
1813! !
1814! !-- Writing out input parameters that are not part of chemistry_parameters
1815! !-- namelist (namelist parameters are anyway read in again in case of restart)
1816!     DO lsp = 1, nvar
1817!        CALL wrd_write_string( 'conc_pr_init_'//chem_species(lsp)%name )
1818!        WRITE ( 14 )  chem_species(lsp)%conc_pr_init
1819!     ENDDO
1820!
1821!
1822!  END SUBROUTINE chem_wrd_global
1823
1824
1825!------------------------------------------------------------------------------!
1826!
1827! Description:
1828! ------------
1829!> Subroutine for emission
1830!------------------------------------------------------------------------------!
1831 SUBROUTINE chem_emissions
1832
1833    USE chem_modules
1834   
1835    USE netcdf_data_input_mod,                                                 &
1836        ONLY:  street_type_f
1837   
1838    USE surface_mod,                                                           &
1839        ONLY:  surf_lsm_h
1840
1841
1842    IMPLICIT NONE
1843
1844    INTEGER(iwp) ::  i    !< running index for grid in x-direction
1845    INTEGER(iwp) ::  j    !< running index for grid in y-direction
1846    INTEGER(iwp) ::  m    !< running index for horizontal surfaces
1847    INTEGER(iwp) ::  lsp  !< running index for chem spcs
1848
1849!
1850!-- Comment??? (todo)
1851    IF ( street_type_f%from_file )  THEN
1852!
1853!--    Streets are lsm surfaces, hence, no usm surface treatment required
1854       DO  m = 1, surf_lsm_h%ns
1855          i = surf_lsm_h%i(m)
1856          j = surf_lsm_h%j(m)
1857         
1858          IF ( street_type_f%var(j,i) >= main_street_id  .AND.                 &
1859               street_type_f%var(j,i) < max_street_id )  THEN
1860             DO  lsp = 1, nvar
1861                surf_lsm_h%cssws(lsp,m) = emiss_factor_main * surface_csflux(lsp)
1862             ENDDO
1863          ELSEIF ( street_type_f%var(j,i) >= side_street_id  .AND.             &
1864                   street_type_f%var(j,i) < main_street_id )  THEN
1865             DO  lsp = 1, nvar
1866                surf_lsm_h%cssws(lsp,m) = emiss_factor_side * surface_csflux(lsp)
1867             ENDDO
1868          ELSE
1869             surf_lsm_h%cssws(:,m) = 0.0_wp
1870          ENDIF
1871       ENDDO
1872       
1873    ENDIF
1874
1875 END SUBROUTINE chem_emissions
1876
1877
1878 END MODULE chemistry_model_mod
1879
Note: See TracBrowser for help on using the repository browser.