PartMC  2.2.1
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     call check_time_multiple("t_max", run_sect_opt%t_max, &
00092          "del_t", run_sect_opt%del_t)
00093     call check_time_multiple("t_output", run_sect_opt%t_output, &
00094          "del_t", run_sect_opt%del_t)
00095     call check_time_multiple("t_progress", run_sect_opt%t_progress, &
00096          "del_t", run_sect_opt%del_t)
00097 
00098     ! g         : spectral mass distribution (mg/cm^3)
00099     ! e         : droplet mass grid (mg)
00100     ! r         : droplet radius grid (um)
00101     ! log_width : constant grid distance of logarithmic grid 
00102 
00103     if (aero_data%n_spec /= 1) then
00104        call die_msg(844211192, &
00105             'run_sect() can only use one aerosol species')
00106     end if
00107 
00108     ! output data structure
00109     call aero_binned_allocate_size(aero_binned, bin_grid%n_bin, &
00110          aero_data%n_spec)
00111     aero_binned%vol_conc = 0d0
00112     call gas_state_allocate_size(gas_state, gas_data%n_spec)
00113     
00114     ! mass and radius grid
00115     do i = 1,bin_grid%n_bin
00116        r(i) = bin_grid%center_radius(i) * 1d6 ! radius in m to um
00117        e(i) = rad2vol(bin_grid%center_radius(i)) &
00118             * aero_data%density(1) * 1d6 ! vol in m^3 to mass in mg
00119     end do
00120     
00121     ! initial mass distribution
00122     call aero_binned_add_aero_dist(aero_binned, bin_grid, aero_data, &
00123          aero_dist)
00124     
00125     call courant(bin_grid%n_bin, bin_grid%log_width, e, ima, c)
00126     
00127     ! initialize time
00128     last_progress_time = 0d0
00129     time = 0d0
00130     i_summary = 1
00131     
00132     ! precompute kernel values for all pairs of bins
00133     call bin_kernel(bin_grid%n_bin, bin_grid%center_radius, aero_data, &
00134          run_sect_opt%coag_kernel_type, env_state, k_bin)
00135     call smooth_bin_kernel(bin_grid%n_bin, k_bin, ck)
00136     do i = 1,bin_grid%n_bin
00137        do j = 1,bin_grid%n_bin
00138           ck(i,j) = ck(i,j) * 1d6  ! m^3/s to cm^3/s
00139        end do
00140     end do
00141     
00142     ! multiply kernel with constant timestep and logarithmic grid distance
00143     do i = 1,bin_grid%n_bin
00144        do j = 1,bin_grid%n_bin
00145           ck(i,j) = ck(i,j) * run_sect_opt%del_t * bin_grid%log_width
00146        end do
00147     end do
00148     
00149     ! initial output
00150     call check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
00151          last_output_time, do_output)
00152     if (do_output) then
00153        call output_sectional(run_sect_opt%prefix, bin_grid, aero_data, &
00154             aero_binned, gas_data, gas_state, env_state, i_summary, &
00155             time, run_sect_opt%t_output, run_sect_opt%uuid)
00156     end if
00157     
00158     ! main time-stepping loop
00159     num_t = nint(run_sect_opt%t_max / run_sect_opt%del_t)
00160     call env_state_allocate(old_env_state)
00161     do i_time = 1, num_t
00162 
00163        if (run_sect_opt%do_coagulation) then
00164           g = aero_binned%vol_conc(:,1) * aero_data%density(1)
00165           call coad(bin_grid%n_bin, run_sect_opt%del_t, taug, taup, taul, &
00166                tauu, prod, ploss, c, ima, g, r, e, ck, ec)
00167           aero_binned%vol_conc(:,1) = g / aero_data%density(1)
00168           aero_binned%num_conc = aero_binned%vol_conc(:,1) &
00169                / rad2vol(bin_grid%center_radius)
00170        end if
00171 
00172        time = run_sect_opt%t_max * real(i_time, kind=dp) 
00173             / real(num_t, kind=dp)
00174 
00175        call env_state_copy(env_state, old_env_state)
00176        call env_data_update_state(env_data, env_state, time, &
00177             update_rel_humid = .true.)
00178        call env_state_update_gas_state(env_state, run_sect_opt%del_t, &
00179             old_env_state, gas_data, gas_state)
00180        call env_state_update_aero_binned(env_state, run_sect_opt%del_t, &
00181             old_env_state, bin_grid, aero_data, aero_binned)
00182        
00183        ! print output
00184        call check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
00185             last_output_time, do_output)
00186        if (do_output) then
00187           i_summary = i_summary + 1
00188           call output_sectional(run_sect_opt%prefix, bin_grid, aero_data, &
00189                aero_binned, gas_data, gas_state, env_state, i_summary, &
00190                time, run_sect_opt%t_output, run_sect_opt%uuid)
00191        end if
00192        
00193        ! print progress to stdout
00194        call check_event(time, run_sect_opt%del_t, run_sect_opt%t_progress, &
00195             last_progress_time, do_progress)
00196        if (do_progress) then
00197           write(*,'(a6,a8)') 'step', 'time'
00198           write(*,'(i6,f8.1)') i_time, time
00199        end if
00200     end do
00201 
00202     call env_state_deallocate(old_env_state)
00203     call aero_binned_deallocate(aero_binned)
00204     call gas_state_deallocate(gas_state)
00205 
00206   end subroutine run_sect
00207   
00208 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00209 
00210   !> Collision subroutine, exponential approach.
00211   subroutine coad(n_bin, dt, taug, taup, taul, tauu, prod, ploss, &
00212        c, ima, g, r, e, ck, ec)
00213     
00214     integer n_bin
00215     real(kind=dp) dt
00216     real(kind=dp) taug(n_bin)
00217     real(kind=dp) taup(n_bin)
00218     real(kind=dp) taul(n_bin)
00219     real(kind=dp) tauu(n_bin)
00220     real(kind=dp) prod(n_bin)
00221     real(kind=dp) ploss(n_bin)
00222     real(kind=dp) c(n_bin,n_bin)
00223     integer ima(n_bin,n_bin)
00224     real(kind=dp) g(n_bin)
00225     real(kind=dp) r(n_bin)
00226     real(kind=dp) e(n_bin)
00227     real(kind=dp) ck(n_bin,n_bin)
00228     real(kind=dp) ec(n_bin,n_bin)
00229     
00230     real(kind=dp), parameter :: gmin = 1d-60
00231     
00232     integer i, i0, i1, j, k, kp
00233     real(kind=dp) x0, gsi, gsj, gsk, gk, x1, flux
00234 
00235     do i = 1,n_bin
00236        prod(i) = 0d0
00237        ploss(i) = 0d0
00238     end do
00239     
00240     ! lower and upper integration limit i0,i1
00241     do i0 = 1,(n_bin - 1)
00242        if (g(i0) .gt. gmin) goto 2000
00243     end do
00244 2000 continue
00245     do i1 = (n_bin - 1),1,-1
00246        if (g(i1) .gt. gmin) goto 2010
00247     end do
00248 2010 continue
00249     
00250     do i = i0,i1
00251        do j = i,i1
00252           k = ima(i,j) ! k = 0 means that i + j goes nowhere
00253           kp = k + 1
00254           
00255           x0 = ck(i,j) * g(i) * g(j)
00256           x0 = min(x0, g(i) * e(j))
00257           
00258           if (j .ne. k) x0 = min(x0, g(j) * e(i))
00259           gsi = x0 / e(j)
00260           gsj = x0 / e(i)
00261           gsk = gsi + gsj
00262              
00263           ! loss from positions i, j
00264           ploss(i) = ploss(i) + gsi
00265           ploss(j) = ploss(j) + gsj
00266           g(i) = g(i) - gsi
00267           g(j) = g(j) - gsj
00268 
00269           if (k > 0) then ! do we have a valid bin for the coagulation result?
00270              gk = g(k) + gsk
00271              
00272              if (gk .gt. gmin) then
00273                 x1 = log(g(kp) / gk + 1d-60)
00274                 flux = gsk / x1 * (exp(0.5d0 * x1) &
00275                      - exp(x1 * (0.5d0 - c(i,j))))
00276                 flux = min(flux, gk)
00277                 g(k) = gk - flux
00278                 g(kp) = g(kp) + flux
00279                 ! gain for positions i, j
00280                 prod(k) =  prod(k) + gsk - flux           
00281                 prod(kp) = prod(kp) + flux
00282              end if
00283           end if
00284        end do
00285     end do
00286     
00287   end subroutine coad
00288   
00289 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00290 
00291   !> Determines the Courant number for each bin pair.
00292   subroutine courant(n_bin, log_width, e, ima, c)
00293 
00294     !> Number of bins.
00295     integer, intent(in) :: n_bin
00296     !> Bin scale factor.
00297     real(kind=dp), intent(in) :: log_width
00298     !> Droplet mass grid (mg).
00299     real(kind=dp), intent(in) :: e(n_bin)
00300     !> i + j goes in bin ima(i,j).
00301     integer, intent(out) :: ima(n_bin,n_bin)
00302     !> Courant number for bin pairs.
00303     real(kind=dp), intent(out) :: c(n_bin,n_bin)
00304     
00305     integer i, j, k, kk
00306     real(kind=dp) x0
00307 
00308     c = 0d0 ! added to avoid uninitialized access errors
00309     ima = 0 ! ima(i,j) = 0 means that particles i + j go nowhere
00310     do i = 1,n_bin
00311        do j = i,n_bin
00312           x0 = e(i) + e(j)
00313           ! this is basically the same as particle_in_bin(), but that
00314           ! is actually slightly different than what was always done
00315           ! here
00316 
00317           ! MW 2011-04-28: I think the above comment no longer
00318           ! applies, and we can make this just
00319           ! bin_grid_particle_in_bin(). FIXME.
00320           k = find_1d(n_bin, e, x0)
00321           if (k < n_bin) then
00322              k = k + 1
00323              if (c(i,j) .lt. 1d0 - 1d-08) then
00324                 kk = k - 1
00325                 c(i,j) = log(x0 / e(k-1)) / (3d0 * log_width)
00326              else
00327                 c(i,j) = 0d0
00328                 kk = k
00329              end if
00330              ima(i,j) = min(n_bin - 1, kk)
00331           end if
00332           c(j,i) = c(i,j)
00333           ima(j,i) = ima(i,j)
00334        end do
00335     end do
00336     
00337   end subroutine courant
00338   
00339 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00340 
00341   !> Smooths kernel values for bin pairs, and halves the self-rate.
00342   subroutine smooth_bin_kernel(n_bin, k, k_smooth)
00343     
00344     !> Number of bins.
00345     integer, intent(in) :: n_bin
00346     !> Kernel values.
00347     real(kind=dp), intent(in) :: k(n_bin,n_bin)
00348     !> Smoothed kernel values.
00349     real(kind=dp), intent(out) :: k_smooth(n_bin,n_bin)
00350     
00351     integer i, j, im, ip, jm, jp
00352     
00353     do i = 1,n_bin
00354        do j = 1,n_bin
00355           jm = max0(j - 1, 1)
00356           im = max0(i - 1, 1)
00357           jp = min0(j + 1, n_bin)
00358           ip = min0(i + 1, n_bin)
00359           k_smooth(i,j) = 0.125d0 * (k(i,jm) + k(im,j) &
00360                + k(ip,j) + k(i,jp)) &
00361                + 0.5d0 * k(i,j)
00362           if (i .eq. j) then
00363              k_smooth(i,j) = 0.5d0 * k_smooth(i,j)
00364           end if
00365        end do
00366     end do
00367     
00368   end subroutine smooth_bin_kernel
00369   
00370 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00371   
00372 end module pmc_run_sect