Source

yt-3.0 / yt / geometry / oct_container.pyx

  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
"""
Oct container

Author: Matthew Turk <matthewturk@gmail.com>
Affiliation: Columbia University
Homepage: http://yt.enzotools.org/
License:
  Copyright (C) 2011 Matthew Turk.  All Rights Reserved.

  This file is part of yt.

  yt is free software; you can redistribute it and/or modify
  it under the terms of the GNU General Public License as published by
  the Free Software Foundation; either version 3 of the License, or
  (at your option) any later version.

  This program is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.

  You should have received a copy of the GNU General Public License
  along with this program.  If not, see <http://www.gnu.org/licenses/>.
"""

from libc.stdlib cimport malloc, free, qsort
cimport numpy as np
import numpy as np
from oct_container cimport Oct, OctAllocationContainer, OctreeContainer
from selection_routines cimport SelectorObject
cimport cython

cdef extern from "stdlib.h":
    # NOTE that size_t might not be int
    void *alloca(int)

cdef OctAllocationContainer *allocate_octs(
        int n_octs, OctAllocationContainer *prev):
    cdef OctAllocationContainer *n_cont
    cdef Oct *oct
    cdef int n, i, j, k
    n_cont = <OctAllocationContainer *> malloc(
        sizeof(OctAllocationContainer))
    if prev == NULL:
        n_cont.offset = 0
    else:
        n_cont.offset = prev.offset + prev.n
    n_cont.my_octs = <Oct *> malloc(sizeof(Oct) * n_octs)
    if n_cont.my_octs == NULL:
        raise MemoryError
    n_cont.n = n_octs
    n_cont.n_assigned = 0
    for n in range(n_octs):
        oct = &n_cont.my_octs[n]
        oct.parent = NULL
        oct.ind = oct.domain = -1
        oct.local_ind = n + n_cont.offset
        oct.level = -1
        for i in range(2):
            for j in range(2):
                for k in range(2):
                    oct.children[i][j][k] = NULL
    if prev != NULL:
        prev.next = n_cont
    n_cont.next = NULL
    return n_cont

cdef void free_octs(
        OctAllocationContainer *first):
    cdef OctAllocationContainer *cur
    while first != NULL:
        cur = first
        free(first.my_octs)
        first = cur.next
        free(cur)

# Here is the strategy for RAMSES containers:
#   * Read each domain individually, creating *all* octs found in that domain
#     file, even if they reside on other CPUs.
#   * Only allocate octs that reside on >= domain
#   * For all octs, insert into tree, which may require traversing existing
#     octs
#   * Note that this doesn ot allow OctAllocationContainer to exactly be a
#     chunk, but it is close.  For IO chunking, we can theoretically examine
#     those octs that live inside a given allocator.

