COOLFluiD  Release kernel
COOLFluiD is a Collaborative Simulation Environment (CSE) focused on complex MultiPhysics simulations.
utest-lss-distributed-matrix.cpp
Go to the documentation of this file.
1 // Copyright (C) 2010-2013 von Karman Institute for Fluid Dynamics, Belgium
2 //
3 // This software is distributed under the terms of the
4 // GNU Lesser General Public License version 3 (LGPLv3).
5 // See doc/lgpl.txt and doc/gpl.txt for the license text.
6 //
7 
8 #define BOOST_TEST_DYN_LINK
9 #define BOOST_TEST_MODULE "Test module for cf3::math::LSS where testing on a distributed matrix."
10 
12 
13 #include <fstream>
14 
15 #include <boost/test/unit_test.hpp>
16 #include <boost/lexical_cast.hpp>
17 
18 #include "common/Log.hpp"
19 #include "common/PE/Comm.hpp"
21 #include "math/LSS/System.hpp"
23 
25 
26 using namespace cf3;
27 using namespace cf3::math;
28 using namespace cf3::math::LSS;
29 
31 
33 {
36  {
37  m_argc = boost::unit_test::framework::master_test_suite().argc;
38  m_argv = boost::unit_test::framework::master_test_suite().argv;
39  }
40 
43  {
44  }
45 
47  int m_argc;
48  char** m_argv;
49 };
50 
52 
53 BOOST_FIXTURE_TEST_SUITE( LSSDistributedMatrixSuite, LSSDistributedMatrixFixture )
54 
55 
57 BOOST_AUTO_TEST_CASE( init_parallel_environment )
58 {
59  common::PE::Comm::instance().init(m_argc,m_argv);
60  BOOST_CHECK_EQUAL(common::PE::Comm::instance().is_active(),true);
61  CFinfo.setFilterRankZero(false);
62 }
63 
65 
66 BOOST_AUTO_TEST_CASE( system_solve )
67 {
68  // test matrix
69  test_matrix m;
70 
71  // commpattern
72  boost::shared_ptr<common::PE::CommPattern> cp_ptr = common::allocate_component<common::PE::CommPattern>("commpattern");
73  common::PE::CommPattern& cp = *cp_ptr;
74  cp.insert("gid",m.global_numbering,1,false);
75  cp.setup(Handle<common::PE::CommWrapper>(cp.get_child("gid")),m.irank_updatable);
76 
77  // system
78  boost::shared_ptr<LSS::System> sys_ptr = common::allocate_component<LSS::System>("system");
79  LSS::System& sys = *sys_ptr;
80  sys.options().option("matrix_builder").change_value(boost::lexical_cast<std::string>(m_argv[1]));
81  sys.create(cp,m.nbeqs,m.column_indices,m.rowstart_positions);
82  sys.reset();
83 
84  // mimicing the assembly for performance measuring
85  /*
86  with az_msrmatrix and native written fill, the suminto part takes the following number of ticks in morpheus (for one step):
87  process | assembly_ticks | matrixfill_ticks | setbc_ticks | solve_ticks
88  -----------------------------------------------------------------------
89  0 | 330499 | 239996 | 24888 | 295230350
90  1 | 337050 | 248826 | 49497 | 306481891
91  2 | 329032 | 242850 | 37709 | 298537437
92  3 | 388584 | 280844 | 13722 | 299270122
93  */
95  ba.resize(3,m.nbeqs);
96  ba.reset(1.);
97  for (int i=0; i<(const int)(m.elem_nodes.size()/3); i++)
98  {
99  ba.indices[0]=m.elem_nodes[3*i+0];
100  ba.indices[1]=m.elem_nodes[3*i+1];
101  ba.indices[2]=m.elem_nodes[3*i+2];
102  sys.add_values(ba);
103  }
104 
105 sys.print("sys_test_assembly_bc_" + boost::lexical_cast<std::string>(m.irank) + ".plt");
106 
107 
108  // fastly filling with pre-boundary condition values, THIS IS NOT THE WAY YOU NORMALLY ASSEMBLE
109  Real* vals=&m.mat_prebc[0];
110  sys.reset();
111  for (int i=0; i<(const int)m.global_numbering.size(); i++)
112  if (cp.isUpdatable()[i])
113  for (int j=0; j<(const int)m.nbeqs; j++)
114  for (int k=m.rowstart_positions[i]; k<(const int)(m.rowstart_positions[i+1]); k++)
115  for (int l=0; l<(const int)m.nbeqs; l++)
116  sys.matrix()->set_value(m.column_indices[k]*m.nbeqs+l,i*m.nbeqs+j,*vals++);
117  vals=&m.rhs_prebc[0];
118  for (int i=0; i<(const int)m.global_numbering.size(); i++)
119  for (int j=0; j<(const int)m.nbeqs; j++)
120  sys.rhs()->set_value(i,j,*vals++);
121  vals=&m.sol_prebc[0];
122  for (int i=0; i<(const int)m.global_numbering.size(); i++)
123  for (int j=0; j<(const int)m.nbeqs; j++)
124  sys.solution()->set_value(i,j,*vals++);
125 
126  // applying boundary conditions
127  for (int i=0; i<(const int)m.bc_node.size(); i++)
128  if (cp.isUpdatable()[m.bc_node[i]])
129  sys.dirichlet(m.bc_node[i],m.bc_eqn[i],m.bc_value[i]);
130  for (int i=0; i<(const int)m.periodic_pairs.size(); i+=2)
131 // if (cp.isUpdatable()[m.periodic_pairs[i]])
132  sys.periodicity(m.periodic_pairs[i+0],m.periodic_pairs[i+1]);
133 
134 
135 //sys.print("sys_prebc_bc_" + boost::lexical_cast<std::string>(m.irank) + ".plt");
136 
137 // // filling the system with prescribed values
138 // sys.reset();
139 // vals=&m.mat_presolve[0];
140 // for (int i=0; i<(const int)m.global_numbering.size(); i++)
141 // if (cp.isUpdatable()[i])
142 // for (int j=0; j<(const int)m.nbeqs; j++)
143 // for (int k=m.rowstart_positions[i]; k<(const int)(m.rowstart_positions[i+1]); k++)
144 // for (int l=0; l<(const int)m.nbeqs; l++)
145 // sys.matrix()->set_value(m.column_indices[k]*m.nbeqs+l,i*m.nbeqs+j,*vals++);
146 // vals=&m.rhs_presolve[0];
147 // for (int i=0; i<(const int)m.global_numbering.size(); i++)
148 // for (int j=0; j<(const int)m.nbeqs; j++)
149 // sys.rhs()->set_value(i,j,*vals++);
150 // vals=&m.sol_presolve[0];
151 // for (int i=0; i<(const int)m.global_numbering.size(); i++)
152 // for (int j=0; j<(const int)m.nbeqs; j++)
153 // sys.solution()->set_value(i,j,*vals++);
154 
155 //sys.print("sys_orig_" + boost::lexical_cast<std::string>(m.irank) + ".plt");
156 
157  // and solve the system
158  sys.solve();
159 
160  // check results
161  std::vector<Real> v;
162  sys.solution()->debug_data(v);
163  for (int i=0; i<(const int)m.global_numbering.size(); i++)
164  if (cp.isUpdatable()[i])
165  for (int j=0; j<(const int)m.nbeqs; j++)
166  BOOST_CHECK_CLOSE(v[i*m.nbeqs+j],m.result[i*m.nbeqs+j],1e-3);
167 
168 // sys.print("distributed_system_" + boost::lexical_cast<std::string>(m.irank) + ".plt");
169 
170  // print a test file
171  std::ofstream fres(std::string("coords_with_sol_" + boost::lexical_cast<std::string>(m.irank) + ".plt").c_str());
172  fres << "VARIABLES=\"X\",\"Y\",\"V0\",\"V1\",\"V2\"\nZONE T=\"Solve results of utest-lss-distributed-matrix.\"" << std::flush;
173  fres.precision(15);
174  sys.solution()->debug_data(v);
175  for (int i=0; i<(const int)m.global_numbering.size(); i++)
176  if (cp.isUpdatable()[i])
177  fres << m.nodal_coordinates[2*i+0] << " " << m.nodal_coordinates[2*i+1] << " " << v[m.nbeqs*i+0] << " " << v[m.nbeqs*i+1] << " " << v[m.nbeqs*i+2] << "\n" << std::flush;
178  fres.close();
179 
180 }
181 
183 
184 BOOST_AUTO_TEST_CASE( finalize_parallel_environment )
185 {
186  CFinfo.setFilterRankZero(true);
188  BOOST_CHECK_EQUAL(common::PE::Comm::instance().is_active(),false);
189 }
190 
192 
193 BOOST_AUTO_TEST_SUITE_END()
194 
195 
std::vector< cf3::Uint > bc_node
local node index of boundaries (without extension by nbeqs sub-matrix)
#define CFinfo
these are always defined
Definition: Log.hpp:104
void resize(Uint numnodes, Uint numeqs)
setting up sizes
Safe pointer to an object. This is the supported method for referring to components.
Definition: Handle.hpp:39
This header collects all the headers needed for the linear system solver, also including configure-ti...
std::vector< cf3::Real > bc_value
what value to fix
Holding all the data and the constructor sets everything up, everything is public.
Basic Classes for Mathematical applications used by COOLFluiD.
Parallel Communication Pattern. This class provides functionality to collect communication. For efficiency it works such a way that you submit your request via the constructor or the add/remove/move magic triangle and then call setup to modify the commpattern. The data needed to be kept synchronous can be registered via the insert function. The word node here means any kind of "point of storage", in this context it is not directly related with the computational mesh.
std::vector< cf3::Uint > rowstart_positions
connectivity structure, CSR-style (without extension by nbeqs sub-matrix)
std::vector< cf3::Uint > global_numbering
global numbering of the nodes (without extension by nbeqs sub-matrix)
void insert(const std::string &name, T *&data, const int size, const unsigned int stride=1, const bool needs_update=true)
Real e()
Definition of the Unit charge [C].
Definition: Consts.hpp:30
LSSDistributedMatrixFixture()
common setup for each test case
std::vector< Uint > indices
local numbering of the unknowns
virtual void change_value(const boost::any &value)
change the value of this option
Definition: Option.cpp:89
std::vector< cf3::Uint > column_indices
connectivity structure, CSR-style (without extension by nbeqs sub-matrix)
std::vector< cf3::Real > mat_prebc
matrix before applying boundary conditions
BOOST_AUTO_TEST_CASE(init_parallel_environment)
void init(int argc=0, char **args=0)
Definition: Comm.cpp:80
Top-level namespace for coolfluid.
Definition: Action.cpp:18
std::vector< cf3::Uint > irank_updatable
rank where the node is updatable (without extension by nbeqs sub-matrix)
std::vector< cf3::Real > sol_prebc
solution before applying boundary conditions
std::vector< cf3::Real > result
reference result after solving the system
void reset(Real reset_to=0.)
reset the values to the value of reset_to
std::vector< cf3::Real > rhs_prebc
right hand side before applying boundary conditions
cf3::Uint nbeqs
number of equations
std::vector< cf3::Uint > periodic_pairs
pairs of periodic local node indices (obviously, size of vector is twice of number of periodic pairs)...
~LSSDistributedMatrixFixture()
common tear-down for each test case
std::vector< cf3::Real > nodal_coordinates
node coordinates for visualization
OptionList & options()
Definition: Component.cpp:856
std::vector< cf3::Uint > bc_eqn
which equation to apply the dirichlet bc to (within nbeqs)
cf3::Uint irank
rank of current process
std::vector< cf3::Uint > elem_nodes
node ids of the elements, for mimicing the assembly
static Comm & instance()
Return a reference to the current PE.
Definition: Comm.cpp:44
const Option & option(const std::string &pname) const
get a constant option from the list
Definition: OptionList.cpp:52
Send comments to:
COOLFluiD Web Admin