-
Notifications
You must be signed in to change notification settings - Fork 24
Improve runtime and peak memory usage of PVEngine.run_full_mode #140
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
base: master
Are you sure you want to change the base?
Changes from 6 commits
38922c8
4450fcc
077a160
e0e1e3e
5268d5f
12f93e5
0d2e6f6
480b350
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
|
|
@@ -2,11 +2,51 @@ | |
| timeseries simulations.""" | ||
|
|
||
| import numpy as np | ||
| from scipy.sparse import csc_matrix | ||
| from scipy.sparse.linalg import spsolve | ||
| from pvfactors.viewfactors import VFCalculator | ||
| from pvfactors.irradiance import HybridPerezOrdered | ||
| from pvfactors.config import DEFAULT_RHO_FRONT, DEFAULT_RHO_BACK | ||
|
|
||
|
|
||
| def _sparse_solve_3D(A, b): | ||
| """ | ||
| Solve the linear system A*x=b for x, where A is sparse (contains mostly | ||
| zeros). | ||
|
|
||
| For large matrices with many zeros, this is faster and uses less memory | ||
| than the dense analog `np.linalg.solve`. | ||
|
|
||
| Parameters | ||
| ---------- | ||
| A : np.ndarray of shape (M, N, N) | ||
| First dimension is timestamps, second and third are surface dimensions | ||
| b : np.ndarray of shape (M, N) | ||
| First dimension is timestamps, second is surfaces | ||
|
|
||
| Returns | ||
| ------- | ||
| x : np.ndarray of shape (M, N) | ||
| First dimension is timestamps, second is surfaces | ||
| """ | ||
| # implementation notes: | ||
| # - csc_matrix seemed to be the fastest option of the various | ||
| # sparse matrix formats in scipy.sparse | ||
| # - unfortunately the sparse matrix formats are 2-D only, so | ||
| # iteration across the time dimension is required | ||
| # - scipy 1.8.0 added new sparse arrays (as opposed to sparse matrices); | ||
| # they are 2-D only at time of writing, but in the future if they | ||
| # become N-D then it may make sense to use them here. | ||
|
||
| xs = [] | ||
| for A_slice, b_slice in zip(A, b): | ||
| A_sparse = csc_matrix(A_slice) | ||
| b_sparse = csc_matrix(b_slice).T | ||
| x = spsolve(A_sparse, b_sparse) | ||
| xs.append(x) | ||
| x = np.stack(xs) | ||
| return x | ||
|
|
||
|
|
||
| class PVEngine(object): | ||
| """Class putting all of the calculations together into simple | ||
| workflows. | ||
|
|
@@ -215,33 +255,36 @@ def run_full_mode(self, fn_build_report=None): | |
| invrho_ts_diag = np.zeros(shape_system) | ||
| for idx_surf in range(shape_system[1]): | ||
| invrho_ts_diag[:, idx_surf, idx_surf] = invrho_mat[idx_surf, :] | ||
| del invrho_mat | ||
| # Subtract matrices: will rely on broadcasting | ||
| # shape = n_timesteps, n_surfaces + 1, n_surfaces + 1 | ||
| a_mat = invrho_ts_diag - ts_vf_matrix_reshaped | ||
| # Calculate inverse, requires specific shape | ||
| # shape = n_timesteps, n_surfaces + 1, n_surfaces + 1 | ||
| inv_a_mat = np.linalg.inv(a_mat) | ||
| # Use einstein sum to get final timeseries radiosities | ||
| # shape = n_surfaces + 1, n_timesteps | ||
| q0 = np.einsum('ijk,ki->ji', inv_a_mat, irradiance_mat) | ||
| del ts_vf_matrix_reshaped | ||
| # solve the linear system a_mat * q0 = irradiance_mat for q0 | ||
| q0 = _sparse_solve_3D(a_mat, irradiance_mat.T).T | ||
| del a_mat | ||
| # Calculate incident irradiance: will rely on broadcasting | ||
| # shape = n_surfaces + 1, n_timesteps | ||
| qinc = np.einsum('ijk,ki->ji', invrho_ts_diag, q0) | ||
| del invrho_ts_diag | ||
|
|
||
| # --- Derive other irradiance terms | ||
| # shape = n_surfaces, n_timesteps | ||
| isotropic_mat = ts_vf_matrix[:-1, -1, :] * irradiance_mat[-1, :] | ||
| del ts_vf_matrix | ||
| reflection_mat = qinc[:-1, :] - irradiance_mat[:-1, :] - isotropic_mat | ||
|
|
||
| del irradiance_mat | ||
| # --- Calculate AOI losses and absorbed irradiance | ||
| # Create tiled reflection matrix of | ||
| # shape = n_surfaces + 1, n_surfaces + 1, n_timestamps | ||
| rho_ts_tiled = np.moveaxis(np.tile(rho_mat.T, (shape_system[1], 1, 1)), | ||
| -1, 0) | ||
| del rho_mat | ||
| # Get vf AOI matrix | ||
| # shape [n_surfaces + 1, n_surfaces + 1, n_timestamps] | ||
| vf_aoi_matrix = (self.vf_calculator | ||
| .build_ts_vf_aoi_matrix(pvarray, rho_ts_tiled)) | ||
| del rho_ts_tiled | ||
| pvarray.ts_vf_aoi_matrix = vf_aoi_matrix | ||
| # Get absorbed irradiance matrix | ||
| # shape [n_surfaces, n_timestamps] | ||
|
|
@@ -250,6 +293,8 @@ def run_full_mode(self, fn_build_report=None): | |
| # Calculate absorbed irradiance | ||
| qabs = (np.einsum('ijk,jk->ik', vf_aoi_matrix, q0)[:-1, :] | ||
| + irradiance_abs_mat) | ||
| del irradiance_abs_mat | ||
| del vf_aoi_matrix | ||
|
|
||
| # --- Update surfaces with values: the lists are ordered by index | ||
| for idx_surf, ts_surface in enumerate(pvarray.all_ts_surfaces): | ||
|
Contributor
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. awesome 👍 it doesn't shock me to have the |
||
|
|
||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -1,4 +1,4 @@ | ||
| from pvfactors.engine import PVEngine | ||
| from pvfactors.engine import PVEngine, _sparse_solve_3D | ||
| from pvfactors.geometry.pvarray import OrderedPVArray | ||
| from pvfactors.irradiance import IsotropicOrdered, HybridPerezOrdered | ||
| from pvfactors.irradiance.utils import breakup_df_inputs | ||
|
|
@@ -709,3 +709,13 @@ def test_engine_variable_albedo(params, df_inputs_clearsky_8760): | |
| expected_bfg_after_aoi = 14.670709 | ||
| np.testing.assert_allclose(bfg, expected_bfg) | ||
| np.testing.assert_allclose(bfg_after_aoi, expected_bfg_after_aoi) | ||
|
|
||
|
|
||
| def test__sparse_solve_3d(): | ||
| """Verify the sparse solver returns same results as np.linalg.solve""" | ||
| A = np.random.rand(5, 3, 3) | ||
| b = np.random.rand(5, 3) | ||
| x_dense = np.linalg.solve(A, b) | ||
| x_sparse = _sparse_solve_3D(A, b) | ||
| assert x_sparse.shape == b.shape | ||
| np.testing.assert_allclose(x_dense, x_sparse) | ||
|
Contributor
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. awesome 👍 |
||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -1,5 +1,6 @@ | ||
| pvlib>=0.9.0,<0.10.0 | ||
| shapely>=1.6.4.post2,<2 | ||
| scipy>=1.2.0 | ||
| matplotlib | ||
| future | ||
| six |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
any chance you could add the arg and returned types as type hints please? (instead of in the docstrings) most of the code doesn't have it but maybe we could start adding it for new implementations
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Done, but I've not used type hints before so please let me know if I've messed it up :)