Changeset 2773 for palm/trunk/SOURCE/pmc_interface_mod.f90
- Timestamp:
- Jan 30, 2018 2:12:54 PM (7 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
palm/trunk/SOURCE/pmc_interface_mod.f90
r2718 r2773 25 25 ! ----------------- 26 26 ! $Id$ 27 ! - Nesting for chemical species 28 ! - Bugfix in setting boundary condition at downward-facing walls for passive 29 ! scalar 30 ! - Some formatting adjustments 31 ! 32 ! 2718 2018-01-02 08:49:38Z maronga 27 33 ! Corrected "Former revisions" section 28 34 ! … … 220 226 221 227 USE control_parameters, & 222 ONLY: cloud_physics, coupling_char, dt_3d, dz, humidity,&228 ONLY: air_chemistry, cloud_physics, coupling_char, dt_3d, dz, humidity,& 223 229 message_string, microphysics_morrison, microphysics_seifert, & 224 230 nest_bound_l, nest_bound_r, nest_bound_s, nest_bound_n, & … … 226 232 simulated_time, topography, volume_flow 227 233 228 USE cpulog, & 234 USE chem_modules, & 235 ONLY: nspec 236 237 USE chemistry_model_mod, & 238 ONLY: chem_species, spec_conc_2 239 240 USE cpulog, & 229 241 ONLY: cpu_log, log_point_s 230 242 231 USE grid_variables, 243 USE grid_variables, & 232 244 ONLY: dx, dy 233 245 234 USE indices, 235 ONLY: nbgp, nx, nxl, nxlg, nxlu, nxr, nxrg, ny, nyn, nyng, nys, nysg, 246 USE indices, & 247 ONLY: nbgp, nx, nxl, nxlg, nxlu, nxr, nxrg, ny, nyn, nyng, nys, nysg, & 236 248 nysv, nz, nzb, nzt, wall_flags_0 237 249 … … 245 257 #endif 246 258 247 USE pegrid, 248 ONLY: collective_wait, comm1dx, comm1dy, comm2d, myid, myidx, myidy, 259 USE pegrid, & 260 ONLY: collective_wait, comm1dx, comm1dy, comm2d, myid, myidx, myidy, & 249 261 numprocs 250 262 251 USE pmc_child, 252 ONLY: pmc_childinit, pmc_c_clear_next_array_list, 253 pmc_c_getnextarray, pmc_c_get_2d_index_list, pmc_c_getbuffer, 254 pmc_c_putbuffer, pmc_c_setind_and_allocmem, 263 USE pmc_child, & 264 ONLY: pmc_childinit, pmc_c_clear_next_array_list, & 265 pmc_c_getnextarray, pmc_c_get_2d_index_list, pmc_c_getbuffer, & 266 pmc_c_putbuffer, pmc_c_setind_and_allocmem, & 255 267 pmc_c_set_dataarray, pmc_set_dataarray_name 256 268 257 USE pmc_general, 269 USE pmc_general, & 258 270 ONLY: da_namelen 259 271 260 USE pmc_handle_communicator, 261 ONLY: pmc_get_model_info, pmc_init_model, pmc_is_rootmodel, 272 USE pmc_handle_communicator, & 273 ONLY: pmc_get_model_info, pmc_init_model, pmc_is_rootmodel, & 262 274 pmc_no_namelist_found, pmc_parent_for_child 263 275 264 USE pmc_mpi_wrapper, 265 ONLY: pmc_bcast, pmc_recv_from_child, pmc_recv_from_parent, 276 USE pmc_mpi_wrapper, & 277 ONLY: pmc_bcast, pmc_recv_from_child, pmc_recv_from_parent, & 266 278 pmc_send_to_child, pmc_send_to_parent 267 279 268 USE pmc_parent, 269 ONLY: pmc_parentinit, pmc_s_clear_next_array_list, pmc_s_fillbuffer, 270 pmc_s_getdata_from_buffer, pmc_s_getnextarray, 271 pmc_s_setind_and_allocmem, pmc_s_set_active_data_array, 280 USE pmc_parent, & 281 ONLY: pmc_parentinit, pmc_s_clear_next_array_list, pmc_s_fillbuffer, & 282 pmc_s_getdata_from_buffer, pmc_s_getnextarray, & 283 pmc_s_setind_and_allocmem, pmc_s_set_active_data_array, & 272 284 pmc_s_set_dataarray, pmc_s_set_2d_index_list 273 285 274 286 #endif 275 287 276 USE surface_mod, 288 USE surface_mod, & 277 289 ONLY: get_topography_top_index_ji, surf_def_h, surf_lsm_h, surf_usm_h 278 290 … … 338 350 REAL(wp), SAVE, DIMENSION(:,:,:), ALLOCATABLE, TARGET :: ncc !: 339 351 REAL(wp), SAVE, DIMENSION(:,:,:), ALLOCATABLE, TARGET :: sc !: 352 353 REAL(wp), SAVE, DIMENSION(:,:,:,:), ALLOCATABLE, TARGET :: chem_spec_c !< child coarse data array for chemical species 354 340 355 ! 341 356 !-- Child interpolation coefficients and child-array indices to be … … 651 666 INTEGER(iwp) :: m !: 652 667 INTEGER(iwp) :: mm !: 668 INTEGER(iwp) :: n = 1 !< running index for chemical species 653 669 INTEGER(iwp) :: nest_overlap !: 654 670 INTEGER(iwp) :: nomatch !: … … 836 852 CALL pmc_s_clear_next_array_list 837 853 DO WHILE ( pmc_s_getnextarray( child_id, myname ) ) 838 CALL pmci_set_array_pointer( myname, child_id = child_id, & 839 nz_cl = nz_cl ) 854 IF ( INDEX( TRIM( myname ), 'chem_' ) /= 0 ) THEN 855 CALL pmci_set_array_pointer( myname, child_id = child_id, & 856 nz_cl = nz_cl, n = n ) 857 n = n + 1 858 ELSE 859 CALL pmci_set_array_pointer( myname, child_id = child_id, & 860 nz_cl = nz_cl ) 861 ENDIF 840 862 ENDDO 841 863 CALL pmc_s_setind_and_allocmem( child_id ) … … 965 987 SUBROUTINE pmci_setup_child 966 988 989 967 990 #if defined( __parallel ) 968 991 IMPLICIT NONE 969 992 970 CHARACTER(LEN=da_namelen) :: myname !: 971 972 INTEGER(iwp) :: i !: 973 INTEGER(iwp) :: ierr !: 974 INTEGER(iwp) :: icl !: 975 INTEGER(iwp) :: icr !: 976 INTEGER(iwp) :: j !: 977 INTEGER(iwp) :: jcn !: 978 INTEGER(iwp) :: jcs !: 993 CHARACTER(LEN=da_namelen) :: myname !< 994 995 INTEGER(iwp) :: i !< 996 INTEGER(iwp) :: ierr !< 997 INTEGER(iwp) :: icl !< 998 INTEGER(iwp) :: icr !< 999 INTEGER(iwp) :: j !< 1000 INTEGER(iwp) :: jcn !< 1001 INTEGER(iwp) :: jcs !< 1002 INTEGER(iwp) :: n !< running index for number of chemical species 979 1003 980 1004 INTEGER(iwp), DIMENSION(5) :: val !: … … 1028 1052 IF ( passive_scalar ) THEN 1029 1053 CALL pmc_set_dataarray_name( 'coarse', 's' ,'fine', 's', ierr ) 1054 ENDIF 1055 1056 IF ( air_chemistry ) THEN 1057 DO n = 1, nspec 1058 CALL pmc_set_dataarray_name( 'coarse', & 1059 'chem_' // & 1060 TRIM( chem_species(n)%name ), & 1061 'fine', & 1062 'chem_' // & 1063 TRIM( chem_species(n)%name ), & 1064 ierr ) 1065 ENDDO 1030 1066 ENDIF 1031 1067 … … 1116 1152 !-- TO_DO: Klaus: better explain the above comment (what is child content?) 1117 1153 CALL pmc_c_clear_next_array_list 1154 1155 n = 1 1118 1156 DO WHILE ( pmc_c_getnextarray( myname ) ) 1119 1157 !-- Note that cg%nz is not th eoriginal nz of parent, but the highest 1120 !-- parent-grid level needed for nesting. 1121 CALL pmci_create_child_arrays ( myname, icl, icr, jcs, jcn, cg%nz ) 1158 !-- parent-grid level needed for nesting. 1159 !-- Please note, in case of chemical species an additional parameter 1160 !-- need to be passed, which is required to set the pointer correctly 1161 !-- to the chemical-species data structure. Hence, first check if current 1162 !-- variable is a chemical species. If so, pass index id of respective 1163 !-- species and increment this subsequently. 1164 IF ( INDEX( TRIM( myname ), 'chem_' ) /= 0 ) THEN 1165 CALL pmci_create_child_arrays ( myname, icl, icr, jcs, jcn, cg%nz, n ) 1166 n = n + 1 1167 ELSE 1168 CALL pmci_create_child_arrays ( myname, icl, icr, jcs, jcn, cg%nz ) 1169 ENDIF 1122 1170 ENDDO 1123 1171 CALL pmc_c_setind_and_allocmem … … 2785 2833 2786 2834 2787 SUBROUTINE pmci_set_array_pointer( name, child_id, nz_cl )2835 SUBROUTINE pmci_set_array_pointer( name, child_id, nz_cl, n ) 2788 2836 2789 2837 IMPLICIT NONE … … 2791 2839 INTEGER, INTENT(IN) :: child_id !: 2792 2840 INTEGER, INTENT(IN) :: nz_cl !: 2841 INTEGER, INTENT(IN),OPTIONAL :: n !< index of chemical species 2842 2793 2843 CHARACTER(LEN=*), INTENT(IN) :: name !: 2794 2844 … … 2819 2869 IF ( TRIM(name) == "nc" ) p_3d => nc 2820 2870 IF ( TRIM(name) == "s" ) p_3d => s 2871 IF ( INDEX( TRIM(name), "chem_" ) /= 0 ) p_3d => chem_species(n)%conc 2872 2821 2873 ! 2822 2874 !-- Next line is just an example for a 2D array (not active for coupling!) … … 2855 2907 IF ( TRIM(name) == "nc" ) p_3d_sec => nc_2 2856 2908 IF ( TRIM(name) == "s" ) p_3d_sec => s_2 2909 IF ( INDEX( TRIM(name), "chem_" ) /= 0 ) p_3d_sec => spec_conc_2(:,:,:,n) 2857 2910 2858 2911 IF ( ASSOCIATED( p_3d ) ) THEN … … 2883 2936 2884 2937 2885 SUBROUTINE pmci_create_child_arrays( name, is, ie, js, je, nzc 2938 SUBROUTINE pmci_create_child_arrays( name, is, ie, js, je, nzc, n ) 2886 2939 2887 2940 IMPLICIT NONE … … 2895 2948 INTEGER(iwp), INTENT(IN) :: nzc !: Note that nzc is cg%nz 2896 2949 2950 INTEGER(iwp), INTENT(IN), OPTIONAL :: n !< number of chemical species 2951 2897 2952 #if defined( __parallel ) 2898 2953 INTEGER(iwp) :: ierr !: … … 2908 2963 !-- List of array names, which can be coupled 2909 2964 IF ( TRIM( name ) == "u" ) THEN 2910 IF ( .NOT. ALLOCATED( uc ) ) ALLOCATE( uc(0:nzc+1, js:je,is:ie) )2965 IF ( .NOT. ALLOCATED( uc ) ) ALLOCATE( uc(0:nzc+1,js:je,is:ie) ) 2911 2966 p_3d => uc 2912 2967 ELSEIF ( TRIM( name ) == "v" ) THEN 2913 IF ( .NOT. ALLOCATED( vc ) ) ALLOCATE( vc(0:nzc+1, js:je,is:ie) )2968 IF ( .NOT. ALLOCATED( vc ) ) ALLOCATE( vc(0:nzc+1,js:je,is:ie) ) 2914 2969 p_3d => vc 2915 2970 ELSEIF ( TRIM( name ) == "w" ) THEN 2916 IF ( .NOT. ALLOCATED( wc ) ) ALLOCATE( wc(0:nzc+1, js:je,is:ie) )2971 IF ( .NOT. ALLOCATED( wc ) ) ALLOCATE( wc(0:nzc+1,js:je,is:ie) ) 2917 2972 p_3d => wc 2918 2973 ELSEIF ( TRIM( name ) == "e" ) THEN 2919 IF ( .NOT. ALLOCATED( ec ) ) ALLOCATE( ec(0:nzc+1, js:je,is:ie) )2974 IF ( .NOT. ALLOCATED( ec ) ) ALLOCATE( ec(0:nzc+1,js:je,is:ie) ) 2920 2975 p_3d => ec 2921 2976 ELSEIF ( TRIM( name ) == "pt") THEN 2922 IF ( .NOT. ALLOCATED( ptc ) ) ALLOCATE( ptc(0:nzc+1, js:je,is:ie) )2977 IF ( .NOT. ALLOCATED( ptc ) ) ALLOCATE( ptc(0:nzc+1,js:je,is:ie) ) 2923 2978 p_3d => ptc 2924 2979 ELSEIF ( TRIM( name ) == "q") THEN 2925 IF ( .NOT. ALLOCATED( q_c ) ) ALLOCATE( q_c(0:nzc+1, js:je,is:ie) )2980 IF ( .NOT. ALLOCATED( q_c ) ) ALLOCATE( q_c(0:nzc+1,js:je,is:ie) ) 2926 2981 p_3d => q_c 2927 2982 ELSEIF ( TRIM( name ) == "qc") THEN 2928 IF ( .NOT. ALLOCATED( qcc ) ) ALLOCATE( qcc(0:nzc+1, js:je,is:ie) )2983 IF ( .NOT. ALLOCATED( qcc ) ) ALLOCATE( qcc(0:nzc+1,js:je,is:ie) ) 2929 2984 p_3d => qcc 2930 2985 ELSEIF ( TRIM( name ) == "qr") THEN 2931 IF ( .NOT. ALLOCATED( qrc ) ) ALLOCATE( qrc(0:nzc+1, js:je,is:ie) )2986 IF ( .NOT. ALLOCATED( qrc ) ) ALLOCATE( qrc(0:nzc+1,js:je,is:ie) ) 2932 2987 p_3d => qrc 2933 2988 ELSEIF ( TRIM( name ) == "nr") THEN 2934 IF ( .NOT. ALLOCATED( nrc ) ) ALLOCATE( nrc(0:nzc+1, js:je,is:ie) )2989 IF ( .NOT. ALLOCATED( nrc ) ) ALLOCATE( nrc(0:nzc+1,js:je,is:ie) ) 2935 2990 p_3d => nrc 2936 2991 ELSEIF ( TRIM( name ) == "nc") THEN 2937 IF ( .NOT. ALLOCATED( ncc ) ) ALLOCATE( ncc(0:nzc+1, js:je,is:ie) )2992 IF ( .NOT. ALLOCATED( ncc ) ) ALLOCATE( ncc(0:nzc+1,js:je,is:ie) ) 2938 2993 p_3d => ncc 2939 2994 ELSEIF ( TRIM( name ) == "s") THEN 2940 IF ( .NOT. ALLOCATED( sc ) ) ALLOCATE( sc(0:nzc+1, js:je,is:ie) )2995 IF ( .NOT. ALLOCATED( sc ) ) ALLOCATE( sc(0:nzc+1,js:je,is:ie) ) 2941 2996 p_3d => sc 2997 ELSEIF ( TRIM( name(1:5) ) == "chem_" ) THEN 2998 IF ( .NOT. ALLOCATED( chem_spec_c ) ) & 2999 ALLOCATE( chem_spec_c(0:nzc+1,js:je,is:ie,1:nspec) ) 3000 p_3d => chem_spec_c(:,:,:,n) 2942 3001 !ELSEIF (trim(name) == "z0") then 2943 3002 !IF (.not.allocated(z0c)) allocate(z0c(js:je, is:ie)) … … 3001 3060 IMPLICIT NONE 3002 3061 3003 INTEGER(iwp) :: i !: 3004 INTEGER(iwp) :: icl !: 3005 INTEGER(iwp) :: icr !: 3006 INTEGER(iwp) :: j !: 3007 INTEGER(iwp) :: jcn !: 3008 INTEGER(iwp) :: jcs !: 3009 INTEGER(iwp) :: k !: 3010 3011 REAL(wp) :: waittime !: 3062 INTEGER(iwp) :: i !< 3063 INTEGER(iwp) :: icl !< 3064 INTEGER(iwp) :: icr !< 3065 INTEGER(iwp) :: j !< 3066 INTEGER(iwp) :: jcn !< 3067 INTEGER(iwp) :: jcs !< 3068 INTEGER(iwp) :: k !< 3069 INTEGER(iwp) :: n !< running index for chemical species 3070 3071 REAL(wp) :: waittime !< 3012 3072 3013 3073 ! … … 3063 3123 CALL pmci_interp_tril_all ( s, sc, ico, jco, kco, r1xo, r2xo, r1yo, & 3064 3124 r2yo, r1zo, r2zo, 's' ) 3125 ENDIF 3126 3127 IF ( air_chemistry ) THEN 3128 DO n = 1, nspec 3129 CALL pmci_interp_tril_all ( chem_species(n)%conc, & 3130 chem_spec_c(:,:,:,n), & 3131 ico, jco, kco, r1xo, r2xo, r1yo, & 3132 r2yo, r1zo, r2zo, 's' ) 3133 ENDDO 3065 3134 ENDIF 3066 3135 … … 3707 3776 3708 3777 IMPLICIT NONE 3778 3779 INTEGER(iwp) :: n !< running index for number of chemical species 3709 3780 3710 3781 ! … … 3787 3858 ENDIF 3788 3859 3860 IF ( air_chemistry ) THEN 3861 DO n = 1, nspec 3862 CALL pmci_interp_tril_lr( chem_species(n)%conc, & 3863 chem_spec_c(:,:,:,n), & 3864 ico, jco, kco, r1xo, r2xo, & 3865 r1yo, r2yo, r1zo, r2zo, & 3866 logc_u_l, logc_ratio_u_l, & 3867 nzt_topo_nestbc_l, 'l', 's' ) 3868 ENDDO 3869 ENDIF 3870 3789 3871 IF ( TRIM( nesting_mode ) == 'one-way' ) THEN 3790 3872 CALL pmci_extrap_ifoutflow_lr( u, 'l', 'u' ) … … 3818 3900 IF ( passive_scalar ) THEN 3819 3901 CALL pmci_extrap_ifoutflow_lr( s, 'l', 's' ) 3902 ENDIF 3903 3904 IF ( air_chemistry ) THEN 3905 DO n = 1, nspec 3906 CALL pmci_extrap_ifoutflow_lr( chem_species(n)%conc, & 3907 'l', 's' ) 3908 ENDDO 3820 3909 ENDIF 3821 3910 … … 3903 3992 nzt_topo_nestbc_r, 'r', 's' ) 3904 3993 3994 IF ( air_chemistry ) THEN 3995 DO n = 1, nspec 3996 CALL pmci_interp_tril_lr( chem_species(n)%conc, & 3997 chem_spec_c(:,:,:,n), & 3998 ico, jco, kco, r1xo, r2xo, & 3999 r1yo, r2yo, r1zo, r2zo, & 4000 logc_u_r, logc_ratio_u_r, & 4001 nzt_topo_nestbc_r, 'r', 's' ) 4002 ENDDO 4003 ENDIF 4004 3905 4005 ENDIF 3906 4006 … … 3932 4032 IF ( passive_scalar ) THEN 3933 4033 CALL pmci_extrap_ifoutflow_lr( s, 'r', 's' ) 4034 ENDIF 4035 4036 IF ( air_chemistry ) THEN 4037 DO n = 1, nspec 4038 CALL pmci_extrap_ifoutflow_lr( chem_species(n)%conc, & 4039 'r', 's' ) 4040 ENDDO 3934 4041 ENDIF 3935 4042 ENDIF … … 4010 4117 ENDIF 4011 4118 4119 IF ( air_chemistry ) THEN 4120 DO n = 1, nspec 4121 CALL pmci_interp_tril_sn( chem_species(n)%conc, & 4122 chem_spec_c(:,:,:,n), & 4123 ico, jco, kco, r1xo, r2xo, & 4124 r1yo, r2yo, r1zo, r2zo, & 4125 logc_u_s, logc_ratio_u_s, & 4126 nzt_topo_nestbc_s, 's', 's' ) 4127 ENDDO 4128 ENDIF 4129 4012 4130 IF ( TRIM( nesting_mode ) == 'one-way' ) THEN 4013 4131 CALL pmci_extrap_ifoutflow_sn( u, 's', 'u' ) … … 4038 4156 IF ( passive_scalar ) THEN 4039 4157 CALL pmci_extrap_ifoutflow_sn( s, 's', 's' ) 4158 ENDIF 4159 4160 IF ( air_chemistry ) THEN 4161 DO n = 1, nspec 4162 CALL pmci_extrap_ifoutflow_sn( chem_species(n)%conc, & 4163 's', 's' ) 4164 ENDDO 4040 4165 ENDIF 4041 4166 … … 4120 4245 ENDIF 4121 4246 4247 IF ( air_chemistry ) THEN 4248 DO n = 1, nspec 4249 CALL pmci_interp_tril_sn( chem_species(n)%conc, & 4250 chem_spec_c(:,:,:,n), & 4251 ico, jco, kco, r1xo, r2xo, & 4252 r1yo, r2yo, r1zo, r2zo, & 4253 logc_u_n, logc_ratio_u_n, & 4254 nzt_topo_nestbc_n, 'n', 's' ) 4255 ENDDO 4256 ENDIF 4257 4122 4258 IF ( TRIM( nesting_mode ) == 'one-way' ) THEN 4123 4259 CALL pmci_extrap_ifoutflow_sn( u, 'n', 'u' ) … … 4147 4283 IF ( passive_scalar ) THEN 4148 4284 CALL pmci_extrap_ifoutflow_sn( s, 'n', 's' ) 4285 ENDIF 4286 4287 IF ( air_chemistry ) THEN 4288 DO n = 1, nspec 4289 CALL pmci_extrap_ifoutflow_sn( chem_species(n)%conc, & 4290 'n', 's' ) 4291 ENDDO 4149 4292 ENDIF 4150 4293 … … 4203 4346 ENDIF 4204 4347 4348 IF ( air_chemistry ) THEN 4349 DO n = 1, nspec 4350 CALL pmci_interp_tril_t( chem_species(n)%conc, & 4351 chem_spec_c(:,:,:,n), & 4352 ico, jco, kco, r1xo, r2xo, & 4353 r1yo, r2yo, r1zo, r2zo, & 4354 's' ) 4355 ENDDO 4356 ENDIF 4357 4205 4358 IF ( TRIM( nesting_mode ) == 'one-way' ) THEN 4206 4359 … … 4235 4388 ENDIF 4236 4389 4390 IF ( air_chemistry ) THEN 4391 DO n = 1, nspec 4392 CALL pmci_extrap_ifoutflow_t( chem_species(n)%conc, 's' ) 4393 ENDDO 4394 ENDIF 4395 4237 4396 ENDIF 4238 4397 … … 4247 4406 !-- Note that TKE is not anterpolated. 4248 4407 IMPLICIT NONE 4408 4409 INTEGER(iwp) :: n !< running index for number of chemical species 4410 4411 4249 4412 4250 4413 CALL pmci_anterp_tophat( u, uc, kctu, iflu, ifuu, jflo, jfuo, kflo, & … … 4290 4453 CALL pmci_anterp_tophat( s, sc, kctu, iflo, ifuo, jflo, jfuo, kflo, & 4291 4454 kfuo, ijfc_s, kfc_s, 's' ) 4455 ENDIF 4456 4457 IF ( air_chemistry ) THEN 4458 DO n = 1, nspec 4459 CALL pmci_anterp_tophat( chem_species(n)%conc, & 4460 chem_spec_c(:,:,:,n), & 4461 kctu, iflo, ifuo, jflo, jfuo, kflo, & 4462 kfuo, ijfc_s, kfc_s, 's' ) 4463 ENDDO 4292 4464 ENDIF 4293 4465 … … 5201 5373 SUBROUTINE pmci_boundary_conds 5202 5374 5375 USE chem_modules, & 5376 ONLY: ibc_cs_b 5377 5203 5378 USE control_parameters, & 5204 5379 ONLY: ibc_pt_b, ibc_q_b, ibc_s_b, ibc_uv_b … … 5213 5388 INTEGER(iwp) :: k !< Index along z-direction 5214 5389 INTEGER(iwp) :: m !< Running index for surface type 5390 INTEGER(iwp) :: n !< running index for number of chemical species 5215 5391 5216 5392 ! … … 5338 5514 j = bc_h(1)%j(m) 5339 5515 k = bc_h(1)%k(m) 5340 s(k -1,j,i) = s(k,j,i)5516 s(k+1,j,i) = s(k,j,i) 5341 5517 ENDDO 5342 5518 ENDIF 5343 5519 ENDIF 5520 ! 5521 !-- Set Neumann boundary conditions for chemical species 5522 IF ( air_chemistry ) THEN 5523 IF ( ibc_cs_b == 1 ) THEN 5524 DO n = 1, nspec 5525 DO m = 1, bc_h(0)%ns 5526 i = bc_h(0)%i(m) 5527 j = bc_h(0)%j(m) 5528 k = bc_h(0)%k(m) 5529 chem_species(n)%conc(k-1,j,i) = chem_species(n)%conc(k,j,i) 5530 ENDDO 5531 DO m = 1, bc_h(1)%ns 5532 i = bc_h(1)%i(m) 5533 j = bc_h(1)%j(m) 5534 k = bc_h(1)%k(m) 5535 chem_species(n)%conc(k+1,j,i) = chem_species(n)%conc(k,j,i) 5536 ENDDO 5537 ENDDO 5538 ENDIF 5539 ENDIF 5344 5540 5345 5541 END SUBROUTINE pmci_boundary_conds
Note: See TracChangeset
for help on using the changeset viewer.