emissionmodel.py 28.7 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
#!/usr/bin/env python
# stilt_tools.py

"""
Author : I. Super

Revision History:
Newly developed code, September 2017

This module holds an emission model that prepares emission files used by the observation operator and
to create pseudo-data

"""

import shutil
import os
import logging
import datetime as dtm
import numpy as np
from numpy import array, logical_and
import da.tools.io4 as io
import math

import da.tools.rc as rc
from da.tools.general import create_dirs, to_datetime
import netCDF4 as nc
identifier = 'EmissionModel ensemble '
version = '1.0'

30
#from da.ccffdas.observationoperator import get_time_indices
31
# Improve: should not be a python file? Defenitely not be stored in this folder!
32
from da.ccffdas import energy_use_country
33
34
energy_use_per_country = energy_use_country.energy_use_per_country 

35
from da.ccffdas import category_info
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
categories = category_info.categories
################### Begin Class Emission model ###################

class EmisModel(object):

    def __init__(self, dacycle=None):
        if dacycle != None:
            self.dacycle = dacycle
        else:
            self.dacycle = {}

    def setup(self, dacycle):
        self.dacycle = dacycle
        self.startdate = self.dacycle['time.fxstart']
        self.enddate = self.dacycle['time.finish']
        self.emisdir = dacycle.dasystem['datadir']
        self.inputdir = dacycle.dasystem['inputdir']
        self.proxyfile = dacycle.dasystem['emis.input.spatial']
        self.tempfileo = dacycle.dasystem['emis.input.tempobs']
        self.tempfilep = dacycle.dasystem['emis.input.tempprior']
        self.btime = int(dacycle.dasystem['run.backtime'])
        self.obsfile = dacycle.dasystem['obs.input.id']
        self.nrspc = int(dacycle.dasystem['obs.spec.nr'])
        self.species = dacycle.dasystem['obs.spec.name'].split(',')
        self.nrcat = int(dacycle.dasystem['obs.cat.nr'])
        self.nparams = int(dacycle.dasystem['nparameters'])
        self.nmembers = int(dacycle['da.optimizer.nmembers'])
        self.pparam = dacycle.dasystem['emis.pparam']
        self.paramfile = dacycle.dasystem['emis.paramfiles']
        self.countries = [country.strip() for country in dacycle.dasystem['countries'].split(';')]
        areafile = dacycle.dasystem['area.file']
        self.area = nc.Dataset(areafile)['Area'][:]

        self.stations2country = rc.read('/home/awoude/ffdas/RINGO/Data/station_countries.rc')
        self.energy_use_per_country = energy_use_per_country
        self.categories = categories
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
        self.paramdict = rc.read(dacycle.dasystem['paramdict'])

        # --- Settings
        year = self.startdate.year
        ndays = (dtm.datetime(year, 12, 31) - dtm.datetime(year - 1, 12, 31)).days
        self.ndays = ndays
        times = np.array([dtm.datetime(year, 1, 1, 0, 0, 0) + dtm.timedelta(hours = i) for i in range(ndays*24)])
        self.times = times
        datapath   = '/projects/0/ctdas/RINGO/EmissionInventories/DynEmMod_TimeProfiles'
        infilename = '{}/RINGO_ECMWF_DailyMeteo{}.nc'.format(datapath, year)
        with nc.Dataset(infilename) as infile:
            self.T2myr    = infile.variables['T2m'     ][:ndays] - 273.15         # T2m           K --> oC
            self.U10      = infile.variables['U10m'    ][:ndays]                  # u             m/s
            self.Rinyr    = infile.variables['Rsin'    ][:ndays] / (24. * 1.0e4)  # Radiation (?) J/m2/day --> J/cm2/hr
            self.T2myr_av = infile.variables['T2m_avg' ][:ndays] - 273.15         #
            self.U10_av   = infile.variables['U10m_avg'][:ndays]
            self.Rinyr_av = infile.variables['Rsin_avg'][:ndays]

