00001 module mrg_x2o_mct
00002
00003 use shr_kind_mod, only: r8 => shr_kind_r8, cl => shr_kind_cl
00004 use mct_mod
00005 use seq_flds_mod
00006 use seq_flds_indices
00007 use seq_comm_mct
00008 use seq_cdata_mod
00009 use seq_infodata_mod
00010
00011 implicit none
00012 save
00013 private
00014
00015
00016
00017
00018
00019 public :: mrg_x2o_init_mct
00020 public :: mrg_x2o_run_mct
00021 public :: mrg_x2o_final_mct
00022
00023
00024
00025
00026
00027
00028 contains
00029
00030
00031 subroutine mrg_x2o_init_mct( cdata_o, a2x_o, i2x_o, r2x_o )
00032
00033
00034 type(seq_cdata), intent(in) :: cdata_o
00035 type(mct_aVect), intent(inout) :: a2x_o
00036 type(mct_aVect), intent(inout) :: i2x_o
00037 type(mct_aVect), intent(inout) :: r2x_o
00038
00039 type(mct_gsMap), pointer :: gsMap_ocn
00040 integer :: mpicom
00041
00042
00043
00044 call seq_cdata_setptrs(cdata_o, gsMap=gsMap_ocn, mpicom=mpicom)
00045
00046
00047
00048 call mct_aVect_init(a2x_o, rList=seq_flds_a2x_fields, &
00049 lsize=MCT_GSMap_lsize(GSMap_ocn, mpicom))
00050 call mct_aVect_zero(a2x_o)
00051
00052
00053
00054 call mct_aVect_init(i2x_o, rList=seq_flds_i2x_fields, &
00055 lsize=mct_GSMap_lsize(GSMap_ocn, mpicom))
00056 call mct_aVect_zero(i2x_o)
00057
00058
00059
00060 call mct_aVect_init(r2x_o, rList=seq_flds_r2x_fields, &
00061 lsize=mct_GSMap_lsize(GSMap_ocn, mpicom))
00062 call mct_aVect_zero(r2x_o)
00063
00064 end subroutine mrg_x2o_init_mct
00065
00066
00067
00068
00069 subroutine mrg_x2o_run_mct( cdata_o, a2x_o, i2x_o, xao_o, fractions_o, x2o_o )
00070
00071
00072
00073
00074
00075
00076 type(seq_cdata), intent(in) :: cdata_o
00077 type(mct_aVect), intent(in) :: a2x_o
00078 type(mct_aVect), intent(in) :: i2x_o
00079
00080 type(mct_aVect), intent(in) :: xao_o
00081 type(mct_aVect), intent(in) :: fractions_o
00082 type(mct_aVect), intent(inout) :: x2o_o
00083
00084
00085
00086 integer :: n, ki, ko, kir, kor
00087 integer :: lsize
00088 real(r8) :: ifrac,ifracr
00089 real(r8) :: afrac,afracr
00090 logical :: usevector
00091 real(r8) :: flux_epbalfact
00092 character(len=cl) :: flux_epbal
00093 type(seq_infodata_type),pointer :: infodata
00094 real(r8) :: frac_sum
00095 real(r8) :: avsdr, anidr, avsdf, anidf
00096 real(r8) :: fswabsv, fswabsi
00097
00098
00099 character(*),parameter :: F01 = "('(mrg_x2o_run_mct) ',a,3e11.3,a,f9.6)"
00100 character(*),parameter :: subName = '(mrg_x2o_run_mct) '
00101
00102
00103
00104
00105 #ifdef CPP_VECTOR
00106 usevector = .true.
00107 #else
00108 usevector = .false.
00109 #endif
00110
00111 call seq_cdata_setptrs(cdata_o, infodata=infodata)
00112 call mct_aVect_zero(x2o_o)
00113
00114 call mct_aVect_copy(aVin=a2x_o, aVout=x2o_o, vector=usevector)
00115 call mct_aVect_copy(aVin=i2x_o, aVout=x2o_o, vector=usevector)
00116
00117
00118 call mct_aVect_copy(aVin=xao_o, aVout=x2o_o, vector=usevector)
00119
00120 call seq_infodata_GetData(infodata,flux_epbal=flux_epbal)
00121 flux_epbalfact = 1.0_r8
00122 if (trim(flux_epbal) == 'ocn') then
00123 call seq_infodata_GetData(infodata, precip_fact = flux_epbalfact)
00124 flux_epbalfact = flux_epbalfact * 1.0e-6_r8
00125 end if
00126 if (flux_epbalfact <= 0.0_R8) then
00127 if (loglevel > 0) write(logunit,F01) 'WARNING: factor from ocn = ',flux_epbalfact
00128 if (loglevel > 0) write(logunit,F01) 'WARNING: resetting flux_epbalfact to 1.0'
00129 flux_epbalfact = 1.0_r8
00130 end if
00131
00132
00133
00134 ki = mct_aVect_indexRa(fractions_o,"ifrac",perrWith=subName)
00135 ko = mct_aVect_indexRa(fractions_o,"ofrac",perrWith=subName)
00136 kir = mct_aVect_indexRa(fractions_o,"ifrad",perrWith=subName)
00137 kor = mct_aVect_indexRa(fractions_o,"ofrad",perrWith=subName)
00138 lsize = mct_aVect_lsize(x2o_o)
00139 do n = 1,lsize
00140
00141 ifrac = fractions_o%rAttr(ki,n)
00142 afrac = fractions_o%rAttr(ko,n)
00143 frac_sum = ifrac + afrac
00144 if ((frac_sum) /= 0._r8) then
00145 ifrac = ifrac / (frac_sum)
00146 afrac = afrac / (frac_sum)
00147 endif
00148
00149 ifracr = fractions_o%rAttr(kir,n)
00150 afracr = fractions_o%rAttr(kor,n)
00151 frac_sum = ifracr + afracr
00152 if ((frac_sum) /= 0._r8) then
00153 ifracr = ifracr / (frac_sum)
00154 afracr = afracr / (frac_sum)
00155 endif
00156
00157 x2o_o%rAttr(index_x2o_Foxx_taux ,n) = xao_o%rAttr(index_xao_Faox_taux ,n) * afrac + &
00158 i2x_o%rAttr(index_i2x_Fioi_taux ,n) * ifrac
00159
00160 x2o_o%rAttr(index_x2o_Foxx_tauy ,n) = xao_o%rAttr(index_xao_Faox_tauy ,n) * afrac + &
00161 i2x_o%rAttr(index_i2x_Fioi_tauy ,n) * ifrac
00162
00163
00164 avsdr = xao_o%rAttr(index_xao_So_avsdr,n)
00165 anidr = xao_o%rAttr(index_xao_So_anidr,n)
00166 avsdf = xao_o%rAttr(index_xao_So_avsdf,n)
00167 anidf = xao_o%rAttr(index_xao_So_anidf,n)
00168 fswabsv = a2x_o%rAttr(index_a2x_Faxa_swvdr,n) * (1.0_R8 - avsdr) &
00169 + a2x_o%rAttr(index_a2x_Faxa_swvdf,n) * (1.0_R8 - avsdf)
00170 fswabsi = a2x_o%rAttr(index_a2x_Faxa_swndr,n) * (1.0_R8 - anidr) &
00171 + a2x_o%rAttr(index_a2x_Faxa_swndf,n) * (1.0_R8 - anidf)
00172
00173 x2o_o%rAttr(index_x2o_Foxx_swnet,n) = (fswabsv + fswabsi) * afracr + &
00174 i2x_o%rAttr(index_i2x_Fioi_swpen,n) * ifrac
00175
00176 x2o_o%rAttr(index_x2o_Foxx_lat ,n) = xao_o%rAttr(index_xao_Faox_lat ,n) * afrac
00177
00178 x2o_o%rAttr(index_x2o_Foxx_sen ,n) = xao_o%rAttr(index_xao_Faox_sen ,n) * afrac
00179
00180 x2o_o%rAttr(index_x2o_Foxx_evap ,n) = xao_o%rAttr(index_xao_Faox_evap ,n) * afrac
00181
00182 x2o_o%rAttr(index_x2o_Foxx_lwup ,n) = xao_o%rAttr(index_xao_Faox_lwup ,n) * afrac
00183
00184 x2o_o%rAttr(index_x2o_Foxx_lwdn ,n) = a2x_o%rAttr(index_a2x_Faxa_lwdn ,n) * afrac
00185
00186 x2o_o%rAttr(index_x2o_Foxx_snow ,n) = a2x_o%rAttr(index_a2x_Faxa_snowc,n) * afrac + &
00187 a2x_o%rAttr(index_a2x_Faxa_snowl,n) * afrac
00188
00189 x2o_o%rAttr(index_x2o_Foxx_rain ,n) = a2x_o%rAttr(index_a2x_Faxa_rainc,n) * afrac + &
00190 a2x_o%rAttr(index_a2x_Faxa_rainl,n) * afrac
00191
00192 x2o_o%rAttr(index_x2o_Foxx_melth,n) = i2x_o%rAttr(index_i2x_Fioi_melth,n) * ifrac
00193
00194 x2o_o%rAttr(index_x2o_Foxx_meltw,n) = i2x_o%rAttr(index_i2x_Fioi_meltw,n) * ifrac
00195
00196 x2o_o%rAttr(index_x2o_Foxx_salt ,n) = i2x_o%rAttr(index_i2x_Fioi_salt ,n) * ifrac
00197
00198
00199
00200
00201 x2o_o%rAttr(index_x2o_Foxx_rain ,n) = x2o_o%rAttr(index_x2o_Foxx_rain ,n) * flux_epbalfact
00202 x2o_o%rAttr(index_x2o_Foxx_snow ,n) = x2o_o%rAttr(index_x2o_Foxx_snow ,n) * flux_epbalfact
00203 x2o_o%rAttr(index_x2o_Forr_roff ,n) = x2o_o%rAttr(index_x2o_Forr_roff ,n) * flux_epbalfact
00204 x2o_o%rAttr(index_x2o_Forr_ioff ,n) = x2o_o%rAttr(index_x2o_Forr_ioff ,n) * flux_epbalfact
00205
00206 x2o_o%rAttr(index_x2o_Foxx_prec ,n) = x2o_o%rAttr(index_x2o_Foxx_rain ,n) + &
00207 x2o_o%rAttr(index_x2o_Foxx_snow ,n)
00208
00209 x2o_o%rAttr(index_x2o_Foxx_bcphidry,n) = a2x_o%rAttr(index_a2x_Faxa_bcphidry,n) * afrac
00210 x2o_o%rAttr(index_x2o_Foxx_bcphodry,n) = a2x_o%rAttr(index_a2x_Faxa_bcphodry,n) * afrac
00211 x2o_o%rAttr(index_x2o_Foxx_bcphiwet,n) = a2x_o%rAttr(index_a2x_Faxa_bcphiwet,n) * afrac
00212 x2o_o%rAttr(index_x2o_Foxx_ocphidry,n) = a2x_o%rAttr(index_a2x_Faxa_ocphidry,n) * afrac
00213 x2o_o%rAttr(index_x2o_Foxx_ocphodry,n) = a2x_o%rAttr(index_a2x_Faxa_ocphodry,n) * afrac
00214 x2o_o%rAttr(index_x2o_Foxx_ocphiwet,n) = a2x_o%rAttr(index_a2x_Faxa_ocphiwet,n) * afrac
00215 x2o_o%rAttr(index_x2o_Foxx_dstwet1,n) = a2x_o%rAttr(index_a2x_Faxa_dstwet1,n) * afrac
00216 x2o_o%rAttr(index_x2o_Foxx_dstwet2,n) = a2x_o%rAttr(index_a2x_Faxa_dstwet2,n) * afrac
00217 x2o_o%rAttr(index_x2o_Foxx_dstwet3,n) = a2x_o%rAttr(index_a2x_Faxa_dstwet3,n) * afrac
00218 x2o_o%rAttr(index_x2o_Foxx_dstwet4,n) = a2x_o%rAttr(index_a2x_Faxa_dstwet4,n) * afrac
00219 x2o_o%rAttr(index_x2o_Foxx_dstdry1,n) = a2x_o%rAttr(index_a2x_Faxa_dstdry1,n) * afrac
00220 x2o_o%rAttr(index_x2o_Foxx_dstdry2,n) = a2x_o%rAttr(index_a2x_Faxa_dstdry2,n) * afrac
00221 x2o_o%rAttr(index_x2o_Foxx_dstdry3,n) = a2x_o%rAttr(index_a2x_Faxa_dstdry3,n) * afrac
00222 x2o_o%rAttr(index_x2o_Foxx_dstdry4,n) = a2x_o%rAttr(index_a2x_Faxa_dstdry4,n) * afrac
00223
00224 end do
00225 end subroutine mrg_x2o_run_mct
00226
00227
00228
00229 subroutine mrg_x2o_final_mct
00230
00231
00232
00233 end subroutine mrg_x2o_final_mct
00234
00235 end module mrg_x2o_mct
00236
00237