PartMC 2.1.2
run_sect.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2005-2011 Nicole Riemer and Matthew West
00002 ! Copyright (C) Andreas Bott
00003 ! Licensed under the GNU General Public License version 2 or (at your
00004 ! option) any later version. See the file COPYING for details.
00005 
00006 !> \file
00007 !> The pmc_run_sect module.
00008 
00009 !> 1D sectional simulation.
00010 !!
00011 !! Sectional code based on \c coad1d.f by Andreas Bott
00012 !!     - http://www.meteo.uni-bonn.de/mitarbeiter/ABott/
00013 !!     - Released under the GPL to Nicole Riemer (personal communication)
00014 !!     - A. Bott, A flux method for the numerical solution of the
00015 !!       stochastic collection equation, J. Atmos. Sci. 55, 2284-2293,
00016 !!       1998.
00017 module pmc_run_sect
00018 
00019   use pmc_bin_grid
00020   use pmc_aero_binned
00021   use pmc_util
00022   use pmc_aero_dist
00023   use pmc_env_data
00024   use pmc_env_state
00025   use pmc_aero_data
00026   use pmc_coag_kernel
00027   use pmc_output
00028   use pmc_gas_data
00029   use pmc_gas_state
00030 
00031   !> Options controlling the operation of run_sect().
00032   type run_sect_opt_t
00033      !> Final time (s).
00034     real(kind=dp) :: t_max
00035     !> Timestep for coagulation (s).
00036     real(kind=dp) :: del_t
00037     !> Output interval (0 disables) (s).
00038     real(kind=dp) :: t_output
00039     !> Progress interval (0 disables) (s).
00040     real(kind=dp) :: t_progress
00041     !> Whether to do coagulation.
00042     logical :: do_coagulation
00043     !> Output prefix.
00044      character(len=300) :: prefix
00045     !> Type of coagulation kernel.
00046     integer :: coag_kernel_type
00047      !> UUID of the simulation.
00048      character(len=PMC_UUID_LEN) :: uuid
00049   end type run_sect_opt_t
00050 
00051 contains
00052 
00053 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00054 
00055   !> Run a sectional simulation.
00056   subroutine run_sect(bin_grid, gas_data, aero_data, aero_dist, &
00057        env_data, env_state, run_sect_opt)
00058   
00059     !> Bin grid.
00060     type(bin_grid_t), intent(in) :: bin_grid
00061     !> Gas data.
00062     type(gas_data_t), intent(in) :: gas_data
00063     !> Aerosol data.
00064     type(aero_data_t), intent(in) :: aero_data
00065     !> Aerosol distribution.
00066     type(aero_dist_t), intent(inout) :: aero_dist
00067     !> Environment data.
00068     type(env_data_t), intent(inout) :: env_data
00069     !> Environment state.
00070     type(env_state_t), intent(inout) :: env_state
00071     !> Options.
00072     type(run_sect_opt_t), intent(in) :: run_sect_opt
00073     
00074     real(kind=dp) c(bin_grid%n_bin,bin_grid%n_bin)
00075     integer ima(bin_grid%n_bin,bin_grid%n_bin)
00076     real(kind=dp) g(bin_grid%n_bin), r(bin_grid%n_bin), e(bin_grid%n_bin)
00077     real(kind=dp) k_bin(bin_grid%n_bin,bin_grid%n_bin)
00078     real(kind=dp) ck(bin_grid%n_bin,bin_grid%n_bin)
00079     real(kind=dp) ec(bin_grid%n_bin,bin_grid%n_bin)
00080     real(kind=dp) taug(bin_grid%n_bin), taup(bin_grid%n_bin)
00081     real(kind=dp) taul(bin_grid%n_bin), tauu(bin_grid%n_bin)
00082     real(kind=dp) prod(bin_grid%n_bin), ploss(bin_grid%n_bin)
00083     real(kind=dp) time, last_output_time, last_progress_time
00084     type(env_state_t) :: old_env_state
00085     type(aero_binned_t) :: aero_binned
00086     type(gas_state_t) :: gas_state
00087     
00088     integer i, j, i_time, num_t, i_summary
00089     logical do_output, do_progress
00090   
00091     ! g         : spectral mass distribution (mg/cm^3)
00092     ! e         : droplet mass grid (mg)
00093     ! r         : droplet radius grid (um)
00094     ! log_width : constant grid distance of logarithmic grid 
00095 
00096     if (aero_data%n_spec /= 1) then
00097        call die_msg(844211192, &
00098             'run_sect() can only use one aerosol species')
00099     end if
00100 
00101     ! output data structure
00102     call aero_binned_allocate_size(aero_binned, bin_grid%n_bin, &
00103          aero_data%n_spec)
00104     aero_binned%vol_conc = 0d0
00105     call gas_state_allocate_size(gas_state, gas_data%n_spec)
00106     
00107     ! mass and radius grid
00108     do i = 1,bin_grid%n_bin
00109        r(i) = bin_grid%center_radius(i) * 1d6 ! radius in m to um
00110        e(i) = rad2vol(bin_grid%center_radius(i)) &
00111             * aero_data%density(1) * 1d6 ! vol in m^3 to mass in mg
00112     end do
00113     
00114     ! initial mass distribution
00115     call aero_binned_add_aero_dist(aero_binned, bin_grid, aero_data, &
00116          aero_dist)
00117     
00118     call courant(bin_grid%n_bin, bin_grid%log_width, e, ima, c)
00119     
00120     ! initialize time
00121     last_progress_time = 0d0
00122     time = 0d0
00123     i_summary = 1
00124     
00125     ! precompute kernel values for all pairs of bins
00126     call bin_kernel(bin_grid%n_bin, bin_grid%center_radius, aero_data, &
00127          run_sect_opt%coag_kernel_type, env_state, k_bin)
00128     call smooth_bin_kernel(bin_grid%n_bin, k_bin, ck)
00129     do i = 1,bin_grid%n_bin
00130        do j = 1,bin_grid%n_bin
00131           ck(i,j) = ck(i,j) * 1d6  ! m^3/s to cm^3/s
00132        end do
00133     end do
00134     
00135     ! multiply kernel with constant timestep and logarithmic grid distance
00136     do i = 1,bin_grid%n_bin
00137        do j = 1,bin_grid%n_bin
00138           ck(i,j) = ck(i,j) * run_sect_opt%del_t * bin_grid%log_width
00139        end do
00140     end do
00141     
00142     ! initial output
00143     call check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
00144          last_output_time, do_output)
00145     if (do_output) then
00146        call output_sectional(run_sect_opt%prefix, bin_grid, aero_data, &
00147             aero_binned, gas_data, gas_state, env_state, i_summary, &
00148             time, run_sect_opt%t_output, run_sect_opt%uuid)
00149     end if
00150     
00151     ! main time-stepping loop
00152     num_t = nint(run_sect_opt%t_max / run_sect_opt%del_t)
00153     call env_state_allocate(old_env_state)
00154     do i_time = 1, num_t
00155 
00156        if (run_sect_opt%do_coagulation) then
00157           g = aero_binned%vol_conc(:,1) * aero_data%density(1)
00158           call coad(bin_grid%n_bin, run_sect_opt%del_t, taug, taup, taul, &
00159                tauu, prod, ploss, c, ima, g, r, e, ck, ec)
00160           aero_binned%vol_conc(:,1) = g / aero_data%density(1)
00161           aero_binned%num_conc = aero_binned%vol_conc(:,1) &
00162                / rad2vol(bin_grid%center_radius)
00163        end if
00164 
00165        time = run_sect_opt%t_max * real(i_time, kind=dp) 
00166             / real(num_t, kind=dp)
00167 
00168        call env_state_copy(env_state, old_env_state)
00169        call env_data_update_state(env_data, env_state, time, &
00170             update_rel_humid = .true.)
00171        call env_state_update_gas_state(env_state, run_sect_opt%del_t, &
00172             old_env_state, gas_data, gas_state)
00173        call env_state_update_aero_binned(env_state, run_sect_opt%del_t, &
00174             old_env_state, bin_grid, aero_data, aero_binned)
00175        
00176        ! print output
00177        call check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
00178             last_output_time, do_output)
00179        if (do_output) then
00180           i_summary = i_summary + 1
00181           call output_sectional(run_sect_opt%prefix, bin_grid, aero_data, &
00182                aero_binned, gas_data, gas_state, env_state, i_summary, &
00183                time, run_sect_opt%t_output, run_sect_opt%uuid)
00184        end if
00185        
00186        ! print progress to stdout
00187        call check_event(time, run_sect_opt%del_t, run_sect_opt%t_progress, &
00188             last_progress_time, do_progress)
00189        if (do_progress) then
00190           write(*,'(a6,a8)') 'step', 'time'
00191           write(*,'(i6,f8.1)') i_time, time
00192        end if
00193     end do
00194 
00195     call env_state_deallocate(old_env_state)
00196     call aero_binned_deallocate(aero_binned)
00197     call gas_state_deallocate(gas_state)
00198 
00199   end subroutine run_sect
00200   
00201 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00202 
00203   !> Collision subroutine, exponential approach.
00204   subroutine coad(n_bin, dt, taug, taup, taul, tauu, prod, ploss, &
00205        c, ima, g, r, e, ck, ec)
00206     
00207     integer n_bin
00208     real(kind=dp) dt
00209     real(kind=dp) taug(n_bin)
00210     real(kind=dp) taup(n_bin)
00211     real(kind=dp) taul(n_bin)
00212     real(kind=dp) tauu(n_bin)
00213     real(kind=dp) prod(n_bin)
00214     real(kind=dp) ploss(n_bin)
00215     real(kind=dp) c(n_bin,n_bin)
00216     integer ima(n_bin,n_bin)
00217     real(kind=dp) g(n_bin)
00218     real(kind=dp) r(n_bin)
00219     real(kind=dp) e(n_bin)
00220     real(kind=dp) ck(n_bin,n_bin)
00221     real(kind=dp) ec(n_bin,n_bin)
00222     
00223     real(kind=dp), parameter :: gmin = 1d-60
00224     
00225     integer i, i0, i1, j, k, kp
00226     real(kind=dp) x0, gsi, gsj, gsk, gk, x1, flux
00227 
00228     do i = 1,n_bin
00229        prod(i) = 0d0
00230        ploss(i) = 0d0
00231     end do
00232     
00233     ! lower and upper integration limit i0,i1
00234     do i0 = 1,(n_bin - 1)
00235        if (g(i0) .gt. gmin) goto 2000
00236     end do
00237 2000 continue
00238     do i1 = (n_bin - 1),1,-1
00239        if (g(i1) .gt. gmin) goto 2010
00240     end do
00241 2010 continue
00242     
00243     do i = i0,i1
00244        do j = i,i1
00245           k = ima(i,j) ! k = 0 means that i + j goes nowhere
00246           kp = k + 1
00247           
00248           x0 = ck(i,j) * g(i) * g(j)
00249           x0 = min(x0, g(i) * e(j))
00250           
00251           if (j .ne. k) x0 = min(x0, g(j) * e(i))
00252           gsi = x0 / e(j)
00253           gsj = x0 / e(i)
00254           gsk = gsi + gsj
00255              
00256           ! loss from positions i, j
00257           ploss(i) = ploss(i) + gsi
00258           ploss(j) = ploss(j) + gsj
00259           g(i) = g(i) - gsi
00260           g(j) = g(j) - gsj
00261 
00262           if (k > 0) then ! do we have a valid bin for the coagulation result?
00263              gk = g(k) + gsk
00264              
00265              if (gk .gt. gmin) then
00266                 x1 = log(g(kp) / gk + 1d-60)
00267                 flux = gsk / x1 * (exp(0.5d0 * x1) &
00268                      - exp(x1 * (0.5d0 - c(i,j))))
00269                 flux = min(flux, gk)
00270                 g(k) = gk - flux
00271                 g(kp) = g(kp) + flux
00272                 ! gain for positions i, j
00273                 prod(k) =  prod(k) + gsk - flux           
00274                 prod(kp) = prod(kp) + flux
00275              end if
00276           end if
00277        end do
00278     end do
00279     
00280   end subroutine coad
00281   
00282 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00283 
00284   !> Determines the Courant number for each bin pair.
00285   subroutine courant(n_bin, log_width, e, ima, c)
00286 
00287     !> Number of bins.
00288     integer, intent(in) :: n_bin
00289     !> Bin scale factor.
00290     real(kind=dp), intent(in) :: log_width
00291     !> Droplet mass grid (mg).
00292     real(kind=dp), intent(in) :: e(n_bin)
00293     !> i + j goes in bin ima(i,j).
00294     integer, intent(out) :: ima(n_bin,n_bin)
00295     !> Courant number for bin pairs.
00296     real(kind=dp), intent(out) :: c(n_bin,n_bin)
00297     
00298     integer i, j, k, kk
00299     real(kind=dp) x0
00300 
00301     c = 0d0 ! added to avoid uninitialized access errors
00302     ima = 0 ! ima(i,j) = 0 means that particles i + j go nowhere
00303     do i = 1,n_bin
00304        do j = i,n_bin
00305           x0 = e(i) + e(j)
00306           ! this is basically the same as particle_in_bin(), but that
00307           ! is actually slightly different than what was always done
00308           ! here
00309           k = find_1d(n_bin, e, x0)
00310           if (k < n_bin) then
00311              k = k + 1
00312              if (c(i,j) .lt. 1d0 - 1d-08) then
00313                 kk = k - 1
00314                 c(i,j) = log(x0 / e(k-1)) / (3d0 * log_width)
00315              else
00316                 c(i,j) = 0d0
00317                 kk = k
00318              end if
00319              ima(i,j) = min(n_bin - 1, kk)
00320           end if
00321           c(j,i) = c(i,j)
00322           ima(j,i) = ima(i,j)
00323        end do
00324     end do
00325     
00326   end subroutine courant
00327   
00328 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00329 
00330   !> Smooths kernel values for bin pairs, and halves the self-rate.
00331   subroutine smooth_bin_kernel(n_bin, k, k_smooth)
00332     
00333     !> Number of bins.
00334     integer, intent(in) :: n_bin
00335     !> Kernel values.
00336     real(kind=dp), intent(in) :: k(n_bin,n_bin)
00337     !> Smoothed kernel values.
00338     real(kind=dp), intent(out) :: k_smooth(n_bin,n_bin)
00339     
00340     integer i, j, im, ip, jm, jp
00341     
00342     do i = 1,n_bin
00343        do j = 1,n_bin
00344           jm = max0(j - 1, 1)
00345           im = max0(i - 1, 1)
00346           jp = min0(j + 1, n_bin)
00347           ip = min0(i + 1, n_bin)
00348           k_smooth(i,j) = 0.125d0 * (k(i,jm) + k(im,j) &
00349                + k(ip,j) + k(i,jp)) &
00350                + 0.5d0 * k(i,j)
00351           if (i .eq. j) then
00352              k_smooth(i,j) = 0.5d0 * k_smooth(i,j)
00353           end if
00354        end do
00355     end do
00356     
00357   end subroutine smooth_bin_kernel
00358   
00359 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00360   
00361 end module pmc_run_sect