90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106

    def find_in_state(self, station, cat, name=None, return_index=False):
        """Function that finds the index in the state vector"""
        if not name:
            key = station + '.' + cat
        else:
            key = station + '.' + cat + '.' + name

        if key in  self.paramdict:
            i_in_state = int(self.paramdict[key])
            if return_index: return i_in_state
            else:
                value = self.prm[i_in_state]
                return value
        elif return_index: return False
        return 1

107
    def get_emis(self, dacycle, samples, indices, do_pseudo):
108
109
110
111
112
113
114
115
116
117
118
119
        """set up emission information for pseudo-obs (do_pseudo=1) and ensemble runs (do_pseudo=0)"""
        
        if do_pseudo==1:
            priorparam=os.path.join(self.emisdir,self.pparam)
            f = io.ct_read(priorparam, 'read')
            self.prm = f.get_variable('prior_values')[:self.nparams]
            f.close()
            self.get_spatial(dacycle, 999, samples, infile=os.path.join(self.emisdir, self.paramfile))
            self.get_temporal(dacycle, 999, samples, do_pseudo=1)
        elif do_pseudo==0:
            self.timestartkey = self.dacycle['time.sample.start']
            self.timefinishkey = self.dacycle['time.sample.end']
120
            time_profiles = self.make_time_profiles(indices=indices)
121
122
123
124
125
126
127
128
129
130
131
132
            for member in range(self.nmembers):
                #first remove old files, then create new emission files per ensemble member
                if self.startdate == self.timestartkey:
                    file = os.path.join(dacycle.dasystem['datadir'],'temporal_data_%03d.nc'%member)
                    try:
                        os.remove(file)
                    except OSError:
                        pass
                prmfile=os.path.join(dacycle['dir.input'],'parameters.%03d.nc'%member)
                f = io.ct_read(prmfile, 'read')
                self.prm = f.get_variable('parametervalues')
                f.close()
133
134
                self.get_yearly_emissions()
                self.get_emissions(dacycle, member, time_profiles, infile=os.path.join(self.emisdir, self.paramfile))
135
136
        logging.debug('Succesfully wrote prior spatial and temporal emission files')
        
137
    def get_yearly_emissions(self):
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
        yremis = np.zeros((len(self.categories), len(self.countries), len(self.species)))
        for i_country, country in enumerate(self.countries):
            for i_cat, (cat, values) in enumerate(self.categories.items()):
                emission_factor = values['emission_factors']

                emission_factor *= self.find_in_state(country, cat, 'emission_factors')
                fraction_of_total = values['fraction_of_total']
                fraction_of_total *= self.find_in_state(country, cat, 'fraction_of_total')
               
                e_use = self.energy_use_per_country[country][values['spatial']] 
                for i_species, specie in enumerate(self.species):
                    emission_ratio = values[specie]
                    uncertainty = values[specie+'.uncertainty']
                    if uncertainty == 'l':
                        emission_ratio = np.exp(emission_ratio)
                        if emission_ratio > 1:
                            logging.debug('{} {} {} {}'.format(country, cat, specie, emission_ratio))
                    emission_ratio *= self.find_in_state(country, cat, specie)
                    emission = e_use * fraction_of_total * emission_factor * emission_ratio
                    yremis[i_cat, i_country, i_species] = emission
        self.yremis = yremis

160
    def get_emissions(self, dacycle, member, time_profiles, infile=None):
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
        """read in proxy data used for spatial distribution of the gridded emissions, disaggregate yearly totals for the area"""
        # Create a recalculation factor from kg/km2/yr to umol/m2/sec
        M_mass = [44e-9, 28e-9]
        sec_year = 24*366*3600. #seconds in a year (leapyear)
        kgperkmperyr2umolperm2pers = np.array(M_mass)[:, None, None] * sec_year * self.area[None, :, :] 
        self.kgperkmperyr2umolperm2pers = kgperkmperyr2umolperm2pers

        #read in proxy data for spatial disaggregation
        infile = os.path.join(self.emisdir, self.proxyfile)
        proxy = io.ct_read(infile, method='read')
        proxy_category_names = proxy['emis_category_names'][:]
        proxy_category_names = [b''.join(category_name).strip().decode() for category_name in proxy_category_names]
         
        proxy_country_names = proxy['country_names'][:]
        proxy_country_names = [b''.join(country_name).strip().decode() for country_name in proxy_country_names]

        spatial_distributions = np.zeros((self.nrcat, len(self.countries), self.area.shape[0], self.area.shape[1]))
