[2562] | 1 | module nucleaco2_mod |
---|
| 2 | implicit none |
---|
| 3 | contains |
---|
| 4 | subroutine nucleaco2(pco2,temp,sat,n_ccn,nucrate,vo2co2, teta) |
---|
| 5 | USE comcstfi_h |
---|
| 6 | |
---|
| 7 | implicit none |
---|
| 8 | !* * |
---|
| 9 | !* This subroutine computes the nucleation rate * |
---|
| 10 | !* as given in Pruppacher & Klett (1978) in the * |
---|
| 11 | !* case of water ice forming on a solid substrate. * |
---|
| 12 | !* Definition refined by Keese (jgr,1989) * |
---|
| 13 | !* Authors: F. Montmessin * |
---|
| 14 | !* Adapted for the LMD/GCM by J.-B. Madeleine * |
---|
| 15 | !* (October 2011) * |
---|
| 16 | !* Optimisation by A. Spiga (February 2012) * |
---|
| 17 | !* CO2 nucleation routine dev. by Constantino * |
---|
| 18 | !* Listowski and Joachim Audouard (2016-2017), * |
---|
| 19 | !* adapted from the water ice nucleation |
---|
| 20 | !* It computes two different nucleation rates : one |
---|
| 21 | !* on the dust CCN distribution and the other one on |
---|
| 22 | !* the water ice particles distribution |
---|
| 23 | !******************************************************* |
---|
| 24 | ! nucrate = output |
---|
| 25 | ! nucrate_h2o en sortie aussi : |
---|
| 26 | !nucleation sur dust et h2o separement ici |
---|
| 27 | |
---|
| 28 | include "microphys.h" |
---|
| 29 | include "callkeys.h" |
---|
| 30 | |
---|
| 31 | ! Inputs |
---|
| 32 | DOUBLE PRECISION, intent(in) :: pco2,sat,vo2co2 |
---|
| 33 | DOUBLE PRECISION, intent(in) :: n_ccn(nbinco2_cld) |
---|
| 34 | REAL, intent(in) :: temp !temperature |
---|
| 35 | REAL, intent(in) :: teta |
---|
| 36 | ! Output |
---|
| 37 | DOUBLE PRECISION, intent(out) :: nucrate(nbinco2_cld) |
---|
| 38 | |
---|
| 39 | ! Local variables |
---|
| 40 | DOUBLE PRECISION nco2 |
---|
| 41 | DOUBLE PRECISION rstar ! Radius of the critical germ (m) |
---|
| 42 | DOUBLE PRECISION gstar ! # of molecules forming a critical embryo |
---|
| 43 | DOUBLE PRECISION fistar ! Activation energy required to form a critical embryo (J) |
---|
| 44 | !DOUBLE PRECISION fshapeco2 ! function defined at the end of the file |
---|
| 45 | DOUBLE PRECISION deltaf |
---|
| 46 | double precision mtetalocal ! local teta in double precision |
---|
| 47 | double precision fshapeco2simple,zefshapeco2 |
---|
| 48 | integer i |
---|
| 49 | !************************************************* |
---|
| 50 | |
---|
| 51 | mtetalocal = dble(teta) |
---|
| 52 | |
---|
| 53 | nco2 = pco2 / kbz / temp |
---|
| 54 | rstar = 2. * sigco2 * vo2co2 / (kbz*temp*log(sat)) |
---|
| 55 | gstar = 4. * pi * (rstar * rstar * rstar) / (3.*vo2co2) |
---|
| 56 | |
---|
| 57 | fshapeco2simple = (2.+mtetalocal)*(1.-mtetalocal)*(1.-mtetalocal) / 4. |
---|
| 58 | |
---|
| 59 | !c Loop over size bins |
---|
| 60 | do i=1,nbinco2_cld |
---|
| 61 | |
---|
| 62 | if ( n_ccn(i) .lt. 1e-10 ) then |
---|
| 63 | !c no dust, no need to compute nucleation! |
---|
| 64 | nucrate(i)=0. |
---|
| 65 | else |
---|
| 66 | if (rad_cldco2(i).gt.3000.*rstar) then |
---|
| 67 | zefshapeco2 = fshapeco2simple |
---|
| 68 | else |
---|
| 69 | zefshapeco2 = fshapeco2(mtetalocal, rad_cldco2(i)/rstar) |
---|
| 70 | endif |
---|
| 71 | |
---|
| 72 | fistar = (4./3.*pi) * sigco2 * (rstar * rstar) * zefshapeco2 |
---|
| 73 | deltaf = (2.*desorpco2-surfdifco2-fistar) / (kbz*temp) |
---|
| 74 | deltaf = min( max(deltaf, -100.d0), 100.d0) |
---|
| 75 | |
---|
| 76 | if (deltaf.eq.-100.) then |
---|
| 77 | nucrate(i) = 0. |
---|
| 78 | else |
---|
| 79 | nucrate(i) = dble(sqrt ( fistar / (3.*pi*kbz*temp*(gstar*gstar)) ) * kbz * temp * rstar * rstar * 4. * pi & |
---|
| 80 | * ( nco2*rad_cldco2(i) ) * ( nco2*rad_cldco2(i) ) / ( zefshapeco2 * nusco2 * m0co2 ) * exp (deltaf)) |
---|
| 81 | endif |
---|
| 82 | endif ! if n_ccn(i) .lt. 1e-10 |
---|
| 83 | |
---|
| 84 | enddo |
---|
| 85 | |
---|
| 86 | end subroutine nucleaco2 |
---|
| 87 | |
---|
| 88 | |
---|
| 89 | !********************************************************* |
---|
| 90 | double precision function fshapeco2(cost,rap) |
---|
| 91 | implicit none |
---|
| 92 | !* function computing the f(m,x) factor * |
---|
| 93 | !* related to energy required to form a critical embryo * |
---|
| 94 | !********************************************************* |
---|
| 95 | |
---|
| 96 | double precision cost,rap |
---|
| 97 | double precision yeah |
---|
| 98 | |
---|
| 99 | !! PHI |
---|
| 100 | yeah = sqrt( 1. - 2.*cost*rap + rap*rap ) |
---|
| 101 | !! FSHAPECO2 = TERM A |
---|
| 102 | fshapeco2 = (1.-cost*rap) / yeah |
---|
| 103 | fshapeco2 = fshapeco2 * fshapeco2 * fshapeco2 |
---|
| 104 | fshapeco2 = 1. + fshapeco2 |
---|
| 105 | !! ... + TERM B |
---|
| 106 | yeah = (rap-cost)/yeah |
---|
| 107 | fshapeco2 = fshapeco2 + rap*rap*rap*(2.-3.*yeah+yeah*yeah*yeah) |
---|
| 108 | !! ... + TERM C |
---|
| 109 | fshapeco2 = fshapeco2 + 3. * cost * rap * rap * (yeah-1.) |
---|
| 110 | !! FACTOR 1/2 |
---|
| 111 | fshapeco2 = 0.5*fshapeco2 |
---|
| 112 | |
---|
| 113 | end function fshapeco2 |
---|
| 114 | |
---|
| 115 | end module nucleaco2_mod |
---|