cdef class OctreeContainer:

    def __init__(self, domain_dimensions, domain_left_edge, domain_right_edge):
        self.nn[0], self.nn[1], self.nn[2] = domain_dimensions
        cdef int i, j, k, p
        self.max_domain = -1
        p = 0
        self.nocts = 0 # Increment when initialized
        self.root_mesh = <Oct****> malloc(sizeof(void*) * self.nn[0])
        for i in range(3):
            self.nn[i] = domain_dimensions[i]
        for i in range(self.nn[0]):
            self.root_mesh[i] = <Oct ***> malloc(sizeof(void*) * self.nn[1])
            for j in range(self.nn[1]):
                self.root_mesh[i][j] = <Oct **> malloc(sizeof(void*) * self.nn[2])
                for k in range(self.nn[2]):
                    self.root_mesh[i][j][k] = NULL
        # We don't initialize the octs yet
        for i in range(3):
            self.DLE[i] = domain_left_edge[i]
            self.DRE[i] = domain_right_edge[i]

    def __dealloc__(self):
        free_octs(self.cont)
        for i in range(self.nn[0]):
            for j in range(self.nn[1]):
                if self.root_mesh[i][j] == NULL: continue
                free(self.root_mesh[i][j])
            if self.root_mesh[i] == NULL: continue
            free(self.root_mesh[i])
        if self.root_mesh == NULL: return
        free(self.root_mesh)

    def __iter__(self):
        cdef OctAllocationContainer *cur = self.cont
        cdef Oct *this
        cdef int i
        while cur != NULL:
            for i in range(cur.n):
                this = &cur.my_octs[i]
                yield (this.ind, this.local_ind, this.domain)
            cur = cur.next

    cdef void oct_bounds(self, Oct *o, np.float64_t *corner, np.float64_t *size):
        cdef np.float64_t base_dx[3]
        cdef int i
        for i in range(3):
            size[i] = (self.DRE[i] - self.DLE[i]) / (self.nn[i] << o.level)
            corner[i] = o.pos[i] * size[i] + self.DLE[i]

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    cdef Oct *get(self, ppos):
        cdef np.int64_t ind[3]
        cdef np.float64_t dds[3], cp[3], pp[3]
        cdef Oct *cur
        cdef int i
        for i in range(3):
            pp[i] = ppos[i] - self.DLE[i]
            dds[i] = (self.DRE[i] - self.DLE[i])/self.nn[i]
            ind[i] = <np.int64_t> (pp[i]/dds[i])
            cp[i] = (ind[i] + 0.5) * dds[i]
        cur = self.root_mesh[ind[0]][ind[1]][ind[2]]
        while cur.children[0][0][0] != NULL:
            for i in range(3):
                dds[i] = dds[i] / 2.0
                if cp[i] > pp[i]:
                    ind[i] = 0
                    cp[i] -= dds[i] / 2.0
                else:
                    ind[i] = 1
                    cp[i] += dds[i]/2.0
            cur = cur.children[ind[0]][ind[1]][ind[2]]
        return cur

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def count_cells(self, SelectorObject selector,
              np.ndarray[np.uint8_t, ndim=2, cast=True] mask):
        cdef int i, j, k, oi
        # pos here is CELL center, not OCT center.
        cdef np.float64_t pos[3]
        cdef int n = mask.shape[0]
        cdef np.float64_t base_dx[3], dx[3]
        cdef int eterm[3]
        cdef np.ndarray[np.int64_t, ndim=1] count
        count = np.zeros(self.max_domain, 'int64')
        for i in range(3):
            # This is the base_dx, but not the base distance from the center
            # position.  Note that the positions will also all be offset by
            # dx/2.0.
            base_dx[i] = (self.DRE[i] - self.DLE[i])/self.nn[i]
        cur = self.cont
        for oi in range(n):
            if oi - cur.offset >= cur.n:
                cur = cur.next
            o = &cur.my_octs[oi - cur.offset]
            for i in range(8):
                count[o.domain - 1] += mask[oi,i]
        return count

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    cdef void neighbors(self, Oct* o, Oct* neighbors[27]):
        cdef np.int64_t curopos[3]
        cdef np.int64_t curnpos[3]
        cdef np.int64_t npos[3]
        cdef int i, j, k, ni, nj, nk, ind[3], nn, dl, skip
        cdef np.float64_t dds[3], cp[3], pp[3]
        cdef Oct* candidate
        for i in range(27): neighbors[i] = NULL
        nn = 0
        for ni in range(3):
            for nj in range(3):
                for nk in range(3):
                    #print "Handling neighbor", nn
                    if ni == nj == nk == 1:
                        neighbors[nn] = o
                        nn += 1
                        continue
                    npos[0] = o.pos[0] + (ni - 1)
                    npos[1] = o.pos[1] + (nj - 1)
                    npos[2] = o.pos[2] + (nk - 1)
                    for i in range(3):
                        # Periodicity
                        if npos[i] == -1:
                            npos[i] = (self.nn[i]  << o.level) - 1
                        elif npos[i] == (self.nn[i] << o.level):
                            npos[i] = 0
                        curopos[i] = o.pos[i]
                        curnpos[i] = npos[i] 
                    # Now we have our neighbor position and a safe place to
                    # keep it.  curnpos will be the root index of the neighbor
                    # at a given level, and npos will be constant.  curopos is
                    # the candidate root at a level.
                    candidate = o
                    while candidate != NULL:
                        if ((curopos[0] == curnpos[0]) and 
                            (curopos[1] == curnpos[1]) and
                            (curopos[2] == curnpos[2])):
                            break
                        # This one doesn't meet it, so we pop up a level.
                        # First we update our positions, then we update our
                        # candidate.
                        for i in range(3):
                            # We strip a digit off the right
                            curopos[i] = (curopos[i] >> 1)
                            curnpos[i] = (curnpos[i] >> 1)
                        # Now we update to the candidate's parent, which should
                        # have a matching position to curopos[]
                        candidate = candidate.parent
                    if candidate == NULL:
                        # Worst case scenario
                        for i in range(3):
                            ind[i] = (npos[i] >> (o.level))
                        #print "NULL", npos[0], npos[1], npos[2], ind[0], ind[1], ind[2], o.level
                        candidate = self.root_mesh[ind[0]][ind[1]][ind[2]]
                    # Now we have the common root, which may be NULL
                    while candidate.level < o.level:
                        dl = o.level - (candidate.level + 1)
                        for i in range(3):
                            ind[i] = (npos[i] >> dl) & 1
                        if candidate.children[0][0][0] == NULL: break
                        candidate = candidate.children[ind[0]][ind[1]][ind[2]]
                    neighbors[nn] = candidate
                    nn += 1

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def get_neighbor_boundaries(self, ppos):
        cdef Oct *main = self.get(ppos)
        cdef Oct* neighbors[27]
        self.neighbors(main, neighbors)
        cdef np.ndarray[np.float64_t, ndim=2] bounds
        cdef np.float64_t corner[3], size[3]
        bounds = np.zeros((27,6), dtype="float64")
        cdef int i, ii
        cdef np.float64_t base_dx[3]
        tnp = 0
        for i in range(27):
            self.oct_bounds(neighbors[i], corner, size)
            for ii in range(3):
                bounds[i, ii] = corner[ii]
                bounds[i, 3+ii] = size[ii]
        return bounds