178
179
180
181
182
        for i, category in enumerate(self.categories):
            spatial_name = self.categories[category]['spatial']
            cat_index = proxy_category_names.index(spatial_name)
            for country in self.countries:
                country_index = self.countries.index(country)
183
184
185
186
187
188
189
190
                # Get the emission distribution for the category
                category_distribution_country = proxy['proxy_maps'][cat_index, country_index, :, :]
                spatial_distributions[i, country_index, :, :] = category_distribution_country
                 
        # Multiply spatial distributions with the yearly emissions in the country to get spatially distributed emissions
        spatial_emissions = spatial_distributions[:, :, None, :, :] * self.yremis[:, :, :, None, None] # cat, country, species, lat, lon
        # Sum over the countries to overlay them.
        spatial_emissions = spatial_emissions.sum(axis=1) # Indices: category, species, lat, lon
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
        self.spatial_emissions = spatial_emissions 

        emissions = []
        for i, category in enumerate(self.categories):
            spatial_name = self.categories[category]['spatial']
            cat_index = proxy_category_names.index(spatial_name)
            temporal_name = self.categories[category]['temporal']
            temporal_profile = time_profiles[temporal_name]
            self.temporal_profile = temporal_profile
            self.emissions = emissions
            emissions.append(spatial_emissions[cat_index, :, :, :] * temporal_profile[None, :, :, :])

        emissions = np.asarray(emissions)
        emissions = np.swapaxes(emissions, 0,1) # Indices: species, category, time, lat, lon
        
206
        # Recalculate spatial emissions to umol/sec/m2
207
        emissions = emissions / kgperkmperyr2umolperm2pers[:, None, None, :, :]
208
209
210
211
212
        ## create output file
        prior_file = os.path.join(self.inputdir, 'prior_spatial_{0:03d}.nc'.format(member))
        f = io.CT_CDF(prior_file, method='create')
        dimid = f.add_dim('ncat', self.nrcat)
        dimid2 = f.add_dim('ops',2 )
213
        dimtime= f.add_dim('time', emissions.shape[2])
214
215
216
217
218
219
220
221
222
        dimlat = f.add_dim('lat', self.area.shape[0])
        dimlon = f.add_dim('lon', self.area.shape[1])
    
        #loop over all tracers
        for i, species in enumerate(self.species): 
            savedict = io.std_savedict.copy() 
            savedict['name'] = species
            savedict['long_name'] = "Spatially distributed emissions"
            savedict['units'] = "micromole/m2/s"
223
224
            savedict['dims'] = dimid + dimtime + dimlat + dimlon
            savedict['values'] = emissions[i,:,:,:]
