Skip to content

Commit 2030c78

Browse files
committed
Merge branch 'develop' into feature-imgui
2 parents 0ae9348 + 11fe938 commit 2030c78

File tree

14 files changed

+329
-247
lines changed

14 files changed

+329
-247
lines changed

core/include/Spirit/System.h

+6-1
Original file line numberDiff line numberDiff line change
@@ -53,8 +53,13 @@ PREFIX float System_Get_Rx(State * state, int idx_image=-1, int idx_chain=-1) SU
5353
// Returns the energy of a spin system.
5454
PREFIX float System_Get_Energy(State * state, int idx_image=-1, int idx_chain=-1) SUFFIX;
5555

56+
// Retrieves the names of the energy contributions, represented as a single string and separated by "|". E.g "Zeeman|Exchange|DMI"
57+
// If 'names' is a nullptr, the required length of the char array is returned.
58+
PREFIX int System_Get_Energy_Array_Names(State * state, char* names, int idx_image=-1, int idx_chain=-1) SUFFIX;
59+
5660
// Retrieves the energy contributions of a spin system.
57-
PREFIX void System_Get_Energy_Array(State * state, float * energies, int idx_image=-1, int idx_chain=-1) SUFFIX;
61+
// If 'energies' is a nullptr, the required length of the energies array is returned.
62+
PREFIX int System_Get_Energy_Array(State * state, float * energies, bool divide_by_nspins=true, int idx_image=-1, int idx_chain=-1) SUFFIX;
5863

5964
// Retrieves the eigenvalues of a spin system
6065
PREFIX void System_Get_Eigenvalues(State * state, float * eigenvalues, int idx_image=-1, int idx_chain=-1) SUFFIX;

core/include/engine/FFT.hpp

+8-7
Original file line numberDiff line numberDiff line change
@@ -3,6 +3,7 @@
33
#define FFT_H
44
#include <Eigen/Core>
55
#include <data/Geometry.hpp>
6+
#include <utility/Logging.hpp>
67
#include <engine/Vectormath_Defines.hpp>
78
#include <iostream>
89
#include <complex>
@@ -25,7 +26,7 @@
2526
#include <cufft.h>
2627
#endif
2728

28-
namespace Engine
29+
namespace Engine
2930
{
3031
namespace FFT
3132
{
@@ -84,7 +85,7 @@ namespace Engine
8485
FFT_cpx_type res;
8586
res[0] = d1[0] * s1[0] + d2[0] * s2[0] + d3[0] * s3[0] - d1[1] * s1[1] - d2[1] * s2[1] - d3[1] * s3[1];
8687
res[1] = d1[0] * s1[1] + d2[0] * s2[1] + d3[0] * s3[1] + d1[1] * s1[0] + d2[1] * s2[0] + d3[1] * s3[0];
87-
return res;
88+
return res;
8889
}
8990

9091
inline void addTo(FFT_cpx_type & a, const FFT_cpx_type & b, bool overwrite)
@@ -112,7 +113,7 @@ namespace Engine
112113
FFT_cpx_type res;
113114
res.x = d1.x * s1.x + d2.x * s2.x + d3.x * s3.x - d1.y * s1.y - d2.y * s2.y - d3.y * s3.y;
114115
res.y = d1.x * s1.y + d2.x * s2.y + d3.x * s3.y + d1.y * s1.x + d2.y * s2.x + d3.y * s3.x;
115-
return res;
116+
return res;
116117
}
117118

118119
inline __device__ void addTo(FFT_cpx_type & a, const FFT_cpx_type & b, bool overwrite = false)
@@ -135,7 +136,7 @@ namespace Engine
135136
fftw_plan_with_nthreads(omp_get_max_threads());
136137
#endif
137138
}
138-
139+
139140
inline void get_strides(field<int*> & strides, const field<int> & maxVal)
140141
{
141142
strides.resize(maxVal.size());
@@ -173,7 +174,7 @@ namespace Engine
173174

174175
//Constructor delegation
175176
FFT_Plan() : FFT_Plan({2,2,2}, true, 1, 8)
176-
{}
177+
{}
177178

178179
FFT_Plan(std::vector<int> dims, bool inverse, int n_transforms, int len) :
179180
dims(dims),
@@ -187,7 +188,7 @@ namespace Engine
187188
}
188189