cdef class RAMSESOctreeContainer(OctreeContainer):

    def allocate_domains(self, domain_counts):
        cdef int count, i
        cdef OctAllocationContainer *cur = self.cont
        assert(cur == NULL)
        self.max_domain = len(domain_counts) # 1-indexed
        self.domains = <OctAllocationContainer **> malloc(
            sizeof(OctAllocationContainer *) * len(domain_counts))
        for i, count in enumerate(domain_counts):
            cur = allocate_octs(count, cur)
            if self.cont == NULL: self.cont = cur
            self.domains[i] = cur
        
    def __dealloc__(self):
        # This gets called BEFORE the superclass deallocation.  But, both get
        # called.
        if self.domains != NULL: free(self.domains)

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def count(self, np.ndarray[np.uint8_t, ndim=1, cast=True] mask,
                     split = False):
        cdef int n = mask.shape[0]
        cdef int i, dom
        cdef OctAllocationContainer *cur
        cdef np.ndarray[np.int64_t, ndim=1] count
        count = np.zeros(self.max_domain, 'int64')
        # This is the idiom for iterating over many containers.
        cur = self.cont
        for i in range(n):
            if i - cur.offset >= cur.n: cur = cur.next
            if mask[i] == 1:
                count[cur.my_octs[i - cur.offset].domain - 1] += 1
        return count

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def add(self, int curdom, int curlevel, int ng,
            np.ndarray[np.float64_t, ndim=2] pos,
            int local_domain, int skip_boundary = 1):
        cdef int level, no, p, i, j, k, ind[3]
        cdef int local = (local_domain == curdom)
        cdef Oct *cur, *next = NULL
        cdef np.float64_t pp[3], cp[3], dds[3]
        no = pos.shape[0] #number of octs
        if curdom > self.max_domain: curdom = local_domain
        cdef OctAllocationContainer *cont = self.domains[curdom - 1]
        cdef int initial = cont.n_assigned
        cdef int in_boundary = 0
        # How do we bootstrap ourselves?
        for p in range(no):
            #for every oct we're trying to add find the 
            #floating point unitary position on this level
            in_boundary = 0
            for i in range(3):
                pp[i] = pos[p, i]
                dds[i] = (self.DRE[i] + self.DLE[i])/self.nn[i]
                ind[i] = <np.int64_t> (pp[i]/dds[i])
                cp[i] = (ind[i] + 0.5) * dds[i]
                if ind[i] < 0 or ind[i] >= self.nn[i]:
                    in_boundary = 1
            if skip_boundary == in_boundary == 1: continue
            cur = self.root_mesh[ind[0]][ind[1]][ind[2]]
            if cur == NULL:
                if curlevel != 0:
                    raise RuntimeError
                cur = &cont.my_octs[cont.n_assigned]
                cur.parent = NULL
                cur.level = 0
                for i in range(3):
                    cur.pos[i] = ind[i]
                cont.n_assigned += 1
                self.nocts += 1
                self.root_mesh[ind[0]][ind[1]][ind[2]] = cur
            # Now we find the location we want
            # Note that RAMSES I think 1-findiceses levels, but we don't.
            for level in range(curlevel):
                # At every level, find the cell this oct
                # lives inside
                for i in range(3):
                    #as we get deeper, oct size halves
                    dds[i] = dds[i] / 2.0
                    if cp[i] > pp[i]: 
                        ind[i] = 0
                        cp[i] -= dds[i]/2.0
                    else:
                        ind[i] = 1
                        cp[i] += dds[i]/2.0
                # Check if it has not been allocated
                next = cur.children[ind[0]][ind[1]][ind[2]]
                if next == NULL:
                    next = &cont.my_octs[cont.n_assigned]
                    cur.children[ind[0]][ind[1]][ind[2]] = next
                    cont.n_assigned += 1
                    next.parent = cur
                    for i in range(3):
                        next.pos[i] = ind[i] + (cur.pos[i] << 1)
                    next.level = level + 1
                    self.nocts += 1
                cur = next
            cur.domain = curdom
            if local == 1:
                cur.ind = p
            cur.level = curlevel
        return cont.n_assigned - initial

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def icoords(self, int domain_id,
                np.ndarray[np.uint8_t, ndim=2, cast=True] mask,
                np.int64_t cell_count,
                np.ndarray[np.int64_t, ndim=1] level_counts):
        # Wham, bam, it's a scam
        cdef np.int64_t i, j, k, oi, ci, n, ii, level
        cdef OctAllocationContainer *cur = self.domains[domain_id - 1]
        cdef Oct *o
        n = mask.shape[0]
        cdef np.ndarray[np.int64_t, ndim=2] coords
        coords = np.empty((cell_count, 3), dtype="int64")
        for oi in range(cur.n):
            o = &cur.my_octs[oi]
            for i in range(2):
                for j in range(2):
                    for k in range(2):
                        ii = ((k*2)+j)*2+i
                        if mask[o.local_ind, ii] == 0: continue
                        ci = level_counts[o.level]
                        coords[ci, 0] = (o.pos[0] << 1) + i
                        coords[ci, 1] = (o.pos[1] << 1) + j
                        coords[ci, 2] = (o.pos[2] << 1) + k
                        level_counts[o.level] += 1
        return coords

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def ires(self, int domain_id,
                np.ndarray[np.uint8_t, ndim=2, cast=True] mask,
                np.int64_t cell_count,
                np.ndarray[np.int64_t, ndim=1] level_counts):
        # Wham, bam, it's a scam
        cdef np.int64_t i, j, k, oi, ci, n
        cdef OctAllocationContainer *cur = self.domains[domain_id - 1]
        cdef Oct *o
        n = mask.shape[0]
        cdef np.ndarray[np.int64_t, ndim=1] levels
        levels = np.empty(cell_count, dtype="int64")
        ci = 0
        for oi in range(cur.n):
            o = &cur.my_octs[oi]
            for i in range(8):
                if mask[oi + cur.offset, i] == 0: continue
                ci = level_counts[o.level]
                levels[ci] = o.level
                level_counts[o.level] += 1
        return levels

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def count_levels(self, int max_level, int domain_id,
                     np.ndarray[np.uint8_t, ndim=2, cast=True] mask):
        cdef np.ndarray[np.int64_t, ndim=1] level_count
        cdef OctAllocationContainer *cur = self.domains[domain_id - 1]
        cdef Oct *o
        cdef int oi, i
        level_count = np.zeros(max_level+1, 'int64')
        for oi in range(cur.n):
            o = &cur.my_octs[oi]
            for i in range(8):
                if mask[o.local_ind, i] == 0: continue
                level_count[o.level] += 1
        return level_count

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def fcoords(self, int domain_id,
                np.ndarray[np.uint8_t, ndim=2, cast=True] mask,
                np.int64_t cell_count,
                np.ndarray[np.int64_t, ndim=1] level_counts):
        # Wham, bam, it's a scam
        cdef np.int64_t i, j, k, oi, ci, n, ii
        cdef OctAllocationContainer *cur = self.domains[domain_id - 1]
        cdef Oct *o
        cdef np.float64_t pos[3]
        cdef np.float64_t base_dx[3], dx[3]
        n = mask.shape[0]
        cdef np.ndarray[np.float64_t, ndim=2] coords
        coords = np.empty((cell_count, 3), dtype="float64")
        for i in range(3):
            # This is the base_dx, but not the base distance from the center
            # position.  Note that the positions will also all be offset by
            # dx/2.0.  This is also for *oct grids*, not cells.
            base_dx[i] = (self.DRE[i] - self.DLE[i])/self.nn[i]
        for oi in range(cur.n):
            o = &cur.my_octs[oi]
            for i in range(3):
                # This gives the *grid* width for this level
                dx[i] = base_dx[i] / (1 << o.level)
                # o.pos is the *grid* index, so pos[i] is the center of the
                # first cell in the grid
                pos[i] = self.DLE[i] + o.pos[i]*dx[i] + dx[i]/4.0
                dx[i] = dx[i] / 2.0 # This is now the *offset* 
            for i in range(2):
                for j in range(2):
                    for k in range(2):
                        ii = ((k*2)+j)*2+i
                        if mask[o.local_ind, ii] == 0: continue
                        ci = level_counts[o.level]
                        coords[ci, 0] = pos[0] + dx[0] * i
                        coords[ci, 1] = pos[1] + dx[1] * j
                        coords[ci, 2] = pos[2] + dx[2] * k
                        level_counts[o.level] += 1
        return coords

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def fill_level(self, int domain, int level, dest_fields, source_fields,
                   np.ndarray[np.uint8_t, ndim=2, cast=True] mask, int offset):
        cdef np.ndarray[np.float64_t, ndim=2] source
        cdef np.ndarray[np.float64_t, ndim=1] dest
        cdef OctAllocationContainer *dom = self.domains[domain - 1]
        cdef Oct *o
        cdef int n
        cdef int i, j, k, ii
        cdef int local_pos, local_filled
        cdef np.float64_t val
        for key in dest_fields:
            local_filled = 0
            dest = dest_fields[key]
            source = source_fields[key]
            for n in range(dom.n):
                o = &dom.my_octs[n]
                if o.level != level: continue
                for i in range(2):
                    for j in range(2):
                        for k in range(2):
                            ii = ((k*2)+j)*2+i
                            if mask[o.local_ind, ii] == 0: continue
                            dest[local_filled + offset] = source[o.ind, ii]
                            local_filled += 1
        return local_filled

