Skip to content

Commit 20c7fc5

Browse files
authored
Test: add unit test for class LCAO_Orbitals (#2037)
* LCAO_Orbitals unit test first commit * LCAO_Orbitals unit test finished
1 parent 47afeb1 commit 20c7fc5

File tree

2 files changed

+361
-16
lines changed

2 files changed

+361
-16
lines changed

source/module_orbital/test/CMakeLists.txt

Lines changed: 10 additions & 16 deletions
Original file line numberDiff line numberDiff line change
@@ -57,24 +57,18 @@ AddTest(
5757
AddTest(
5858
TARGET ORB_atomic_lm_test
5959
SOURCES ORB_atomic_lm_test.cpp
60-
${depend_files}
61-
../../module_base/parallel_reduce.cpp
62-
../../module_base/parallel_global.cpp
63-
../../module_base/parallel_common.cpp
64-
../../module_cell/parallel_kpoints.cpp
65-
LIBS ${math_libs} device
60+
../ORB_atomic_lm.cpp
61+
LIBS ${math_libs} device base
6662
)
6763

68-
#AddTest(
69-
# TARGET ORB_read_test
70-
# SOURCES ORB_read_test.cpp
71-
# ${depend_files}
72-
# ../../module_base/parallel_reduce.cpp
73-
# ../../module_base/parallel_global.cpp
74-
# ../../module_base/parallel_common.cpp
75-
# ../../module_cell/parallel_kpoints.cpp
76-
# LIBS ${math_libs} device
77-
#)
64+
AddTest(
65+
TARGET ORB_read_test
66+
SOURCES ORB_read_test.cpp
67+
../ORB_read.cpp
68+
../ORB_atomic.cpp
69+
../ORB_atomic_lm.cpp
70+
LIBS ${math_libs} device base
71+
)
7872

7973
install(DIRECTORY lcao_H2O DESTINATION ${CMAKE_CURRENT_BINARY_DIR})
8074
install(DIRECTORY lcao_H2O DESTINATION ${CMAKE_CURRENT_BINARY_DIR}/../../../tests)
Lines changed: 351 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,351 @@
1+
#include <fstream>
2+
#include "gtest/gtest.h"
3+
#include "module_base/global_variable.h"
4+
#include "module_orbital/ORB_atomic.h"
5+
#include "module_orbital/ORB_atomic_lm.h"
6+
7+
#define private public
8+
#include "module_orbital/ORB_read.h"
9+
#undef private
10+
11+
#ifdef __MPI
12+
#include <mpi.h>
13+
#endif
14+
15+
/***********************************************************
16+
* unit test of class "LCAO_Orbitals"
17+
***********************************************************/
18+
19+
/**
20+
* Tested functions:
21+
*
22+
* - Read_Orbitals
23+
* read orbital & (optionally) descriptor files and pour the data
24+
* into Numerical_Orbital (and its member Numerical_Orbital_Lm) objects
25+
*
26+
* - bcast_files
27+
* broadcast the orbital_file from rank-0 to all other processes
28+
*
29+
* - all "getters"
30+
* get access to class members
31+
*
32+
***********************************************************/
33+
34+
class LcaoOrbitalsTest : public ::testing::Test
35+
{
36+
protected:
37+
38+
void SetUp();
39+
void TearDown();
40+
41+
// object under unit test
42+
LCAO_Orbitals lcao_;
43+
44+
// initialize lcao_ with parameters below & call Read_Orbitals
45+
void lcao_read();
46+
47+
// parameters to pass to lcao_
48+
std::ofstream ofs_log_;
49+
int ntype_;
50+
int lmax_;
51+
int out_mat_r_; // unused variable
52+
bool force_flag_;
53+
int my_rank_;
54+
bool deepks_setorb_;
55+
bool read_in_flag_;
56+
std::string descriptor_file_;
57+
std::vector<std::string> orbital_file_;
58+
double ecutwfc_;
59+
double dk_;
60+
double dR_;
61+
double Rmax_;
62+
};
63+
64+
65+
void LcaoOrbitalsTest::SetUp() {
66+
ofs_log_.open("ORB_read_test.log");
67+
ntype_ = 2;
68+
lmax_ = 2;
69+
out_mat_r_ = 0; // unused variable
70+
force_flag_ = true;
71+
my_rank_ = 0;
72+
deepks_setorb_ = true;
73+
74+
read_in_flag_ = true;
75+
descriptor_file_ = "./lcao_H2O/jle.orb";
76+
orbital_file_.push_back("./lcao_H2O/H_gga_8au_60Ry_2s1p.orb");
77+
orbital_file_.push_back("./lcao_H2O/O_gga_7au_60Ry_2s2p1d.orb");
78+
79+
//if (GlobalV::MY_RANK != 0) {
80+
// std::cout << "rank=" << GlobalV::MY_RANK
81+
// << " orbital_file_0 = " << orbital_file_[0] << std::endl;
82+
//}
83+
84+
// below we mimic ORB_control::read_orb_first
85+
ecutwfc_ = 123.0;
86+
dk_ = 0.01;
87+
dR_ = 0.01;
88+
Rmax_ = 20;
89+
}
90+
91+
92+
void LcaoOrbitalsTest::lcao_read() {
93+
94+
// see UnitCell::read_atom_species in module_cell/read_atoms.cpp
95+
lcao_.read_in_flag = read_in_flag_;
96+
lcao_.descriptor_file = descriptor_file_;
97+
lcao_.orbital_file = orbital_file_;
98+
#ifdef __MPI
99+
lcao_.bcast_files(ntype_, GlobalV::MY_RANK);
100+
#endif
101+
102+
// see ORB_control::read_orb_first
103+
lcao_.ecutwfc = ecutwfc_;
104+
lcao_.dk = dk_;
105+
lcao_.dR = dR_;
106+
lcao_.Rmax = Rmax_;
107+
108+
lcao_.Read_Orbitals(ofs_log_, ntype_, lmax_, deepks_setorb_, out_mat_r_,
109+
force_flag_, my_rank_);
110+
}
111+
112+
113+
void LcaoOrbitalsTest::TearDown() {
114+
}
115+
116+
117+
TEST_F(LcaoOrbitalsTest, ReadInFlag) {
118+
read_in_flag_ = false;
119+
EXPECT_EXIT(this->lcao_read(), testing::ExitedWithCode(0), "");
120+
}
121+
122+
123+
TEST_F(LcaoOrbitalsTest, WrongOrbFile) {
124+
orbital_file_[0] = "./lcao_H2O/H_gga_8au_60Ry_2s1.orb";
125+
EXPECT_EXIT(this->lcao_read(), testing::ExitedWithCode(0), "");
126+
}
127+
128+
129+
TEST_F(LcaoOrbitalsTest, WrongDescFile) {
130+
descriptor_file_ = "./lcao_H2O/jl.orb";
131+
EXPECT_EXIT(this->lcao_read(), testing::ExitedWithCode(0), "");
132+
}
133+
134+
135+
TEST_F(LcaoOrbitalsTest, BcastFiles) {
136+
#ifdef __MPI
137+
if ( GlobalV::MY_RANK == 0 ) {
138+
lcao_.orbital_file = orbital_file_;
139+
}
140+
141+
if ( GlobalV::MY_RANK != 0) {
142+
EXPECT_EQ(lcao_.orbital_file, std::vector<std::string>{});
143+
}
144+
145+
lcao_.read_in_flag = read_in_flag_;
146+
lcao_.bcast_files(2, GlobalV::MY_RANK);
147+
148+
EXPECT_EQ(lcao_.orbital_file, orbital_file_);
149+
#endif
150+
}
151+
152+
153+
TEST_F(LcaoOrbitalsTest, ReadOrbitals) {
154+
155+
this->lcao_read();
156+
157+
// This test checks whether Read_Orbitals behaves as expected.
158+
159+
// alias
160+
Numerical_Orbital& ao0 = lcao_.Phi[0];
161+
Numerical_Orbital& ao1 = lcao_.Phi[1];
162+
Numerical_Orbital& aod = lcao_.Alpha[0];
163+
164+
// we now check whether the contents of the above
165+
// Numerical_Orbital objects are correct
166+
167+
// maximum tolerance for double-type comparison
168+
double max_tol = 1e-12;
169+
170+
// H
171+
EXPECT_EQ(ao0.getType(), 0);
172+
EXPECT_EQ(ao0.getLabel(), "H");
173+
EXPECT_EQ(ao0.getLmax(), 1);
174+
EXPECT_EQ(ao0.getNchi(0), 2);
175+
EXPECT_EQ(ao0.getNchi(1), 1);
176+
ASSERT_EQ(ao0.getTotal_nchi(), 3);
177+
178+
std::vector<int> L0_list{0,0,1};
179+
std::vector<int> N0_list{0,1,0};
180+
181+
for (size_t i = 0; i != 3; ++i) {
182+
int L = L0_list[i], N = N0_list[i];
183+
EXPECT_EQ(ao0.PhiLN(L,N).getLabel(), "H");
184+
EXPECT_EQ(ao0.PhiLN(L,N).getType(), 0);
185+
EXPECT_EQ(ao0.PhiLN(L,N).getL(), L);
186+
EXPECT_EQ(ao0.PhiLN(L,N).getChi(), N);
187+
EXPECT_EQ(ao0.PhiLN(L,N).getNr(), 801);
188+
EXPECT_EQ(ao0.PhiLN(L,N).getNk(), lcao_.kmesh);
189+
EXPECT_EQ(ao0.PhiLN(L,N).getDk(), lcao_.dk);
190+
EXPECT_EQ(ao0.PhiLN(L,N).getDruniform(), lcao_.dr_uniform);
191+
192+
for (int ir = 0; ir != 801; ++ir) {
193+
EXPECT_DOUBLE_EQ(ao0.PhiLN(L,N).getRab(ir), 0.01);
194+
EXPECT_DOUBLE_EQ(ao0.PhiLN(L,N).getRadial(ir), 0.01*ir);
195+
}
196+
}
197+
198+
EXPECT_NEAR(ao0.PhiLN(0,0).getPsi(0 ), 1.837183001954e+00, max_tol);
199+
EXPECT_NEAR(ao0.PhiLN(0,0).getPsi(1 ), 1.836944589913e+00, max_tol);
200+
EXPECT_NEAR(ao0.PhiLN(0,0).getPsi(4 ), 1.833374417163e+00, max_tol);
201+
EXPECT_NEAR(ao0.PhiLN(0,0).getPsi(799), 3.037233152557e-07, max_tol);
202+
EXPECT_NEAR(ao0.PhiLN(0,0).getPsi(800), 0.000000000000e+00, max_tol);
203+
204+
EXPECT_NEAR(ao0.PhiLN(0,1).getPsi(0 ), -2.482045090982e+00, max_tol);
205+
EXPECT_NEAR(ao0.PhiLN(0,1).getPsi(1 ), -2.481575045574e+00, max_tol);
206+
EXPECT_NEAR(ao0.PhiLN(0,1).getPsi(4 ), -2.474535579529e+00, max_tol);
207+
EXPECT_NEAR(ao0.PhiLN(0,1).getPsi(799), 1.115867959482e-06, max_tol);
208+
EXPECT_NEAR(ao0.PhiLN(0,1).getPsi(800), 0.000000000000e+00, max_tol);
209+
210+
EXPECT_NEAR(ao0.PhiLN(1,0).getPsi(0 ), 0.000000000000e+00, max_tol);
211+
EXPECT_NEAR(ao0.PhiLN(1,0).getPsi(1 ), -2.619148756396e-02, max_tol);
212+
EXPECT_NEAR(ao0.PhiLN(1,0).getPsi(4 ), -1.045849793771e-01, max_tol);
213+
EXPECT_NEAR(ao0.PhiLN(1,0).getPsi(799), 3.217573100688e-06, max_tol);
214+
EXPECT_NEAR(ao0.PhiLN(1,0).getPsi(800), 0.000000000000e+00, max_tol);
215+
216+
217+
// O
218+
EXPECT_EQ(ao1.getType(), 1);
219+
EXPECT_EQ(ao1.getLabel(), "O");
220+
EXPECT_EQ(ao1.getLmax(), 2);
221+
EXPECT_EQ(ao1.getNchi(0), 2);
222+
EXPECT_EQ(ao1.getNchi(1), 2);
223+
EXPECT_EQ(ao1.getNchi(2), 1);
224+
ASSERT_EQ(ao1.getTotal_nchi(), 5);
225+
226+
std::vector<int> L1_list{0,0,1,1,2};
227+
std::vector<int> N1_list{0,1,0,1,0};
228+
229+
for (size_t i = 0; i != 5; ++i) {
230+
int L = L1_list[i], N = N1_list[i];
231+
EXPECT_EQ(ao1.PhiLN(L,N).getLabel(), "O");
232+
EXPECT_EQ(ao1.PhiLN(L,N).getType(), 1);
233+
EXPECT_EQ(ao1.PhiLN(L,N).getL(), L);
234+
EXPECT_EQ(ao1.PhiLN(L,N).getChi(), N);
235+
EXPECT_EQ(ao1.PhiLN(L,N).getNr(), 701);
236+
EXPECT_EQ(ao1.PhiLN(L,N).getNk(), lcao_.kmesh);
237+
EXPECT_EQ(ao1.PhiLN(L,N).getDk(), lcao_.dk);
238+
EXPECT_EQ(ao1.PhiLN(L,N).getDruniform(), lcao_.dr_uniform);
239+
240+
for (int ir = 0; ir != 701; ++ir) {
241+
EXPECT_DOUBLE_EQ(ao1.PhiLN(L,N).getRab(ir), 0.01);
242+
EXPECT_DOUBLE_EQ(ao1.PhiLN(L,N).getRadial(ir), 0.01*ir);
243+
}
244+
}
245+
246+
EXPECT_NEAR(ao1.PhiLN(0,0).getPsi(0), 1.208504975904e+00, max_tol);
247+
EXPECT_NEAR(ao1.PhiLN(0,0).getPsi(1), 1.208605373194e+00, max_tol);
248+
EXPECT_NEAR(ao1.PhiLN(0,0).getPsi(4), 1.210103935461e+00, max_tol);
249+
EXPECT_NEAR(ao1.PhiLN(0,0).getPsi(699), 4.465396560257e-08, max_tol);
250+
EXPECT_NEAR(ao1.PhiLN(0,0).getPsi(700), 0.0, max_tol);
251+
252+
EXPECT_NEAR(ao1.PhiLN(0,1).getPsi(0), 7.254873428942e-01, max_tol);
253+
EXPECT_NEAR(ao1.PhiLN(0,1).getPsi(1), 7.256666701836e-01, max_tol);
254+
EXPECT_NEAR(ao1.PhiLN(0,1).getPsi(4), 7.283448557011e-01, max_tol);
255+
EXPECT_NEAR(ao1.PhiLN(0,1).getPsi(699), -1.916246212603e-06, max_tol);
256+
EXPECT_NEAR(ao1.PhiLN(0,1).getPsi(700), 0.0, max_tol);
257+
258+
EXPECT_NEAR(ao1.PhiLN(1,0).getPsi(0), 0.0, max_tol);
259+
EXPECT_NEAR(ao1.PhiLN(1,0).getPsi(1), 4.626669306440e-02, max_tol);
260+
EXPECT_NEAR(ao1.PhiLN(1,0).getPsi(4), 1.845014292772e-01, max_tol);
261+
EXPECT_NEAR(ao1.PhiLN(1,0).getPsi(699), 2.870401658966e-07, max_tol);
262+
EXPECT_NEAR(ao1.PhiLN(1,0).getPsi(700), 0.0, max_tol);
263+
264+
EXPECT_NEAR(ao1.PhiLN(1,1).getPsi(0), 0.0, max_tol);
265+
EXPECT_NEAR(ao1.PhiLN(1,1).getPsi(1), 3.375340101333e-02, max_tol);
266+
EXPECT_NEAR(ao1.PhiLN(1,1).getPsi(4), 1.346256082234e-01, max_tol);
267+
EXPECT_NEAR(ao1.PhiLN(1,1).getPsi(699), -2.771091616120e-06, max_tol);
268+
EXPECT_NEAR(ao1.PhiLN(1,1).getPsi(700), 0.0, max_tol);
269+
270+
EXPECT_NEAR(ao1.PhiLN(2,0).getPsi(0), 0.0, max_tol);
271+
EXPECT_NEAR(ao1.PhiLN(2,0).getPsi(1), -3.343626342662e-04, max_tol);
272+
EXPECT_NEAR(ao1.PhiLN(2,0).getPsi(4), -5.337546547975e-03, max_tol);
273+
EXPECT_NEAR(ao1.PhiLN(2,0).getPsi(699), 1.396308876444e-06, max_tol);
274+
EXPECT_NEAR(ao1.PhiLN(2,0).getPsi(700), 0.0, max_tol);
275+
276+
277+
// Descriptor
278+
279+
EXPECT_EQ(aod.getType(), 0);
280+
EXPECT_EQ(aod.getLabel(), "");
281+
EXPECT_EQ(aod.getLmax(), 2);
282+
EXPECT_EQ(aod.getNchi(0), 2);
283+
EXPECT_EQ(aod.getNchi(1), 2);
284+
EXPECT_EQ(aod.getNchi(2), 2);
285+
ASSERT_EQ(aod.getTotal_nchi(), 6);
286+
287+
std::vector<int> Ld_list{0,0,1,1,2,2};
288+
std::vector<int> Nd_list{0,1,0,1,0,1};
289+
290+
for (size_t i = 0; i != 6; ++i) {
291+
int L = Ld_list[i], N = Nd_list[i];
292+
EXPECT_EQ(aod.PhiLN(L,N).getLabel(), "");
293+
EXPECT_EQ(aod.PhiLN(L,N).getType(), 0);
294+
EXPECT_EQ(aod.PhiLN(L,N).getL(), L);
295+
EXPECT_EQ(aod.PhiLN(L,N).getChi(), N);
296+
EXPECT_EQ(aod.PhiLN(L,N).getNr(), 205);
297+
EXPECT_EQ(aod.PhiLN(L,N).getNk(), lcao_.kmesh);
298+
EXPECT_EQ(aod.PhiLN(L,N).getDk(), lcao_.dk);
299+
EXPECT_EQ(aod.PhiLN(L,N).getDruniform(), lcao_.dr_uniform);
300+
301+
for (int ir = 0; ir != 205; ++ir) {
302+
EXPECT_DOUBLE_EQ(aod.PhiLN(L,N).getRab(ir), 0.01);
303+
EXPECT_DOUBLE_EQ(aod.PhiLN(L,N).getRadial(ir), 0.01*ir);
304+
}
305+
}
306+
307+
// TODO chi value check is skipped for now
308+
// orbitals in jle.orb are not normalized
309+
// getPsi() does not gives the numbers in jle.orb
310+
}
311+
312+
313+
TEST_F(LcaoOrbitalsTest, Getters) {
314+
315+
this->lcao_read();
316+
317+
EXPECT_EQ(lcao_.get_ecutwfc(), lcao_.ecutwfc);
318+
EXPECT_EQ(lcao_.get_kmesh(), lcao_.kmesh);
319+
EXPECT_EQ(lcao_.get_dk(), lcao_.dk);
320+
EXPECT_EQ(lcao_.get_dR(), lcao_.dR);
321+
EXPECT_EQ(lcao_.get_Rmax(), lcao_.Rmax);
322+
EXPECT_EQ(lcao_.get_lmax(), lcao_.lmax);
323+
EXPECT_EQ(lcao_.get_lmax_d(), lcao_.lmax_d);
324+
EXPECT_EQ(lcao_.get_nchimax(), lcao_.nchimax);
325+
EXPECT_EQ(lcao_.get_nchimax_d(), lcao_.nchimax_d);
326+
EXPECT_EQ(lcao_.get_ntype(), lcao_.ntype);
327+
EXPECT_EQ(lcao_.get_dr_uniform(), lcao_.dr_uniform);
328+
EXPECT_EQ(lcao_.get_rcutmax_Phi(), lcao_.rcutmax_Phi);
329+
}
330+
331+
332+
int main(int argc, char **argv)
333+
{
334+
335+
#ifdef __MPI
336+
MPI_Init(&argc, &argv);
337+
MPI_Comm_size(MPI_COMM_WORLD,&GlobalV::NPROC);
338+
MPI_Comm_rank(MPI_COMM_WORLD,&GlobalV::MY_RANK);
339+
#endif
340+
341+
testing::InitGoogleTest(&argc, argv);
342+
int result = RUN_ALL_TESTS();
343+
344+
#ifdef __MPI
345+
MPI_Finalize();
346+
#endif
347+
348+
return result;
349+
}
350+
351+

0 commit comments

Comments
 (0)