189190
//copy constructor
190-
FFT_Plan(FFT_Plan const & other)
191+
FFT_Plan(FFT_Plan const & other)
191192
{
192193
this->dims = other.dims;
193194
this->inverse = other.inverse;
@@ -219,7 +220,7 @@ namespace Engine
219220
}
220221
return *this;
221222
}
222-
223+
223224
//move assignment operator
224225
FFT_Plan& operator=(FFT_Plan const && other)
225226
{

core/include/engine/Sparse_HTST.hpp

+1-1
Original file line numberDiff line numberDiff line change
@@ -52,7 +52,7 @@ namespace Engine
5252
SpMatrixX & hessian_geodesic_3N, SpMatrixX & hessian_geodesic_2N, SpMatrixX & tangent_basis, VectorX & eigenvalues, MatrixX & eigenvectors);
5353

5454
void sparse_hessian_bordered_3N(const vectorfield & image, const vectorfield & gradient, const SpMatrixX & hessian, SpMatrixX & hessian_out);
55-
55+
5656
};
5757
}
5858

core/include/engine/Vectormath.hpp

+2-3
Original file line numberDiff line numberDiff line change
@@ -272,7 +272,7 @@ namespace Engine
272272
// Translations (cell) of spin i
273273
int nic = ispin / (N*Na*Nb);
274274
int nib = (ispin - nic*N*Na*Nb) / (N*Na);
275-
int nia = ispin - nic*N*Na*Nb - nib*N*Na;
275+
int nia = (ispin - nic*N*Na*Nb - nib*N*Na) / N;
276276

277277
// Translations (cell) of spin j (possibly outside of non-periodical domain)
278278
int pm = 1;
@@ -336,8 +336,7 @@ namespace Engine
336336
// Calculate the index of spin j according to it's translations
337337
int jspin = pair.j + (nja)*N + (njb)*N*Na + (njc)*N*Na*Nb;
338338

339-
// Invalid index if atom type of spin j is not correct
340-
if ( pair.j != jspin%N || !cu_check_atom_type(atom_types[jspin]) )
339+
if ( !cu_check_atom_type(atom_types[jspin]) )
341340
return -1;
342341

343342
// Return a valid index

core/python/spirit/system.py

+31-10
Original file line numberDiff line numberDiff line change
@@ -98,16 +98,37 @@ def get_eigenvalues(p_state, idx_image=-1, idx_chain=-1):
9898
_Get_Eigenvalues(ctypes.c_void_p(p_state), eigenvalues, ctypes.c_int(idx_image), ctypes.c_int(idx_chain))
9999
return eigenvalues
100100

101-
# NOTE: excluded since there is no clean way to get the C++ pairs
102-
### Get Energy array
103-
# _Get_Energy_Array = _spirit.System_Get_Energy_Array
104-
# _Get_Energy_Array.argtypes = [ctypes.c_void_p, ctypes.POINTER(ctypes.c_float),
105-
# ctypes.c_int, ctypes.c_int]
106-
# _Get_Energy_Array.restype = None
107-
# def Get_Energy_Array(p_state, idx_image=-1, idx_chain=-1):
108-
# Energies
109-
# _Get_Energy_Array(ctypes.c_void_p(p_state), energies,
110-
# ctypes.c_int(idx_image), ctypes.c_int(idx_chain))
101+
### Get Energy Contributions
102+
### The result is a dictionary with strings as keys and floats as values
103+
### The keys are the names of the energy contributions, the values the energy_contribution in meV
104+
_Get_Energy_Array = _spirit.System_Get_Energy_Array
105+
_Get_Energy_Array.argtypes = [ctypes.c_void_p, ctypes.POINTER(ctypes.c_float), ctypes.c_bool,
106+
ctypes.c_int, ctypes.c_int]
107+
_Get_Energy_Array.restype = None
108+
109+
_Get_Energy_Array_Names = _spirit.System_Get_Energy_Array_Names
110+
_Get_Energy_Array_Names.argtypes = [ctypes.c_void_p, ctypes.POINTER(ctypes.c_char),
111+
ctypes.c_int, ctypes.c_int]
112+
_Get_Energy_Array_Names.restype = ctypes.c_int
113+
def get_energy_contributions(p_state, divide_by_nspins = True, idx_image=-1, idx_chain=-1):
114+
NULL = ctypes.POINTER(ctypes.c_char)()
115+
116+
n_char_array = _Get_Energy_Array_Names(ctypes.c_void_p(p_state), NULL,
117+
ctypes.c_int(idx_image), ctypes.c_int(idx_chain))
118+
119+
energy_array_names = (n_char_array*ctypes.c_char)()
120+
121+
_Get_Energy_Array_Names(ctypes.c_void_p(p_state), energy_array_names,
122+
ctypes.c_int(idx_image), ctypes.c_int(idx_chain))
123+
124+
contrib_names = str(energy_array_names[:].decode("utf-8")).split("|")
125+
n_contribs = len(contrib_names)
126+
energies = (n_contribs*ctypes.c_float)()
127+
128+
_Get_Energy_Array(ctypes.c_void_p(p_state), energies, divide_by_nspins,
129+
ctypes.c_int(idx_image), ctypes.c_int(idx_chain))
130+
131+
return dict(zip(contrib_names, energies))
111132

