v0.15.0
Loading...
Searching...
No Matches
Functions | Variables
mesh_data_transfer.cpp File Reference

Transfer data from source mesh to target mesh. More...

#include <MoFEM.hpp>
#include <mbcoupler/Coupler.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 

Detailed Description

Transfer data from source mesh to target mesh.

The data is transferred from tags on 3D elements of the source mesh to tags on 3D elements in the target mesh. This code uses MBCoupler functionality of MOAB and is based on mbcoupler_test.cpp from MOAB.

For this tool to work, source and target files should be partitioned using mbpart tool from MOAB, e.g. mbpart -z PHG -t 4 sslv116_hex.cub sslv116_hex_4mp.h5m

Definition in file mesh_data_transfer.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 21 of file mesh_data_transfer.cpp.

21 {
22 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
23
24 try {
25
26 moab::Core mb_instance;
27 moab::Interface &moab = mb_instance;
28
29 // Create MoFEM database
30 MoFEM::Core core(moab);
31 MoFEM::Interface &m_field = core;
32
33 char mesh_source_file[255] = "source.h5m";
34 char mesh_target_file[255] = "target.h5m";
35 char mesh_out_file[255] = "out.h5m";
36 char iterp_tag_name[255] = "INTERNAL_STRESS";
37
38 int interp_order = 1;
39
40 PetscBool use_target_verts = PETSC_TRUE;
41
42 int atom_test = 0;
43
44 double toler = 5.e-10;
45
46 PetscOptionsBegin(m_field.get_comm(), "", "mesh data transfer tool",
47 "none");
48 CHKERR PetscOptionsString("-source_file", "source mesh file name", "",
49 "source.h5m", mesh_source_file, 255, PETSC_NULLPTR);
50 CHKERR PetscOptionsString("-target_file", "target mesh file name", "",
51 "target.h5m", mesh_target_file, 255, PETSC_NULLPTR);
52 CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
53 "out.h5m", mesh_out_file, 255, PETSC_NULLPTR);
54 CHKERR PetscOptionsString("-interp_tag", "interpolation tag name", "",
55 "INTERNAL_STRESS", iterp_tag_name, 255,
56 PETSC_NULLPTR);
57 CHKERR PetscOptionsInt("-interp_order", "interpolation order", "", 0,
58 &interp_order, PETSC_NULLPTR);
59 CHKERR PetscOptionsBool("-use_target_verts",
60 "use target vertices for interpolation", "",
61 use_target_verts, &use_target_verts, PETSC_NULLPTR);
62 CHKERR PetscOptionsInt("-atom_test", "atom test number", "", 0, &atom_test,
63 PETSC_NULLPTR);
64
65 PetscOptionsEnd();
66
67 MOFEM_LOG_CHANNEL("WORLD");
68 MOFEM_LOG_TAG("WORLD", "mesh_data_transfer");
69
70 Coupler::Method method;
71 switch (interp_order) {
72 case 0:
73 method = Coupler::CONSTANT;
74 break;
75 case 1:
76 method = Coupler::LINEAR_FE;
77 break;
78 default:
79 SETERRQ(PETSC_COMM_WORLD, MOFEM_NOT_IMPLEMENTED,
80 "Unsupported interpolation order");
81 }
82
83 std::vector<std::string> mesh_files(2);
84 mesh_files[0] = string(mesh_source_file);
85 mesh_files[1] = string(mesh_target_file);
86
87 int nprocs, rank;
88 ierr = MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
89 CHKERRQ(ierr);
90 ierr = MPI_Comm_rank(MPI_COMM_WORLD, &rank);
91 CHKERRQ(ierr);
92
93 std::string read_opts, write_opts;
94 read_opts = "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARTITION_"
95 "DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS";
96 if (nprocs > 1)
97 read_opts += ";PARALLEL_GHOSTS=3.0.1";
98 write_opts = (nprocs > 1) ? "PARALLEL=WRITE_PART" : "";
99
100 std::vector<ParallelComm *> pcs(mesh_files.size());
101 std::vector<EntityHandle> roots(mesh_files.size());
102
103 for (unsigned int i = 0; i < mesh_files.size(); i++) {
104 pcs[i] = new ParallelComm(&moab, MPI_COMM_WORLD);
105 int index = pcs[i]->get_id();
106 std::string newread_opts;
107 std::ostringstream extra_opt;
108 extra_opt << ";PARALLEL_COMM=" << index;
109 newread_opts = read_opts + extra_opt.str();
110
111 CHKERR moab.create_meshset(MESHSET_SET, roots[i]);
112
113 CHKERR moab.load_file(mesh_files[i].c_str(), &roots[i],
114 newread_opts.c_str());
115 }
116
117 Range src_elems, targ_elems, targ_verts;
118 CHKERR pcs[0]->get_part_entities(src_elems, 3);
119
120 Tag interp_tag;
121 CHKERR moab.tag_get_handle(iterp_tag_name, interp_tag);
122
123 int interp_tag_len;
124 CHKERR moab.tag_get_length(interp_tag, interp_tag_len);
125
126 if (interp_tag_len != 1 && interp_tag_len != 3 && interp_tag_len != 9) {
127 SETERRQ(PETSC_COMM_WORLD, MOFEM_NOT_IMPLEMENTED,
128 "Unsupported interpolation tag length: %d", interp_tag_len);
129 }
130
131 Coupler mbc(&moab, pcs[0], src_elems, 0, true);
132
133 std::vector<double> vpos; // the positions we are interested in
134 int num_pts = 0;
135
136 Range tmp_verts;
137
138 // First get all vertices adj to partition entities in target mesh
139 CHKERR pcs[1]->get_part_entities(targ_elems, 3);
140 if (!use_target_verts) {
141 targ_verts = targ_elems;
142 } else {
143 CHKERR moab.get_adjacencies(targ_elems, 0, false, targ_verts,
144 moab::Interface::UNION);
145 }
146
147 // Then get non-owned verts and subtract
148 CHKERR pcs[1]->get_pstatus_entities(0, PSTATUS_NOT_OWNED, tmp_verts);
149 targ_verts = subtract(targ_verts, tmp_verts);
150
151 // get position of these entities; these are the target points
152 num_pts = (int)targ_verts.size();
153 vpos.resize(3 * targ_verts.size());
154 CHKERR moab.get_coords(targ_verts, &vpos[0]);
155
156 // Locate those points in the source mesh
157 CHKERR mbc.locate_points(&vpos[0], num_pts, 0, toler);
158
159 std::vector<double> source_data(interp_tag_len * src_elems.size(), 0.0);
160 std::vector<double> target_data(interp_tag_len * num_pts, 0.0);
161
162 CHKERR moab.tag_get_data(interp_tag, src_elems, &source_data[0]);
163
164 Tag scalar_tag, adj_count_tag;
165 double def_scl = 0;
166 string scalar_tag_name = string(iterp_tag_name) + "_COMP";
167 CHKERR moab.tag_get_handle(scalar_tag_name.c_str(), 1, MB_TYPE_DOUBLE,
168 scalar_tag, MB_TAG_CREAT | MB_TAG_DENSE,
169 &def_scl);
170
171 string adj_count_tag_name = "ADJ_COUNT";
172 double def_adj = 0;
173 CHKERR moab.tag_get_handle(adj_count_tag_name.c_str(), 1, MB_TYPE_DOUBLE,
174 adj_count_tag, MB_TAG_CREAT | MB_TAG_DENSE,
175 &def_adj);
176
177 std::vector<double> source_data_scalar(src_elems.size());
178
179 // MBCoupler functionality supports only scalar meshes. For the case of
180 // vector or tensor tags we need to save each component as a scalar tag
181 for (int itag = 0; itag < interp_tag_len; itag++) {
182
183 for (int ielem = 0; ielem < src_elems.size(); ielem++) {
184 source_data_scalar[ielem] = source_data[itag + ielem * interp_tag_len];
185 }
186
187 CHKERR moab.tag_set_data(scalar_tag, src_elems, &source_data_scalar[0]);
188
189 if (interp_order == 1) {
190 // In case of linear interpolation, data needs to be saved on source
191 // mesh vertices. Here we compute average value of the data on the
192 // vertices from all adjacent 3D elements
193
194 Range src_verts;
195 CHKERR moab.get_connectivity(src_elems, src_verts, true);
196
197 CHKERR moab.tag_clear_data(scalar_tag, src_verts, &def_scl);
198 CHKERR moab.tag_clear_data(adj_count_tag, src_verts, &def_adj);
199
200 for (auto &tet : src_elems) {
201 double tet_data = 0;
202 CHKERR moab.tag_get_data(scalar_tag, &tet, 1, &tet_data);
203
204 Range adj_verts;
205 CHKERR moab.get_connectivity(&tet, 1, adj_verts, true);
206
207 std::vector<double> adj_vert_data(adj_verts.size(), 0.0);
208 std::vector<double> adj_vert_count(adj_verts.size(), 0.0);
209
210 CHKERR moab.tag_get_data(scalar_tag, adj_verts, &adj_vert_data[0]);
211 CHKERR moab.tag_get_data(adj_count_tag, adj_verts,
212 &adj_vert_count[0]);
213
214 for (int ivert = 0; ivert < adj_verts.size(); ivert++) {
215 adj_vert_data[ivert] += tet_data;
216 adj_vert_count[ivert] += 1;
217 }
218
219 CHKERR moab.tag_set_data(scalar_tag, adj_verts, &adj_vert_data[0]);
220 CHKERR moab.tag_set_data(adj_count_tag, adj_verts,
221 &adj_vert_count[0]);
222 }
223
224 // We need to reduce tags for the parallel case
225 std::vector<Tag> tags;
226 tags.push_back(scalar_tag);
227 tags.push_back(adj_count_tag);
228 pcs[0]->reduce_tags(tags, tags, MPI_SUM, src_verts);
229
230 std::vector<double> src_vert_data(src_verts.size(), 0.0);
231 std::vector<double> src_vert_adj_count(src_verts.size(), 0.0);
232
233 CHKERR moab.tag_get_data(scalar_tag, src_verts, &src_vert_data[0]);
234 CHKERR moab.tag_get_data(adj_count_tag, src_verts,
235 &src_vert_adj_count[0]);
236
237 for (int ivert = 0; ivert < src_verts.size(); ivert++) {
238 src_vert_data[ivert] /= src_vert_adj_count[ivert];
239 }
240 CHKERR moab.tag_set_data(scalar_tag, src_verts, &src_vert_data[0]);
241 }
242
243 std::vector<double> target_data_scalar(num_pts, 0.0);
244 CHKERR mbc.interpolate(method, scalar_tag_name, &target_data_scalar[0]);
245
246 for (int ielem = 0; ielem < num_pts; ielem++) {
247 target_data[itag + ielem * interp_tag_len] = target_data_scalar[ielem];
248 }
249 }
250
251 CHKERR moab.tag_delete(scalar_tag);
252 CHKERR moab.tag_delete(adj_count_tag);
253
254 // Use original tag
255 CHKERR moab.tag_set_data(interp_tag, targ_verts, &target_data[0]);
256
257 if (use_target_verts) {
258 // If the data was interpolated to target mesh vertices, we need to
259 // average the vertex values for each 3D element and save on its tag
260
261 // We need to reduce tags for the parallel case
262 std::vector<Tag> tags;
263 tags.push_back(interp_tag);
264 pcs[1]->reduce_tags(tags, tags, MPI_SUM, targ_verts);
265
266 for (auto &tet : targ_elems) {
267 Range adj_verts;
268 CHKERR moab.get_connectivity(&tet, 1, adj_verts, true);
269
270 std::vector<double> adj_vert_data(adj_verts.size() * interp_tag_len,
271 0.0);
272 CHKERR moab.tag_get_data(interp_tag, adj_verts, &adj_vert_data[0]);
273
274 std::vector<double> tet_data(interp_tag_len, 0.0);
275 for (int itag = 0; itag < interp_tag_len; itag++) {
276 for (int i = 0; i < adj_verts.size(); i++) {
277 tet_data[itag] += adj_vert_data[i * interp_tag_len + itag];
278 }
279 tet_data[itag] /= adj_verts.size();
280 }
281 CHKERR moab.tag_set_data(interp_tag, &tet, 1, &tet_data[0]);
282 }
283 }
284
285 if (atom_test == 1 || atom_test == 2) {
286 // Check the integrated stress for the SSLV116 test
287
288 std::vector<double> data_integ(interp_tag_len, 0.0);
289 for (auto &tet : targ_elems) {
290
291 std::vector<double> tet_data(interp_tag_len, 0.0);
292 CHKERR moab.tag_get_data(interp_tag, &tet, 1, &tet_data[0]);
293
294 const EntityHandle *vert_conn;
295 int vert_num;
296 CHKERR moab.get_connectivity(tet, vert_conn, vert_num, true);
297 std::vector<double> vpos(3 * vert_num);
298 CHKERR moab.get_coords(vert_conn, vert_num, &vpos[0]);
299 double vol = Tools::tetVolume(vpos.data());
300
301 for (int itag = 0; itag < interp_tag_len; itag++) {
302 data_integ[itag] += tet_data[itag] * vol;
303 }
304 }
305
306 std::vector<double> global_data_integ(interp_tag_len, 0.0);
307 MPI_Allreduce(&data_integ[0], &global_data_integ[0], interp_tag_len,
308 MPI_DOUBLE, MPI_SUM, m_field.get_comm());
309
310 std::string temp = "Integrated stress for sslv116 test: ";
311 for (int i = 0; i < 9; ++i) {
312 std::stringstream ss;
313 ss << std::scientific << std::setprecision(3) << global_data_integ[i];
314 temp += ss.str() + " ";
315 }
316 MOFEM_LOG("WORLD", Sev::inform) << temp;
317
318 double non_zero_val = 1.655e12;
319 double non_zero_tol = 5.1e-2;
320 if (interp_order == 1) {
321 non_zero_tol = 1e-3;
322 }
323 std::vector<int> non_zero_inds(3);
324 std::vector<int> zero_inds(6);
325 if (atom_test == 1) {
326 non_zero_inds = {0, 4, 8};
327 zero_inds = {1, 2, 3, 5, 6, 7};
328 } else {
329 non_zero_inds = {0, 1, 2};
330 zero_inds = {3, 4, 5, 6, 7, 8};
331 }
332 bool non_zero_check =
333 all_of(begin(non_zero_inds), end(non_zero_inds), [&](int i) {
334 return abs(global_data_integ[i] - non_zero_val) / non_zero_val <
335 non_zero_tol;
336 });
337 bool zero_check = all_of(begin(zero_inds), end(zero_inds), [&](int i) {
338 return abs(global_data_integ[i]) < 1e-12;
339 });
340 if (!non_zero_check || !zero_check) {
341 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
342 "Wrong value of the integrated stress");
343 }
344 } else if (atom_test) {
345 SETERRQ(PETSC_COMM_WORLD, MOFEM_NOT_IMPLEMENTED,
346 "Unknown test number: %d", atom_test);
347 }
348
349 Range part_sets;
350 // Only save the target mesh
351 part_sets.insert(roots[1]);
352 std::string new_write_opts;
353 std::ostringstream extra_opt;
354 if (nprocs > 1) {
355 int index = pcs[1]->get_id();
356 extra_opt << ";PARALLEL_COMM=" << index;
357 }
358 new_write_opts = write_opts + extra_opt.str();
359 CHKERR moab.write_file(mesh_out_file, NULL, new_write_opts.c_str(),
360 part_sets);
361 if (0 == rank) {
362 MOFEM_LOG("WORLD", Sev::inform) << "Wrote file " << mesh_out_file;
363 }
364
365 for (unsigned int i = 0; i < mesh_files.size(); i++)
366 delete pcs[i];
367 }
369
371
372 return 0;
373}
#define CATCH_ERRORS
Catch errors.
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define CHKERR
Inline error check.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
static char help[]
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
int atom_test
Atom test.
Definition plastic.cpp:121
void temp(int x, int y=10)
Definition simple.cpp:4
virtual MPI_Comm & get_comm() const =0
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
static double tetVolume(const double *coords)
Calculate volume of tetrahedron.
Definition Tools.cpp:30

Variable Documentation

◆ help

char help[] = "...\n\n"
static

Definition at line 19 of file mesh_data_transfer.cpp.