-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmesh.hh
1697 lines (1374 loc) · 59.9 KB
/
mesh.hh
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
mesh.hh - This file is part of MUSIC -
a code to generate multi-scale initial conditions
for cosmological simulations
Copyright (C) 2010 Oliver Hahn
*/
#ifndef __MESH_HH
#define __MESH_HH
#include <iomanip>
#include <iostream>
#include <stdexcept>
#include <vector>
#include <math.h>
#include "config_file.hh"
#include "log.hh"
#include "region_generator.hh"
class refinement_mask {
protected:
std::vector<short> mask_;
size_t nx_, ny_, nz_;
public:
refinement_mask(void) : nx_(0), ny_(0), nz_(0) {}
refinement_mask(size_t nx, size_t ny, size_t nz, short value = 0.) : nx_(nx), ny_(ny), nz_(nz) {
mask_.assign(nx_ * ny_ * nz_, value);
}
refinement_mask(const refinement_mask &r) {
nx_ = r.nx_;
ny_ = r.ny_;
nz_ = r.nz_;
mask_ = r.mask_;
}
refinement_mask &operator=(const refinement_mask &r) {
nx_ = r.nx_;
ny_ = r.ny_;
nz_ = r.nz_;
mask_ = r.mask_;
return *this;
}
void init(size_t nx, size_t ny, size_t nz, short value = 0.) {
nx_ = nx;
ny_ = ny;
nz_ = nz;
mask_.assign(nx_ * ny_ * nz_, value);
}
const short &operator()(size_t i, size_t j, size_t k) const { return mask_[(i * ny_ + j) * nz_ + k]; }
short &operator()(size_t i, size_t j, size_t k) { return mask_[(i * ny_ + j) * nz_ + k]; }
size_t count_flagged(void) {
size_t count = 0;
for (size_t i = 0; i < mask_.size(); ++i)
if (mask_[i])
++count;
return count;
}
size_t count_notflagged(void) {
size_t count = 0;
for (size_t i = 0; i < mask_.size(); ++i)
if (!mask_[i])
++count;
return count;
}
};
//! base class for all things that have rectangular mesh structure
template <typename T> class Meshvar {
public:
typedef T real_t;
size_t m_nx, //!< x-extent of the rectangular mesh
m_ny, //!< y-extent of the rectangular mesh
m_nz; //!< z-extent of the rectangular mesh
int m_offx, //!< x-offset of the grid (just as a helper, not used inside the class)
m_offy, //!< y-offset of the grid (just as a helper, not used inside the class)
m_offz; //!< z-offset of the grid (just as a helper, not used inside the class)
real_t *m_pdata; //!< pointer to the dynamic data array
//! constructor for cubic mesh
explicit Meshvar(size_t n, int offx, int offy, int offz)
: m_nx(n), m_ny(n), m_nz(n), m_offx(offx), m_offy(offy), m_offz(offz) {
m_pdata = new real_t[m_nx * m_ny * m_nz];
}
//! constructor for rectangular mesh
Meshvar(size_t nx, size_t ny, size_t nz, int offx, int offy, int offz)
: m_nx(nx), m_ny(ny), m_nz(nz), m_offx(offx), m_offy(offy), m_offz(offz) {
m_pdata = new real_t[m_nx * m_ny * m_nz];
}
//! variant copy constructor with optional copying of the actual data
Meshvar(const Meshvar<real_t> &m, bool copy_over = true) {
m_nx = m.m_nx;
m_ny = m.m_ny;
m_nz = m.m_nz;
m_offx = m.m_offx;
m_offy = m.m_offy;
m_offz = m.m_offz;
m_pdata = new real_t[m_nx * m_ny * m_nz];
if (copy_over)
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] = m.m_pdata[i];
}
//! standard copy constructor
explicit Meshvar(const Meshvar<real_t> &m) {
m_nx = m.m_nx;
m_ny = m.m_ny;
m_nz = m.m_nz;
m_offx = m.m_offx;
m_offy = m.m_offy;
m_offz = m.m_offz;
m_pdata = new real_t[m_nx * m_ny * m_nz];
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] = m.m_pdata[i];
}
//! destructor
~Meshvar() {
if (m_pdata != NULL)
delete[] m_pdata;
}
//! deallocate the data, but keep the structure
inline void deallocate(void) {
if (m_pdata != NULL)
delete[] m_pdata;
m_pdata = NULL;
}
//! get extent of the mesh along a specified dimension (const)
inline size_t size(unsigned dim) const {
if (dim == 0)
return m_nx;
if (dim == 1)
return m_ny;
return m_nz;
}
//! get offset of the mesh along a specified dimension (const)
inline int offset(unsigned dim) const {
if (dim == 0)
return m_offx;
if (dim == 1)
return m_offy;
return m_offz;
}
//! get extent of the mesh along a specified dimension
inline int &offset(unsigned dim) {
if (dim == 0)
return m_offx;
if (dim == 1)
return m_offy;
return m_offz;
}
//! set all the data to zero values
void zero(void) {
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] = 0.0;
}
//! direct array random acces to the data block
inline real_t *operator[](const size_t i) { return &m_pdata[i]; }
//! direct array random acces to the data block (const)
inline const real_t *operator[](const size_t i) const { return &m_pdata[i]; }
//! 3D random access to the data block via index 3-tuples
inline real_t &operator()(const int ix, const int iy, const int iz) {
#ifdef DEBUG
if (ix < 0 || ix >= (int)m_nx || iy < 0 || iy >= (int)m_ny || iz < 0 || iz >= (int)m_nz)
LOGERR("Array index (%d,%d,%d) out of bounds", ix, iy, iz);
#endif
return m_pdata[((size_t)ix * m_ny + (size_t)iy) * m_nz + (size_t)iz];
}
//! 3D random access to the data block via index 3-tuples (const)
inline const real_t &operator()(const int ix, const int iy, const int iz) const {
#ifdef DEBUG
if (ix < 0 || ix >= (int)m_nx || iy < 0 || iy >= (int)m_ny || iz < 0 || iz >= (int)m_nz)
LOGERR("Array index (%d,%d,%d) out of bounds", ix, iy, iz);
#endif
return m_pdata[((size_t)ix * m_ny + (size_t)iy) * m_nz + (size_t)iz];
}
//! direct multiplication of the whole data block with a number
Meshvar<real_t> &operator*=(real_t x) {
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] *= x;
return *this;
}
//! direct addition of a number to the whole data block
Meshvar<real_t> &operator+=(real_t x) {
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] += x;
return *this;
}
//! direct element-wise division of the whole data block by a number
Meshvar<real_t> &operator/=(real_t x) {
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] /= x;
return *this;
}
//! direct subtraction of a number from the whole data block
Meshvar<real_t> &operator-=(real_t x) {
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] -= x;
return *this;
}
//! direct element-wise multiplication with another compatible mesh
Meshvar<real_t> &operator*=(const Meshvar<real_t> &v) {
if (v.m_nx * v.m_ny * v.m_nz != m_nx * m_ny * m_nz) {
LOGERR("Meshvar::operator*= : attempt to operate on incompatible data");
throw std::runtime_error("Meshvar::operator*= : attempt to operate on incompatible data");
}
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] *= v.m_pdata[i];
return *this;
}
//! direct element-wise division with another compatible mesh
Meshvar<real_t> &operator/=(const Meshvar<real_t> &v) {
if (v.m_nx * v.m_ny * v.m_nz != m_nx * m_ny * m_nz) {
LOGERR("Meshvar::operator/= : attempt to operate on incompatible data");
throw std::runtime_error("Meshvar::operator/= : attempt to operate on incompatible data");
}
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] /= v.m_pdata[i];
return *this;
}
//! direct element-wise addition of another compatible mesh
Meshvar<real_t> &operator+=(const Meshvar<real_t> &v) {
if (v.m_nx * v.m_ny * v.m_nz != m_nx * m_ny * m_nz) {
LOGERR("Meshvar::operator+= : attempt to operate on incompatible data");
throw std::runtime_error("Meshvar::operator+= : attempt to operate on incompatible data");
}
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] += v.m_pdata[i];
return *this;
}
//! direct element-wise subtraction of another compatible mesh
Meshvar<real_t> &operator-=(const Meshvar<real_t> &v) {
if (v.m_nx * v.m_ny * v.m_nz != m_nx * m_ny * m_nz) {
LOGERR("Meshvar::operator-= : attempt to operate on incompatible data");
throw std::runtime_error("Meshvar::operator-= : attempt to operate on incompatible data");
}
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] -= v.m_pdata[i];
return *this;
}
//! assignment operator for rectangular meshes
Meshvar<real_t> &operator=(const Meshvar<real_t> &m) {
m_nx = m.m_nx;
m_ny = m.m_ny;
m_nz = m.m_nz;
m_offx = m.m_offx;
m_offy = m.m_offy;
m_offz = m.m_offz;
if (m_pdata != NULL)
delete m_pdata;
m_pdata = new real_t[m_nx * m_ny * m_nz];
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
m_pdata[i] = m.m_pdata[i];
return *this;
}
real_t *get_ptr(void) { return m_pdata; }
};
//! MeshvarBnd derived class adding boundary ghost cell functionality
template <typename T> class MeshvarBnd : public Meshvar<T> {
using Meshvar<T>::m_nx;
using Meshvar<T>::m_ny;
using Meshvar<T>::m_nz;
using Meshvar<T>::m_pdata;
public:
typedef T real_t;
//! number of boundary (ghost) cells
int m_nbnd;
//! most general constructor
MeshvarBnd(int nbnd, size_t nx, size_t ny, size_t nz, size_t xoff, size_t yoff, size_t zoff)
: Meshvar<real_t>(nx + 2 * nbnd, ny + 2 * nbnd, nz + 2 * nbnd, xoff, yoff, zoff), m_nbnd(nbnd) {}
//! zero-offset constructor
MeshvarBnd(size_t nbnd, size_t nx, size_t ny, size_t nz)
: Meshvar<real_t>(nx + 2 * nbnd, ny + 2 * nbnd, nz + 2 * nbnd, 0, 0, 0), m_nbnd(nbnd) {}
//! constructor for cubic meshes
MeshvarBnd(size_t nbnd, size_t n, size_t xoff, size_t yoff, size_t zoff)
: Meshvar<real_t>(n + 2 * nbnd, xoff, yoff, zoff), m_nbnd(nbnd) {}
//! constructor for cubic meshes with zero offset
MeshvarBnd(size_t nbnd, size_t n) : Meshvar<real_t>(n + 2 * nbnd, 0, 0, 0), m_nbnd(nbnd) {}
//! modified copy constructor, allows to avoid copying actual data
MeshvarBnd(const MeshvarBnd<real_t> &v, bool copyover) : Meshvar<real_t>(v, copyover), m_nbnd(v.m_nbnd) {}
//! copy constructor
explicit MeshvarBnd(const MeshvarBnd<real_t> &v) : Meshvar<real_t>(v, true), m_nbnd(v.m_nbnd) {}
//! get extent of the mesh along a specified dimension
inline size_t size(unsigned dim = 0) const {
if (dim == 0)
return m_nx - 2 * m_nbnd;
if (dim == 1)
return m_ny - 2 * m_nbnd;
return m_nz - 2 * m_nbnd;
}
//! 3D random access to the data block via index 3-tuples
inline real_t &operator()(const int ix, const int iy, const int iz) {
size_t iix(ix + m_nbnd), iiy(iy + m_nbnd), iiz(iz + m_nbnd);
return m_pdata[(iix * m_ny + iiy) * m_nz + iiz];
}
//! 3D random access to the data block via index 3-tuples (const)
inline const real_t &operator()(const int ix, const int iy, const int iz) const {
size_t iix(ix + m_nbnd), iiy(iy + m_nbnd), iiz(iz + m_nbnd);
return m_pdata[(iix * m_ny + iiy) * m_nz + iiz];
}
//! assignment operator for rectangular meshes with ghost zones
MeshvarBnd<real_t> &operator=(const MeshvarBnd<real_t> &m) {
if (this->m_nx != m.m_nx || this->m_ny != m.m_ny || this->m_nz != m.m_nz) {
this->m_nx = m.m_nx;
this->m_ny = m.m_ny;
this->m_nz = m.m_nz;
if (m_pdata != NULL)
delete[] m_pdata;
m_pdata = new real_t[m_nx * m_ny * m_nz];
}
for (size_t i = 0; i < m_nx * m_ny * m_nz; ++i)
this->m_pdata[i] = m.m_pdata[i];
return *this;
}
//! sets the value of all ghost zones to zero
void zero_bnd(void) {
int nx, ny, nz;
nx = this->size(0);
ny = this->size(1);
nz = this->size(2);
for (int j = -m_nbnd; j < ny + m_nbnd; ++j)
for (int k = -m_nbnd; k < nz + m_nbnd; ++k) {
for (int i = -m_nbnd; i < 0; ++i) {
(*this)(i, j, k) = 0.0;
(*this)(nx - 1 - i, j, k) = 0.0;
}
}
for (int i = -m_nbnd; i < nx + m_nbnd; ++i)
for (int k = -m_nbnd; k < nz + m_nbnd; ++k) {
for (int j = -m_nbnd; j < 0; ++j) {
(*this)(i, j, k) = 0.0;
(*this)(i, ny - j - 1, k) = 0.0;
}
}
for (int i = -m_nbnd; i < nx + m_nbnd; ++i)
for (int j = -m_nbnd; j < ny + m_nbnd; ++j) {
for (int k = -m_nbnd; k < 0; ++k) {
(*this)(i, j, k) = 0.0;
(*this)(i, j, nz - k - 1) = 0.0;
}
}
}
//! outputs the data, for debugging only, not practical for large datasets
void print(void) const {
int nbnd = m_nbnd;
std::cout << "size is [" << this->size(0) << ", " << this->size(1) << ", " << this->size(2) << "]\n";
std::cout << "ghost region has length of " << nbnd << std::endl;
std::cout.precision(3);
for (int i = -nbnd; i < (int)this->size(0) + nbnd; ++i) {
std::cout << "ix = " << i << ": \n";
for (int j = -nbnd; j < (int)this->size(1) + nbnd; ++j) {
for (int k = -nbnd; k < (int)this->size(2) + nbnd; ++k) {
if (i < 0 || i >= this->size(0) || j < 0 || j >= this->size(1) || k < 0 || k >= this->size(2))
std::cout << "[" << std::setw(6) << (*this)(i, j, k) << "] ";
else
std::cout << std::setw(8) << (*this)(i, j, k) << " ";
}
std::cout << std::endl;
}
std::cout << std::endl;
}
}
};
//! class that subsumes a nested grid collection
template <typename T> class GridHierarchy {
public:
//! number of ghost cells on boundary
size_t m_nbnd;
//! highest level without adaptive refinement
unsigned m_levelmin;
//! vector of pointers to the underlying rectangular mesh data for each level
std::vector<MeshvarBnd<T> *> m_pgrids;
std::vector<int> m_xoffabs, //!< vector of x-offsets of a level mesh relative to the coarser level
m_yoffabs, //!< vector of x-offsets of a level mesh relative to the coarser level
m_zoffabs; //!< vector of x-offsets of a level mesh relative to the coarser level
std::vector<refinement_mask *> m_ref_masks;
bool bhave_refmask;
protected:
//! check whether a given grid has identical hierarchy, dimensions to this
bool is_consistent(const GridHierarchy<T> &gh) {
if (gh.levelmax() != levelmax())
return false;
if (gh.levelmin() != levelmin())
return false;
for (unsigned i = levelmin(); i <= levelmax(); ++i)
for (int j = 0; j < 3; ++j) {
if (size(i, j) != gh.size(i, j))
return false;
if (offset(i, j) != gh.offset(i, j))
return false;
}
return true;
}
public:
//! return a pointer to the MeshvarBnd object representing data for one level
MeshvarBnd<T> *get_grid(unsigned ilevel) {
if (ilevel >= m_pgrids.size()) {
LOGERR("Attempt to access level %d but maxlevel = %d", ilevel, m_pgrids.size() - 1);
throw std::runtime_error("Fatal: attempt to access non-existent grid");
}
return m_pgrids[ilevel];
}
//! return a pointer to the MeshvarBnd object representing data for one level (const)
const MeshvarBnd<T> *get_grid(unsigned ilevel) const {
if (ilevel >= m_pgrids.size()) {
LOGERR("Attempt to access level %d but maxlevel = %d", ilevel, m_pgrids.size() - 1);
throw std::runtime_error("Fatal: attempt to access non-existent grid");
}
return m_pgrids[ilevel];
}
//! constructor for a collection of rectangular grids representing a multi-level hierarchy
/*! creates an empty hierarchy, levelmin is initially zero, no grids are stored
* @param nbnd number of ghost zones added at the boundary
*/
explicit GridHierarchy(size_t nbnd) : m_nbnd(nbnd), m_levelmin(0), bhave_refmask(false) { m_pgrids.clear(); }
//! copy constructor
explicit GridHierarchy(const GridHierarchy<T> &gh) {
for (unsigned i = 0; i <= gh.levelmax(); ++i)
m_pgrids.push_back(new MeshvarBnd<T>(*gh.get_grid(i)));
m_nbnd = gh.m_nbnd;
m_levelmin = gh.m_levelmin;
m_xoffabs = gh.m_xoffabs;
m_yoffabs = gh.m_yoffabs;
m_zoffabs = gh.m_zoffabs;
// ref_mask = gh.ref_mask;
bhave_refmask = gh.bhave_refmask;
if (bhave_refmask) {
for (size_t i = 0; i < gh.m_ref_masks.size(); ++i)
m_ref_masks.push_back(new refinement_mask(*(gh.m_ref_masks[i])));
}
}
//! destructor
~GridHierarchy() { this->deallocate(); }
//! free all memory occupied by the grid hierarchy
void deallocate() {
for (unsigned i = 0; i < m_pgrids.size(); ++i)
delete m_pgrids[i];
m_pgrids.clear();
std::vector<MeshvarBnd<T> *>().swap(m_pgrids);
m_xoffabs.clear();
m_yoffabs.clear();
m_zoffabs.clear();
m_levelmin = 0;
for (size_t i = 0; i < m_ref_masks.size(); ++i)
delete m_ref_masks[i];
m_ref_masks.clear();
}
// meaning of the mask:
// -1 = outside of mask
// 0.5 = in mask and refined (i.e. cell exists also on finer level)
// 1 = in mask and not refined (i.e. cell exists only on this level)
void add_refinement_mask(const double *shift) {
bhave_refmask = false;
//! generate a mask
if (m_levelmin != levelmax()) {
for (int ilevel = (int)levelmax(); ilevel >= (int)levelmin(); --ilevel) {
double xq[3], dx = 1.0 / (1ul << ilevel);
m_ref_masks[ilevel]->init(size(ilevel, 0), size(ilevel, 1), size(ilevel, 2), 0);
for (size_t i = 0; i < size(ilevel, 0); i += 2) {
xq[0] = (offset_abs(ilevel, 0) + i) * dx + 0.5 * dx + shift[0];
for (size_t j = 0; j < size(ilevel, 1); j += 2) {
xq[1] = (offset_abs(ilevel, 1) + j) * dx + 0.5 * dx + shift[1];
for (size_t k = 0; k < size(ilevel, 2); k += 2) {
xq[2] = (offset_abs(ilevel, 2) + k) * dx + 0.5 * dx + shift[2];
short mask_val = -1; // outside mask
if (the_region_generator->query_point(xq, ilevel) || ilevel == (int)levelmin())
mask_val = 1; // inside mask
(*m_ref_masks[ilevel])(i + 0, j + 0, k + 0) = mask_val;
(*m_ref_masks[ilevel])(i + 0, j + 0, k + 1) = mask_val;
(*m_ref_masks[ilevel])(i + 0, j + 1, k + 0) = mask_val;
(*m_ref_masks[ilevel])(i + 0, j + 1, k + 1) = mask_val;
(*m_ref_masks[ilevel])(i + 1, j + 0, k + 0) = mask_val;
(*m_ref_masks[ilevel])(i + 1, j + 0, k + 1) = mask_val;
(*m_ref_masks[ilevel])(i + 1, j + 1, k + 0) = mask_val;
(*m_ref_masks[ilevel])(i + 1, j + 1, k + 1) = mask_val;
}
}
}
}
bhave_refmask = true;
for (int ilevel = (int)levelmin(); ilevel < (int)levelmax(); ++ilevel) {
for (size_t i = 0; i < size(ilevel, 0); i++)
for (size_t j = 0; j < size(ilevel, 1); j++)
for (size_t k = 0; k < size(ilevel, 2); k++) {
bool fine_is_flagged = false;
int ifine[] = {
2 * (int)i - 2 * (int)offset(ilevel + 1, 0), 2 * (int)j - 2 * (int)offset(ilevel + 1, 1),
2 * (int)k - 2 * (int)offset(ilevel + 1, 2),
};
if (ifine[0] >= 0 && ifine[0] < (int)size(ilevel + 1, 0) && ifine[1] >= 0 &&
ifine[1] < (int)size(ilevel + 1, 1) && ifine[2] >= 0 && ifine[2] < (int)size(ilevel + 1, 2)) {
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 0, ifine[2] + 0) > 0;
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 0, ifine[2] + 1) > 0;
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 1, ifine[2] + 0) > 0;
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 1, ifine[2] + 1) > 0;
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 0, ifine[2] + 0) > 0;
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 0, ifine[2] + 1) > 0;
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 1, ifine[2] + 0) > 0;
fine_is_flagged |= (*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 1, ifine[2] + 1) > 0;
if (fine_is_flagged) {
(*m_ref_masks[ilevel])(i, j, k) = 2; // cell is refined
(*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 0, ifine[2] + 0) = 1;
(*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 0, ifine[2] + 1) = 1;
(*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 1, ifine[2] + 0) = 1;
(*m_ref_masks[ilevel + 1])(ifine[0] + 0, ifine[1] + 1, ifine[2] + 1) = 1;
(*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 0, ifine[2] + 0) = 1;
(*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 0, ifine[2] + 1) = 1;
(*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 1, ifine[2] + 0) = 1;
(*m_ref_masks[ilevel + 1])(ifine[0] + 1, ifine[1] + 1, ifine[2] + 1) = 1;
}
}
}
}
}
}
//! get offset of a grid at specified refinement level
/*! the offset describes the shift of a refinement grid with respect to its coarser parent grid
* @param ilevel the level for which the offset is to be determined
* @param idim the dimension along which the offset is to be determined
* @return integer value denoting the offset in units of coarse grid cells
* @sa offset_abs
*/
int offset(int ilevel, int idim) const { return m_pgrids[ilevel]->offset(idim); }
//! get size of a grid at specified refinement level
/*! the size describes the number of cells along one dimension of a grid
* @param ilevel the level for which the size is to be determined
* @param idim the dimension along which the size is to be determined
* @return integer value denoting the size of refinement grid at level ilevel along dimension idim
*/
size_t size(int ilevel, int idim) const { return m_pgrids[ilevel]->size(idim); }
//! get the absolute offset of a grid at specified refinement level
/*! the absolute offset describes the shift of a refinement grid with respect to the simulation volume
* @param ilevel the level for which the offset is to be determined
* @param idim the dimension along which the offset is to be determined
* @return integer value denoting the absolute offset in units of fine grid cells
* @sa offset
*/
int offset_abs(int ilevel, int idim) const {
if (idim == 0)
return m_xoffabs[ilevel];
if (idim == 1)
return m_yoffabs[ilevel];
return m_zoffabs[ilevel];
}
//! get the coordinate posisition of a grid cell
/*! returns the position of a grid cell at specified level relative to the simulation volume
* @param ilevel the refinement level of the grid cell
* @param i the x-index of the cell in the level grid
* @param j the y-index of the cell in the level grid
* @param k the z-index of the cell in the level grid
* @param ppos pointer to a double[3] array to which the coordinates are written
* @return none
*/
void cell_pos(unsigned ilevel, int i, int j, int k, double *ppos) const {
double h = 1.0 / (1 << ilevel); //, htop = h*2.0;
ppos[0] = h * ((double)offset_abs(ilevel, 0) + (double)i + 0.5);
ppos[1] = h * ((double)offset_abs(ilevel, 1) + (double)j + 0.5);
ppos[2] = h * ((double)offset_abs(ilevel, 2) + (double)k + 0.5);
if (ppos[0] >= 1.0 || ppos[1] >= 1.0 || ppos[2] >= 1.0)
std::cerr << " - Cell seems outside domain! : (" << ppos[0] << ", " << ppos[1] << ", " << ppos[2] << "\n";
}
//! get the bounding box of a grid in code units
/*! returns the bounding box of a grid at specified level in code units
* @param ilevel the refinement level of the grid
* @param left pointer to a double[3] array to which the left corner is written
* @param right pointer to a double[3] array to which the right corner is written
* @return none
*/
void grid_bbox(unsigned ilevel, double *left, double *right) const {
double h = 1.0 / (1 << ilevel); //, htop = h*2.0;
left[0] = h * ((double)offset_abs(ilevel, 0));
left[1] = h * ((double)offset_abs(ilevel, 1));
left[2] = h * ((double)offset_abs(ilevel, 2));
right[0] = left[0] + h * ((double)size(ilevel, 0));
right[1] = left[1] + h * ((double)size(ilevel, 1));
right[2] = left[2] + h * ((double)size(ilevel, 2));
}
//! checks whether a given grid cell is refined
/*! a grid cell counts as refined if it is divided into 8 cells at the next higher level
* @param ilevel the refinement level of the grid cell
* @param i the x-index of the cell in the level grid
* @param j the y-index of the cell in the level grid
* @param k the z-index of the cell in the level grid
* @return true if cell is refined, false otherwise
*/
bool is_refined(unsigned ilevel, int i, int j, int k) const {
// meaning of the mask:
// -1 = outside of mask
// 2 = in mask and refined (i.e. cell exists also on finer level)
// 1 = in mask and not refined (i.e. cell exists only on this level)
if (bhave_refmask) {
return (*m_ref_masks[ilevel])(i, j, k) == 2;
}
if (!bhave_refmask && ilevel == levelmax())
return false;
if (i < offset(ilevel + 1, 0) || i >= offset(ilevel + 1, 0) + (int)size(ilevel + 1, 0) / 2 ||
j < offset(ilevel + 1, 1) || j >= offset(ilevel + 1, 1) + (int)size(ilevel + 1, 1) / 2 ||
k < offset(ilevel + 1, 2) || k >= offset(ilevel + 1, 2) + (int)size(ilevel + 1, 2) / 2)
return false;
return true;
}
bool is_in_mask(unsigned ilevel, int i, int j, int k) const {
// meaning of the mask:
// -1 = outside of mask
// 2 = in mask and refined (i.e. cell exists also on finer level)
// 1 = in mask and not refined (i.e. cell exists only on this level)
if (bhave_refmask) {
return ((*m_ref_masks[ilevel])(i, j, k) >= 0);
}
return true;
}
//! sets the values of all grids on all levels to zero
void zero(void) {
for (unsigned i = 0; i < m_pgrids.size(); ++i)
m_pgrids[i]->zero();
}
//! count the number of cells that are not further refined (=leafs)
/*! for allocation purposes it is useful to query the number of cells to be expected
* @param lmin the minimum refinement level to consider
* @param lmax the maximum refinement level to consider
* @return the integer number of cells between lmin and lmax that are not further refined
*/
size_t count_leaf_cells(unsigned lmin, unsigned lmax) const {
size_t npcount = 0;
for (int ilevel = lmax; ilevel >= (int)lmin; --ilevel)
for (unsigned i = 0; i < get_grid(ilevel)->size(0); ++i)
for (unsigned j = 0; j < get_grid(ilevel)->size(1); ++j)
for (unsigned k = 0; k < get_grid(ilevel)->size(2); ++k)
if (is_in_mask(ilevel, i, j, k) && !is_refined(ilevel, i, j, k))
++npcount;
return npcount;
}
//! count the number of cells that are not further refined (=leafs)
/*! for allocation purposes it is useful to query the number of cells to be expected
* @return the integer number of cells in the hierarchy that are not further refined
*/
size_t count_leaf_cells(void) const { return count_leaf_cells(levelmin(), levelmax()); }
//! creates a hierarchy of coextensive grids, refined by factors of 2
/*! creates a hierarchy of lmax grids, each extending over the whole simulation volume with
* grid length 2^n for level 0<=n<=lmax
* @param lmax the maximum refinement level to be added (sets the resolution to 2^lmax for each dim)
* @return none
*/
void create_base_hierarchy(unsigned lmax) {
size_t n = 1;
this->deallocate();
m_pgrids.clear();
m_xoffabs.clear();
m_yoffabs.clear();
m_zoffabs.clear();
for (unsigned i = 0; i <= lmax; ++i) {
// std::cout << "....adding level " << i << " (" << n << ", " << n << ", " << n << ")" << std::endl;
m_pgrids.push_back(new MeshvarBnd<T>(m_nbnd, n, n, n, 0, 0, 0));
m_pgrids[i]->zero();
n *= 2;
m_xoffabs.push_back(0);
m_yoffabs.push_back(0);
m_zoffabs.push_back(0);
}
m_levelmin = lmax;
for (unsigned i = 0; i <= lmax; ++i)
m_ref_masks.push_back(new refinement_mask(size(i, 0), size(i, 1), size(i, 2), (short)(i != lmax)));
}
//! multiply entire grid hierarchy by a constant
GridHierarchy<T> &operator*=(T x) {
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) *= x;
return *this;
}
//! divide entire grid hierarchy by a constant
GridHierarchy<T> &operator/=(T x) {
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) /= x;
return *this;
}
//! add a constant to the entire grid hierarchy
GridHierarchy<T> &operator+=(T x) {
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) += x;
return *this;
}
//! subtract a constant from the entire grid hierarchy
GridHierarchy<T> &operator-=(T x) {
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) -= x;
return *this;
}
//! multiply (element-wise) two grid hierarchies
GridHierarchy<T> &operator*=(const GridHierarchy<T> &gh) {
if (!is_consistent(gh)) {
LOGERR("GridHierarchy::operator*= : attempt to operate on incompatible data");
throw std::runtime_error("GridHierarchy::operator*= : attempt to operate on incompatible data");
}
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) *= *gh.get_grid(i);
return *this;
}
//! divide (element-wise) two grid hierarchies
GridHierarchy<T> &operator/=(const GridHierarchy<T> &gh) {
if (!is_consistent(gh)) {
LOGERR("GridHierarchy::operator/= : attempt to operate on incompatible data");
throw std::runtime_error("GridHierarchy::operator/= : attempt to operate on incompatible data");
}
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) /= *gh.get_grid(i);
return *this;
}
//! add (element-wise) two grid hierarchies
GridHierarchy<T> &operator+=(const GridHierarchy<T> &gh) {
if (!is_consistent(gh))
throw std::runtime_error("GridHierarchy::operator+= : attempt to operate on incompatible data");
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) += *gh.get_grid(i);
return *this;
}
//! subtract (element-wise) two grid hierarchies
GridHierarchy<T> &operator-=(const GridHierarchy<T> &gh) {
if (!is_consistent(gh)) {
LOGERR("GridHierarchy::operator-= : attempt to operate on incompatible data");
throw std::runtime_error("GridHierarchy::operator-= : attempt to operate on incompatible data");
}
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) -= *gh.get_grid(i);
return *this;
}
//! assign (element-wise) two grid hierarchies
GridHierarchy<T> &operator=(const GridHierarchy<T> &gh) {
bhave_refmask = gh.bhave_refmask;
if (bhave_refmask) {
for (unsigned i = 0; i <= gh.levelmax(); ++i)
m_ref_masks.push_back(new refinement_mask(*(gh.m_ref_masks[i])));
}
if (!is_consistent(gh)) {
for (unsigned i = 0; i < m_pgrids.size(); ++i)
delete m_pgrids[i];
m_pgrids.clear();
for (unsigned i = 0; i <= gh.levelmax(); ++i)
m_pgrids.push_back(new MeshvarBnd<T>(*gh.get_grid(i)));
m_levelmin = gh.levelmin();
m_nbnd = gh.m_nbnd;
m_xoffabs = gh.m_xoffabs;
m_yoffabs = gh.m_yoffabs;
m_zoffabs = gh.m_zoffabs;
return *this;
} // throw std::runtime_error("GridHierarchy::operator= : attempt to operate on incompatible data");
for (unsigned i = 0; i < m_pgrids.size(); ++i)
(*m_pgrids[i]) = *gh.get_grid(i);
return *this;
}
/*
//! assignment operator
GridHierarchy& operator=( const GridHierarchy<T>& gh )
{
for( unsigned i=0; i<m_pgrids.size(); ++i )
delete m_pgrids[i];
m_pgrids.clear();
for( unsigned i=0; i<=gh.levelmax(); ++i )
m_pgrids.push_back( new MeshvarBnd<T>( *gh.get_grid(i) ) );
m_levelmin = gh.levelmin();
m_nbnd = gh.m_nbnd;
m_xoffabs = gh.m_xoffabs;
m_yoffabs = gh.m_yoffabs;
m_zoffabs = gh.m_zoffabs;
return *this;
}
*/
/*! add a new refinement patch to the so-far finest level
* @param xoff x-offset in units of the coarser grid (finest grid before adding new patch)
* @param yoff y-offset in units of the coarser grid (finest grid before adding new patch)
* @param zoff z-offset in units of the coarser grid (finest grid before adding new patch)
* @param nx x-extent in fine grid cells
* @param ny y-extent in fine grid cells
* @param nz z-extent in fine grid cells
*/
void add_patch(unsigned xoff, unsigned yoff, unsigned zoff, unsigned nx, unsigned ny, unsigned nz) {
m_pgrids.push_back(new MeshvarBnd<T>(m_nbnd, nx, ny, nz, xoff, yoff, zoff));
m_pgrids.back()->zero();
//.. add absolute offsets (in units of current level grid cells)
m_xoffabs.push_back(2 * (m_xoffabs.back() + xoff));
m_yoffabs.push_back(2 * (m_yoffabs.back() + yoff));
m_zoffabs.push_back(2 * (m_zoffabs.back() + zoff));
m_ref_masks.push_back(new refinement_mask(nx, ny, nz, 0));
}
/*! cut a refinement patch to the specified size
* @param ilevel grid level on which to perform the size adjustment
* @param xoff new x-offset in units of the coarser grid (finest grid before adding new patch)
* @param yoff new y-offset in units of the coarser grid (finest grid before adding new patch)
* @param zoff new z-offset in units of the coarser grid (finest grid before adding new patch)
* @param nx new x-extent in fine grid cells
* @param ny new y-extent in fine grid cells
* @param nz new z-extent in fine grid cells
*/
void cut_patch(unsigned ilevel, unsigned xoff, unsigned yoff, unsigned zoff, unsigned nx, unsigned ny, unsigned nz) {
unsigned dx, dy, dz, dxtop, dytop, dztop;
dx = xoff - m_xoffabs[ilevel];
dy = yoff - m_yoffabs[ilevel];
dz = zoff - m_zoffabs[ilevel];
assert(dx % 2 == 0 && dy % 2 == 0 && dz % 2 == 0);
dxtop = m_pgrids[ilevel]->offset(0) + dx / 2;
dytop = m_pgrids[ilevel]->offset(1) + dy / 2;
dztop = m_pgrids[ilevel]->offset(2) + dz / 2;
MeshvarBnd<T> *mnew = new MeshvarBnd<T>(m_nbnd, nx, ny, nz, dxtop, dytop, dztop);
//... copy data
for (unsigned i = 0; i < nx; ++i)
for (unsigned j = 0; j < ny; ++j)
for (unsigned k = 0; k < nz; ++k)
(*mnew)(i, j, k) = (*m_pgrids[ilevel])(i + dx, j + dy, k + dz);
//... replace in hierarchy
delete m_pgrids[ilevel];
m_pgrids[ilevel] = mnew;
//... update offsets
m_xoffabs[ilevel] += dx;
m_yoffabs[ilevel] += dy;
m_zoffabs[ilevel] += dz;
if (ilevel < levelmax()) {
m_pgrids[ilevel + 1]->offset(0) -= dx;
m_pgrids[ilevel + 1]->offset(1) -= dy;
m_pgrids[ilevel + 1]->offset(2) -= dz;
}