225
226
            f.add_data(savedict)
        f.close()
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
            
    def make_time_profiles(self, indices):
        """Function that calculates the time profiles based on pre-specified
        monthly, daily and hourly profiles. Temperature and radiation affect 
        household heating and powerplant time profiles.
        Input: 
            year: int: The year for which the profiles should be calculated
        Returns: 
            dict of np.arrays:
                The temporal profiles (one for each hour) for each gridcel and timestep """
    
        ndays, nlat, nlon = self.T2myr.shape
        # --- calculate degree day sum and corresponding time profiles
        fr_cons = 0.2  #    constant emission  for consumers (cooking, warm water)
        fr_gls  = 0.   # no constant emission  for glasshouses
        fr_coal = 0.8  #    constant emission  for coal-fired power plants
        fr_gas  = 0.1  #    constant emission  for gas-fired power plants
    
        T0_cons = 18.  # temperature threshold for consumer heating
        T0_gls  = 15.  # temperature threshold for glasshouse heating
        T0_coal = 25.  # temperature threshold for coal-fired power plants
        U0_gas  = 10.  # wind speed  threshold for gas-fired power plants
        R0_gas  = 50.  # radiation   threshold for gas-fired power plants
    
        HDC_cons    = np.empty((ndays, nlat, nlon)) # Heating Demand Category 1 (Household heating)
        HDC_gls     = np.empty((ndays, nlat, nlon)) # Heating Demand Category 2 (Glasshouse heating)
        HDC_coal    = np.empty((ndays, nlat, nlon)) # Heating Demand Category 3 (Power plants)
        HDC_gas     = np.empty((ndays, nlat, nlon)) # Heating Demand Category 4 (Renewable activity)
    
        for i in range(ndays):
            HDC_cons[i] = np.fmax(T0_cons - self.T2myr[i, :, :], 0) # Daily demand for consumers / household heating
            HDC_gls [i] = np.fmax(T0_gls  - self.T2myr[i, :, :], 0) # Daily demand for glasshouse heating
            HDC_coal[i] = np.fmax(T0_coal - self.T2myr[i, :, :], 0) # Daily demand for coal-fired powerplant productvity
            dum1        = np.fmax(U0_gas  - self.U10[i, :, :], 0) # Wind energy
            dum2        = np.fmax(R0_gas  - self.Rinyr[i, :, :], 0) # Solar energy
            HDC_gas [i] = dum1 * dum2
    
        HC_cons         = HDC_cons.mean(axis=0)
        HC_gls          = HDC_gls.mean(axis=0)
        HC_coal         = HDC_coal.mean(axis=0)
        HC_gas          = HDC_gas.mean(axis=0)
    
        t_consd  = (HDC_cons + fr_cons * HC_cons) / ((1 + fr_cons) * HC_cons) # daily time profile for consumer/household heating
        t_glsd   = (HDC_gls  + fr_gls  * HC_gls) / ((1 + fr_gls ) * HC_gls)   #                        glasshouse
        t_coald  = (HDC_coal + fr_coal * HC_coal) / ((1 + fr_coal) * HC_coal) #                        coal-fired powerplant
        t_gasd   = (HDC_gas  + fr_gas  * HC_gas) / ((1 + fr_gas ) * HC_gas)   #                        gas-fired  powerplant
    
        # Hourly time profiles for energy use and consumers
        engy_hr  = np.array([0.79, 0.72, 0.72, 0.71, 0.74, 0.80,
                             0.92, 1.08, 1.19, 1.22, 1.21, 1.21,
                             1.17, 1.15, 1.14, 1.13, 1.10, 1.07, 
                             1.04, 1.02, 1.02, 1.01, 0.96, 0.88])
        cons_hr  = np.array([0.38, 0.36, 0.36, 0.36, 0.37, 0.50, 
                             1.19, 1.53, 1.57, 1.56, 1.35, 1.16,
                             1.07, 1.06, 1.00, 0.98, 0.99, 1.12, 
                             1.41, 1.52, 1.39, 1.35, 1.00, 0.42])
    
        engy_hrm = np.tile(engy_hr, (nlat, nlon, 1)).transpose(2, 0, 1) #Repeat the daily time profile for all days
        cons_hrm = np.tile(cons_hr, (nlat, nlon, 1)).transpose(2, 0, 1)
    
        nt       = len(self.times)
        t_gas    = np.empty((nt, nlat, nlon))
        t_coal   = np.empty((nt, nlat, nlon))
        t_cons   = np.empty((nt, nlat, nlon))
        t_gls    = np.empty((nt, nlat, nlon))
    
        #Repeat the daily time profile for all days
        for iday in range(ndays):
            t_gas [iday * 24:(iday + 1) * 24, :, :] = np.tile(t_gasd [iday], (24, 1, 1)) * engy_hrm # hourly time profile = daily * hourly
            t_coal[iday * 24:(iday + 1) * 24, :, :] = np.tile(t_coald[iday], (24, 1, 1)) * engy_hrm
            t_cons[iday * 24:(iday + 1) * 24, :, :] = np.tile(t_consd[iday], (24, 1, 1)) * cons_hrm
            t_gls [iday * 24:(iday + 1) * 24, :, :] = np.tile(t_glsd [iday], (24, 1, 1)) * cons_hrm
    
        #######
        # road traffic time profiles
        #######
        carhw_mnt = [1.04, 1.05, 0.98, 1.00, 1.00, 1.04, 1.01, 0.93, 0.86, 1.01, 1.07, 1.02] # monthly time profiles car highway
        carmr_mnt = [1.05, 1.06, 0.93, 1.01, 0.99, 1.01, 1.00, 1.00, 0.86, 1.02, 1.08, 1.03] #                       car middle road
        carur_mnt = [1.06, 1.10, 0.97, 1.00, 1.01, 1.00, 1.00, 0.98, 0.85, 1.01, 1.07, 1.01] #                       car urban road
        hdvhw_mnt = [0.92, 0.89, 1.03, 1.02, 1.07, 1.05, 1.03, 1.02, 0.90, 0.97, 0.98, 0.94]
        hdvmr_mnt = [1.06, 1.09, 0.93, 0.98, 0.99, 1.01, 1.02, 1.02, 0.78, 1.09, 1.08, 1.03]
        hdvur_mnt = [1.10, 1.17, 0.91, 0.94, 0.76, 1.02, 1.04, 1.04, 0.82, 1.13, 1.14, 1.03]
        carhw_wk  = [1.00, 1.09, 1.10, 1.10, 1.11, 0.87, 0.73] # weekly time profile (mon-sun)
        carmr_wk  = [0.99, 1.07, 1.08, 1.07, 1.09, 0.93, 0.76] 
        carur_wk  = [0.98, 1.06, 1.07, 1.07, 1.09, 0.92, 0.80] 
        hdvhw_wk  = [1.23, 1.30, 1.27, 1.27, 1.33, 0.42, 0.19] 
        hdvmr_wk  = [1.27, 1.25, 1.26, 1.29, 1.29, 0.45, 0.19] 
        hdvur_wk  = [1.07, 1.24, 1.23, 1.21, 1.23, 0.63, 0.38] 
        carhw_hrw = [0.17, 0.09, 0.05, 0.05, 0.09, 0.42, 1.38, 1.87, # hourly time profile weekdays
                     1.81, 1.31, 1.11, 1.13, 1.26, 1.30, 1.39, 1.63, 
                     1.95, 1.97, 1.46, 1.05, 0.78, 0.69, 0.62, 0.42] 
        carmr_hrw = [0.20, 0.09, 0.05, 0.05, 0.11, 0.36, 1.04, 1.42, 
                     1.49, 1.26, 1.20, 1.25, 1.37, 1.42, 1.49, 1.62, 
                     1.83, 1.92, 1.56, 1.23, 0.92, 0.81, 0.77, 0.56]
        carur_hrw = [0.30, 0.20, 0.12, 0.09, 0.10, 0.25, 0.63, 1.22, 
                     1.53, 1.17, 1.10, 1.21, 1.37, 1.41, 1.52, 1.71, 
                     1.94, 2.03, 1.59, 1.18, 0.98, 0.91, 0.83, 0.61]
        hdvhw_hrw = [0.13, 0.12, 0.12, 0.15, 0.28, 0.68, 1.59, 1.50, 
                     1.46, 1.67, 1.72, 1.77, 1.70, 1.77, 1.82, 1.85, 
                     1.60, 1.18, 0.97, 0.73, 0.49, 0.33, 0.22, 0.16]
        hdvmr_hrw = [0.13, 0.09, 0.11, 0.10, 0.21, 0.52, 1.58, 1.75, 
                     1.61, 1.74, 1.70, 1.69, 1.58, 1.70, 1.79, 1.96, 
                     1.82, 1.23, 0.93, 0.61, 0.43, 0.30, 0.23, 0.19]
        hdvur_hrw = [0.14, 0.08, 0.06, 0.08, 0.15, 0.36, 0.95, 1.61, 
                     1.86, 1.77, 1.76, 1.79, 1.74, 1.75, 1.83, 2.16, 
                     1.76, 1.27, 0.92, 0.65, 0.47, 0.37, 0.29, 0.21]
        carhw_hrd = [0.57, 0.34, 0.20, 0.13, 0.11, 0.17, 0.29, 0.38, # hourly time profile weekends
                     0.65, 1.03, 1.30, 1.54, 1.75, 1.95, 1.94, 1.83, 
                     1.81, 1.72, 1.32, 1.27, 1.18, 1.00, 0.84, 0.68] 
        carmr_hrd = [0.58, 0.35, 0.22, 0.14, 0.14, 0.25, 0.32, 0.38, 
                     0.61, 1.01, 1.29, 1.50, 1.73, 1.90, 1.96, 1.88, 
                     1.84, 1.76, 1.33, 1.18, 1.10, 0.93, 0.88, 0.72]
        carur_hrd = [0.79, 0.58, 0.45, 0.33, 0.28, 0.27, 0.28, 0.32, 
                     0.48, 0.79, 1.04, 1.29, 1.60, 1.77, 1.85, 1.87, 
                     1.86, 1.76, 1.38, 1.18, 1.10, 0.99, 0.93, 0.81]
        hdvhw_hrd = [0.44, 0.38, 0.36, 0.38, 0.42, 0.67, 0.99, 1.19, 
                     1.51, 1.69, 1.77, 1.68, 1.63, 1.53, 1.39, 1.34, 
                     1.33, 1.29, 1.12, 0.86, 0.69, 0.55, 0.44, 0.35]
        hdvmr_hrd = [0.49, 0.33, 0.35, 0.23, 0.36, 0.56, 0.93, 1.21, 
                     1.31, 1.47, 1.58, 1.76, 1.63, 1.52, 1.52, 1.28, 
                     1.24, 1.09, 0.89, 0.67, 0.61, 0.58, 0.55, 0.52]
        hdvur_hrd = [0.47, 0.41, 0.33, 0.33, 0.34, 0.42, 0.61, 0.85, 
                     1.04, 1.27, 1.43, 1.55, 1.72, 1.72, 1.69, 1.69, 
                     1.66, 1.56, 1.34, 1.03, 0.83, 0.68, 0.58, 0.46]
    
        t_carhw   = np.empty((nt, nlat, nlon))
        t_carmr   = np.empty((nt, nlat, nlon))
        t_carur   = np.empty((nt, nlat, nlon))
        t_hdvhw   = np.empty((nt, nlat, nlon))
        t_hdvmr   = np.empty((nt, nlat, nlon))
        t_hdvur   = np.empty((nt, nlat, nlon))
        t_ship    = np.empty((nt, nlat, nlon))
    
        # Multiply the hourly, daily and monthly time profiles to create new profile
        for i, t in enumerate(self.times):
            month_idx = t.month - 1
            day_idx   = t.weekday()
            hour_idx  = t.hour
            if day_idx < 5: # weekdays
                t_carhw[i, None, None] = carhw_mnt[month_idx] * carhw_wk[day_idx] * carhw_hrw[hour_idx]
                t_carmr[i, None, None] = carmr_mnt[month_idx] * carmr_wk[day_idx] * carmr_hrw[hour_idx]
                t_carur[i, None, None] = carur_mnt[month_idx] * carur_wk[day_idx] * carur_hrw[hour_idx]
                t_hdvhw[i, None, None] = hdvhw_mnt[month_idx] * hdvhw_wk[day_idx] * hdvhw_hrw[hour_idx]
                t_hdvmr[i, None, None] = hdvmr_mnt[month_idx] * hdvmr_wk[day_idx] * hdvmr_hrw[hour_idx]
                t_hdvur[i, None, None] = hdvur_mnt[month_idx] * hdvur_wk[day_idx] * hdvur_hrw[hour_idx]
                t_ship [i, None, None] = hdvhw_mnt[month_idx] * hdvhw_wk[day_idx] * hdvhw_hrw[hour_idx]
            else: # weekends
                t_carhw[i, None, None] = carhw_mnt[month_idx] * carhw_wk[day_idx] * carhw_hrd[hour_idx]
                t_carmr[i, None, None] = carmr_mnt[month_idx] * carmr_wk[day_idx] * carmr_hrd[hour_idx]
                t_carur[i, None, None] = carur_mnt[month_idx] * carur_wk[day_idx] * carur_hrd[hour_idx]
                t_hdvhw[i, None, None] = hdvhw_mnt[month_idx] * hdvhw_wk[day_idx] * hdvhw_hrd[hour_idx]
                t_hdvmr[i, None, None] = hdvmr_mnt[month_idx] * hdvmr_wk[day_idx] * hdvmr_hrd[hour_idx]
                t_hdvur[i, None, None] = hdvur_mnt[month_idx] * hdvur_wk[day_idx] * hdvur_hrd[hour_idx]
                t_ship [i, None, None] = hdvhw_mnt[month_idx] * hdvhw_wk[day_idx] * hdvhw_hrd[hour_idx]
        # --- normalize time profiles
        t_gas   /= nt
        t_coal  /= nt
        t_cons  /= nt
        t_gls   /= nt
        t_carhw /= nt
        t_carmr /= nt
        t_carur /= nt
        t_hdvhw /= nt
        t_hdvmr /= nt
        t_hdvur /= nt
        t_ship  /= nt
   
        datepoint = self.startdate
        enddate = self.enddate
        #assert False
        #time_indices = get_time_indices(datepoint, startdate=None, backtimes=24)
 
        t_ind = np.ones((nt, nlat, nlon))
        # Make container (dict) containing all time profiles
        time_profiles = {'t_gas': t_gas[indices], 
                         't_coal': t_coal[indices],
                         't_ind': t_ind[indices],
                         't_cons': t_cons[indices], 
                         't_gls': t_gls[indices],
                         't_carhw': t_carhw[indices],
                         't_carmr': t_carmr[indices], 
                         't_carur': t_carur[indices], 
                         't_hdvhw': t_hdvhw[indices], 
                         't_hdvmr': t_hdvmr[indices], 
                         't_hdvur': t_hdvur[indices], 
                         't_ship': t_ship[indices]}
        logging.debug('Time profiles created')
        return time_profiles

