! source file: /sfs/fs1/work-geomar6/smomw113/UVic_ESCM/2.9_dk_0814/updates/07_aouassim/source/common/co2data.F subroutine co2emitdata return end subroutine co2ccndata return end subroutine satdata return end subroutine co2ccnforcdata !======================================================================= ! routine to read and interpolate one dimensional forcing data !======================================================================= ! ! I will use this file only to develop the subroutine co2ccnforcdata. ! Later this code shall be moved into the file co2data.F. ! ! wkoeve@geomar.de, 7.09.15 ! ! In order to keep this as simple as possible, I will not consider co-operability ! with several other options, namely ! O_carbon_co2_2d, O_co2emit_track_co2, defined O_co2emit_track_sat ! ! I need to add a respective break statement early in this routine, or in UVic_ESCM, ! to make sure there is no conflict with these options. See lines 698ff for a similar example. ! ! As I have coded it right now, one should set both 1 and 1. ! If only 1 is set the value of co2_yr (=pyear) will be used continuously. ! ! In this version, I have changed variable name to fit the std of UVic better. ! implicit none character(120) :: fname, name, new_file_name, text integer iou, n, ln, ib(10), ic(10) logical inqvardef, exists, track real avg_co2, dat(3), data_time, pk, tim(3), wt1, wt3, fa real, allocatable :: data(:), time(:) save dat, data, ln, tim, time include "size.h" include "param.h" include "pconst.h" include "stdunits.h" include "calendar.h" include "cembm.h" include "switch.h" include "tmngr.h" include "atm.h" ! WK: The following block is for the call from *subroutine setembm*, I think. if (.not. allocated (time)) then print*, "####################################################" print*, "WK: loading data from A_co2_forc.nc for forced mode." name = "A_co2_forc.nc" fname = new_file_name (name) print*,"load from the following file: ",trim(fname) inquire (file=trim(fname), exist=exists) if (.not. exists) then print*, "==> Error: ", trim(fname), " does not exist." stop '=>co2ccnforc' else call openfile (fname, iou) call getdimlen ('time', iou, ln) allocate ( time(ln) ) allocate ( data(ln) ) ib(:) = 1 ic(:) = ln call getvara ('time', iou, ln, ib, ic, time, c1, c0) text = 'years' call getatttext (iou, 'time', 'units', text) if (trim(text) .eq. "days since 1-1-1") & time(:) = time(:)/yrlen - 1. if (trim(text) .eq. "days since 0-1-1") & time(:) = time(:)/yrlen if (trim(text) .eq. "years since 1-1-1") & time(:) = time(:) - 1. exists = inqvardef('A_co2forc', iou) if (.not. exists) then print*, "==> Warning: A_co2forc data does not exist." else call getvara ('A_co2forc', iou, ln, ib, ic, data, c1, c0) endif endif tim(:) = time(1) dat(:) = data(1) endif ! WK: Here the interpolation starts at any call of this routine from *subroutine gasbc* ! WK: Initially, co2_yr = pyear by default. data_time = year0 + accel_yr0 + (relyr - accel_yr0)*accel tim(2) = min(time(ln), max(time(1), data_time)) if (tim(2) .le. time(1)) then dat(2) = data(1) elseif (tim(2) .ge. time(ln)) then dat(2) = data(ln) else if (tim(2) .gt. tim(3)) then do n=2,ln if (time(n-1) .le. tim(2) .and. time(n) .ge. tim(2)) then tim(1) = time(n-1) dat(1) = data(n-1) tim(3) = time(n) dat(3) = data(n) endif enddo endif wt1 = 1. if (tim(3) .ne. tim(1)) wt1 = (tim(3)-tim(2))/(tim(3)-tim(1)) wt1 = max(0., min(1., wt1)) wt3 = 1. - wt1 dat(2) = dat(1)*wt1 + dat(3)*wt3 endif co2ccnforc = dat(2) if (tsits) then print*,"co2ccnforc: ",co2ccnforc endif return end