112133
### Get Chain number of images
113134
_Update_Data = _spirit.System_Update_Data

core/python/test/system.py

+11-3
Original file line numberDiff line numberDiff line change
@@ -5,7 +5,7 @@
55
spirit_py_dir = os.path.abspath(os.path.join(os.path.dirname( __file__ ), ".."))
66
sys.path.insert(0, spirit_py_dir)
77

8-
from spirit import state, system, configuration
8+
from spirit import state, system, configuration, hamiltonian
99

1010
import unittest
1111

@@ -42,8 +42,16 @@ def test_get_spin_directions(self):
4242
def test_get_energy(self):
4343
# NOTE: that test is trivial
4444
E = system.get_energy(self.p_state)
45-
46-
45+
46+
def test_get_energy_contributions(self):
47+
configuration.plus_z(self.p_state)
48+
configuration.domain(self.p_state, [0,0,-1], border_cylindrical=2)
49+
system.update_data(self.p_state)
50+
E_contribs = system.get_energy_contributions(self.p_state, divide_by_nspins=False)
51+
E = system.get_energy(self.p_state)
52+
system.print_energy_array(p_state)
53+
self.assertEqual( len(E_contribs.values()), 3 ) # There should be 3 contributions
54+
self.assertAlmostEqual( sum(E_contribs.values()), E, places=5 ) #TODO: Apparently we can not go higher with the number of decimal places, because the order of summation differs. This Should be invesitgated.
4755
# NOTE: there is no way to test the system.Update_Data() and system.Print_Energy_Array()
4856

4957
#########

core/src/Spirit/System.cpp

+51-2
Original file line numberDiff line numberDiff line change
@@ -133,7 +133,7 @@ catch( ... )
133133
return 0;
134134
}
135135

136-
void System_Get_Energy_Array(State * state, float * energies, int idx_image, int idx_chain) noexcept
136+
int System_Get_Energy_Array_Names(State * state, char* names, int idx_image, int idx_chain) noexcept
137137
try
138138
{
139139
std::shared_ptr<Data::Spin_System> image;
@@ -142,14 +142,63 @@ try
142142
// Fetch correct indices and pointers
143143
from_indices( state, idx_image, idx_chain, image, chain );
144144

145+
int n_char_array = -1; // Start of with offset -1, because the last contributions gets no "|" delimiter
145146
for (unsigned int i=0; i<image->E_array.size(); ++i)
146147
{
147-
energies[i] = (float)image->E_array[i].second;
148+
n_char_array += image->E_array[i].first.size() + 1; // Add +1 because we separate the contribution names with the character "|"
149+
}
150+
151+
// If 'names' is a nullptr, we return the required length of the names array
152+
if(names==nullptr)
153+
{
154+
return n_char_array;
155+
} else { // Else we try to fill the provided char array
156+
int idx=0;
157+
for (unsigned int i=0; i<image->E_array.size(); ++i)
158+
{
159+
for(const char & cur_char : (image->E_array[i]).first)
160+
{
161+
names[idx++] = cur_char;
162+
}
163+
if(i != image->E_array.size()-1)
164+
names[idx++] = '|';
165+
}
166+
return -1;
148167
}
149168
}
150169
catch( ... )
151170
{
152171
spirit_handle_exception_api(idx_image, idx_chain);
172+
return -1;
173+
}
174+
175+
176+
int System_Get_Energy_Array(State * state, float * energies, bool divide_by_nspins, int idx_image, int idx_chain) noexcept
177+
try
178+
{
179+
std::shared_ptr<Data::Spin_System> image;
180+
std::shared_ptr<Data::Spin_System_Chain> chain;
181+
182+
// Fetch correct indices and pointers
183+
from_indices( state, idx_image, idx_chain, image, chain );
184+
185+
scalar nd = divide_by_nspins ? 1/(scalar)image->nos : 1;
186+
187+
if(energies == nullptr)
188+
{
189+
return image->E_array.size();
190+
} else {
191+
for (unsigned int i=0; i<image->E_array.size(); ++i)
192+
{
193+
energies[i] = nd * (float)image->E_array[i].second;
194+
}
195+
return -1;
196+
}
197+
}
198+
catch( ... )
199+
{
200+
spirit_handle_exception_api(idx_image, idx_chain);
201+
return -1;
153202
}
154203