416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
    def get_temporal(self, dacycle, member, samples, do_pseudo):
        """read in time profiles used for temporal distribution of the emissions"""
        # First, get the station names from the smaples. For these stations, the time profiles will be written.
        codes = samples.getvalues('code')
        self.codes = codes
        ## SI: Get the country names
        stationnames = []
        for code in codes:
            two_names = any(x in code for x in ['UW', 'DW'])
            stationnames.append('_'.join(code.split('_')[1:2 + two_names]))
        stationnames = list(set(stationnames))
        # For pseudo-observation (do_pseudo==1) or when no time profiles need to be optimized the profiles are simply read from the
        # input file and copied to another file. Otherwise create a new file per ensemble member at t=0 and update the profiles for each time step
        # Check if the ensemble file exists. Otherwise create.
        ## SI: tmeporal data should include only countries.
        ensfile = os.path.join(self.emisdir, 'temporal_data_%03d.nc'%member)
        if not os.path.exists(ensfile):
            dumfile = os.path.join(self.emisdir, self.tempfilep)
            shutil.copy2(dumfile,ensfile)
        time_profiles_ds = nc.Dataset(ensfile)
        times = time_profiles_ds['Times'][:]
        times = np.array([dtm.datetime.strptime(time, '%Y-%m-%d %H:%M:%S') for time in np.array(times)])
        self.times = times
        subselect = logical_and(times >= self.timestartkey , times < self.timefinishkey).nonzero()[0]
        date_selection = times.take(subselect, axis=0)
        # The time profiles should always cover at least one full year
        start_date = dtm.datetime(self.timestartkey.year,1,1,0,0) #first time included
        end_date = dtm.datetime(self.timestartkey.year,12,31,23,0) #last time included
        dt = dtm.timedelta(0,3600)
        starttime_index = np.where(times==self.timestartkey)[0][0]
        startdate_index = np.where(times==self.startdate)[0][0]
        end_index = np.where(times==self.timefinishkey)[0][0]
        
        """ Time profiles should, for a full year, always have an average value of 1.0. Therefore, a new method has been developed
        to optimize time profiles such that we comply with this and the time profiles do not affect the total emissions.
        For this purpose we apply the scaling factor (statevector) to the period covered in this cycle. The time profile for all dates 
        outside this period are scaled equally such that the time profile remains its average value of 1.0. Except previously updated
        dates (from previous cycles) are maintained (they are already optimized!)."""

        unselected_day = np.where((times<self.startdate) |  (times>self.timefinishkey))[0]
        category_names = list(time_profiles_ds['category_name'][:]) 
        self.category_names = category_names
        station_names_ds = list(time_profiles_ds['station_names'][:])
        profiles = np.zeros(time_profiles_ds['time_profile'][:].shape)
        for category, values in self.categories.items():
            cat_index = category_names.index(values['temporal'])
            for station in stationnames:
            ## SI: for country in countries:
                paramvalue = self.find_in_state(station, category, values['temporal'])
                if paramvalue != 1:
                    station_index = station_names_ds.index(station)
                    original_profile = time_profiles_ds['time_profile'][station_index, cat_index, :]
                    selected_profile = time_profiles_ds['time_profile'][station_index, cat_index, :].take(subselect)
                    new_profile      = selected_profile[:] * paramvalue
                    daily_sum = np.array(original_profile[unselected_day]).sum()

                    original_profile[:startdate_index] = original_profile[:startdate_index] - (original_profile[:startdate_index] / daily_sum) * (new_profile.sum() - selected_profile.sum())
                    original_profile[end_index:]       = original_profile[end_index:]       - (original_profile[end_index:]       / daily_sum) * (new_profile.sum() - selected_profile.sum())
                    original_profile[starttime_index:end_index] = new_profile
                    profiles[station_index, cat_index, :] = original_profile
        time_profiles_ds.close()
        # Now, write the output
        tofile = nc.Dataset(ensfile, 'r+')

        for category, values in self.categories.items():
            cat_index = category_names.index(values['temporal'])
            for station in stationnames:
            ## SI: for country in countries:
                if self.find_in_state(station, category, 'time_profile') != 1:
                    station_index = station_names_ds.index(station)
                    tofile['time_profile'][station_index, cat_index, :] = profiles[station_index, cat_index, :]
        tofile.close()
        # Now read in the correct profiles, select the correct time period and write the profiles into one file per ensemble member
        time_profiles_ds = nc.Dataset(ensfile)
        subselect = logical_and(times >= times[0] , times <= times[-1]).nonzero()[0]
        date_selection = times.take(subselect, axis=0)

        
        for station in stationnames:
        ## SI: for country in countries:
            station_index = station_names_ds.index(station)
            prior_file = os.path.join(self.inputdir, 'prior_temporal_{0}_{1:03}.nc'.format(station, member))
            f = io.CT_CDF(prior_file, method='create')
            dimtime = f.add_dim('Times', len(date_selection))
                
            cat_names_done = []
            for category, values in self.categories.items():
                cat_name = values['temporal']
                cat_index = category_names.index(cat_name)
                if not cat_name in cat_names_done:
                    profile = np.array(time_profiles_ds['time_profile'][station_index, cat_index, :].take(subselect))
                    savedict = io.std_savedict.copy()
                    savedict['name'] = cat_name
                    savedict['long_name'] = "Temporal distribution"
                    savedict['units'] = ""
                    savedict['dims'] = dimtime
                    savedict['values'] = profile
                    f.add_data(savedict)
                    cat_names_done.append(cat_name)
            f.close()

        time_profiles_ds.close() 
################### End Class Emission model ###################


if __name__ == "__main__":
    pass