cdef public int compare_octs(void *vo1, void *vo2) nogil:
    cdef Oct *o1 = (<Oct**> vo1)[0]
    cdef Oct *o2 = (<Oct**> vo2)[0]
    if o1.domain < o2.domain: return -1
    elif o1.domain == o2.domain:
        if o1.level < o2.level: return -1
        if o1.level > o2.level: return 1
        else: return 0
    elif o1.domain > o2.domain: return 1

cdef class ParticleOctreeContainer(OctreeContainer):
    cdef ParticleArrays *first_sd
    cdef ParticleArrays *last_sd
    cdef Oct** oct_list
    cdef np.int64_t *dom_offsets

    def allocate_root(self):
        cdef int i, j, k
        cdef Oct *cur
        for i in range(self.nn[0]):
            for j in range(self.nn[1]):
                for k in range(self.nn[2]):
                    cur = self.allocate_oct()
                    cur.level = 0
                    cur.pos[0] = i
                    cur.pos[1] = j
                    cur.pos[2] = k
                    cur.parent = NULL
                    self.root_mesh[i][j][k] = cur

    def __dealloc__(self):
        cdef i, j, k
        for i in range(self.nn[0]):
            for j in range(self.nn[1]):
                for k in range(self.nn[2]):
                    self.visit_free(self.root_mesh[i][j][k])
        free(self.oct_list)
        free(self.dom_offsets)

    cdef void visit_free(self, Oct *o):
        cdef int i, j, k
        for i in range(2):
            for j in range(2):
                for k in range(2):
                    if o.children[i][j][k] == NULL: continue
                    self.visit_free(o.children[i][j][k])
        if o.sd.np >= 0:
            if o.sd.pos != NULL:
                for i in range(3):
                    free(o.sd.pos[i])
                free(o.sd.pos)
        free(o)

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def icoords(self, np.ndarray[np.uint8_t, ndim=2, cast=True] mask,
                np.int64_t cell_count):
        cdef np.ndarray[np.int64_t, ndim=2] coords
        coords = np.empty((cell_count, 3), dtype="int64")
        cdef int oi, i, ci, ii
        ci = 0
        for oi in range(self.nocts):
            o = self.oct_list[oi]
            for i in range(2):
                for j in range(2):
                    for k in range(2):
                        ii = ((k*2)+j)*2+i
                        if mask[oi, ii] == 1:
                            coords[ci, 0] = (o.pos[0] << 1) + i
                            coords[ci, 1] = (o.pos[1] << 1) + j
                            coords[ci, 2] = (o.pos[2] << 1) + k
                            ci += 1
        return coords

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def ires(self, np.ndarray[np.uint8_t, ndim=2, cast=True] mask,
                np.int64_t cell_count):
        cdef np.ndarray[np.int64_t, ndim=1] res
        res = np.empty(cell_count, dtype="int64")
        cdef int oi, i, ci
        ci = 0
        for oi in range(self.nocts):
            o = self.oct_list[oi]
            for i in range(8):
                if mask[oi, i] == 1:
                    res[ci] = o.level
                    ci += 1
        return res

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def fcoords(self, np.ndarray[np.uint8_t, ndim=2, cast=True] mask,
                np.int64_t cell_count):
        cdef np.ndarray[np.float64_t, ndim=2] coords
        coords = np.empty((cell_count, 3), dtype="float64")
        cdef int oi, i, ci
        cdef np.float64_t base_dx[3], dx[3], pos[3]
        for i in range(3):
            # This is the base_dx, but not the base distance from the center
            # position.  Note that the positions will also all be offset by
            # dx/2.0.  This is also for *oct grids*, not cells.
            base_dx[i] = (self.DRE[i] - self.DLE[i])/self.nn[i]
        ci = 0
        cdef int proc
        for oi in range(self.nocts):
            proc = 0
            for i in range(8):
                if mask[oi, i] == 1:
                    proc = 1
                    break
            if proc == 0: continue
            o = self.oct_list[oi]
            for i in range(3):
                # This gives the *grid* width for this level
                dx[i] = base_dx[i] / (1 << o.level)
                # o.pos is the *grid* index, so pos[i] is the center of the
                # first cell in the grid
                pos[i] = self.DLE[i] + o.pos[i]*dx[i] + dx[i]/4.0
                dx[i] = dx[i] / 2.0 # This is now the *offset* 
            for i in range(2):
                for j in range(2):
                    for k in range(2):
                        ii = ((k*2)+j)*2+i
                        if mask[oi, ii] == 0: continue
                        coords[ci, 0] = pos[0] + dx[0] * i
                        coords[ci, 1] = pos[1] + dx[1] * j
                        coords[ci, 2] = pos[2] + dx[2] * k
                        ci += 1
        return coords

    def allocate_domains(self, domain_counts):
        pass

    def finalize(self):
        self.oct_list = <Oct**> malloc(sizeof(Oct*)*self.nocts)
        cdef np.int64_t i = 0
        cdef ParticleArrays *c = self.first_sd
        while c != NULL:
            self.oct_list[i] = c.oct
            if c.oct.domain == -1:
                assert(c.oct.children[0][0][0] == NULL)
                assert(c.np == 0)
            if c.np >= 0:
                for j in range(3):
                    free(c.pos[j])
                free(c.pos)
                c.pos = NULL
                # We should also include a shortening of the domain IDs here
            c = c.next
            i += 1
        assert(i == self.nocts)
        qsort(self.oct_list, self.nocts, sizeof(Oct*), &compare_octs)
        cdef int cur_dom = -1
        self.dom_offsets = <np.int64_t *>malloc(sizeof(np.int64_t) *
                                                self.max_domain + 1)
        for i in range(self.nocts):
            if self.oct_list[i].domain > cur_dom:
                cur_dom = self.oct_list[i].domain
                self.dom_offsets[cur_dom] = i

    cdef Oct* allocate_oct(self):
        self.nocts += 1
        cdef Oct *my_oct = <Oct*> malloc(sizeof(Oct))
        cdef ParticleArrays *sd = <ParticleArrays*> \
            malloc(sizeof(ParticleArrays))
        cdef int i, j, k
        my_oct.ind = my_oct.domain = -1
        my_oct.domain = -1
        my_oct.local_ind = self.nocts - 1
        my_oct.pos[0] = my_oct.pos[1] = my_oct.pos[2] = -1
        my_oct.level = -1
        my_oct.sd = sd
        for i in range(2):
            for j in range(2):
                for k in range(2):
                    my_oct.children[i][j][k] = NULL
        my_oct.parent = NULL
        if self.first_sd == NULL:
            self.first_sd = sd
        if self.last_sd != NULL:
            self.last_sd.next = sd
        self.last_sd = sd
        sd.oct = my_oct
        sd.next = NULL
        sd.pos = <np.float64_t **> malloc(sizeof(np.float64_t*) * 3)
        for i in range(3):
            sd.pos[i] = <np.float64_t *> malloc(sizeof(np.float64_t) * 32)
        for i in range(32):
            sd.pos[0][i] = sd.pos[1][i] = sd.pos[2][i] = 0.0
        sd.np = 0
        return my_oct

    def linearly_count(self):
        cdef np.int64_t total = 0
        cdef ParticleArrays *c = self.first_sd
        while c != NULL:
            total += 1
            c = c.next
        return total

    def add(self, np.ndarray[np.float64_t, ndim=2] pos, np.int64_t domain_id):
        cdef int no = pos.shape[0]
        cdef int p, i, level
        cdef np.float64_t dds[3], cp[3], pp[3]
        cdef int ind[3]
        self.max_domain = max(self.max_domain, domain_id)
        cdef int mid, mad
        if self.root_mesh[0][0][0] == NULL: self.allocate_root()
        for p in range(no):
            level = 0
            for i in range(3):
                pp[i] = pos[p, i]
                dds[i] = (self.DRE[i] + self.DLE[i])/self.nn[i]
                ind[i] = <np.int64_t> (pp[i]/dds[i])
                cp[i] = (ind[i] + 0.5) * dds[i]
            cur = self.root_mesh[ind[0]][ind[1]][ind[2]]
            if cur == NULL:
                raise RuntimeError
            if self._check_refine(cur, cp, domain_id) == 1:
                self.refine_oct(cur, cp)
            while cur.sd.np < 0:
                for i in range(3):
                    dds[i] = dds[i] / 2.0
                    if cp[i] > pp[i]:
                        ind[i] = 0
                        cp[i] -= dds[i] / 2.0
                    else:
                        ind[i] = 1
                        cp[i] += dds[i]/2.0
                cur = cur.children[ind[0]][ind[1]][ind[2]]
                level += 1
                if self._check_refine(cur, cp, domain_id) == 1:
                    self.refine_oct(cur, cp)
            # Now we copy in our particle 
            pi = cur.sd.np
            cur.level = level
            for i in range(3):
                cur.sd.pos[i][pi] = pp[i]
            cur.domain = domain_id
            cur.sd.np += 1

    cdef int _check_refine(self, Oct *cur, np.float64_t cp[3], int domain_id):
        if cur.children[0][0][0] != NULL:
            return 0
        elif cur.sd.np == 32:
            return 1
        elif cur.domain >= 0 and cur.domain != domain_id:
            return 1
        return 0

    cdef void refine_oct(self, Oct *o, np.float64_t pos[3]):
        cdef int i, j, k, m, ind[3]
        cdef Oct *noct
        for i in range(2):
            for j in range(2):
                for k in range(2):
                    noct = self.allocate_oct()
                    noct.level = o.level + 1
                    noct.pos[0] = (o.pos[0] << 1) + i
                    noct.pos[1] = (o.pos[1] << 1) + j
                    noct.pos[2] = (o.pos[2] << 1) + k
                    noct.parent = o
                    o.children[i][j][k] = noct
        for m in range(o.sd.np):
            for i in range(3):
                if o.sd.pos[i][m] < pos[i]:
                    ind[i] = 0
                else:
                    ind[i] = 1
            noct = o.children[ind[0]][ind[1]][ind[2]]
            k = noct.sd.np
            for i in range(3):
                noct.sd.pos[i][k] = o.sd.pos[i][m]
            noct.domain = o.domain
            noct.sd.np += 1
        o.sd.np = -1
        for i in range(3):
            free(o.sd.pos[i])
        free(o.sd.pos)

    def recursively_count(self):
        cdef int i, j, k
        cdef np.int64_t counts[128]
        for i in range(128): counts[i] = 0
        for i in range(self.nn[0]):
            for j in range(self.nn[1]):
                for k in range(self.nn[2]):
                    if self.root_mesh[i][j][k] != NULL:
                        self.visit(self.root_mesh[i][j][k], counts)
        level_counts = {}
        for i in range(128):
            if counts[i] == 0: break
            level_counts[i] = counts[i]
        return level_counts
        
    cdef visit(self, Oct *o, np.int64_t *counts, level = 0):
        cdef int i, j, k
        counts[level] += 1
        for i in range(2):
            for j in range(2):
                for k in range(2):
                    if o.children[i][j][k] != NULL:
                        self.visit(o.children[i][j][k], counts, level + 1)
        return

    def domain_identify(self, np.ndarray[np.uint8_t, ndim=2, cast=True] mask):
        cdef int i, oi, m
        cdef Oct *o
        cdef np.ndarray[np.uint8_t, ndim=1, cast=True] dmask
        dmask = np.zeros(self.max_domain+1, dtype='uint8')
        for oi in range(self.nocts):
            m = 0
            o = self.oct_list[oi]
            if o.sd.np <= 0 or o.domain == -1: continue
            for i in range(8):
                if mask[oi, i] == 1:
                    m = 1
                    break
            if m == 0: continue
            dmask[o.domain] = 1
        return dmask.astype("bool")

    @cython.boundscheck(False)
    @cython.wraparound(False)
    @cython.cdivision(True)
    def count_neighbor_particles(self, ppos):
        cdef Oct *main = self.get(ppos)
        cdef Oct* neighbors[27]
        self.neighbors(main, neighbors)
        cdef int i, ni, dl, tnp
        tnp = 0
        for i in range(27):
            if neighbors[i].sd != NULL:
                tnp += neighbors[i].sd.np
        return tnp