MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
particleset.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "particleset.hpp"
13
14#if defined(MFEM_USE_MPI) && defined(MFEM_USE_GSLIB)
15
16// Ignore warnings from the gslib header (GCC version)
17#ifdef MFEM_HAVE_GCC_PRAGMA_DIAGNOSTIC
18#pragma GCC diagnostic push
19#pragma GCC diagnostic ignored "-Wunused-function"
20#endif
21
22namespace gslib
23{
24extern "C"
25{
26#include <gslib.h>
27} // extern C
28} // namespace gslib
29
30#ifdef MFEM_HAVE_GCC_PRAGMA_DIAGNOSTIC
31#pragma GCC diagnostic pop
32#endif
33
34#endif // MFEM_USE_MPI && MFEM_USE_GSLIB
35
36
37namespace mfem
38{
39
40Particle::Particle(int dim, const Array<int> &field_vdims, int num_tags)
41 : coords(dim), fields(), tags()
42{
43 coords = 0.0;
44
45 fields.reserve(field_vdims.Size());
46 for (int f = 0; f < field_vdims.Size(); f++)
47 {
48 fields.emplace_back(field_vdims[f]);
49 fields.back() = 0.0;
50 }
51
52 tags.reserve(num_tags);
53 for (int t = 0; t < num_tags; t++)
54 {
55 tags.emplace_back(1);
56 tags.back()[0] = 0;
57 }
58}
59
60void Particle::SetTagRef(int t, int *tag_data)
61{
62 MFEM_ASSERT(t >= 0 &&
63 static_cast<size_t>(t) < tags.size(), "Invalid tag index");
64 tags[t].MakeRef(tag_data, 1);
65}
66
67void Particle::SetFieldRef(int f, real_t *field_data)
68{
69 MFEM_ASSERT(f >= 0 &&
70 static_cast<size_t>(f) < fields.size(), "Invalid field "
71 "index");
72 Vector temp(field_data, fields[f].Size());
73 fields[f].MakeRef(temp, 0, fields[f].Size());
74}
75
76bool Particle::operator==(const Particle &rhs) const
77{
78 // Compare coordinate size and values
79 if (coords.Size() != rhs.coords.Size())
80 {
81 return false;
82 }
83 for (int d = 0; d < coords.Size(); d++)
84 {
85 if (coords[d] != rhs.coords[d])
86 {
87 return false;
88 }
89 }
90 // Compare fields vdim and values
91 if (fields.size() != rhs.fields.size())
92 {
93 return false;
94 }
95 for (size_t f = 0; f < fields.size(); f++)
96 {
97 if (fields[f].Size() != rhs.fields[f].Size())
98 {
99 return false;
100 }
101 for (int c = 0; c < fields[f].Size(); c++)
102 {
103 if (fields[f][c] != rhs.fields[f][c])
104 {
105 return false;
106 }
107 }
108 }
109 // Compare tags size and values
110 if (tags.size() != rhs.tags.size())
111 {
112 return false;
113 }
114 for (size_t t = 0; t < tags.size(); t++)
115 {
116 if (tags[t][0] != rhs.tags[t][0])
117 {
118 return false;
119 }
120 }
121
122 return true;
123}
124
125void Particle::Print(std::ostream &os) const
126{
127 os << "Coords: (";
128 for (int d = 0; d < coords.Size(); d++)
129 {
130 os << coords[d] << ( (d+1 < coords.Size()) ? "," : ")\n");
131 }
132 for (size_t f = 0; f < fields.size(); f++)
133 {
134 os << "Field " << f << ": (";
135 for (int c = 0; c < fields[f].Size(); c++)
136 {
137 os << fields[f][c] << ( (c+1 < fields[f].Size()) ? "," : ")\n");
138 }
139 }
140 for (size_t t = 0; t < tags.size(); t++)
141 {
142 os << "Tag " << t << ": " << tags[t][0] << "\n";
143 }
144}
145
146Array<Ordering::Type> ParticleSet::GetOrderingArray(Ordering::Type o, int N)
147{
148 Array<Ordering::Type> ordering_arr(N);
149 ordering_arr = o;
150 return ordering_arr;
151}
152std::string ParticleSet::GetDefaultFieldName(int i)
153{
154 return "Field_" + std::to_string(i);
155}
156
157std::string ParticleSet::GetDefaultTagName(int i)
158{
159 return "Tag_" + std::to_string(i);
160}
161
162Array<const char*> ParticleSet::GetEmptyNameArray(int N)
163{
164 Array<const char*> names(N);
165 names = nullptr;
166 return names;
167}
168
169#ifdef MFEM_USE_MPI
170int ParticleSet::GetRank(MPI_Comm comm_)
171{
172 int r; MPI_Comm_rank(comm_, &r);
173 return r;
174}
175int ParticleSet::GetSize(MPI_Comm comm_)
176{
177 int s; MPI_Comm_size(comm_, &s);
178 return s;
179}
180#endif // MFEM_USE_MPI
181
183{
184 ids.Reserve(res);
185
186 // Reserve fields
187 for (int f = -1; f < GetNFields(); f++)
188 {
189 ParticleVector &pv = (f == -1 ? coords : *fields[f]);
190 pv.Reserve(res*pv.GetVDim());
191 }
192
193 // Reserve tags
194 for (int t = 0; t < GetNTags(); t++)
195 {
196 tags[t]->Reserve(res);
197 }
198
199}
200
202{
203 Array<int> field_vdims(GetNFields());
204 for (int f = 0; f < GetNFields(); f++)
205 {
206 field_vdims[f] = Field(f).GetVDim();
207 }
208 return field_vdims;
209}
210
212 Array<int> *new_indices)
213{
214 int num_add = new_ids.Size();
215 int old_np = GetNParticles();
216 int new_np = old_np + num_add;
217
218 // Set indices of new particles
219 if (new_indices)
220 {
221 new_indices->SetSize(num_add);
222 for (int i = 0; i < num_add; i++)
223 {
224 (*new_indices)[i] = ids.Size() + i;
225 }
226 }
227 // Add new ids
228 ids.Append(new_ids);
229
230 // Update data
231 for (int f = -1; f < GetNFields(); f++)
232 {
233 ParticleVector &pv = (f == -1 ? coords : *fields[f]);
234 pv.SetNumParticles(new_np); // does not delete existing data
235 }
236
237 // Update tags
238 for (int t = 0; t < GetNTags(); t++)
239 {
240 tags[t]->SetSize(new_np);
241 }
242}
243
244#if defined(MFEM_USE_MPI) && defined(MFEM_USE_GSLIB)
245
246/// \cond DO_NOT_DOCUMENT
247template<size_t NBytes>
248void ParticleSet::TransferParticlesImpl(ParticleSet &pset,
249 const Array<int> &send_idxs,
250 const Array<unsigned int> &send_ranks)
251{
252 struct pdata_t
253 {
254 alignas(real_t) std::array<std::byte, NBytes> data;
255 IDType id;
256 };
257
258 int nreals = pset.GetFieldVDims().Sum() + pset.Coords().GetVDim();
259 int ntags = pset.GetNTags();
260 size_t nbytes = nreals*sizeof(real_t) + ntags*sizeof(int);
261 MFEM_VERIFY(nbytes <= NBytes, "More data than can be packed.");
262
263 using parr_t = pdata_t;
264 gslib::array gsl_arr;
265 parr_t *pdata_arr;
266 array_init(parr_t, &gsl_arr, send_idxs.Size());
267 pdata_arr = (parr_t*) gsl_arr.ptr;
268
269 gsl_arr.n = send_idxs.Size();
270 for (int i = 0; i < send_idxs.Size(); i++)
271 {
272 parr_t &pdata = pdata_arr[i];
273 pdata.id = pset.GetIDs()[send_idxs[i]];
274
275 // Copy particle data directly into pdata
276 size_t counter = 0;
277 for (int f = -1; f < pset.GetNFields(); f++)
278 {
279 ParticleVector &pv = (f == -1 ? pset.Coords() : pset.Field(f));
280 for (int c = 0; c < pv.GetVDim(); c++)
281 {
282 std::memcpy(pdata.data.data() + counter, &pv(send_idxs[i], c),
283 sizeof(real_t));
284 counter += sizeof(real_t);
285 }
286 }
287
288 // Copy tags
289 for (int t = 0; t < pset.GetNTags(); t++)
290 {
291 Array<int> &tag_arr = pset.Tag(t);
292 std::memcpy(pdata.data.data() + counter, &tag_arr[send_idxs[i]],
293 sizeof(int));
294 counter += sizeof(int);
295 }
296 }
297
298 int nparticles = pset.GetNParticles();
299 int nsend = send_idxs.Size();
300
301 // Transfer particles
302 sarray_transfer_ext(parr_t, &gsl_arr, send_ranks.GetData(),
303 sizeof(unsigned int), pset.cr);
304
305 // Make sure we have enough space for received particles
306 int nrecv = (int) gsl_arr.n;
307 int ndelete = nsend - nrecv;
308 if (ndelete > 0)
309 {
310 // Remove unneeded particles
311 auto datap = const_cast<int*>(send_idxs.GetData());
312 Array<int> delete_idxs(datap + nrecv, ndelete);
313 pset.RemoveParticles(delete_idxs);
314 }
315 else
316 {
317 pset.Reserve(nparticles-ndelete);
318 }
319
320 pdata_arr = (parr_t*) gsl_arr.ptr;
321
322 // Add newly-recvd data directly to active state
323 for (int i = 0; i < nrecv; i++)
324 {
325 parr_t &pdata = pdata_arr[i];
326 IDType id = pdata.id;
327
328 int new_loc_idx;
329 if (i < nsend) // update existing particle
330 {
331 new_loc_idx = send_idxs[i];
332 pset.UpdateID(new_loc_idx, id);
333 }
334 else
335 {
336 // add new particle
337 Array<int> idx_temp;
338 pset.AddParticles(Array<IDType>({id}), &idx_temp);
339 new_loc_idx = idx_temp[0]; // Get index of newly-added particle
340 }
341
342 size_t counter = 0;
343 for (int f = -1; f < pset.GetNFields(); f++)
344 {
345 ParticleVector &pv = (f == -1 ? pset.Coords() : pset.Field(f));
346 for (int c = 0; c < pv.GetVDim(); c++)
347 {
348 real_t& val = pv(new_loc_idx, c);
349 std::memcpy(&val, pdata.data.data() + counter, sizeof(real_t));
350 counter += sizeof(real_t);
351 }
352 }
353
354 for (int t = 0; t < pset.GetNTags(); t++)
355 {
356 Array<int> &tag_arr = pset.Tag(t);
357 std::memcpy(&tag_arr[new_loc_idx],
358 pdata.data.data() + counter, sizeof(int));
359 counter += sizeof(int);
360 }
361 }
362 array_free(&gsl_arr);
363}
364
365template<size_t NBytes>
366ParticleSet::TransferParticlesType ParticleSet::TransferParticles::Kernel()
367{
368 return &ParticleSet::TransferParticlesImpl<NBytes>;
369}
370
371ParticleSet::Kernels::Kernels()
372{
373 constexpr size_t sizd = sizeof(real_t);
374 TransferParticles::Specialization<2*sizd>::Add();
375 TransferParticles::Specialization<3*sizd>::Add();
376 TransferParticles::Specialization<4*sizd>::Add();
377 TransferParticles::Specialization<8*sizd>::Add();
378 TransferParticles::Specialization<12*sizd>::Add();
379 TransferParticles::Specialization<16*sizd>::Add();
380 TransferParticles::Specialization<20*sizd>::Add();
381 TransferParticles::Specialization<24*sizd>::Add();
382 TransferParticles::Specialization<28*sizd>::Add();
383 TransferParticles::Specialization<32*sizd>::Add();
384 TransferParticles::Specialization<36*sizd>::Add();
385 TransferParticles::Specialization<40*sizd>::Add();
386}
387
388auto ParticleSet::TransferParticles::Fallback(size_t bufsize)
389-> ParticleSet::TransferParticlesType
390{
391 constexpr size_t sizd = sizeof(real_t);
392 if (bufsize < 4*sizd)
393 {
394 return &ParticleSet::TransferParticlesImpl<4*sizd>;
395 }
396 else if (bufsize < 8*sizd)
397 {
398 return &ParticleSet::TransferParticlesImpl<8*sizd>;
399 }
400 else if (bufsize < 12*sizd)
401 {
402 return &ParticleSet::TransferParticlesImpl<12*sizd>;
403 }
404 else if (bufsize < 16*sizd)
405 {
406 return &ParticleSet::TransferParticlesImpl<16*sizd>;
407 }
408 else if (bufsize < 20*sizd)
409 {
410 return &ParticleSet::TransferParticlesImpl<20*sizd>;
411 }
412 else if (bufsize < 24*sizd)
413 {
414 return &ParticleSet::TransferParticlesImpl<24*sizd>;
415 }
416 else if (bufsize < 28*sizd)
417 {
418 return &ParticleSet::TransferParticlesImpl<28*sizd>;
419 }
420 else if (bufsize < 32*sizd)
421 {
422 return &ParticleSet::TransferParticlesImpl<32*sizd>;
423 }
424 else if (bufsize < 36*sizd)
425 {
426 return &ParticleSet::TransferParticlesImpl<36*sizd>;
427 }
428 else if (bufsize < 40*sizd)
429 {
430 return &ParticleSet::TransferParticlesImpl<40*sizd>;
431 }
432 return &ParticleSet::TransferParticlesImpl<60*sizd>;
433}
434/// \endcond DO_NOT_DOCUMENT
435
437{
438 MFEM_ASSERT(rank_list.Size() == GetNParticles(),
439 "rank_list must be of size GetNParticles().");
440
441 int rank = GetRank(comm);
442
443 // Get particles to be transferred
444 // (Avoid unnecessary copies of particle data into and out of buffers)
445 Array<int> send_idxs;
446 Array<unsigned int> send_ranks;
447 send_idxs.Reserve(rank_list.Size());
448 send_ranks.Reserve(rank_list.Size());
449 for (int i = 0; i < rank_list.Size(); i++)
450 {
451 if (rank != static_cast<int>(rank_list[i]))
452 {
453 send_idxs.Append(i);
454 send_ranks.Append(rank_list[i]);
455 }
456 }
457
458 // Compute number of bytes of a single particle
459 int nreals = GetFieldVDims().Sum() + coords.GetVDim();
460 int ntags = GetNTags();
461 size_t nbytes = nreals*sizeof(real_t) + ntags*sizeof(int);
462
463 // Dispatch to appropriate redistribution function for this size
464 TransferParticles::Run(nbytes, *this, send_idxs, send_ranks);
465}
466
467#endif // MFEM_USE_MPI && MFEM_USE_GSLIB
468
473
474void ParticleSet::WriteToFile(const char *fname,
475 const std::stringstream &ss_header, const std::stringstream &ss_data)
476{
477
478#ifdef MFEM_USE_MPI
479 // Parallel:
480 int rank = GetRank(comm);
481
482 MPI_File_delete(fname, MPI_INFO_NULL); // delete old file if it exists
483 MPI_File file;
484 int mpi_err = MPI_File_open(comm, fname, MPI_MODE_CREATE | MPI_MODE_WRONLY,
485 MPI_INFO_NULL, &file);
486 MFEM_VERIFY(mpi_err == MPI_SUCCESS, "MPI_File_open failed.");
487
488 // Print header
489 if (rank == 0)
490 {
491 MPI_File_write_at(file, 0, ss_header.str().data(), ss_header.str().size(),
492 MPI_CHAR, MPI_STATUS_IGNORE);
493 }
494
495 // Compute the data size in bytes
496 MPI_Offset data_size = ss_data.str().size();
497 MPI_Offset offset;
498
499 // Compute the offsets using an exclusive scan
500 MPI_Exscan(&data_size, &offset, 1, MPI_OFFSET, MPI_SUM, comm);
501 if (rank == 0)
502 {
503 offset = 0;
504 }
505
506 // Add offset from the header
507 offset += ss_header.str().size();
508
509 // Write data collectively
510 MPI_File_write_at_all(file, offset, ss_data.str().data(),
511 data_size, MPI_BYTE, MPI_STATUS_IGNORE);
512
513 // Close file
514 MPI_File_close(&file);
515#else
516 // Serial:
517 std::ofstream ofs(fname);
518 MFEM_VERIFY(ofs.is_open() && !ofs.fail(),
519 "Error: Could not open file " << fname << " for writing.");
520 ofs << ss_header.str() << ss_data.str();
521 ofs.close();
522#endif // MFEM_USE_MPI
523}
524
525ParticleSet::ParticleSet(int id_stride_, IDType id_counter_, int num_particles,
526 int dim, Ordering::Type coords_ordering, const Array<int> &field_vdims,
527 const Array<Ordering::Type> &field_orderings,
528 const Array<const char*> &field_names_, int num_tags,
529 const Array<const char*> &tag_names_)
530 : id_stride(id_stride_),
531 id_counter(id_counter_),
532 coords(dim, coords_ordering)
533{
534 // Initialize fields
535 for (int f = 0; f < field_vdims.Size(); f++)
536 {
537 AddField(field_vdims[f], field_orderings[f], field_names_[f]);
538 }
539
540 // Initialize tags
541 for (int t = 0; t < num_tags; t++)
542 {
543 AddTag(tag_names_[t]);
544 }
545
546 // Add num_particles
547 Array<IDType> init_ids(num_particles);
548 for (int i = 0; i < num_particles; i++)
549 {
550 init_ids[i] = id_counter;
552 }
553 AddParticles(init_ids);
554}
555
557{
558 if (p.GetDim() != GetDim())
559 {
560 return false;
561 }
562 if (p.GetNFields() != GetNFields())
563 {
564 return false;
565 }
566 for (int f = 0; f < GetNFields(); f++)
567 {
568 if (p.GetFieldVDim(f) != Field(f).GetVDim())
569 {
570 return false;
571 }
572 }
573 if (p.GetNTags() != GetNTags())
574 {
575 return false;
576 }
577
578 return true;
579
580}
581
582ParticleSet::ParticleSet(int num_particles, int dim,
583 Ordering::Type coords_ordering)
584 : ParticleSet(1, 0, num_particles, dim, coords_ordering, Array<int>(),
585 Array<Ordering::Type>(), Array<const char*>(), 0,
586 Array<const char*>())
587{
588
589}
590
591ParticleSet::ParticleSet(int num_particles, int dim,
592 const Array<int> &field_vdims, int num_tags,
593 Ordering::Type all_ordering)
594 : ParticleSet(1, 0, num_particles, dim, all_ordering, field_vdims,
595 GetOrderingArray(all_ordering, field_vdims.Size()),
596 GetEmptyNameArray(field_vdims.Size()), num_tags,
597 GetEmptyNameArray(num_tags))
598{
599}
600
601ParticleSet::ParticleSet(int num_particles, int dim,
602 const Array<int> &field_vdims, const Array<const
603 char*> &field_names_, int num_tags,
604 const Array<const char*> &tag_names_,
605 Ordering::Type all_ordering)
606 : ParticleSet(1, 0, num_particles, dim, all_ordering, field_vdims,
607 GetOrderingArray(all_ordering, field_vdims.Size()),
608 field_names_, num_tags,
609 tag_names_)
610{
611
612}
613
614ParticleSet::ParticleSet(int num_particles, int dim,
615 Ordering::Type coords_ordering,
616 const Array<int> &field_vdims,
617 const Array<Ordering::Type> &field_orderings,
618 const Array<const char*> &field_names_, int num_tags,
619 const Array<const char*> &tag_names_)
620 : ParticleSet(1, 0, num_particles, dim, coords_ordering, field_vdims,
621 field_orderings, field_names_, num_tags, tag_names_)
622{
623
624}
625
626
627
628#ifdef MFEM_USE_MPI
629ParticleSet::ParticleSet(MPI_Comm comm_, int rank_num_particles, int dim,
630 Ordering::Type coords_ordering)
631 : ParticleSet(comm_, rank_num_particles, dim, coords_ordering, Array<int>(),
632 Array<Ordering::Type>(), Array<const char*>(), 0,
633 Array<const char*>())
634{
635
636};
637
638ParticleSet::ParticleSet(MPI_Comm comm_, int rank_num_particles, int dim,
639 const Array<int> &field_vdims, int num_tags,
640 Ordering::Type all_ordering)
641 : ParticleSet(comm_, rank_num_particles, dim, all_ordering, field_vdims,
642 GetOrderingArray(all_ordering, field_vdims.Size()),
643 GetEmptyNameArray(field_vdims.Size()), num_tags,
644 GetEmptyNameArray(num_tags))
645{
646
647}
648
649ParticleSet::ParticleSet(MPI_Comm comm_, int rank_num_particles, int dim,
650 const Array<int> &field_vdims, const Array<const
651 char*> &field_names_,
652 int num_tags, const Array<const char*> &tag_names_,
653 Ordering::Type all_ordering)
654 : ParticleSet(comm_, rank_num_particles, dim, all_ordering, field_vdims,
655 GetOrderingArray(all_ordering, field_vdims.Size()),
656 field_names_, num_tags,
657 tag_names_)
658{
659
660}
661
662ParticleSet::ParticleSet(MPI_Comm comm_, int rank_num_particles, int dim,
663 Ordering::Type coords_ordering,
664 const Array<int> &field_vdims,
665 const Array<Ordering::Type> &field_orderings,
666 const Array<const char*> &field_names_, int num_tags,
667 const Array<const char*> &tag_names_)
668 : ParticleSet(GetSize(comm_), (IDType)GetRank(comm_),
669 rank_num_particles,
670 dim,
671 coords_ordering,
672 field_vdims,
673 field_orderings,
674 field_names_,
675 num_tags,
676 tag_names_)
677{
678 comm = comm_;
679#ifdef MFEM_USE_GSLIB
680 gsl_comm = new gslib::comm;
681 cr = new gslib::crystal;
682 comm_init(gsl_comm, comm);
683 crystal_init(cr, gsl_comm);
684#endif // MFEM_USE_GSLIB
685}
686#endif // MFEM_USE_MPI
687
689{
690 IDType total = (IDType)GetNParticles();
691#ifdef MFEM_USE_MPI
692 MPI_Allreduce(MPI_IN_PLACE, &total, 1, MPI_UNSIGNED_LONG_LONG,
693 MPI_SUM, comm);
694#endif // MFEM_USE_MPI
695 return total;
696}
697
698int ParticleSet::AddField(int vdim, Ordering::Type field_ordering,
699 const char* field_name)
700{
701 std::string field_name_str(field_name ? field_name : "");
702 if (!field_name)
703 {
704 field_name_str = GetDefaultFieldName(field_names.size());
705 }
706 fields.emplace_back(std::make_unique<ParticleVector>(vdim, field_ordering,
707 GetNParticles()));
708 field_names.emplace_back(field_name_str);
709
710 return GetNFields() - 1;
711}
712
713int ParticleSet::AddTag(const char* tag_name)
714{
715 std::string tag_name_str(tag_name ? tag_name : "");
716 if (!tag_name)
717 {
718 tag_name_str = GetDefaultTagName(tag_names.size());
719 }
720 tags.emplace_back(std::make_unique<Array<int>>(GetNParticles()));
721 tag_names.emplace_back(tag_name_str);
722
723 return GetNTags() - 1;
724}
725
727{
728 MFEM_ASSERT(IsValidParticle(p),
729 "Particle is incompatible with ParticleSet.");
730
731 // Add the particle
732 Array<int> idxs;
735
736 // Set the new particle data
737 int idx = idxs[0];
738 SetParticle(idx, p);
739}
740
741void ParticleSet::AddParticles(int num_particles, Array<int> *new_indices)
742{
743 Array<IDType> add_ids(num_particles);
744 for (int i = 0; i < num_particles; i++)
745 {
746 add_ids[i] = id_counter;
748 }
749
750 AddParticles(add_ids, new_indices);
751}
752
754{
755 // Delete IDs
756 ids.DeleteAt(list);
757
758 // Delete data
759 for (int f = -1; f < GetNFields(); f++)
760 {
761 ParticleVector &pv = (f == -1 ? coords : *fields[f]);
762 pv.DeleteParticles(list);
763 }
764
765 // Delete tags
766 for (int t = 0; t < GetNTags(); t++)
767 {
768 tags[t]->DeleteAt(list);
769 }
770}
771
773{
775
776 Coords().GetValues(i, p.Coords());
777
778 for (int f = 0; f < GetNFields(); f++)
779 {
780 Field(f).GetValues(i, p.Field(f));
781 }
782
783 for (int t = 0; t < GetNTags(); t++)
784 {
785 p.Tag(t) = Tag(t)[i];
786 }
787
788 return p;
789}
790
792{
794 {
795 return false;
796 }
797 for (int f = 0; f < GetNFields(); f++)
798 {
799 if (fields[f]->GetOrdering() == Ordering::byNODES)
800 {
801 return false;
802 }
803 }
804 return true;
805}
806
808{
810
811 Coords().GetValuesRef(i, p.Coords());
812
813 for (int f = 0; f < GetNFields(); f++)
814 {
815 MFEM_ASSERT(Field(f).GetOrdering() == Ordering::byVDIM,
816 "GetParticleRef only valid when all fields ordered byVDIM.");
817 p.SetFieldRef(f, Field(f).GetData() + i*Field(f).GetVDim());
818 }
819
820 for (int t = 0; t < GetNTags(); t++)
821 {
822 p.SetTagRef(t, &(*tags[t])[i]);
823 }
824
825 return p;
826}
827
829{
830 MFEM_ASSERT(IsValidParticle(p),
831 "Particle is incompatible with ParticleSet.");
832
833 Coords().SetValues(i, p.Coords());
834
835 for (int f = 0; f < GetNFields(); f++)
836 {
837 Field(f).SetValues(i, p.Field(f));
838 }
839
840 for (int t = 0; t < GetNTags(); t++)
841 {
842 Tag(t)[i] = p.Tag(t);
843 }
844}
845
846void ParticleSet::PrintCSV(const char *fname, int precision)
847{
848 Array<int> all_field_idxs(GetNFields()), all_tag_idxs(GetNTags());
849
850 for (int f = 0; f < GetNFields(); f++)
851 {
852 all_field_idxs[f] = f;
853 }
854
855 for (int t = 0; t < GetNTags(); t++)
856 {
857 all_tag_idxs[t] = t;
858 }
859
860 PrintCSV(fname, all_field_idxs, all_tag_idxs, precision);
861}
862
863void ParticleSet::PrintCSV(const char *fname, const Array<int> &field_idxs,
864 const Array<int> &tag_idxs, int precision)
865{
866 std::stringstream ss_header;
867
868 // Configure header:
869 ss_header << "id";
870
871#ifdef MFEM_USE_MPI
872 ss_header << ",rank";
873#endif // MFEM_USE_MPI
874
875 std::array<char, 3> ax = {'X', 'Y', 'Z'};
876 for (int c = 0; c < coords.GetVDim(); c++)
877 {
878 ss_header << "," << ax[c];
879 }
880
881 for (int f = 0; f < field_idxs.Size(); f++)
882 {
883 ParticleVector &pv = *fields[field_idxs[f]];
884 for (int c = 0; c < pv.GetVDim(); c++)
885 {
886 ss_header << "," << field_names[field_idxs[f]] <<
887 (pv.GetVDim() > 1 ? "_" + std::to_string(c) : "");
888 }
889 }
890
891 for (int t = 0; t < tag_idxs.Size(); t++)
892 {
893 ss_header << "," << tag_names[tag_idxs[t]];
894 }
895 ss_header << "\n";
896
897 // Configure data
898 std::stringstream ss_data;
899 ss_data.precision(precision);
900#ifdef MFEM_USE_MPI
901 int rank = GetRank(comm);
902#endif // MFEM_USE_MPI
903 for (int i = 0; i < GetNParticles(); i++)
904 {
905 ss_data << ids[i];
906#ifdef MFEM_USE_MPI
907 ss_data << "," << rank;
908#endif // MFEM_USE_MPI
909
910 for (int c = 0; c < coords.GetVDim(); c++)
911 {
912 ss_data << "," << coords(i, c);
913 }
914 for (int f = 0; f < field_idxs.Size(); f++)
915 {
916 ParticleVector &pv = *fields[field_idxs[f]];
917 for (int c = 0; c < pv.GetVDim(); c++)
918 {
919 ss_data << "," << pv(i, c);
920 }
921 }
922 for (int t = 0; t < tag_idxs.Size(); t++)
923 {
924 ss_data << "," << (*tags[tag_idxs[t]])[i];
925 }
926 ss_data << "\n";
927 }
928
929 // Write
930 WriteToFile(fname, ss_header, ss_data);
931}
932
934{
935#if defined(MFEM_USE_MPI) && defined(MFEM_USE_GSLIB)
936 if (gsl_comm)
937 {
938 if (!Mpi::IsFinalized()) // currently segfaults inside gslib otherwise
939 {
940 crystal_free(cr);
941 comm_free(gsl_comm);
942 delete gsl_comm;
943 delete cr;
944 }
945 }
946#endif // MFEM_USE_MPI && MFEM_USE_GSLIB
947}
948
949
950} // namespace mfem
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
Definition array.hpp:184
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:840
int Size() const
Return the logical size of the array.
Definition array.hpp:166
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:912
T * GetData()
Returns the data.
Definition array.hpp:140
void DeleteAt(const Array< int > &indices)
Delete entries at indices, and resize.
Definition array.hpp:1007
T Sum() const
Return the sum of all the array entries using the '+'' operator for class 'T'.
Definition array.cpp:129
static bool IsFinalized()
Return true if MPI has been finalized.
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition ordering.hpp:13
Type
Ordering methods:
Definition ordering.hpp:17
ParticleSet initializes and manages data associated with particles.
Array< IDType > ids
Global unique IDs of particles owned by this rank.
ParticleVector coords
Spatial coordinates of particles owned by this rank.
~ParticleSet()
Destructor.
unsigned long long IDType
void Redistribute(const Array< unsigned int > &rank_list)
Redistribute particle data to rank_list.
int GetNTags() const
Get the number of tags registered to particles.
void AddParticle(const Particle &p)
Add a particle using Particle .
bool IsValidParticle(const Particle &p) const
Check if a particle could belong in this ParticleSet by comparing field and tag dimension.
Particle GetParticleRef(int i)
Get Particle object whose members reference the actual data associated with particle i in this Partic...
std::vector< std::string > tag_names
Tag names, to be written when PrintCSV() is called.
std::vector< std::unique_ptr< ParticleVector > > fields
All particle fields for particles owned by this rank.
Particle GetParticle(int i) const
Get new Particle object with copy of data associated with particle i .
const Array< int > GetFieldVDims() const
Get an Array<int> of the field vector-dimensions registered to particles.
ParticleVector & Coords()
Get a reference to the coordinates ParticleVector.
const Array< IDType > & GetIDs() const
Get the global IDs of the active particles owned by this ParticleSet.
int GetNFields() const
Get the number of fields registered to particles.
int AddTag(const char *tag_name=nullptr)
Add a tag to the ParticleSet.
const int id_stride
Stride for IDs (used internally when new particles are added).
void UpdateID(int local_idx, IDType new_global_id)
Update global ID of a particle.
bool IsParticleRefValid() const
Determine if GetParticleRef is valid.
ParticleVector & Field(int f)
Get a reference to field f 's ParticleVector.
int GetDim() const
Get the spatial dimension.
int GetNParticles() const
Get the number of active particles currently held by this ParticleSet.
void WriteToFile(const char *fname, const std::stringstream &ss_header, const std::stringstream &ss_data)
Write string in ss_header , followed by ss_data , to a single file; compatible in parallel.
ParticleSet(int id_stride_, IDType id_counter_, int num_particles, int dim, Ordering::Type coords_ordering, const Array< int > &field_vdims, const Array< Ordering::Type > &field_orderings, const Array< const char * > &field_names_, int num_tags, const Array< const char * > &tag_names_)
Hidden main constructor of ParticleSet.
Array< int > & Tag(int t)
Get a reference to tag t 's Array<int>.
int AddField(int vdim, Ordering::Type field_ordering=Ordering::byVDIM, const char *field_name=nullptr)
Add a field to the ParticleSet.
struct gslib::crystal * cr
IDType GetGlobalNParticles() const
Get the global number of active particles across all ranks.
void PrintCSV(const char *fname, int precision=16)
Print all particle data to a comma-delimited CSV file.
void Reserve(int res)
Reserve memory for res particles.
std::vector< std::unique_ptr< Array< int > > > tags
All particle tags for particles owned by this rank.
IDType id_counter
Current globally unique ID to be assigned to the next particle added.
void RemoveParticles(const Array< int > &list)
Remove particle data specified by list of particle indices.
struct gslib::comm * gsl_comm
void AddParticles(const Array< IDType > &new_ids, Array< int > *new_indices=nullptr)
Add particles with global identifiers new_ids and optionally get the local indices of new particles i...
std::vector< std::string > field_names
Field names, to be written when PrintCSV() is called.
Particle CreateParticle() const
Create a Particle object with the same spatial dimension, number of fields and field vdims,...
void SetParticle(int i, const Particle &p)
Set data for particle at index i with data from provided particle p.
ParticleVector carries vector data (of a given vector dimension) for an arbitrary number of particles...
void SetValues(int i, const Vector &nvals)
Set particle i 's data to nvals .
void GetValues(int i, Vector &nvals) const
Get a copy of particle i 's data.
void DeleteParticles(const Array< int > &indices)
Remove particle data at indices.
void SetNumParticles(int num_vectors, bool keep_data=true)
Set the number of particle Vector data to be held by the ParticleVector, keeping existing data.
int GetVDim() const
Get the Vector dimension of the ParticleVector.
void GetValuesRef(int i, Vector &nref)
For GetOrdering == Ordering::byVDIM, set nref to refer to particle i 's data.
Ordering::Type GetOrdering() const
Get the ordering of data in the ParticleVector.
Container for data associated with a single particle.
std::vector< Array< int > > tags
A std::vector of Array<int> where each Array<int> holds data for a given tag.
bool operator==(const Particle &rhs) const
Particle equality operator.
std::vector< Vector > fields
A std::vector of Vector where each Vector holds data for a given field (e.g., mass,...
Vector coords
Spatial coordinates.
Particle(int dim, const Array< int > &field_vdims, int num_tags)
Construct a Particle instance.
void Print(std::ostream &os=mfem::out) const
Print all particle data to os.
void SetTagRef(int t, int *tag_data)
Set tag t to reference external data.
void SetFieldRef(int f, real_t *field_data)
Set field f to reference external data.
Vector data type.
Definition vector.hpp:82
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
void Reserve(int res)
Update Capacity() to res (if less than current), keeping existing entries.
Definition vector.hpp:633
int dim
Definition ex24.cpp:53
mfem::real_t real_t
float real_t
Definition config.hpp:46
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
real_t p(const Vector &x, real_t t)