Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
particle_accessor.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_particles_particle_accessor_h
17#define dealii_particles_particle_accessor_h
18
19#include <deal.II/base/config.h>
20
22
23#include <deal.II/grid/tria.h>
24
26
28
29namespace Particles
30{
31 // Forward declarations
32#ifndef DOXYGEN
33 template <int, int>
34 class ParticleIterator;
35 template <int, int>
36 class ParticleHandler;
37#endif
38
42 template <int dim, int spacedim = dim>
44 {
45 public:
49 void *
51
52
56 const void *
58
76 void
77 set_location(const Point<spacedim> &new_location);
78
84 const Point<spacedim> &
85 get_location() const;
86
104 void
105 set_reference_location(const Point<dim> &new_reference_location);
106
110 const Point<dim> &
112
117 get_id() const;
118
126 void
128
133 bool
135
153 void
154 set_properties(const std::vector<double> &new_properties);
155
173 void
175
183
191
197 std::size_t
199
209
214 template <class Archive>
215 void
216 serialize(Archive &ar, const unsigned int version);
217
221 void
223
227 void
229
233 bool
235
239 bool
241
242 private:
247
254 const std::multimap<internal::LevelInd, Particle<dim, spacedim>> &map,
255 const typename std::multimap<internal::LevelInd,
256 Particle<dim, spacedim>>::iterator
257 &particle);
258
259 private:
264 std::multimap<internal::LevelInd, Particle<dim, spacedim>> *map;
265
270 typename std::multimap<internal::LevelInd,
272
273 // Make ParticleIterator a friend to allow it constructing
274 // ParticleAccessors.
275 template <int, int>
276 friend class ParticleIterator;
277 template <int, int>
278 friend class ParticleHandler;
279 };
280
281
282
283 template <int dim, int spacedim>
284 template <class Archive>
285 void
287 const unsigned int version)
288 {
289 return particle->second.serialize(ar, version);
290 }
291
292
293 // ------------------------- inline functions ------------------------------
294
295 template <int dim, int spacedim>
297 : map(nullptr)
298 , particle()
299 {}
300
301
302
303 template <int dim, int spacedim>
305 const std::multimap<internal::LevelInd, Particle<dim, spacedim>> &map,
306 const typename std::multimap<internal::LevelInd,
307 Particle<dim, spacedim>>::iterator & particle)
308 : map(const_cast<
309 std::multimap<internal::LevelInd, Particle<dim, spacedim>> *>(&map))
310 , particle(particle)
311 {}
312
313
314
315 template <int dim, int spacedim>
316 inline const void *
318 const void *data)
319 {
320 Assert(particle != map->end(), ExcInternalError());
321
322 return particle->second.read_particle_data_from_memory(data);
323 }
324
325
326
327 template <int dim, int spacedim>
328 inline void *
330 void *data) const
331 {
332 Assert(particle != map->end(), ExcInternalError());
333
334 return particle->second.write_particle_data_to_memory(data);
335 }
336
337
338
339 template <int dim, int spacedim>
340 inline void
342 {
343 Assert(particle != map->end(), ExcInternalError());
344
345 particle->second.set_location(new_loc);
346 }
347
348
349
350 template <int dim, int spacedim>
351 inline const Point<spacedim> &
353 {
354 Assert(particle != map->end(), ExcInternalError());
355
356 return particle->second.get_location();
357 }
358
359
360
361 template <int dim, int spacedim>
362 inline void
364 const Point<dim> &new_loc)
365 {
366 Assert(particle != map->end(), ExcInternalError());
367
368 particle->second.set_reference_location(new_loc);
369 }
370
371
372
373 template <int dim, int spacedim>
374 inline const Point<dim> &
376 {
377 Assert(particle != map->end(), ExcInternalError());
378
379 return particle->second.get_reference_location();
380 }
381
382
383
384 template <int dim, int spacedim>
387 {
388 Assert(particle != map->end(), ExcInternalError());
389
390 return particle->second.get_id();
391 }
392
393
394
395 template <int dim, int spacedim>
396 inline void
398 PropertyPool<dim, spacedim> &new_property_pool)
399 {
400 Assert(particle != map->end(), ExcInternalError());
401
402 particle->second.set_property_pool(new_property_pool);
403 }
404
405
406
407 template <int dim, int spacedim>
408 inline bool
410 {
411 Assert(particle != map->end(), ExcInternalError());
412
413 return particle->second.has_properties();
414 }
415
416
417
418 template <int dim, int spacedim>
419 inline void
421 const std::vector<double> &new_properties)
422 {
423 Assert(particle != map->end(), ExcInternalError());
424
425 particle->second.set_properties(new_properties);
426 }
427
428
429
430 template <int dim, int spacedim>
431 inline void
433 const ArrayView<const double> &new_properties)
434 {
435 Assert(particle != map->end(), ExcInternalError());
436
437 particle->second.set_properties(new_properties);
438 }
439
440
441
442 template <int dim, int spacedim>
443 inline const ArrayView<const double>
445 {
446 Assert(particle != map->end(), ExcInternalError());
447
448 return particle->second.get_properties();
449 }
450
451
452
453 template <int dim, int spacedim>
457 {
458 Assert(particle != map->end(), ExcInternalError());
459
461 &triangulation, particle->first.first, particle->first.second);
462 return cell;
463 }
464
465
466
467 template <int dim, int spacedim>
468 inline const ArrayView<double>
470 {
471 Assert(particle != map->end(), ExcInternalError());
472
473 return particle->second.get_properties();
474 }
475
476
477
478 template <int dim, int spacedim>
479 inline std::size_t
481 {
482 Assert(particle != map->end(), ExcInternalError());
483
484 return particle->second.serialized_size_in_bytes();
485 }
486
487
488
489 template <int dim, int spacedim>
490 inline void
492 {
493 Assert(particle != map->end(), ExcInternalError());
494 ++particle;
495 }
496
497
498
499 template <int dim, int spacedim>
500 inline void
502 {
503 Assert(particle != map->begin(), ExcInternalError());
504 --particle;
505 }
506
507
508
509 template <int dim, int spacedim>
510 inline bool
513 {
514 return (map != other.map) || (particle != other.particle);
515 }
516
517
518
519 template <int dim, int spacedim>
520 inline bool
523 {
524 return (map == other.map) && (particle == other.particle);
525 }
526
527
528} // namespace Particles
529
531
532namespace boost
533{
534 namespace geometry
535 {
536 namespace index
537 {
538 // Forward declaration of bgi::indexable
539 template <class T>
540 struct indexable;
541
546 template <int dim, int spacedim>
547 struct indexable<::Particles::ParticleAccessor<dim, spacedim>>
548 {
553 using result_type = const ::Point<spacedim> &;
554
556 operator()(const ::Particles::ParticleAccessor<dim, spacedim>
557 &accessor) const
558 {
559 return accessor.get_location();
560 }
561 };
562 } // namespace index
563 } // namespace geometry
564} // namespace boost
565
566#endif
std::multimap< internal::LevelInd, Particle< dim, spacedim > > * map
const Point< spacedim > & get_location() const
types::particle_index get_id() const
const void * read_particle_data_from_memory(const void *data)
Triangulation< dim, spacedim >::cell_iterator get_surrounding_cell(const Triangulation< dim, spacedim > &triangulation) const
bool operator!=(const ParticleAccessor< dim, spacedim > &other) const
const ArrayView< double > get_properties()
void set_location(const Point< spacedim > &new_location)
void serialize(Archive &ar, const unsigned int version)
std::multimap< internal::LevelInd, Particle< dim, spacedim > >::iterator particle
ParticleAccessor(const std::multimap< internal::LevelInd, Particle< dim, spacedim > > &map, const typename std::multimap< internal::LevelInd, Particle< dim, spacedim > >::iterator &particle)
const ArrayView< const double > get_properties() const
void set_property_pool(PropertyPool< dim, spacedim > &property_pool)
void set_reference_location(const Point< dim > &new_reference_location)
const Point< dim > & get_reference_location() const
std::size_t serialized_size_in_bytes() const
void set_properties(const std::vector< double > &new_properties)
void * write_particle_data_to_memory(void *data) const
void set_properties(const ArrayView< const double > &new_properties)
bool operator==(const ParticleAccessor< dim, spacedim > &other) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcInternalError()
std::pair< int, int > LevelInd
Definition: particle.h:42
STL namespace.
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
result_type operator()(const ::Particles::ParticleAccessor< dim, spacedim > &accessor) const