Skip to content

OpenMP test branch #349

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

Draft
wants to merge 2 commits into
base: master
Choose a base branch
from
Draft
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
3 changes: 3 additions & 0 deletions .gitattributes
Original file line number Diff line number Diff line change
@@ -1,2 +1,5 @@
*.ipynb filter=nbstripout
*.ipynb diff=ipynb

examples/zundel_i-PI.ipynb filter= diff=

2 changes: 1 addition & 1 deletion CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -81,7 +81,7 @@ SET(TYPE_ARCHITECTURE "native" CACHE STRING
)

########## COMPILATION FLAGS ##########
set(CMAKE_CXX_FLAGS "${CMAKE_CXX_FLAGS} -Wall -Wextra -Weffc++ -Wno-non-virtual-dtor")
set(CMAKE_CXX_FLAGS "${CMAKE_CXX_FLAGS} -Wall -Wextra -Weffc++ -Wno-non-virtual-dtor -fopenmp")
set(CMAKE_C_FLAGS "${CMAKE_C_FLAGS} -Wall -Wextra")

if ("${CMAKE_CXX_COMPILER_ID}" STREQUAL "GNU" OR
Expand Down
2 changes: 2 additions & 0 deletions src/rascal/models/sparse_kernel_predict.hh
Original file line number Diff line number Diff line change
Expand Up @@ -127,6 +127,8 @@ namespace rascal {

math::Vector_t weights_scaled(weights.size());
size_t i_row{0};

#pragma omp for
for (auto center : manager) {
const int a_sp{center.get_atom_type()};
// compute contraction of the model weights with the gradient of the
Expand Down
4 changes: 4 additions & 0 deletions src/rascal/models/sparse_kernels.hh
Original file line number Diff line number Diff line change
Expand Up @@ -178,6 +178,7 @@ namespace rascal {
for (auto & manager : managers) {
auto && propA{*manager->template get_property<Property_t>(
representation_name, true)};
#pragma omp for
for (auto center : manager) {
int sp = center.get_atom_type();
KNM.row(ii_A) =
Expand Down Expand Up @@ -321,6 +322,7 @@ namespace rascal {
if (zeta > 1) {
dkdX_missing_T.setZero();
// zeta * (X * T)**(zeta-1)
#pragma omp for
for (auto center : manager) {
int a_species{center.get_atom_type()};
dkdX_missing_T[center] =
Expand Down Expand Up @@ -362,6 +364,7 @@ namespace rascal {
if (do_block_by_key_dot) {
size_t idx_row{0};
auto repr_grads = prop_repr_grad.get_raw_data_view();
#pragma omp for
for (auto center : manager) {
int a_species{center.get_atom_type()};
Eigen::Vector3d r_i = center.get_position();
Expand Down Expand Up @@ -447,6 +450,7 @@ namespace rascal {
idx_row += nb_rows;
} // center
} else {
#pragma omp for
for (auto center : manager) {
int a_species{center.get_atom_type()};
Eigen::Vector3d r_i = center.get_position();
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -1791,6 +1791,7 @@ namespace rascal {
// coeff C^{ij}_{nlm}
auto c_ij_nlm = math::Matrix_t(n_row, n_col);

#pragma omp for
for (auto center : manager) {
// c^{i}
auto & coefficients_center = expansions_coefficients[center];
Expand Down
10 changes: 9 additions & 1 deletion src/rascal/representations/calculator_spherical_invariants.hh
Original file line number Diff line number Diff line change
Expand Up @@ -605,6 +605,7 @@ namespace rascal {
Eigen::Vector3d soap_vector_dot_gradient{};

// compute the dot product and update the gradients to be normalized
#pragma omp for
for (auto center : manager) {
for (auto neigh : center.pairs_with_self_pair()) {
const auto & soap_vector = soap_vectors[center];
Expand Down Expand Up @@ -744,6 +745,7 @@ namespace rascal {
*manager, "power spectrums inverse norms", true};
soap_vector_norm_inv.resize();

#pragma omp for
for (auto center : manager) {
auto & coefficients{expansions_coefficients[center]};
auto & soap_vector{soap_vectors[center]};
Expand Down Expand Up @@ -976,6 +978,7 @@ namespace rascal {
*manager, "radial spectrum inverse norms", ExcludeGhosts};
soap_vector_norm_inv.resize();

#pragma omp for
for (auto center : manager) {
const auto & coefficients{expansions_coefficients[center]};
auto & soap_vector{soap_vectors[center]};
Expand Down Expand Up @@ -1048,7 +1051,12 @@ namespace rascal {
double mult{1.0};
Key_t trip_type{0, 0, 0};
internal::SortedKey<Key_t> triplet_type{trip_type};
for (auto center : manager) {

#pragma omp for
for (auto iter = manager->begin(); iter<manager->end(); ++iter) {
std::cerr<<"computing "<<iter.get_index()<< " out of "<< manager->end().get_index()<<"\n";
auto center = *iter;
// for (auto center : manager) {
auto & coefficients{expansions_coefficients[center]};
auto & soap_vector{soap_vectors[center]};

Expand Down
30 changes: 30 additions & 0 deletions src/rascal/structure_managers/structure_manager.hh
Original file line number Diff line number Diff line change
Expand Up @@ -1387,6 +1387,20 @@ namespace rascal {
return *this;
}

//! int increment
Iterator & operator+=(size_t increment) {
this->index += increment;
return *this;
}

//! int decrement
Iterator & operator-=(size_t increment) {
this->index -= increment;
return *this;
}
//! subtraction between iterators
auto operator-(const Iterator & other) { return this->index - other.index; }

value_type operator*() {
auto & cluster_indices_properties = std::get<Order - 1>(
this->get_manager().get_cluster_indices_container());
Expand Down Expand Up @@ -1420,6 +1434,22 @@ namespace rascal {
return not(*this == other);
}

//! comparison operators
bool operator<(const Iterator & other) { return this->index < other.index; }

bool operator>(const Iterator & other) { return this->index > other.index; }

bool operator<=(const Iterator & other) {
return this->index <= other.index;
}

bool operator>=(const Iterator & other) {
return this->index >= other.index;
}

size_t get_index() {
return this->index;
}
/**
* const access to container
*/
Expand Down