COOLFluiD  Release kernel
COOLFluiD is a Collaborative Simulation Environment (CSE) focused on complex MultiPhysics simulations.
utest-rdm-merge.cpp
Go to the documentation of this file.
1 
2 // Copyright (C) 2010-2011 von Karman Institute for Fluid Dynamics, Belgium
3 //
4 // This software is distributed under the terms of the
5 // GNU Lesser General Public License version 3 (LGPLv3).
6 // See doc/lgpl.txt and doc/gpl.txt for the license text.
7 
8 #define BOOST_TEST_DYN_LINK
9 #define BOOST_TEST_MODULE "Test module for heat-conduction related proto operations"
10 
11 #include <boost/test/unit_test.hpp>
12 
13 #define BOOST_PROTO_MAX_ARITY 10
14 #ifdef BOOST_MPL_LIMIT_METAFUNCTION_ARITY
15  #undef BOOST_MPL_LIMIT_METAFUNCTION_ARITY
16 #endif
17 #define BOOST_MPL_LIMIT_METAFUNCTION_ARITY 10
18 
19 #include "common/Core.hpp"
20 #include "common/Environment.hpp"
21 
22 #include "math/LSS/System.hpp"
23 
24 #include "mesh/Domain.hpp"
25 
27 #include "solver/ModelUnsteady.hpp"
28 #include "solver/Time.hpp"
29 
30 #include "math/LSS/SolveLSS.hpp"
32 #include "solver/CriterionTime.hpp"
34 
37 
39 
41 #include "UFEM/Solver.hpp"
42 #include "UFEM/Tags.hpp"
43 
44 using namespace cf3;
45 using namespace cf3::solver;
46 using namespace cf3::solver::actions;
47 using namespace cf3::solver::actions::Proto;
48 using namespace cf3::common;
49 using namespace cf3::math::Consts;
50 using namespace cf3::mesh;
51 
52 using namespace boost;
53 
54 
55 typedef std::vector<std::string> StringsT;
56 typedef std::vector<Uint> SizesT;
57 
59 inline void
60 check_close(const Real a, const Real b, const Real threshold)
61 {
62  BOOST_CHECK_CLOSE(a, b, threshold);
63 }
64 
65 static boost::proto::terminal< void(*)(Real, Real, Real) >::type const _check_close = {&check_close};
66 
68 {
70  root( Core::instance().root() )
71  {
72  }
73 
75 
76 };
77 
78 BOOST_FIXTURE_TEST_SUITE( RDMMergeSuite, RDMMergeFixture )
79 
81 {
82  common::PE::Comm::instance().init(boost::unit_test::framework::master_test_suite().argc, boost::unit_test::framework::master_test_suite().argv);
83  BOOST_CHECK_EQUAL(common::PE::Comm::instance().size(), 1);
84 }
85 
86 BOOST_AUTO_TEST_CASE( Heat1DComponent )
87 {
88  Core::instance().environment().options().set("log_level", 4u);
89 
90  // Parameters
91  Real length = 1.;
92  const Uint nb_segments = 25 ;
93 
94  // Setup a model
95  ModelUnsteady& model = *root.create_component<ModelUnsteady>("Model");
96  Domain& domain = model.create_domain("Domain");
97  UFEM::Solver& solver = *model.create_component<UFEM::Solver>("Solver");
98  Handle<UFEM::LSSActionUnsteady> lss_action(solver.add_unsteady_solver("cf3.UFEM.LSSActionUnsteady"));
99  Handle<common::ActionDirector> ic(solver.get_child("InitialConditions"));
100 
101 
102 
103  boost::shared_ptr<solver::actions::Iterate> time_loop = allocate_component<solver::actions::Iterate>("TimeLoop");
104  time_loop->create_component<solver::CriterionTime>("CriterionTime");
105 
106  // Proto placeholders
107  FieldVariable<0, ScalarField> fi("FI", UFEM::Tags::solution());
108 
109  // Allowed elements (reducing this list improves compile times)
110  boost::mpl::vector1<mesh::LagrangeP1::Triag2D> allowed_elements;
111 
112  // BCs
113  boost::shared_ptr<UFEM::BoundaryConditions> bc = allocate_component<UFEM::BoundaryConditions>("BoundaryConditions");
114 
115  FieldVariable<1, VectorField> u_adv("AdvectionSpeed", "linearized_velocity");
116  RealVector u_ref(2); u_ref << 1.,0.;
117 
118  // add the top-level actions (assembly, BC and solve)
119  *ic
121  (
122  "SetInitial",
124  (
125  fi = 5.
126  )
127  )
129  (
130  "InitAdvectionSpeed",
132  (
133  u_adv = coordinates[1] * u_ref
134  )
135  );
136  *lss_action
138  (
139  "Assembly",
141  (
142  allowed_elements,
143  group
144  (
145  _A = _0,
147  (
148  _A(fi,fi) += 0.0025 * transpose(nabla(fi)) * nabla(fi),
149  _A(fi,fi) += transpose(N(fi) /*+ 0.*/ ) * u_adv*nabla(fi)
150  ),
151  lss_action->system_matrix += _A,
152  lss_action->system_rhs += -_A * _x
153  )
154  )
155  )
156  << bc
157  << allocate_component<math::LSS::SolveLSS>("SolveLSS")
158  << create_proto_action("Increment", nodes_expression(fi += lss_action->solution(fi)));
159 
160  // Setup physics
161  model.create_physics("cf3.physics.DynamicModel");
162 
163  // Setup mesh
164  Mesh& mesh = *domain.create_component<Mesh>("Mesh");
165  Tools::MeshGeneration::create_rectangle_tris(mesh, length, length, nb_segments, nb_segments);
166 
167  lss_action->options().set("regions", std::vector<URI>(1, mesh.topology().uri()));
168  ic->get_child("SetInitial")->options().set("regions", std::vector<URI>(1, mesh.topology().uri()));
169  ic->get_child("InitAdvectionSpeed")->options().set("regions", std::vector<URI>(1, mesh.topology().uri()));
170 
171  // Set boundary conditions
172  bc->add_constant_bc("top", "FI", 8.);
173  bc->add_constant_bc("bottom", "FI", 2.);
174  bc->add_constant_bc("left", "FI", 8.);
175  bc->add_constant_bc("right", "FI", 2.);
176 
177  model.create_time();
178  model.time().options().set("end_time", 10.);
179  model.time().options().set("time_step",1.);
180 
181  // Run the solver
182  model.simulate();
183 
184  domain.write_mesh(URI("rdmmerge.plt"));
185 }
186 
188 
189 BOOST_AUTO_TEST_SUITE_END()
190 
191 
boost::proto::terminal< SFOp< ShapeFunctionOp > >::type const N
Time & create_time(const std::string &name="Time")
Create Time component.
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
external boost library namespace
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 ...
void check_close(const Real a, const Real b, const Real threshold)
Check close, for testing purposes.
This header collects all the headers needed for the linear system solver, also including configure-ti...
Handle< common::Action > add_unsteady_solver(const std::string &builder_name)
Definition: Solver.cpp:138
static boost::proto::terminal< void(*)(Real, Real, Real) >::type const _check_close
Basic Classes for Solver applications used by CF.
Definition: Action.cpp:29
URI uri() const
Construct the full path.
Definition: Component.cpp:248
static boost::proto::terminal< ElementQuadratureTag >::type element_quadrature
Use element_quadrature(expr1, expr2, ..., exprN) to evaluate a group of expressions.
void create_rectangle_tris(Mesh &mesh, const Real x_len, const Real y_len, const Uint x_segments, const Uint y_segments)
Create a triangulated version of a 2D rectangular grid.
tuple root
Definition: coolfluid.py:24
boost::shared_ptr< ElementsExpression< ExprT, ElementTypes > > elements_expression(ElementTypes, const ExprT &expr)
Definition: Expression.hpp:292
Manage a collection of UFEM solvers.
Definition: Solver.hpp:31
static boost::proto::terminal< ElementRHS >::type const _x
Terminal for the element RHS vector ("b")
Static functions for mathematical constants.
Definition: Consts.hpp:25
boost::proto::terminal< TransposeFunction >::type const transpose
boost::proto::terminal< SFOp< NablaOp > >::type const nabla
BOOST_AUTO_TEST_CASE(InitMPI)
Basic Classes for Mesh applications used by COOLFluiD.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealVector
Dynamic sized column vector.
Definition: MatrixTypes.hpp:25
tuple model
Global confifuration.
Handle< Component > get_child(const std::string &name)
Definition: Component.cpp:441
void init(int argc=0, char **args=0)
Definition: Comm.cpp:80
Top-level namespace for coolfluid.
Definition: Action.cpp:18
virtual void simulate()
Simulates this model.
boost::proto::terminal< SFOp< CoordinatesOp > >::type const coordinates
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
std::vector< Uint > SizesT
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.
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
OptionList & options()
Definition: Component.cpp:856
Time & time()
Reference to the time.
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
Most basic kernel library.
Definition: Action.cpp:19
std::vector< std::string > StringsT
Send comments to:
COOLFluiD Web Admin