155204
void System_Get_Eigenvalues(State * state, float * eigenvalues, int idx_image, int idx_chain) noexcept

core/src/engine/FFT.cu

+26-5
Original file line numberDiff line numberDiff line change
@@ -13,20 +13,29 @@ namespace Engine
1313
{
1414
std::cerr << "NOT IMPLEMENTED FOR cuFFT" << std::endl;
1515
}
16+
1617
void iFour_3D(FFT_cfg cfg, FFT_cpx_type * in, FFT_real_type * out)
1718
{
1819
std::cerr << "NOT IMPLEMENTED FOR cuFFT" << std::endl;
1920
}
2021

2122
void batch_Four_3D(FFT_Plan & plan)
2223
{
23-
cufftExecR2C(plan.cfg, plan.real_ptr.data(), plan.cpx_ptr.data());
24+
auto res = cufftExecR2C(plan.cfg, plan.real_ptr.data(), plan.cpx_ptr.data());
25+
if(res != CUFFT_SUCCESS)
26+
{
27+
Log(Utility::Log_Level::Error, Utility::Log_Sender::All, fmt::format("cufftExecR2C failed with error: {}", res));
28+
}
2429
cudaDeviceSynchronize();
2530
}
2631

2732
void batch_iFour_3D(FFT_Plan & plan)
2833
{
29-
cufftExecC2R(plan.cfg, plan.cpx_ptr.data(), plan.real_ptr.data());
34+
auto res = cufftExecC2R(plan.cfg, plan.cpx_ptr.data(), plan.real_ptr.data());
35+
if(res != CUFFT_SUCCESS)
36+
{
37+
Log(Utility::Log_Level::Error, Utility::Log_Sender::All, fmt::format("cufftExecC2R failed with error: {}", res));
38+
}
3039
cudaDeviceSynchronize();
3140
}
3241

@@ -47,16 +56,28 @@ namespace Engine
4756

4857
if(this->inverse == false)
4958
{
50-
cufftPlanMany(&this->cfg, rank, n, inembed, istride, idist, onembed, ostride, odist, CUFFT_R2C, n_transforms);
59+
auto res = cufftPlanMany(&this->cfg, rank, n, inembed, istride, idist, onembed, ostride, odist, CUFFT_R2C, n_transforms);
60+
if(res != CUFFT_SUCCESS)
61+
{
62+
Log(Utility::Log_Level::Error, Utility::Log_Sender::All, fmt::format("cufftPlanMany failed with error: {}", res));
63+
}
5164
} else
5265
{
53-
cufftPlanMany(&this->cfg, rank, n, inembed, istride, idist, onembed, ostride, odist, CUFFT_C2R, n_transforms);
66+
auto res = cufftPlanMany(&this->cfg, rank, n, inembed, istride, idist, onembed, ostride, odist, CUFFT_C2R, n_transforms);
67+
if(res != CUFFT_SUCCESS)
68+
{
69+
Log(Utility::Log_Level::Error, Utility::Log_Sender::All, fmt::format("cufftPlanMany failed with error: {}", res));
70+
}
5471
}
5572
}
5673

5774
void FFT_Plan::Free_Configuration()
5875
{
59-
cufftDestroy(this->cfg);
76+
auto res = cufftDestroy(this->cfg);
77+
if(res != CUFFT_SUCCESS)
78+
{
79+
Log(Utility::Log_Level::Error, Utility::Log_Sender::All, fmt::format("cufftDestroy failed with error: {}", res));
80+
}
6081
}
6182
}
6283
}

0 commit comments

Comments
 (0)