Skip to content
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

Remove redundant variable and optimize the copy assignment for the CUDA matrix #636

Merged
merged 5 commits into from
Aug 29, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 5 additions & 3 deletions include/micm/cuda/util/cuda_dense_matrix.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -129,10 +129,12 @@ namespace micm
CudaDenseMatrix& operator=(const CudaDenseMatrix& other)
{
VectorMatrix<T, L>::operator=(other);
if (this->param_.d_data_ != nullptr)
if (this->param_.number_of_elements_ != other.param_.number_of_elements_)
{
CHECK_CUDA_ERROR(micm::cuda::FreeVector(this->param_), "cudaFree");
this->param_ = other.param_;
CHECK_CUDA_ERROR(micm::cuda::MallocVector<T>(this->param_, this->param_.number_of_elements_), "cudaMalloc");
this->param_ = other.param_;
CHECK_CUDA_ERROR(micm::cuda::MallocVector<T>(this->param_, this->param_.number_of_elements_), "cudaMalloc");
}
CHECK_CUDA_ERROR(micm::cuda::CopyToDeviceFromDevice<T>(this->param_, other.param_), "cudaMemcpyDeviceToDevice");
return *this;
}
Expand Down
20 changes: 11 additions & 9 deletions include/micm/solver/rosenbrock.inl
Original file line number Diff line number Diff line change
Expand Up @@ -16,9 +16,8 @@ namespace micm
MatrixPolicy Y = state.variables_;
std::size_t num_rows = Y.NumRows();
std::size_t num_cols = Y.NumColumns();
MatrixPolicy Ynew(num_rows, num_cols, 0.0);
MatrixPolicy initial_forcing(num_rows, num_cols, 0.0);
MatrixPolicy forcing(num_rows, num_cols, 0.0);
MatrixPolicy Ynew(num_rows, num_cols);
MatrixPolicy initial_forcing(num_rows, num_cols);
std::vector<MatrixPolicy> K{};
const double h_max = parameters_.h_max_ == 0.0 ? time_step : std::min(time_step, parameters_.h_max_);
const double h_start =
Expand All @@ -28,7 +27,7 @@ namespace micm

K.reserve(parameters_.stages_);
for (std::size_t i = 0; i < parameters_.stages_; ++i)
K.emplace_back(num_rows, num_cols, 0.0);
K.emplace_back(num_rows, num_cols);

double present_time = 0.0;
double H = std::min(std::max(std::abs(parameters_.h_min_), std::abs(h_start)), std::abs(h_max));
Expand Down Expand Up @@ -59,7 +58,7 @@ namespace micm
// Limit H if necessary to avoid going beyond the specified chemistry time step
H = std::min(H, std::abs(time_step - present_time));

// compute the forcing at the beginning of the current time
// compute the initial forcing at the beginning of the current time
initial_forcing.Fill(0.0);
rates_.AddForcingTerms(state.rate_constants_, Y, initial_forcing);
stats.function_calls_ += 1;
Expand Down Expand Up @@ -88,7 +87,7 @@ namespace micm
double stage_combinations = ((stage + 1) - 1) * ((stage + 1) - 2) / 2;
if (stage == 0)
{
forcing = initial_forcing;
K[stage].Copy(initial_forcing);
}
else
{
Expand All @@ -99,12 +98,15 @@ namespace micm
{
Ynew.Axpy(parameters_.a_[stage_combinations + j], K[j]);
}
forcing.Fill(0.0);
rates_.AddForcingTerms(state.rate_constants_, Ynew, forcing);
K[stage].Fill(0.0);
rates_.AddForcingTerms(state.rate_constants_, Ynew, K[stage]);
stats.function_calls_ += 1;
}
}
K[stage].Copy(forcing);
if (stage + 1 < parameters_.stages_ && !parameters_.new_function_evaluation_[stage + 1])
{
K[stage + 1].Copy(K[stage]);
}
for (uint64_t j = 0; j < stage; ++j)
{
K[stage].Axpy(parameters_.c_[stage_combinations + j] / H, K[j]);
Expand Down
Loading