00001 module mrg_x2a_mct
00002
00003 use shr_kind_mod, only: r8 => shr_kind_r8
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
00010 implicit none
00011 save
00012 private
00013
00014
00015
00016
00017
00018 public :: mrg_x2a_init_mct
00019 public :: mrg_x2a_run_mct
00020 public :: mrg_x2a_final_mct
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 contains
00032
00033
00034 subroutine mrg_x2a_init_mct( cdata_a, l2x_a, o2x_a, i2x_a, xao_a )
00035
00036
00037
00038
00039
00040 type(seq_cdata), intent(in) :: cdata_a
00041 type(mct_aVect), intent(inout) :: l2x_a
00042 type(mct_aVect), intent(inout) :: o2x_a
00043 type(mct_aVect), intent(inout) :: i2x_a
00044 type(mct_aVect), intent(inout) :: xao_a
00045
00046
00047
00048 type(mct_gsmap), pointer :: gsMap_a
00049 integer :: lsize
00050 integer :: mpicom
00051
00052
00053
00054
00055 call seq_cdata_setptrs(cdata_a, gsMap=gsMap_a, mpicom=mpicom)
00056 lsize = mct_GSMap_lsize(GSMap_a, mpicom)
00057
00058
00059
00060 call mct_aVect_init(l2x_a, rList=seq_flds_l2x_fields, lsize=lsize)
00061 call mct_aVect_zero(l2x_a)
00062
00063
00064
00065 call mct_aVect_init(o2x_a, rList=seq_flds_o2x_fields, lsize=lsize)
00066 call mct_aVect_zero(o2x_a)
00067
00068
00069
00070 call mct_aVect_init(i2x_a, rList=seq_flds_i2x_fields, lsize=lsize)
00071 call mct_aVect_zero(i2x_a)
00072
00073
00074
00075 call mct_aVect_init(xao_a, rList=seq_flds_xao_fields, lsize=lsize)
00076 call mct_aVect_zero(xao_a)
00077
00078 end subroutine mrg_x2a_init_mct
00079
00080
00081
00082 subroutine mrg_x2a_run_mct( cdata_a, l2x_a, o2x_a, xao_a, i2x_a, fractions_a, x2a_a )
00083
00084
00085
00086
00087
00088 type(seq_cdata), intent(in) :: cdata_a
00089 type(mct_aVect), intent(in) :: l2x_a
00090 type(mct_aVect), intent(in) :: o2x_a
00091 type(mct_aVect), intent(in) :: xao_a
00092 type(mct_aVect), intent(in) :: i2x_a
00093 type(mct_aVect), intent(in) :: fractions_a
00094 type(mct_aVect), intent(inout) :: x2a_a
00095
00096
00097
00098
00099 logical :: usevector
00100 integer :: n, ki, kl, ko
00101 real(r8) :: frac
00102 integer :: lsize
00103
00104
00105
00106
00107 call mct_avect_zero(x2a_a)
00108
00109
00110
00111
00112 ki=mct_aVect_indexRA(fractions_a,"ifrac")
00113 kl=mct_aVect_indexRA(fractions_a,"lfrac")
00114 ko=mct_aVect_indexRA(fractions_a,"ofrac")
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126 #ifdef CPP_VECTOR
00127 usevector = .true.
00128 #else
00129 usevector = .false.
00130 #endif
00131 call mct_aVect_copy(aVin=l2x_a, aVout=x2a_a, vector=usevector)
00132 call mct_aVect_copy(aVin=o2x_a, aVout=x2a_a, vector=usevector)
00133 call mct_aVect_copy(aVin=i2x_a, aVout=x2a_a, vector=usevector)
00134 call mct_aVect_copy(aVin=xao_a, aVout=x2a_a, vector=usevector)
00135
00136
00137
00138 lsize = mct_avect_lsize(x2a_a)
00139 do n = 1,lsize
00140
00141 x2a_a%rAttr(index_x2a_Sx_lfrac,n) = fractions_a%Rattr(kl,n)
00142 x2a_a%rAttr(index_x2a_Sx_ifrac,n) = fractions_a%Rattr(ki,n)
00143 x2a_a%rAttr(index_x2a_Sx_ofrac,n) = fractions_a%Rattr(ko,n)
00144
00145 frac = x2a_a%rAttr(index_x2a_Sx_lfrac,n)
00146 if (frac > 0._r8) then
00147 x2a_a%rAttr(index_x2a_Sx_avsdr,n) = x2a_a%rAttr(index_x2a_Sx_avsdr,n) + &
00148 l2x_a%rAttr(index_l2x_Sl_avsdr,n) * frac
00149 x2a_a%rAttr(index_x2a_Sx_avsdf,n) = x2a_a%rAttr(index_x2a_Sx_avsdf,n) + &
00150 l2x_a%rAttr(index_l2x_Sl_avsdf,n) * frac
00151 x2a_a%rAttr(index_x2a_Sx_anidr,n) = x2a_a%rAttr(index_x2a_Sx_anidr,n) + &
00152 l2x_a%rAttr(index_l2x_Sl_anidr,n) * frac
00153 x2a_a%rAttr(index_x2a_Sx_anidf,n) = x2a_a%rAttr(index_x2a_Sx_anidf,n) + &
00154 l2x_a%rAttr(index_l2x_Sl_anidf,n) * frac
00155 x2a_a%rAttr(index_x2a_Sx_tref,n) = x2a_a%rAttr(index_x2a_Sx_tref,n) + &
00156 l2x_a%rAttr(index_l2x_Sl_tref,n) * frac
00157 x2a_a%rAttr(index_x2a_Sx_qref,n) = x2a_a%rAttr(index_x2a_Sx_qref,n) + &
00158 l2x_a%rAttr(index_l2x_Sl_qref,n) * frac
00159 x2a_a%rAttr(index_x2a_Faxx_evap,n) = x2a_a%rAttr(index_x2a_Faxx_evap,n) + &
00160 l2x_a%rAttr(index_l2x_Fall_evap,n) * frac
00161 x2a_a%rAttr(index_x2a_Faxx_lwup,n) = x2a_a%rAttr(index_x2a_Faxx_lwup,n) + &
00162 l2x_a%rAttr(index_l2x_Fall_lwup,n) * frac
00163 x2a_a%rAttr(index_x2a_Faxx_sen,n) = x2a_a%rAttr(index_x2a_Faxx_sen,n) + &
00164 l2x_a%rAttr(index_l2x_Fall_sen,n) * frac
00165 x2a_a%rAttr(index_x2a_Faxx_lat,n) = x2a_a%rAttr(index_x2a_Faxx_lat,n) + &
00166 l2x_a%rAttr(index_l2x_Fall_lat,n) * frac
00167 x2a_a%rAttr(index_x2a_Faxx_taux,n) = x2a_a%rAttr(index_x2a_Faxx_taux,n) + &
00168 l2x_a%rAttr(index_l2x_Fall_taux,n) * frac
00169 x2a_a%rAttr(index_x2a_Faxx_tauy,n) = x2a_a%rAttr(index_x2a_Faxx_tauy,n) + &
00170 l2x_a%rAttr(index_l2x_Fall_tauy,n) * frac
00171 x2a_a%rAttr(index_x2a_Sx_t,n) = x2a_a%rAttr(index_x2a_Sx_t,n) + &
00172 l2x_a%rAttr(index_l2x_Sl_t,n) * frac
00173
00174 if ( (index_x2a_Faxx_fco2_lnd /=0) .and. (index_l2x_Fall_nee /=0)) then
00175 x2a_a%rAttr(index_x2a_Faxx_fco2_lnd,n) = l2x_a%rAttr(index_l2x_Fall_nee,n) * frac
00176 end if
00177
00178 if ( index_x2a_Faxx_flxvoc1 /= 0 ) &
00179 x2a_a%rAttr(index_x2a_Faxx_flxvoc1,n) = l2x_a%rAttr(index_l2x_Fall_flxvoc1,n)* frac
00180 if ( index_x2a_Faxx_flxvoc2 /= 0 ) &
00181 x2a_a%rAttr(index_x2a_Faxx_flxvoc2,n) = l2x_a%rAttr(index_l2x_Fall_flxvoc2,n)* frac
00182 if ( index_x2a_Faxx_flxvoc3 /= 0 ) &
00183 x2a_a%rAttr(index_x2a_Faxx_flxvoc3,n) = l2x_a%rAttr(index_l2x_Fall_flxvoc3,n)* frac
00184 if ( index_x2a_Faxx_flxvoc4 /= 0 ) &
00185 x2a_a%rAttr(index_x2a_Faxx_flxvoc4,n) = l2x_a%rAttr(index_l2x_Fall_flxvoc4,n)* frac
00186 if ( index_x2a_Faxx_flxvoc5 /= 0 ) &
00187 x2a_a%rAttr(index_x2a_Faxx_flxvoc5,n) = l2x_a%rAttr(index_l2x_Fall_flxvoc5,n)* frac
00188 end if
00189
00190
00191 frac = x2a_a%rAttr(index_x2a_Sx_ifrac,n)
00192 if (frac > 0._r8) then
00193 x2a_a%rAttr(index_x2a_Sx_avsdr,n) = x2a_a%rAttr(index_x2a_Sx_avsdr,n) + &
00194 i2x_a%rAttr(index_i2x_Si_avsdr,n) * frac
00195 x2a_a%rAttr(index_x2a_Sx_avsdf,n) = x2a_a%rAttr(index_x2a_Sx_avsdf,n) + &
00196 i2x_a%rAttr(index_i2x_Si_avsdf,n) * frac
00197 x2a_a%rAttr(index_x2a_Sx_anidr,n) = x2a_a%rAttr(index_x2a_Sx_anidr,n) + &
00198 i2x_a%rAttr(index_i2x_Si_anidr,n) * frac
00199 x2a_a%rAttr(index_x2a_Sx_anidf,n) = x2a_a%rAttr(index_x2a_Sx_anidf,n) + &
00200 i2x_a%rAttr(index_i2x_Si_anidf,n) * frac
00201 x2a_a%rAttr(index_x2a_Sx_tref,n) = x2a_a%rAttr(index_x2a_Sx_tref,n) + &
00202 i2x_a%rAttr(index_i2x_Si_tref,n) * frac
00203 x2a_a%rAttr(index_x2a_Sx_qref,n) = x2a_a%rAttr(index_x2a_Sx_qref,n) + &
00204 i2x_a%rAttr(index_i2x_Si_qref,n) * frac
00205 x2a_a%rAttr(index_x2a_Faxx_evap,n) = x2a_a%rAttr(index_x2a_Faxx_evap,n) + &
00206 i2x_a%rAttr(index_i2x_Faii_evap,n) * frac
00207 x2a_a%rAttr(index_x2a_Faxx_lwup,n) = x2a_a%rAttr(index_x2a_Faxx_lwup,n) + &
00208 i2x_a%rAttr(index_i2x_Faii_lwup,n) * frac
00209 x2a_a%rAttr(index_x2a_Faxx_sen,n) = x2a_a%rAttr(index_x2a_Faxx_sen,n) + &
00210 i2x_a%rAttr(index_i2x_Faii_sen,n) * frac
00211 x2a_a%rAttr(index_x2a_Faxx_lat,n) = x2a_a%rAttr(index_x2a_Faxx_lat,n) + &
00212 i2x_a%rAttr(index_i2x_Faii_lat,n) * frac
00213 x2a_a%rAttr(index_x2a_Faxx_taux,n) = x2a_a%rAttr(index_x2a_Faxx_taux,n) + &
00214 i2x_a%rAttr(index_i2x_Faii_taux,n) * frac
00215 x2a_a%rAttr(index_x2a_Faxx_tauy,n) = x2a_a%rAttr(index_x2a_Faxx_tauy,n) + &
00216 i2x_a%rAttr(index_i2x_Faii_tauy,n) * frac
00217 x2a_a%rAttr(index_x2a_Sx_t,n) = x2a_a%rAttr(index_x2a_Sx_t,n) + &
00218 i2x_a%rAttr(index_i2x_Si_t,n) * frac
00219 end if
00220
00221 frac = x2a_a%rAttr(index_x2a_Sx_ofrac,n)
00222 if (frac > 0._r8) then
00223 x2a_a%rAttr(index_x2a_Sx_avsdr,n) = x2a_a%rAttr(index_x2a_Sx_avsdr,n) + &
00224 xao_a%rAttr(index_xao_So_avsdr,n) * frac
00225 x2a_a%rAttr(index_x2a_Sx_avsdf,n) = x2a_a%rAttr(index_x2a_Sx_avsdf,n) + &
00226 xao_a%rAttr(index_xao_So_avsdf,n) * frac
00227 x2a_a%rAttr(index_x2a_Sx_anidr,n) = x2a_a%rAttr(index_x2a_Sx_anidr,n) + &
00228 xao_a%rAttr(index_xao_So_anidr,n) * frac
00229 x2a_a%rAttr(index_x2a_Sx_anidf,n) = x2a_a%rAttr(index_x2a_Sx_anidf,n) + &
00230 xao_a%rAttr(index_xao_So_anidf,n) * frac
00231 x2a_a%rAttr(index_x2a_Sx_tref,n) = x2a_a%rAttr(index_x2a_Sx_tref,n) + &
00232 xao_a%rAttr(index_xao_So_tref,n) * frac
00233 x2a_a%rAttr(index_x2a_Sx_qref,n) = x2a_a%rAttr(index_x2a_Sx_qref,n) + &
00234 xao_a%rAttr(index_xao_So_qref,n) * frac
00235 x2a_a%rAttr(index_x2a_Faxx_evap,n) = x2a_a%rAttr(index_x2a_Faxx_evap,n) + &
00236 xao_a%rAttr(index_xao_Faox_evap,n) * frac
00237 x2a_a%rAttr(index_x2a_Faxx_lwup,n) = x2a_a%rAttr(index_x2a_Faxx_lwup,n) + &
00238 xao_a%rAttr(index_xao_Faox_lwup,n) * frac
00239 x2a_a%rAttr(index_x2a_Faxx_sen,n) = x2a_a%rAttr(index_x2a_Faxx_sen,n) + &
00240 xao_a%rAttr(index_xao_Faox_sen,n) * frac
00241 x2a_a%rAttr(index_x2a_Faxx_lat,n) = x2a_a%rAttr(index_x2a_Faxx_lat,n) + &
00242 xao_a%rAttr(index_xao_Faox_lat,n) * frac
00243 x2a_a%rAttr(index_x2a_Faxx_taux,n) = x2a_a%rAttr(index_x2a_Faxx_taux,n) + &
00244 xao_a%rAttr(index_xao_Faox_taux,n) * frac
00245 x2a_a%rAttr(index_x2a_Faxx_tauy,n) = x2a_a%rAttr(index_x2a_Faxx_tauy,n) + &
00246 xao_a%rAttr(index_xao_Faox_tauy,n) * frac
00247 x2a_a%rAttr(index_x2a_Sx_t,n) = x2a_a%rAttr(index_x2a_Sx_t,n) + &
00248 o2x_a%rAttr(index_o2x_So_t,n) * frac
00249 end if
00250
00251 if ( (index_x2a_Faxx_fco2_ocn /=0) .and. (index_o2x_Faoo_fco2 /=0)) then
00252 x2a_a%rAttr(index_x2a_Faxx_fco2_ocn,n) = o2x_a%rAttr(index_o2x_Faoo_fco2,n)
00253 end if
00254
00255
00256 if ( (index_x2a_Faxx_fdms /=0) .and. (index_o2x_Faoo_fdms /=0)) then
00257 x2a_a%rAttr(index_x2a_Faxx_fdms,n) = o2x_a%rAttr(index_o2x_Faoo_fdms,n)
00258 end if
00259
00260 end do
00261
00262 end subroutine mrg_x2a_run_mct
00263
00264
00265
00266 subroutine mrg_x2a_final_mct
00267
00268
00269
00270 end subroutine mrg_x2a_final_mct
00271
00272 end module mrg_x2a_mct
00273