Rheolef  7.2
an efficient C++ finite element environment
Loading...
Searching...
No Matches
geo_mpi_put.cc
Go to the documentation of this file.
1
21#include "rheolef/config.h"
22#ifdef _RHEOLEF_HAVE_MPI
23#include "rheolef/geo.h"
24#include "rheolef/geo_domain.h"
25#include "rheolef/dis_macros.h"
26#include "rheolef/rheostream.h"
27#include "rheolef/iorheo.h"
28#include "rheolef/index_set.h"
29
30namespace rheolef {
31
32// ----------------------------------------------------------------------------
33// output
34// ----------------------------------------------------------------------------
36struct permutation_proxy {
37 typedef geo_element::size_type size_type;
38 permutation_proxy (const hack_array<geo_element_hack>& elts, size_type first_dis_v = 0)
39 : _elts(elts), _first_dis_v(first_dis_v) {}
40 size_type size() const { return _elts.size(); }
41 size_type operator[] (size_type ie) const {
42 const geo_element& K = _elts [ie];
43 return K.ios_dis_ie() - _first_dis_v;
44 }
45 const permutation_proxy& data() const { return *this; }
46// data:
47 const hack_array<geo_element_hack>& _elts;
48 size_type _first_dis_v;
49};
50template <class T>
53{
54 using namespace std;
56 check_macro (format[iorheo::rheo], "geo distributed: unsupported geo output file format=" << format);
57
60 size_type my_proc = comm.rank();
61 size_type nproc = comm.size();
62 size_type dis_nnod = base::_node.dis_size ();
64 //
65 // 1) put header
66 //
67 ops << setprecision(numeric_limits<Float>::digits10)
68 << "#!geo" << endl
69 << endl
70 << "mesh" << endl
72
73 geo_header h;
74 h.dimension = base::_dimension;
75 h.sys_coord = base::_sys_coord;
76 h.order = base::order();
77 h.dis_size_by_variant [0] = base::_node.dis_size();
79 h.dis_size_by_variant [variant] = base::_geo_element [variant].dis_size();
80 }
81 ops << endl << h << endl;
82 //
83 // 2) put nodes
84 //
85 if (base::_dimension > 0) {
86 T rounding_prec = iorheo::getrounding_precision(ops.os());
87 if (rounding_prec == 0) {
88 base::_node.permuted_put_values (ops,
91 } else {
92 base::_node.permuted_put_values (ops,
95 }
96 ops << endl;
97 }
98 //
99 // 3) build a node permutation disarray on io_proc
100 //
101 // elements are permuted back to original order and may
102 // refer to original node numbering
103 // TODO: disarray node_perm = merge_array_on_proc (_inod2ios_dis_inod , io_proc);
104 std::vector<size_type> node_perm ((my_proc == io_proc) ? dis_nnod : 0);
105 size_type tag_gather = distributor::get_new_tag();
106 if (my_proc == io_proc) {
107 size_type i_start = _inod2ios_dis_inod.ownership().first_index(my_proc);
108 size_type i_size = _inod2ios_dis_inod.ownership().size (my_proc);
109 for (size_type i = 0; i < i_size; i++) {
110 node_perm [i_start + i] = _inod2ios_dis_inod [i];
111 }
112 for (size_type iproc = 0; iproc < nproc; iproc++) {
113 if (iproc == my_proc) continue;
114 size_type i_start = _inod2ios_dis_inod.ownership().first_index(iproc);
115 size_type i_size = _inod2ios_dis_inod.ownership().size (iproc);
116 comm.recv (iproc, tag_gather, node_perm.begin().operator->() + i_start, i_size);
117 }
118 } else {
119 comm.send (0, tag_gather, _inod2ios_dis_inod.begin().operator->(), _inod2ios_dis_inod.size());
120 }
121 //
122 // 4) put elements, faces & edges
123 //
124 geo_element_permuted_put put_element (node_perm);
125 for (size_type dim = base::_gs._map_dimension; dim > 0; dim--) {
126 size_type first_dis_v = 0;
129
130 if (base::_gs.ownership_by_variant[variant].dis_size() == 0) continue;
131 permutation_proxy perm (base::_geo_element [variant], first_dis_v);
132 base::_geo_element [variant].permuted_put_values (ops, perm, put_element);
133 first_dis_v += base::_gs.ownership_by_variant[variant].dis_size();
134 }
135 ops << endl;
136 }
137 //
138 // 6) put domains
139 //
141 iter = base::_domains.begin(),
142 last = base::_domains.end();
143 iter != last; ++iter) {
144 ops << endl;
145 (*iter).put (ops, *this);
146 }
147 return ops;
148}
149template <class T>
150void
152{
153 fatal_macro ("dump: not yet");
154#ifdef TODO
155 base::_node.dump (name + "-node");
156 base::_geo_element[variant].dump(name + "-elem-"+ref::name(variant));
157#endif // TODO
158}
159// ----------------------------------------------------------------------------
160// instanciation in library
161// ----------------------------------------------------------------------------
162template class geo_rep<Float,distributed>;
163
164} // namespace rheolef
165#endif // _RHEOLEF_HAVE_MPI
see the communicator page for the full documentation
static tag_type get_new_tag()
returns a new tag
the finite element boundary domain
std::vector< domain_indirect_basic< distributed > > _domains
Definition geo.h:680
disarray< node_type, distributed > _node
Definition geo.h:684
std::array< hack_array< geo_element_hack, distributed >, reference_element::max_variant > _geo_element
Definition geo.h:678
const communicator & comm() const
Definition geo.h:637
size_type dis_size(size_type dim) const
reference_element::size_type size_type
const_iterator end(size_type dim) const
Definition geo.h:1004
base::size_type size_type
Definition geo.h:934
disarray< size_type > _inod2ios_dis_inod
Definition geo.h:1064
base::const_iterator const_iterator
Definition geo.h:941
const_iterator begin(size_type dim) const
Definition geo.h:1003
sequential mesh representation
Definition geo.h:778
std::bitset< last > flag_type
Definition iorheo.h:440
static flag_type flags(std::ios &s)
Definition iorheo.cc:244
static flag_type format_field
Definition iorheo.h:445
odiststream: see the diststream page for the full documentation
Definition diststream.h:137
std::ostream & os()
Definition diststream.h:247
static size_type io_proc()
Definition diststream.cc:79
static const variant_type max_variant
static variant_type last_variant_by_dimension(size_type dim)
static variant_type first_variant_by_dimension(size_type dim)
static const variant_type p
#define fatal_macro(message)
Definition dis_macros.h:33
Expr1::float_type T
Definition field_expr.h:230
check_macro(expr1.have_homogeneous_space(Xh1), "dual(expr1,expr2); expr1 should have homogeneous space. HINT: use dual(interpolate(Xh, expr1),expr2)")
This file is part of Rheolef.
STL namespace.
point output helper
Definition geo.h:163
point output helper, with rounding feature
Definition geo.h:171