COOLFluiD  Release kernel
COOLFluiD is a Collaborative Simulation Environment (CSE) focused on complex MultiPhysics simulations.
utest-proto-heat-parallel.cpp
Go to the documentation of this file.
1 // Copyright (C) 2010-2011 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 heat-conduction related proto operations"
10 
11 #include <boost/assign.hpp>
12 #include <boost/lexical_cast.hpp>
13 #include <boost/test/unit_test.hpp>
14 
15 #define BOOST_PROTO_MAX_ARITY 10
16 #ifdef BOOST_MPL_LIMIT_METAFUNCTION_ARITY
17  #undef BOOST_MPL_LIMIT_METAFUNCTION_ARITY
18  #define BOOST_MPL_LIMIT_METAFUNCTION_ARITY 10
19 #endif
20 
21 #include "common/Core.hpp"
22 #include "common/Environment.hpp"
23 
24 #include "common/PE/debug.hpp"
25 
26 #include "math/LSS/System.hpp"
27 
28 #include "mesh/Domain.hpp"
29 
31 #include "solver/Model.hpp"
32 
33 #include "math/LSS/SolveLSS.hpp"
36 
38 
39 #include "UFEM/LSSAction.hpp"
40 #include "UFEM/Solver.hpp"
41 #include "UFEM/Tags.hpp"
42 
43 using namespace cf3;
44 using namespace cf3::solver;
45 using namespace cf3::solver::actions;
46 using namespace cf3::solver::actions::Proto;
47 using namespace cf3::common;
48 using namespace cf3::math::Consts;
49 using namespace cf3::mesh;
50 
51 using namespace boost::assign;
52 
53 
54 typedef std::vector<std::string> StringsT;
55 typedef std::vector<Uint> SizesT;
56 
58 inline void
59 check_close(const Real a, const Real b, const Real threshold)
60 {
61  BOOST_CHECK_CLOSE(a, b, threshold);
62 }
63 
64 static boost::proto::terminal< void(*)(Real, Real, Real) >::type const _check_close = {&check_close};
65 
67 {
69  root( Core::instance().root() )
70  {
71  }
72 
74 
75 };
76 
77 BOOST_FIXTURE_TEST_SUITE( ProtoHeatSuite, ProtoHeatFixture )
78 
80 {
81  common::PE::Comm::instance().init(boost::unit_test::framework::master_test_suite().argc, boost::unit_test::framework::master_test_suite().argv);
82  //PE::wait_for_debugger(0);
83 }
84 
85 BOOST_AUTO_TEST_CASE( Heat2DParallel)
86 {
87  Core::instance().environment().options().set("log_level", 4u);
88 
89  // Parameters
90  Real length = 5.;
91  const Uint nb_segments = 16 ;
92 
93  // Setup a model
94  Model& model = *root.create_component<Model>("Model");
95  Domain& domain = model.create_domain("Domain");
96  UFEM::Solver& solver = *model.create_component<UFEM::Solver>("Solver");
97 
98  Handle<UFEM::LSSAction> lss_action(solver.add_direct_solver("cf3.UFEM.LSSAction"));
99 
100  // Proto placeholders
101  FieldVariable<0, ScalarField> temperature("Temperature", UFEM::Tags::solution());
102 
103  // Allowed elements (reducing this list improves compile times)
104  boost::mpl::vector1<mesh::LagrangeP1::Quad2D> allowed_elements;
105 
106  // BCs
107  boost::shared_ptr<UFEM::BoundaryConditions> bc = allocate_component<UFEM::BoundaryConditions>("BoundaryConditions");
108 
109  // add the top-level actions (assembly, BC and solve)
110  *lss_action
112  (
113  "Assembly",
115  (
116  allowed_elements,
117  group
118  (
119  _A = _0,
120  element_quadrature( _A(temperature) += transpose(nabla(temperature)) * nabla(temperature) ),
121  lss_action->system_matrix += _A
122  )
123  )
124  )
125  << bc
126  << allocate_component<math::LSS::SolveLSS>("SolveLSS")
127  << create_proto_action("Increment", nodes_expression(temperature += lss_action->solution(temperature)))
128  << create_proto_action("CheckResult", nodes_expression(_check_close(temperature, 10. + 25.*(coordinates(0,0) / length), 1e-6)));
129 
130  // Setup physics
131  model.create_physics("cf3.physics.DynamicModel");
132 
133  // Setup mesh
134  Mesh& mesh = *domain.create_component<Mesh>("Mesh");
135  BlockMesh::BlockArrays& blocks = *domain.create_component<BlockMesh::BlockArrays>("blocks");
136 
137  *blocks.create_points(2, 4) << 0. << 0. << length << 0. << length << length << 0. << length;
138  *blocks.create_blocks(1) << 0 << 1 << 2 << 3;
139  *blocks.create_block_subdivisions() << nb_segments << nb_segments;
140  *blocks.create_block_gradings() << 1. << 1. << 1. << 1.;
141 
142  *blocks.create_patch("bottom", 1) << 0 << 1;
143  *blocks.create_patch("right", 1) << 1 << 2;
144  *blocks.create_patch("top", 1) << 2 << 3;
145  *blocks.create_patch("left", 1) << 3 << 0;
146 
147  blocks.partition_blocks(PE::Comm::instance().size(), YY);
148  blocks.create_mesh(mesh);
149 
150  lss_action->options().set("regions", std::vector<URI>(1, mesh.topology().uri()));
151 
152  // Set boundary conditions
153  bc->add_constant_bc("left", "Temperature", 10.);
154  bc->add_constant_bc("right", "Temperature", 35.);
155 
156  // Run the solver
157  model.simulate();
158 
159  // Write data pertaining to communication
160  Field& comm_field = mesh.geometry_fields().create_field("comm", "ghosts, updatable, rank");
161  const Uint nb_nodes = mesh.geometry_fields().size();
162  for(Uint node = 0; node != nb_nodes; ++node)
163  {
164  comm_field[node][0] = mesh.geometry_fields().is_ghost(node) ? 1. : 0.;
165  comm_field[node][1] = mesh.geometry_fields().comm_pattern().isUpdatable()[node];
166  comm_field[node][2] = mesh.geometry_fields().rank()[node];
167  }
168 
169  // Save
170  model.domain().create_component("writer", "cf3.mesh.VTKXML.Writer");
171  model.domain().write_mesh(URI("utest-proto-heat-parallel_output.pvtu", cf3::common::URI::Scheme::FILE));
172 // lss.matrix()->print("utest-proto-heat-parallel_matrix-" + boost::lexical_cast<std::string>(common::PE::Comm::instance().rank()) + ".plt");
173 // lss.rhs()->print("utest-proto-heat-parallel_rhs-" + boost::lexical_cast<std::string>(common::PE::Comm::instance().rank()) + ".plt");
174 // lss.solution()->print("utest-proto-heat-parallel_solution-" + boost::lexical_cast<std::string>(common::PE::Comm::instance().rank()) + ".plt");
175 }
176 
178 
179 BOOST_AUTO_TEST_SUITE_END()
180 
181 
Field & create_field(const std::string &name, const Uint cols)
Create a new field in this group.
Definition: Dictionary.cpp:178
static boost::proto::terminal< ZeroTag >::type _0
Placeholder for the zero matrix.
Definition: Terminals.hpp:209
virtual mesh::Domain & create_domain(const std::string &name)
creates a domain in this model
Definition: Model.cpp:201
Safe pointer to an object. This is the supported method for referring to components.
Definition: Handle.hpp:39
virtual physics::PhysModel & create_physics(const std::string &builder)
Definition: Model.cpp:182
static boost::proto::terminal< ElementSystemMatrix< boost::mpl::int_< 0 > > >::type const _A
Some predefined element matrices (more can be user-defined, but you have to change the number in the ...
Handle< common::Action > add_direct_solver(const std::string &builder_name)
Definition: Solver.cpp:128
This header collects all the headers needed for the linear system solver, also including configure-ti...
Basic Classes for Solver applications used by CF.
Definition: Action.cpp:29
URI uri() const
Construct the full path.
Definition: Component.cpp:248
common::List< Uint > & rank()
Return the rank of every field row.
Definition: Dictionary.hpp:90
static boost::proto::terminal< ElementQuadratureTag >::type element_quadrature
Use element_quadrature(expr1, expr2, ..., exprN) to evaluate a group of expressions.
tuple root
Definition: coolfluid.py:24
boost::shared_ptr< ElementsExpression< ExprT, ElementTypes > > elements_expression(ElementTypes, const ExprT &expr)
Definition: Expression.hpp:292
void check_close(const Real a, const Real b, const Real threshold)
Check close, for testing purposes.
Manage a collection of UFEM solvers.
Definition: Solver.hpp:31
std::vector< Uint > SizesT
Real e()
Definition of the Unit charge [C].
Definition: Consts.hpp:30
Static functions for mathematical constants.
Definition: Consts.hpp:25
boost::proto::terminal< TransposeFunction >::type const transpose
boost::proto::terminal< SFOp< NablaOp > >::type const nabla
virtual mesh::Domain & domain()
gets the domain from this model
Definition: Model.cpp:161
static boost::proto::terminal< void(*)(Real, Real, Real) >::type const _check_close
std::vector< std::string > StringsT
void write_mesh(const common::URI &file)
write the active mesh
Definition: Domain.cpp:179
Basic Classes for Mesh applications used by COOLFluiD.
tuple model
Global confifuration.
void init(int argc=0, char **args=0)
Definition: Comm.cpp:80
Top-level namespace for coolfluid.
Definition: Action.cpp:18
boost::proto::terminal< SFOp< CoordinatesOp > >::type const coordinates
BOOST_AUTO_TEST_CASE(InitMPI)
static boost::proto::terminal< ExpressionGroupTag >::type group
Use group(expr1, expr2, ..., exprN) to evaluate a group of expressions.
common::Environment & environment() const
Definition: Core.cpp:168
unsigned int Uint
typedef for unsigned int
Definition: CF.hpp:90
virtual void simulate()
Simulates this model.
Definition: Model.cpp:132
boost::shared_ptr< ProtoAction > create_proto_action(const std::string &name, const boost::shared_ptr< Expression > &expression)
Create a new ProtoAction, immediatly setting the expression.
bool is_ghost(const Uint idx) const
Check if a field row is owned by this rank.
Definition: Dictionary.cpp:151
Region & topology() const
Definition: Mesh.hpp:51
static Core & instance()
Definition: Core.cpp:37
boost::shared_ptr< NodesExpression< ExprT, boost::mpl::range_c< Uint, 1, 4 > > > nodes_expression(const ExprT &expr)
Definition: Expression.hpp:308
Definition: Defs.hpp:17
OptionList & options()
Definition: Component.cpp:856
Dictionary & geometry_fields() const
Definition: Mesh.cpp:339
static Comm & instance()
Return a reference to the current PE.
Definition: Comm.cpp:44
Base class for defining CF components.
Definition: Component.hpp:82
void set(const std::string &pname, const boost::any &val)
Definition: OptionList.cpp:132
Handle< Component > create_component(const std::string &name, const std::string &builder)
Build a (sub)component of this component using the extended type_name of the component.
Definition: Component.cpp:568
Uint size() const
Number of rows of contained fields.
Definition: Dictionary.cpp:99
Most basic kernel library.
Definition: Action.cpp:19
std::vector< bool > & isUpdatable()
common::PE::CommPattern & comm_pattern()
Return the comm pattern valid for this field group. Created based on the glb_idx and rank if it didn'...
Definition: Dictionary.cpp:136
Handle< common::Table< Real > > create_points(const Uint dimensions, const Uint nb_points)
Create the table that holds the points for the blocks.
Definition: BlockData.cpp:1287
Send comments to:
COOLFluiD Web Admin