/****************************************************************************** * @brief Read atmospheric forcing data. *****************************************************************************/ void vic_force(void) { extern size_t NF; extern size_t NR; extern size_t current; extern force_data_struct *force; extern dmy_struct *dmy; extern domain_struct global_domain; extern domain_struct local_domain; extern filenames_struct filenames; extern global_param_struct global_param; extern option_struct options; extern soil_con_struct *soil_con; extern veg_con_map_struct *veg_con_map; extern veg_con_struct **veg_con; extern veg_hist_struct **veg_hist; extern parameters_struct param; extern param_set_struct param_set; double *t_offset = NULL; double *dvar = NULL; size_t i; size_t j; size_t v; size_t band; int vidx; size_t d3count[3]; size_t d3start[3]; size_t d4count[4]; size_t d4start[4]; double *Tfactor; // allocate memory for variables to be read dvar = malloc(local_domain.ncells_active * sizeof(*dvar)); check_alloc_status(dvar, "Memory allocation error."); // for now forcing file is determined by the year sprintf(filenames.forcing[0], "%s%4d.nc", filenames.f_path_pfx[0], dmy[current].year); // global_param.forceoffset[0] resets every year since the met file restarts // every year if (current > 1 && (dmy[current].year != dmy[current - 1].year)) { global_param.forceoffset[0] = 0; } // only the time slice changes for the met file reads. The rest is constant d3start[1] = 0; d3start[2] = 0; d3count[0] = 1; d3count[1] = global_domain.n_ny; d3count[2] = global_domain.n_nx; // Air temperature: tas for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[AIR_TEMP].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].air_temp[j] = (double) dvar[i]; } } // Precipitation: prcp for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[PREC].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].prec[j] = (double) dvar[i]; } } // Downward solar radiation: dswrf for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[SWDOWN].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].shortwave[j] = (double) dvar[i]; } } // Downward longwave radiation: dlwrf for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[LWDOWN].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].longwave[j] = (double) dvar[i]; } } // Wind speed: wind for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[WIND].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].wind[j] = (double) dvar[i]; } } // vapor pressure: vp for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[VP].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].vp[j] = (double) dvar[i]; } } // Pressure: pressure for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[PRESSURE].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].pressure[j] = (double) dvar[i]; } } // Optional inputs if (options.LAKES) { // Channel inflow to lake d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[CHANNEL_IN].varname, d3start, d3count, dvar); for (j = 0; j < NF; j++) { for (i = 0; i < local_domain.ncells_active; i++) { force[i].channel_in[j] = (double) dvar[i]; } } } if (options.CARBON) { // Atmospheric CO2 mixing ratio for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[CATM].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].Catm[j] = (double) dvar[i]; } } // Cosine of solar zenith angle for (j = 0; j < NF; j++) { for (i = 0; i < local_domain.ncells_active; i++) { force[i].coszen[j] = compute_coszen( local_domain.locations[i].latitude, local_domain.locations[i].longitude, soil_con[i].time_zone_lng, dmy[current].day_in_year, dmy[current].dayseconds); } } // Fraction of shortwave that is direct for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[FDIR].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].fdir[j] = (double) dvar[i]; } } // Photosynthetically active radiation for (j = 0; j < NF; j++) { d3start[0] = global_param.forceskip[0] + global_param.forceoffset[0] + j; get_scatter_nc_field_double(filenames.forcing[0], param_set.TYPE[PAR].varname, d3start, d3count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { force[i].par[j] = (double) dvar[i]; } } } // Update the offset counter global_param.forceoffset[0] += NF; // Initialize the veg_hist structure with the current climatological // vegetation parameters. This may be overwritten with the historical // forcing time series. for (i = 0; i < local_domain.ncells_active; i++) { for (v = 0; v < options.NVEGTYPES; v++) { vidx = veg_con_map[i].vidx[v]; if (vidx != NODATA_VEG) { for (j = 0; j < NF; j++) { veg_hist[i][vidx].albedo[j] = veg_con[i][vidx].albedo[dmy[current].month - 1]; veg_hist[i][vidx].displacement[j] = veg_con[i][vidx].displacement[dmy[current].month - 1]; veg_hist[i][vidx].fcanopy[j] = veg_con[i][vidx].fcanopy[dmy[current].month - 1]; veg_hist[i][vidx].LAI[j] = veg_con[i][vidx].LAI[dmy[current].month - 1]; veg_hist[i][vidx].roughness[j] = veg_con[i][vidx].roughness[dmy[current].month - 1]; } } } } // Read veg_hist file if (options.LAI_SRC == FROM_VEGHIST || options.FCAN_SRC == FROM_VEGHIST || options.ALB_SRC == FROM_VEGHIST) { // for now forcing file is determined by the year sprintf(filenames.forcing[1], "%s%4d.nc", filenames.f_path_pfx[1], dmy[current].year); // global_param.forceoffset[1] resets every year since the met file restarts // every year if (current > 1 && (dmy[current].year != dmy[current - 1].year)) { global_param.forceoffset[1] = 0; } // only the time slice changes for the met file reads. The rest is constant d4start[2] = 0; d4start[3] = 0; d4count[0] = 1; d4count[1] = 1; d4count[2] = global_domain.n_ny; d4count[3] = global_domain.n_nx; // Leaf Area Index: lai if (options.LAI_SRC == FROM_VEGHIST) { for (j = 0; j < NF; j++) { d4start[0] = global_param.forceskip[1] + global_param.forceoffset[1] + j; for (v = 0; v < options.NVEGTYPES; v++) { d4start[1] = v; get_scatter_nc_field_double(filenames.forcing[1], "lai", d4start, d4count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { vidx = veg_con_map[i].vidx[v]; if (vidx != NODATA_VEG) { veg_hist[i][vidx].LAI[j] = (double) dvar[i]; } } } } } // Partial veg cover fraction: fcov if (options.FCAN_SRC == FROM_VEGHIST) { for (j = 0; j < NF; j++) { d4start[0] = global_param.forceskip[1] + global_param.forceoffset[1] + j; for (v = 0; v < options.NVEGTYPES; v++) { d4start[1] = v; get_scatter_nc_field_double(filenames.forcing[1], "fcov", d4start, d4count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { vidx = veg_con_map[i].vidx[v]; if (vidx != NODATA_VEG) { veg_hist[i][vidx].fcanopy[j] = (double) dvar[i]; } } } } } // Albedo: alb if (options.ALB_SRC == FROM_VEGHIST) { for (j = 0; j < NF; j++) { d4start[0] = global_param.forceskip[1] + global_param.forceoffset[1] + j; for (v = 0; v < options.NVEGTYPES; v++) { d4start[1] = v; get_scatter_nc_field_double(filenames.forcing[1], "alb", d4start, d4count, dvar); for (i = 0; i < local_domain.ncells_active; i++) { vidx = veg_con_map[i].vidx[v]; if (vidx != NODATA_VEG) { veg_hist[i][vidx].albedo[j] = (double) dvar[i]; } } } } } // Update the offset counter global_param.forceoffset[1] += NF; } // allocate memory for t_offset t_offset = malloc(local_domain.ncells_active * sizeof(*t_offset)); check_alloc_status(t_offset, "Memory allocation error."); for (i = 0; i < local_domain.ncells_active; i++) { if (options.SNOW_BAND > 1) { Tfactor = soil_con[i].Tfactor; t_offset[i] = Tfactor[0]; for (band = 1; band < options.SNOW_BAND; band++) { if (Tfactor[band] < t_offset[i]) { t_offset[i] = Tfactor[band]; } } } else { t_offset[i] = 0; } } // Convert forcings into what we need and calculate missing ones for (i = 0; i < local_domain.ncells_active; i++) { for (j = 0; j < NF; j++) { // pressure in Pa force[i].pressure[j] *= PA_PER_KPA; // vapor pressure in Pa force[i].vp[j] *= PA_PER_KPA; // vapor pressure deficit in Pa force[i].vpd[j] = svp(force[i].air_temp[j]) - force[i].vp[j]; if (force[i].vpd[j] < 0) { force[i].vpd[j] = 0; force[i].vp[j] = svp(force[i].air_temp[j]); } // air density in kg/m3 force[i].density[j] = air_density(force[i].air_temp[j], force[i].pressure[j]); // snow flag force[i].snowflag[j] = will_it_snow(&(force[i].air_temp[j]), t_offset[i], param.SNOW_MAX_SNOW_TEMP, &(force[i].prec[j]), 1); } // Check on fcanopy for (v = 0; v < options.NVEGTYPES; v++) { vidx = veg_con_map[i].vidx[v]; if (vidx != NODATA_VEG) { for (j = 0; j < NF; j++) { if ((veg_hist[i][vidx].fcanopy[j] < MIN_FCANOPY) && ((current == 0) || (options.FCAN_SRC == FROM_VEGHIST))) { // Only issue this warning once if not using veg hist fractions log_warn( "cell %zu, veg` %d substep %zu fcanopy %f < minimum of %f; setting = %f", i, vidx, j, veg_hist[i][vidx].fcanopy[j], MIN_FCANOPY, MIN_FCANOPY); veg_hist[i][vidx].fcanopy[j] = MIN_FCANOPY; } } } } } // Put average value in NR field for (i = 0; i < local_domain.ncells_active; i++) { force[i].air_temp[NR] = average(force[i].air_temp, NF); // For precipitation put total force[i].prec[NR] = average(force[i].prec, NF) * NF; force[i].shortwave[NR] = average(force[i].shortwave, NF); force[i].longwave[NR] = average(force[i].longwave, NF); force[i].pressure[NR] = average(force[i].pressure, NF); force[i].wind[NR] = average(force[i].wind, NF); force[i].vp[NR] = average(force[i].vp, NF); force[i].vpd[NR] = (svp(force[i].air_temp[NR]) - force[i].vp[NR]); force[i].density[NR] = air_density(force[i].air_temp[NR], force[i].pressure[NR]); force[i].snowflag[NR] = will_it_snow(force[i].air_temp, t_offset[i], param.SNOW_MAX_SNOW_TEMP, force[i].prec, NF); for (v = 0; v < options.NVEGTYPES; v++) { vidx = veg_con_map[i].vidx[v]; if (vidx != NODATA_VEG) { // not the correct way to calculate average albedo in general, // but leave for now (it's correct if albedo is constant over // the model step) veg_hist[i][vidx].albedo[NR] = average(veg_hist[i][vidx].albedo, NF); veg_hist[i][vidx].displacement[NR] = average( veg_hist[i][vidx].displacement, NF); veg_hist[i][vidx].fcanopy[NR] = average( veg_hist[i][vidx].fcanopy, NF); veg_hist[i][vidx].LAI[NR] = average(veg_hist[i][vidx].LAI, NF); veg_hist[i][vidx].roughness[NR] = average( veg_hist[i][vidx].roughness, NF); } } // Optional inputs if (options.LAKES) { force[i].channel_in[NR] = average(force[i].channel_in, NF) * NF; } if (options.CARBON) { force[i].Catm[NR] = average(force[i].Catm, NF); force[i].fdir[NR] = average(force[i].fdir, NF); force[i].par[NR] = average(force[i].par, NF); // for coszen, use value at noon force[i].coszen[NR] = compute_coszen( local_domain.locations[i].latitude, local_domain.locations[i].longitude, soil_con[i].time_zone_lng, dmy[current].day_in_year, SEC_PER_DAY / 2); } } // cleanup free(dvar); }
/****************************************************************************** * @brief Initialize atmospheric variables for the model and snow time steps. *****************************************************************************/ void vic_force(atmos_data_struct *atmos, dmy_struct *dmy, FILE **infile, veg_con_struct *veg_con, veg_hist_struct **veg_hist, soil_con_struct *soil_con) { extern option_struct options; extern param_set_struct param_set; extern global_param_struct global_param; extern parameters_struct param; extern size_t NR, NF; size_t i; size_t j; size_t v; size_t rec; size_t uidx; double t_offset; double **forcing_data; double ***veg_hist_data; double avgJulyAirTemp; double *Tfactor; bool *AboveTreeLine; /******************************* Check that required inputs were supplied *******************************/ if (!param_set.TYPE[AIR_TEMP].SUPPLIED) { log_err("Air temperature must be supplied as a forcing"); } if (!param_set.TYPE[PREC].SUPPLIED) { log_err("Precipitation must be supplied as a forcing"); } if (!param_set.TYPE[SWDOWN].SUPPLIED) { log_err("Downward shortwave radiation must be supplied as a forcing"); } if (!param_set.TYPE[LWDOWN].SUPPLIED) { log_err("Downward longwave radiation must be supplied as a forcing"); } if (!param_set.TYPE[PRESSURE].SUPPLIED) { log_err("Atmospheric pressure must be supplied as a forcing"); } if (!param_set.TYPE[VP].SUPPLIED) { log_err("Vapor ressure must be supplied as a forcing"); } if (!param_set.TYPE[WIND].SUPPLIED) { log_err("Wind speed must be supplied as a forcing"); } /******************************* Miscellaneous initialization *******************************/ /* Assign local copies of some variables */ avgJulyAirTemp = soil_con->avgJulyAirTemp; Tfactor = soil_con->Tfactor; AboveTreeLine = soil_con->AboveTreeLine; /* Assign N_ELEM for veg-dependent forcings */ if (param_set.TYPE[ALBEDO].SUPPLIED) { param_set.TYPE[ALBEDO].N_ELEM = veg_con[0].vegetat_type_num; } if (param_set.TYPE[LAI_IN].SUPPLIED) { param_set.TYPE[LAI_IN].N_ELEM = veg_con[0].vegetat_type_num; } if (param_set.TYPE[FCANOPY].SUPPLIED) { param_set.TYPE[FCANOPY].N_ELEM = veg_con[0].vegetat_type_num; } /******************************* read in meteorological data *******************************/ forcing_data = read_forcing_data(infile, global_param, &veg_hist_data); log_info("Read meteorological forcing file"); /**************************************************** Variables in the atmos_data structure ****************************************************/ t_offset = Tfactor[0]; for (i = 1; i < options.SNOW_BAND; i++) { if (Tfactor[i] < t_offset) { t_offset = Tfactor[i]; } } for (rec = 0; rec < global_param.nrecs; rec++) { for (i = 0; i < NF; i++) { uidx = rec * NF + i; // temperature in Celsius atmos[rec].air_temp[i] = forcing_data[AIR_TEMP][uidx]; // precipitation in mm/period atmos[rec].prec[i] = forcing_data[PREC][uidx]; // downward shortwave in W/m2 atmos[rec].shortwave[i] = forcing_data[SWDOWN][uidx]; // downward longwave in W/m2 atmos[rec].longwave[i] = forcing_data[LWDOWN][uidx]; // pressure in Pa atmos[rec].pressure[i] = forcing_data[PRESSURE][uidx] * PA_PER_KPA; // vapor pressure in Pa atmos[rec].vp[i] = forcing_data[VP][uidx] * PA_PER_KPA; // vapor pressure deficit in Pa atmos[rec].vpd[i] = svp(atmos[rec].air_temp[i]) - atmos[rec].vp[i]; // air density in kg/m3 atmos[rec].density[i] = air_density(atmos[rec].air_temp[i], atmos[rec].pressure[i]); // wind speed in m/s atmos[rec].wind[i] = forcing_data[WIND][uidx]; // snow flag atmos[rec].snowflag[i] = will_it_snow(&(atmos[rec].air_temp[i]), t_offset, param.SNOW_MAX_SNOW_TEMP, &(atmos[rec].prec[i]), 1); // Optional inputs if (options.LAKES) { // Channel inflow from upstream (into lake) if (param_set.TYPE[CHANNEL_IN].SUPPLIED) { atmos[rec].channel_in[i] = forcing_data[CHANNEL_IN][uidx]; } else { atmos[rec].channel_in[i] = 0; } } if (options.CARBON) { // Atmospheric CO2 concentration atmos[rec].Catm[i] = forcing_data[CATM][uidx]; // Fraction of shortwave that is direct atmos[rec].fdir[i] = forcing_data[FDIR][uidx]; // photosynthetically active radiation atmos[rec].par[i] = forcing_data[PAR][uidx]; // Cosine of solar zenith angle atmos[rec].coszen[i] = compute_coszen(soil_con->lat, soil_con->lng, soil_con->time_zone_lng, dmy[rec].day_in_year, dmy[rec].dayseconds); } } if (NF > 1) { atmos[rec].air_temp[NR] = average(atmos[rec].air_temp, NF); // For precipitation put total atmos[rec].prec[NR] = average(atmos[rec].prec, NF) * NF; atmos[rec].shortwave[NR] = average(atmos[rec].shortwave, NF); atmos[rec].longwave[NR] = average(atmos[rec].longwave, NF); atmos[rec].pressure[NR] = average(atmos[rec].pressure, NF); atmos[rec].vp[NR] = average(atmos[rec].vp, NF); atmos[rec].vpd[NR] = average(atmos[rec].vpd, NF); atmos[rec].density[NR] = average(atmos[rec].density, NF); atmos[rec].wind[NR] = average(atmos[rec].wind, NF); atmos[rec].snowflag[NR] = false; for (i = 0; i < NF; i++) { if (atmos[rec].snowflag[i] == true) { atmos[rec].snowflag[NR] = true; } } if (options.LAKES) { atmos[rec].channel_in[NR] = average(atmos[rec].channel_in, NF) * NF; } if (options.CARBON) { atmos[rec].Catm[NR] = average(atmos[rec].Catm, NF); atmos[rec].fdir[NR] = average(atmos[rec].fdir, NF); atmos[rec].par[NR] = average(atmos[rec].par, NF); // for coszen, use value at noon atmos[rec].coszen[NR] = compute_coszen(soil_con->lat, soil_con->lng, soil_con->time_zone_lng, dmy[rec].day_in_year, SEC_PER_DAY / 2); } } } /**************************************************** Variables in the veg_hist structure ****************************************************/ /* First, assign default climatology */ for (rec = 0; rec < global_param.nrecs; rec++) { for (v = 0; v <= veg_con[0].vegetat_type_num; v++) { for (i = 0; i < NF; i++) { veg_hist[rec][v].albedo[i] = veg_con[v].albedo[dmy[rec].month - 1]; veg_hist[rec][v].displacement[i] = veg_con[v].displacement[dmy[rec].month - 1]; veg_hist[rec][v].fcanopy[i] = veg_con[v].fcanopy[dmy[rec].month - 1]; veg_hist[rec][v].LAI[i] = veg_con[v].LAI[dmy[rec].month - 1]; veg_hist[rec][v].roughness[i] = veg_con[v].roughness[dmy[rec].month - 1]; } } } /* Next, overwrite with veg_hist values, validate, and average */ for (rec = 0; rec < global_param.nrecs; rec++) { for (v = 0; v <= veg_con[0].vegetat_type_num; v++) { for (i = 0; i < NF; i++) { uidx = rec * NF + i; if (param_set.TYPE[ALBEDO].SUPPLIED && options.ALB_SRC == FROM_VEGHIST) { if (veg_hist_data[ALBEDO][v][uidx] != NODATA_VH) { veg_hist[rec][v].albedo[i] = veg_hist_data[ALBEDO][v][uidx]; } } if (param_set.TYPE[LAI_IN].SUPPLIED && options.LAI_SRC == FROM_VEGHIST) { if (veg_hist_data[LAI_IN][v][uidx] != NODATA_VH) { veg_hist[rec][v].LAI[i] = veg_hist_data[LAI_IN][v][uidx]; } } if (param_set.TYPE[FCANOPY].SUPPLIED && options.FCAN_SRC == FROM_VEGHIST) { if (veg_hist_data[FCANOPY][v][uidx] != NODATA_VH) { veg_hist[rec][v].fcanopy[i] = veg_hist_data[FCANOPY][v][uidx]; } } // Check on fcanopy if (veg_hist[rec][v].fcanopy[i] < MIN_FCANOPY) { log_warn( "rec %zu, veg %zu substep %zu fcanopy %f < minimum of %f; setting = %f\n", rec, v, i, veg_hist[rec][v].fcanopy[i], MIN_FCANOPY, MIN_FCANOPY); veg_hist[rec][v].fcanopy[i] = MIN_FCANOPY; } } if (NF > 1) { veg_hist[rec][v].albedo[NR] = average(veg_hist[rec][v].albedo, NF); veg_hist[rec][v].displacement[NR] = average( veg_hist[rec][v].displacement, NF); veg_hist[rec][v].fcanopy[NR] = average( veg_hist[rec][v].fcanopy, NF); veg_hist[rec][v].LAI[NR] = average(veg_hist[rec][v].LAI, NF); veg_hist[rec][v].roughness[NR] = average( veg_hist[rec][v].roughness, NF); } } } /**************************************************** Free forcing_data and veg_hist_data structures ****************************************************/ for (i = 0; i < N_FORCING_TYPES; i++) { if (param_set.TYPE[i].SUPPLIED) { if (i != ALBEDO && i != LAI_IN && i != FCANOPY) { free(forcing_data[i]); } else { for (j = 0; j < param_set.TYPE[i].N_ELEM; j++) { free(veg_hist_data[i][j]); } free(veg_hist_data[i]); } } } free(forcing_data); free(veg_hist_data); /**************************************************** Compute treeline based on July average temperature ****************************************************/ if (options.COMPUTE_TREELINE) { if (!(options.JULY_TAVG_SUPPLIED && avgJulyAirTemp == -999)) { compute_treeline(atmos, dmy, avgJulyAirTemp, Tfactor, AboveTreeLine); } } }