PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
probes.cpp
Go to the documentation of this file.
1
8#include <numeric>
9#include <cstring>
10#include <algorithm>
11
12#include <petscdmcomposite.h>
13#include <petscviewerhdf5.h>
14
15#include <petibm/probes.h>
16
17namespace petibm
18{
19
20namespace misc
21{
22// Factory function to create a probe to monitor the solution.
23PetscErrorCode createProbe(const MPI_Comm &comm,
24 const YAML::Node &node,
25 const type::Mesh &mesh,
26 type::Probe &probe)
27{
28 type::ProbeType probeType;
29
30 PetscFunctionBeginUser;
31
32 probeType = type::str2ProbeType[node["type"].as<std::string>("VOLUME")];
33
34 // either create a probe to monitor a sub-volume or to monitor a single point
35 switch (probeType)
36 {
38 probe = std::make_shared<ProbeVolume>(comm, node, mesh);
39 break;
41 probe = std::make_shared<ProbePoint>(comm, node, mesh);
42 break;
43 default:
44 SETERRQ(comm, PETSC_ERR_ARG_UNKNOWN_TYPE,
45 "Unknown type of probe. Accepted values are: "
46 "VOLUME and POINT");
47 }
48
49 PetscFunctionReturn(0);
50} // createProbe
51
52//***************************************************************************//
53//************************* ProbeBase *************************//
54//***************************************************************************//
55
56// Initialize the probe.
57PetscErrorCode ProbeBase::init(const MPI_Comm &inComm,
58 const YAML::Node &node,
59 const type::Mesh &mesh)
60{
61 PetscErrorCode ierr;
62
63 PetscFunctionBeginUser;
64
65 // store information about the communicator
66 comm = inComm;
67 ierr = MPI_Comm_size(comm, &commSize); CHKERRQ(ierr);
68 ierr = MPI_Comm_rank(comm, &commRank); CHKERRQ(ierr);
69
70 // store information about what to monitor and when to monitor
71 name = node["name"].as<std::string>("unnamed");
72 field = type::str2fd[node["field"].as<std::string>()];
73 path = node["path"].as<std::string>();
74 n_monitor = node["n_monitor"].as<PetscInt>(1);
75 t_start = node["t_start"].as<PetscReal>(0.0);
76 t_end = node["t_end"].as<PetscReal>(1e12);
77
78 viewer = PETSC_NULL;
79
80 PetscFunctionReturn(0);
81} // ProbeBase::init
82
83// Destructor.
85{
86 PetscErrorCode ierr;
87 PetscBool finalized;
88
89 PetscFunctionBeginUser;
90
91 ierr = PetscFinalized(&finalized); CHKERRV(ierr);
92 if (finalized) return;
93
94 ierr = destroy(); CHKERRV(ierr);
95} // ProbeBase::~ProbeBase
96
97// Manually destroy data in the object.
98PetscErrorCode ProbeBase::destroy()
99{
100 PetscErrorCode ierr;
101
102 PetscFunctionBeginUser;
103
104 name = "";
105 path = "";
106 if (viewer != PETSC_NULL) {ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);}
107 comm = MPI_COMM_NULL;
108 commSize = commRank = 0;
109
110 PetscFunctionReturn(0);
111} // ProbeBase::destroy
112
113// Monitor the field solution at points in the index set.
114PetscErrorCode ProbeBase::monitor(const type::Solution &solution,
115 const type::Mesh &mesh,
116 const PetscInt &n,
117 const PetscReal &t)
118{
119 PetscErrorCode ierr;
120
121 PetscFunctionBeginUser;
122
123 // monitor only if appropriate time-step index and appropriate time
124 if (n % n_monitor == 0 && t >= t_start && t <= t_end)
125 {
126 if (field < mesh->dim) // monitor a component of the velocity field
127 {
128 std::vector<Vec> vel(mesh->dim);
129 ierr = DMCompositeGetAccessArray(mesh->UPack, solution->UGlobal,
130 mesh->dim, nullptr,
131 vel.data()); CHKERRQ(ierr);
132 ierr = monitorVec(mesh->da[field], vel[field], n, t); CHKERRQ(ierr);
133 ierr = DMCompositeRestoreAccessArray(mesh->UPack, solution->UGlobal,
134 mesh->dim, nullptr,
135 vel.data()); CHKERRQ(ierr);
136 }
137 else if (field == 3) // monitor the pressure field
138 {
139 ierr = monitorVec(mesh->da[3], solution->pGlobal , n, t); CHKERRQ(ierr);
140 }
141 else
142 SETERRQ(comm, PETSC_ERR_SUP,
143 "Unsupported field. Supported fields are:\n"
144 "\t u (0), v (1), w (2), and p (3).");
145 }
146
147 PetscFunctionReturn(0);
148} // ProbeBase::monitor
149
150//***************************************************************************//
151//************************* ProbeVolume *************************//
152//***************************************************************************//
153
154// Constructor. Initialize the probe.
155ProbeVolume::ProbeVolume(const MPI_Comm &comm,
156 const YAML::Node &node,
157 const type::Mesh &mesh)
158{
159 init(comm, node, mesh);
160} // ProbeVolume::ProbeVolume
161
162// Initialize the probe.
163PetscErrorCode ProbeVolume::init(const MPI_Comm &comm,
164 const YAML::Node &node,
165 const type::Mesh &mesh)
166{
167 PetscErrorCode ierr;
168
169 PetscFunctionBeginUser;
170
171 ierr = ProbeBase::init(comm, node, mesh); CHKERRQ(ierr);
172
173 // store information about the type of PETSc Viewer object to use
174 std::string vtype_str = node["viewer"].as<std::string>("ascii");
175 if (vtype_str == "ascii")
176 viewerType = PETSCVIEWERASCII;
177 else if (vtype_str == "hdf5")
178 viewerType = PETSCVIEWERHDF5;
179 // tolerance to define if a point belong to the sub-volume
180 atol = node["atol"].as<PetscReal>(1e-6);
181 // number of the time-steps over which we accumulate the data
182 // data are added together and we write the time averaged data
183 n_sum = node["n_sum"].as<PetscInt>(0);
184
185 isPetsc = PETSC_NULL;
186 isNatural = PETSC_NULL;
187 dvec = PETSC_NULL;
188
189 // store information about the sub-volume to monitor
190 box = type::RealVec2D(3, type::RealVec1D(2, 0.0));
191 for (auto item : node["box"])
192 {
193 type::Dir dir = type::str2dir[item.first.as<std::string>()];
194 box[dir][0] = item.second[0].as<PetscReal>();
195 box[dir][1] = item.second[1].as<PetscReal>();
196 }
197
198 nPtsDir.resize(3, 1);
199 startIdxDir.resize(3, 0);
200
201 // get information about the location of the sub-mesh
202 ierr = getInfo(mesh, box); CHKERRQ(ierr);
203 // create gridline coordinates for the sub-mesh
204 ierr = createGrid(mesh); CHKERRQ(ierr);
205 // write the sub-mesh to file
206 ierr = writeGrid(path); CHKERRQ(ierr);
207
208 // create a PETSc Viewer to output the data
209 ierr = PetscViewerCreate(comm, &viewer); CHKERRQ(ierr);
210 ierr = PetscViewerSetType(viewer, viewerType); CHKERRQ(ierr);
211 // Note: we set the "append" mode as the output file already exists
212 // (it was created when the sub-mesh was written into it)
213 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_APPEND); CHKERRQ(ierr);
214 ierr = PetscViewerFileSetName(viewer, path.c_str()); CHKERRQ(ierr);
215
216 // create a PETSc Index Set objects to the sub-vector of interest
217 ierr = createIS(mesh); CHKERRQ(ierr);
218 // write Index Set (natural ordering) to file
219 ierr = writeIS(path); CHKERRQ(ierr);
220
221 PetscFunctionReturn(0);
222} // ProbeVolume::init
223
224// Manually destroy data in the object.
225PetscErrorCode ProbeVolume::destroy()
226{
227 PetscErrorCode ierr;
228
229 PetscFunctionBeginUser;
230
231 ierr = ProbeBase::destroy(); CHKERRQ(ierr);
232 if (isPetsc != PETSC_NULL) {ierr = ISDestroy(&isPetsc); CHKERRQ(ierr);}
233 if (isNatural != PETSC_NULL) {ierr = ISDestroy(&isNatural); CHKERRQ(ierr);}
234 if (dvec != PETSC_NULL) {ierr = VecDestroy(&dvec); CHKERRQ(ierr);}
235
236 PetscFunctionReturn(0);
237} // ProbeVolume::destroy
238
239// Get information about the sub-mesh area to monitor.
240PetscErrorCode ProbeVolume::getInfo(const type::Mesh &mesh,
241 const type::RealVec2D &box)
242{
243 std::vector<PetscReal>::iterator low, up;
244
245 PetscFunctionBeginUser;
246
247 // get the starting index along a gridline and the number of points
248 // for each direction of the sub-mesh
249 for (PetscInt d = 0; d < mesh->dim; ++d)
250 {
251 type::RealVec1D line(mesh->coord[field][d],
252 mesh->coord[field][d] + mesh->n[field][d]);
253 low = std::lower_bound(line.begin(), line.end(), box[d][0] - atol);
254 up = std::upper_bound(line.begin(), line.end(), box[d][1] + atol);
255 startIdxDir[d] = low - line.begin();
256 nPtsDir[d] = up - line.begin() - startIdxDir[d];
257 }
258
259 // get the number of points in the sub-volume
260 nPts = std::accumulate(nPtsDir.begin(), nPtsDir.end(),
261 1, std::multiplies<PetscInt>());
262
263 PetscFunctionReturn(0);
264} // ProbeVolume::getInfo
265
266// Create the index sets (PETSc and natural ordering) of the points to monitor.
267PetscErrorCode ProbeVolume::createIS(const type::Mesh &mesh)
268{
269 PetscErrorCode ierr;
270
271 PetscFunctionBeginUser;
272
273 DMDALocalInfo info;
274 ierr = DMDAGetLocalInfo(mesh->da[field], &info); CHKERRQ(ierr);
275 std::vector<PetscInt> indices(nPts);
276 PetscInt count = 0;
277 for (PetscInt k = info.zs; k < info.zs + info.zm; ++k)
278 {
279 for (PetscInt j = info.ys; j < info.ys + info.ym; ++j)
280 {
281 for (PetscInt i = info.xs; i < info.xs + info.xm; ++i)
282 {
283 if (i >= startIdxDir[0] && i < startIdxDir[0] + nPtsDir[0] &&
284 j >= startIdxDir[1] && j < startIdxDir[1] + nPtsDir[1] &&
285 k >= startIdxDir[2] && k < startIdxDir[2] + nPtsDir[2])
286 {
287 indices[count] = k * (info.my * info.mx) + j * info.mx + i;
288 count++;
289 }
290 }
291 }
292 }
293 indices.resize(count);
294
295 // AO will map `indices` from natural ordering to PETSc ordering.
296 // However, we need to keep the vector with natural-ordering indices
297 // so we can post-process the solution in the sub-volume.
298 std::vector<PetscInt> indices2(indices);
299
300 AO ao;
301 ierr = DMDAGetAO(mesh->da[field], &ao); CHKERRQ(ierr);
302 ierr = AOApplicationToPetsc(ao, count, &indices[0]); CHKERRQ(ierr);
303 ierr = ISCreateGeneral(comm, count, &indices[0],
304 PETSC_COPY_VALUES, &isPetsc); CHKERRQ(ierr);
305 // Create IS containing the (natural )index of the points in the sub-volume.
306 ierr = ISCreateGeneral(comm, count, &indices2[0],
307 PETSC_COPY_VALUES, &isNatural); CHKERRQ(ierr);
308
309 PetscFunctionReturn(0);
310} // ProbeVolume::createIS
311
312// Get the sub-mesh area to monitor.
313PetscErrorCode ProbeVolume::createGrid(const type::Mesh &mesh)
314{
315 PetscFunctionBeginUser;
316
317 // store the gridline coordinates of the sub-mesh
318 coord = type::RealVec2D(mesh->dim, type::RealVec1D());
319 for (PetscInt d = 0; d < mesh->dim; ++d)
320 {
321 coord[d].reserve(nPtsDir[d]);
322 PetscInt first = startIdxDir[d];
323 PetscInt last = first + nPtsDir[d];
324 for (PetscInt i = first; i < last; ++i)
325 coord[d].push_back(mesh->coord[field][d][i]);
326 }
327
328 PetscFunctionReturn(0);
329} // ProbeVolume::createGrid
330
331// Write the sub-mesh grid points into a file.
332PetscErrorCode ProbeVolume::writeGrid(const std::string &filePath)
333{
334 PetscErrorCode ierr;
335
336 PetscFunctionBeginUser;
337
338 // only ASCII and HDF5 formats are currently supported
339 if (std::strcmp(viewerType, "ascii") == 0)
340 {
341 ierr = writeGrid_ASCII(filePath); CHKERRQ(ierr);
342 }
343 else if (std::strcmp(viewerType, "hdf5") == 0)
344 {
345 ierr = writeGrid_HDF5(filePath); CHKERRQ(ierr);
346 }
347 else
348 SETERRQ(comm, PETSC_ERR_ARG_UNKNOWN_TYPE,
349 "Unsupported PetscViewerType. Supported types are:\n"
350 "\t PETSCVIEWERASCII and PETSCVIEWERHDF5");
351
352 PetscFunctionReturn(0);
353} // ProbeVolume::writeGrid
354
355// Write the sub mesh into a HDF5 file.
356PetscErrorCode ProbeVolume::writeGrid_HDF5(const std::string &filePath)
357{
358 PetscErrorCode ierr;
359
360 PetscFunctionBeginUser;
361
362 // only the first process in the communicator write the sub-mesh into a file
363 if (commRank == 0)
364 {
365 // because only one process is involved in writing the sub-mesh,
366 // we need to create a temporary viewer
367 PetscViewer viewer2;
368 ierr = PetscViewerCreate(PETSC_COMM_SELF, &viewer2); CHKERRQ(ierr);
369 ierr = PetscViewerSetType(viewer2, PETSCVIEWERHDF5); CHKERRQ(ierr);
370 ierr = PetscViewerFileSetMode(viewer2, FILE_MODE_WRITE); CHKERRQ(ierr);
371 ierr = PetscViewerFileSetName(
372 viewer2, filePath.c_str()); CHKERRQ(ierr);
373 ierr = PetscViewerHDF5PushGroup(viewer2, "mesh"); CHKERRQ(ierr);
374 std::vector<std::string> dirs{"x", "y", "z"};
375 for (unsigned int d = 0; d < coord.size(); ++d)
376 {
377 Vec tmp;
378 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, 1, nPtsDir[d],
379 &coord[d][0], &tmp); CHKERRQ(ierr);
380 ierr = PetscObjectSetName((PetscObject) tmp,
381 dirs[d].c_str()); CHKERRQ(ierr);
382 ierr = VecView(tmp, viewer2); CHKERRQ(ierr);
383 ierr = VecDestroy(&tmp); CHKERRQ(ierr);
384 }
385 ierr = PetscViewerDestroy(&viewer2); CHKERRQ(ierr);
386 }
387
388 PetscFunctionReturn(0);
389} // ProbeVolume::writeGrid_HDF5
390
391// Write the sub mesh into an ASCII file.
392PetscErrorCode ProbeVolume::writeGrid_ASCII(const std::string &filePath)
393{
394 PetscErrorCode ierr;
395
396 PetscFunctionBeginUser;
397
398 // only the first process in the communicator write the sub-mesh into a file
399 if (commRank == 0)
400 {
401 // because only one process is involved in writing the sub-mesh,
402 // we need to create a temporary viewer
403 PetscViewer viewer2;
404 ierr = PetscViewerCreate(PETSC_COMM_SELF, &viewer2); CHKERRQ(ierr);
405 ierr = PetscViewerSetType(viewer2, PETSCVIEWERASCII); CHKERRQ(ierr);
406 ierr = PetscViewerPushFormat(
407 viewer2, PETSC_VIEWER_ASCII_SYMMODU); CHKERRQ(ierr);
408 ierr = PetscViewerFileSetMode(viewer2, FILE_MODE_WRITE); CHKERRQ(ierr);
409 ierr = PetscViewerFileSetName(
410 viewer2, filePath.c_str()); CHKERRQ(ierr);
411 std::vector<std::string> dirs{"x", "y", "z"};
412 for (unsigned int d = 0; d < coord.size(); ++d)
413 {
414 Vec tmp;
415 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, 1, nPtsDir[d],
416 &coord[d][0], &tmp); CHKERRQ(ierr);
417 ierr = PetscObjectSetName((PetscObject) tmp,
418 dirs[d].c_str()); CHKERRQ(ierr);
419 ierr = VecView(tmp, viewer2); CHKERRQ(ierr);
420 ierr = VecDestroy(&tmp); CHKERRQ(ierr);
421 }
422 ierr = PetscViewerPopFormat(viewer2); CHKERRQ(ierr);
423 ierr = PetscViewerDestroy(&viewer2); CHKERRQ(ierr);
424 }
425 ierr = MPI_Barrier(comm); CHKERRQ(ierr);
426
427 PetscFunctionReturn(0);
428} // ProbeVolume::writeGrid_ASCII
429
430// Write the index set (natural ordering) into a file.
431PetscErrorCode ProbeVolume::writeIS(const std::string &filePath)
432{
433 PetscErrorCode ierr;
434
435 PetscFunctionBeginUser;
436
437 // only ASCII and HDF5 formats are currently supported
438 if (std::strcmp(viewerType, "ascii") == 0)
439 {
440 ierr = writeIS_ASCII(filePath); CHKERRQ(ierr);
441 }
442 else if (std::strcmp(viewerType, "hdf5") == 0)
443 {
444 ierr = writeIS_HDF5(filePath); CHKERRQ(ierr);
445 }
446 else
447 SETERRQ(comm, PETSC_ERR_ARG_UNKNOWN_TYPE,
448 "Unsupported PetscViewerType. Supported types are:\n"
449 "\t PETSCVIEWERASCII and PETSCVIEWERHDF5");
450
451 PetscFunctionReturn(0);
452} // ProbeVolume::writeIS
453
454PetscErrorCode ProbeVolume::writeIS_HDF5(const std::string &filePath)
455{
456 PetscErrorCode ierr;
457
458 PetscFunctionBeginUser;
459
460 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_APPEND); CHKERRQ(ierr);
461 ierr = PetscViewerFileSetName(viewer, filePath.c_str()); CHKERRQ(ierr);
462
463 ierr = PetscViewerHDF5PushGroup(viewer, "mesh"); CHKERRQ(ierr);
464 ierr = PetscObjectSetName((PetscObject) isNatural, "IS"); CHKERRQ(ierr);
465 ierr = ISView(isNatural, viewer); CHKERRQ(ierr);
466
467 PetscFunctionReturn(0);
468} // ProbeVolume::writeIS_HDF5
469
470// Write index set (natural ordering) into HDF5 file.
471PetscErrorCode ProbeVolume::writeIS_ASCII(const std::string &filePath)
472{
473 PetscErrorCode ierr;
474
475 PetscFunctionBeginUser;
476
477 ierr = PetscViewerPushFormat(
478 viewer, PETSC_VIEWER_ASCII_SYMMODU); CHKERRQ(ierr);
479 ierr = PetscObjectSetName((PetscObject) isNatural, "IS"); CHKERRQ(ierr);
480 ierr = ISView(isNatural, viewer); CHKERRQ(ierr);
481 ierr = PetscViewerPopFormat(viewer); CHKERRQ(ierr);
482
483 PetscFunctionReturn(0);
484} // ProbeVolume::writeIS_ASCII
485
486// Monitor a sub-region of the full-domain PETSc Vec object.
487PetscErrorCode ProbeVolume::monitorVec(const DM &da,
488 const Vec &fvec,
489 const PetscInt &n,
490 const PetscReal &t)
491{
492 PetscErrorCode ierr;
493
494 PetscFunctionBeginUser;
495
496 // grab the part of the vector that corresponds to the sub-volume
497 Vec svec;
498 ierr = VecGetSubVector(fvec, isPetsc, &svec); CHKERRQ(ierr);
499 if (n_sum != 0) // we accumulate the data over the time-steps
500 {
501 if (dvec == PETSC_NULL)
502 {
503 ierr = VecDuplicate(svec, &dvec); CHKERRQ(ierr);
504 ierr = VecSet(dvec, 0.0); CHKERRQ(ierr);
505 }
506 ierr = VecAXPY(dvec, 1.0, svec); CHKERRQ(ierr);
507 count++;
508 if (count % n_sum == 0) // output the time-averaged data
509 {
510 ierr = VecScale(dvec, 1.0 / count); CHKERRQ(ierr);
511 ierr = writeVec(dvec, t); CHKERRQ(ierr);
512 // reset for the next accumulation cycle
513 ierr = VecSet(dvec, 0.0); CHKERRQ(ierr);
514 count = 0;
515 }
516 }
517 else // output the time-step sub-volume data to file
518 {
519 ierr = writeVec(svec, t); CHKERRQ(ierr);
520 }
521 ierr = VecRestoreSubVector(fvec, isPetsc, &svec); CHKERRQ(ierr);
522
523 PetscFunctionReturn(0);
524} // ProbeVolume::monitorVec
525
526// Write vector data into a file.
527PetscErrorCode ProbeVolume::writeVec(const Vec &vec, const PetscReal &t)
528{
529 PetscErrorCode ierr;
530
531 PetscFunctionBeginUser;
532
533 // only ASCII and HDF5 formats are currently supported
534 if (std::strcmp(viewerType, "ascii") == 0)
535 {
536 ierr = writeVec_ASCII(vec, t); CHKERRQ(ierr);
537 }
538 else if (std::strcmp(viewerType, "hdf5") == 0)
539 {
540 ierr = writeVec_HDF5(vec, t); CHKERRQ(ierr);
541 }
542 else
543 SETERRQ(comm, PETSC_ERR_ARG_UNKNOWN_TYPE,
544 "Unsupported PetscViewerType. Supported types are:\n"
545 "\t PETSCVIEWERASCII and PETSCVIEWERHDF5");
546
547 PetscFunctionReturn(0);
548} // ProbeVolume::writeVec
549
550// Write vector data data into a HDF5 file.
551PetscErrorCode ProbeVolume::writeVec_HDF5(const Vec &vec, const PetscReal &t)
552{
553 PetscErrorCode ierr;
554
555 PetscFunctionBeginUser;
556
557 std::string field_str = type::fd2str[field];
558 ierr = PetscViewerHDF5PushGroup(viewer, field_str.c_str()); CHKERRQ(ierr);
559 ierr = PetscObjectSetName((PetscObject) vec,
560 std::to_string(t).c_str()); CHKERRQ(ierr);
561 ierr = VecView(vec, viewer); CHKERRQ(ierr);
562 if (count != 0) // add number of time-steps accumulated as attribute
563 {
564 ierr = PetscViewerHDF5WriteAttribute(
565 viewer, std::to_string(t).c_str(), "count", PETSC_INT, &count);
566 CHKERRQ(ierr);
567 }
568
569 PetscFunctionReturn(0);
570} // ProbeVolume::writeVec_HDF5
571
572// Write vector data into an ASCII file.
573PetscErrorCode ProbeVolume::writeVec_ASCII(const Vec &vec, const PetscReal &t)
574{
575 PetscErrorCode ierr;
576
577 PetscFunctionBeginUser;
578
579 ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_SYMMODU); CHKERRQ(ierr);
580 ierr = PetscViewerASCIIPrintf(viewer, "\nt = %e\n", t); CHKERRQ(ierr);
581 if (count != 0) // write number of time-steps accumulated
582 {
583 ierr = PetscViewerASCIIPrintf(viewer, "count = %d\n", count); CHKERRQ(ierr);
584 }
585 ierr = VecView(vec, viewer); CHKERRQ(ierr);
586 ierr = PetscViewerPopFormat(viewer); CHKERRQ(ierr);
587
588 PetscFunctionReturn(0);
589} // ProbeVolume::writeVec_ASCII
590
591//***************************************************************************//
592//************************* ProbePoint *************************//
593//***************************************************************************//
594
595// Constructor. Initialize the probe.
596ProbePoint::ProbePoint(const MPI_Comm &comm,
597 const YAML::Node &node,
598 const type::Mesh &mesh)
599{
600 init(comm, node, mesh);
601} // ProbePoint::ProbePoint
602
603// Initialize the probe.
604PetscErrorCode ProbePoint::init(const MPI_Comm &comm,
605 const YAML::Node &node,
606 const type::Mesh &mesh)
607{
608 PetscErrorCode ierr;
609
610 PetscFunctionBeginUser;
611
612 ierr = ProbeBase::init(comm, node, mesh); CHKERRQ(ierr);
613
614 // only ASCII output format is supported
615 viewerType = PETSCVIEWERASCII;
616
617 // get location of the point to interpolate
618 loc = type::RealVec1D(3, 0.0);
619 for (PetscInt d = 0; d < mesh->dim; ++d)
620 loc[d] = node["loc"][d].as<PetscReal>();
621
622 // is the target point in on the current process sub-domain?
623 pointOnLocalProc = mesh->isPointOnLocalProc(loc, field);
624
625 // create a PETSc Viewer to output the data
626 // Note: as the target point belongs to only one sub-domain,
627 // we create a Viewer that involves only one process
628 ierr = PetscViewerCreate(PETSC_COMM_SELF, &viewer); CHKERRQ(ierr);
629 ierr = PetscViewerSetType(viewer, viewerType); CHKERRQ(ierr);
630 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_WRITE); CHKERRQ(ierr);
631 ierr = PetscViewerFileSetName(viewer, path.c_str()); CHKERRQ(ierr);
633 {
634 // create the interpolation object (2D: bi-linear, 3D: tri-linear)
635 ierr = createLinInterp(PETSC_COMM_SELF,
636 loc, mesh, field, interp); CHKERRQ(ierr);
637 }
638 // create a local vector that will have ghost-point values
639 ierr = DMCreateLocalVector(mesh->da[field], &svec); CHKERRQ(ierr);
640
641 PetscFunctionReturn(0);
642} // ProbePoint::init
643
644// Manually destroy data in the object.
645PetscErrorCode ProbePoint::destroy()
646{
647 PetscErrorCode ierr;
648
649 PetscFunctionBeginUser;
650
651 ierr = ProbeBase::destroy(); CHKERRQ(ierr);
653 {
654 ierr = interp->destroy(); CHKERRQ(ierr);
655 }
656 ierr = VecDestroy(&svec); CHKERRQ(ierr);
657
658 PetscFunctionReturn(0);
659} // ProbePoint::destroy
660
661// Monitor a sub-region of the full-domain PETSc Vec object.
662PetscErrorCode ProbePoint::monitorVec(const DM &da,
663 const Vec &fvec,
664 const PetscInt &n,
665 const PetscReal &t)
666{
667 PetscErrorCode ierr;
668
669 PetscFunctionBeginUser;
670
671 // scatter values to local vector
672 ierr = DMGlobalToLocalBegin(da, fvec, INSERT_VALUES, svec); CHKERRQ(ierr);
673 ierr = DMGlobalToLocalEnd(da, fvec, INSERT_VALUES, svec); CHKERRQ(ierr);
674
676 {
677 // interpolate the value at the target point and write it into file
678 ierr = interp->interpolate(da, svec, value); CHKERRQ(ierr);
679 ierr = PetscViewerASCIIPrintf(viewer, "%10.8e\t%10.8e\n", t, value); CHKERRQ(ierr);
680 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_APPEND); CHKERRQ(ierr);
681 }
682
683 PetscFunctionReturn(0);
684} // ProbePoint::monitorVec
685
686} // end of namespace misc
687
688} // end of namespace petibm
virtual PetscErrorCode destroy()
Manually destroy data in the object.
Definition: probes.cpp:98
PetscReal t_end
Monitoring ending time.
Definition: probes.h:81
std::string path
Path of the file to output the solution.
Definition: probes.h:69
PetscErrorCode monitor(const type::Solution &solution, const type::Mesh &mesh, const PetscInt &n, const PetscReal &t)
Monitor the field solution and output data to file.
Definition: probes.cpp:114
PetscReal t_start
Monitoring starting time.
Definition: probes.h:78
PetscMPIInt commRank
Rank of the local process in the MPI communicator.
Definition: probes.h:96
type::Field field
Type of the field to monitor.
Definition: probes.h:72
virtual PetscErrorCode init(const MPI_Comm &comm, const YAML::Node &node, const type::Mesh &mesh)=0
Initialize the probe.
Definition: probes.cpp:57
PetscViewerType viewerType
Type of the PETSc viewer to use.
Definition: probes.h:87
PetscMPIInt commSize
Number of processes in the MPI communicator.
Definition: probes.h:93
PetscViewer viewer
PETSc viewer to output the solution.
Definition: probes.h:84
virtual ~ProbeBase()
Destructor.
Definition: probes.cpp:84
PetscInt n_monitor
Frequency of monitoring the solution (number of time steps).
Definition: probes.h:75
virtual PetscErrorCode monitorVec(const DM &da, const Vec &fvec, const PetscInt &n, const PetscReal &t)=0
Monitor a sub-region of a vector and possibly output data to file.
MPI_Comm comm
MPI communicator.
Definition: probes.h:90
std::string name
Name of the probe as a string.
Definition: probes.h:66
PetscReal value
Interpolated value.
Definition: probes.h:316
PetscErrorCode destroy()
Manually destroy the data.
Definition: probes.cpp:645
Vec svec
Local (ghosted) PETSc Vec object with neighboring values.
Definition: probes.h:313
type::LinInterp interp
Interpolating object to monitor at a single point.
Definition: probes.h:319
PetscErrorCode monitorVec(const DM &da, const Vec &fvec, const PetscInt &n, const PetscReal &t)
Monitor a sub-region of a vector and possibly output data to file.
Definition: probes.cpp:662
PetscErrorCode init(const MPI_Comm &comm, const YAML::Node &node, const type::Mesh &mesh)
Initialize the probe.
Definition: probes.cpp:604
ProbePoint(const MPI_Comm &comm, const YAML::Node &node, const type::Mesh &mesh)
Default constructor.
Definition: probes.cpp:596
type::RealVec1D loc
Coordinates of the point to monitor around.
Definition: probes.h:307
PetscBool pointOnLocalProc
True if target point located on local sub-domain.
Definition: probes.h:310
PetscErrorCode init(const MPI_Comm &comm, const YAML::Node &node, const type::Mesh &mesh)
Initialize the probe.
Definition: probes.cpp:163
PetscErrorCode writeGrid_HDF5(const std::string &filePath)
Write the sub mesh into a HDF5 file.
Definition: probes.cpp:356
type::RealVec2D coord
Grid point coordinates in the volume.
Definition: probes.h:157
ProbeVolume(const MPI_Comm &comm, const YAML::Node &node, const type::Mesh &mesh)
Default constructor.
Definition: probes.cpp:155
PetscErrorCode writeIS_ASCII(const std::string &filePath)
Write index set (natural ordering) into a HDF5 file.
Definition: probes.cpp:471
PetscErrorCode writeGrid(const std::string &filePath)
Write the sub mesh grid points into a file.
Definition: probes.cpp:332
PetscErrorCode getInfo(const type::Mesh &mesh, const type::RealVec2D &box)
Get information about the sub-mesh area to monitor.
Definition: probes.cpp:240
type::IntVec1D nPtsDir
Number of grid points along each direction in the volume.
Definition: probes.h:160
PetscInt n_sum
Frequency of saving the data to file.
Definition: probes.h:172
PetscErrorCode writeIS_HDF5(const std::string &filePath)
Write index set (natural ordering) into a HDF5 file.
Definition: probes.cpp:454
PetscErrorCode monitorVec(const DM &da, const Vec &fvec, const PetscInt &n, const PetscReal &t)
Monitor a sub-region of a vector and possibly output data to file.
Definition: probes.cpp:487
Vec dvec
Sub-vector of the region to monitor.
Definition: probes.h:154
PetscErrorCode destroy()
Manually destroy the data.
Definition: probes.cpp:225
PetscErrorCode writeGrid_ASCII(const std::string &filePath)
Write the sub mesh into an ASCII file.
Definition: probes.cpp:392
PetscInt nPts
Number of grid points in the volume.
Definition: probes.h:163
PetscErrorCode createIS(const type::Mesh &mesh)
Create the index set for the points to monitor.
Definition: probes.cpp:267
PetscErrorCode writeIS(const std::string &filePath)
Write index set (natural ordering) into a file.
Definition: probes.cpp:431
IS isPetsc
Index set for the grid points to monitor (PETSc ordering).
Definition: probes.h:148
PetscErrorCode writeVec_ASCII(const Vec &vec, const PetscReal &t)
Output a PETSc Vec object to an ASCII file.
Definition: probes.cpp:573
type::RealVec2D box
Limits of the volume.
Definition: probes.h:145
PetscReal atol
Absolute tolerance criterion when comparing values.
Definition: probes.h:169
PetscErrorCode createGrid(const type::Mesh &mesh)
Get the sub-mesh area to monitor.
Definition: probes.cpp:313
PetscInt count
Counter to know when to flush to the data to file.
Definition: probes.h:175
PetscErrorCode writeVec(const Vec &vec, const PetscReal &t)
Output a PETSc Vec object to file.
Definition: probes.cpp:527
type::IntVec1D startIdxDir
Index of the first point in the volume in each direction.
Definition: probes.h:166
PetscErrorCode writeVec_HDF5(const Vec &vec, const PetscReal &t)
Output a PETSc Vec object to a HDF5 file.
Definition: probes.cpp:551
IS isNatural
Index set for the grid points to monitor (Natural ordering).
Definition: probes.h:151
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
std::shared_ptr< misc::ProbeBase > Probe
Type definition of Probe.
Definition: probes.h:357
PetscErrorCode createLinInterp(const MPI_Comm &comm, const type::RealVec1D &point, const type::Mesh &mesh, const type::Field &field, type::LinInterp &interp)
Factory function to create a linear interpolation object.
Definition: lininterp.cpp:16
PetscErrorCode createProbe(const MPI_Comm &comm, const YAML::Node &node, const type::Mesh &mesh, type::Probe &probe)
Factory function to create a probe to monitor the solution.
Definition: probes.cpp:23
std::shared_ptr< solution::SolutionBase > Solution
Type definition of solution object.
Definition: solution.h:210
Dir
Legal physical directions.
Definition: type.h:68
std::map< std::string, ProbeType > str2ProbeType
Mapping between std::string and ProbeType.
Definition: type.cpp:53
ProbeType
Type of probe for monitoring solution.
Definition: type.h:123
std::map< Field, std::string > fd2str
Mapping between Field and std::string.
Definition: type.cpp:24
std::vector< RealVec1D > RealVec2D
2D std::vector holding PetscReal.
Definition: type.h:142
std::map< std::string, Dir > str2dir
Mapping between std::string and Dir.
Definition: type.cpp:16
std::vector< PetscReal > RealVec1D
1D std::vector holding PetscReal.
Definition: type.h:140
std::map< std::string, Field > str2fd
Mapping between std::string and Field.
Definition: type.cpp:21
@ VOLUME
Definition: type.h:127
A toolbox for building flow solvers.
Definition: bodypack.h:52
Prototype of probe classes, type::Probe, and factory function.