aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/dh.hh
blob: 70ffe38032904f9c9075faacad1ae05c40a87328 (plain)
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
/***************************************************************************
                          dh.hh  -  Data Hierarchy
			  A grid hierarchy plus ghost zones
                             -------------------
    begin                : Sun Jun 11 2000
    copyright            : (C) 2000 by Erik Schnetter
    email                : schnetter@astro.psu.edu

    $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/dh.hh,v 1.7 2001/04/23 08:10:15 schnetter Exp $

 ***************************************************************************/

/***************************************************************************
 *                                                                         *
 *   This program 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 2 of the License, or     *
 *   (at your option) any later version.                                   *
 *                                                                         *
 ***************************************************************************/

#ifndef DH_HH
#define DH_HH

#include <assert.h>

#include <iostream>
#include <list>
#include <string>
#include <vector>

#include "bbox.hh"
#include "bboxset.hh"
#include "defs.hh"
#include "gh.hh"
#include "vect.hh"

using namespace std;



// Forward declaration
template<int D> class dh;
template<int D> class generic_gf;

// Output
template<int D>
ostream& operator<< (ostream& os, const dh<D>& d);



// A data hierarchy (grid hierarchy plus ghost zones)
template<int D>
class dh {
  
  // Types
  typedef vect<int,D>    ivect;
  typedef bbox<int,D>    ibbox;
  typedef bboxset<int,D> ibset;
  typedef list<ibbox>    iblist;
  typedef vector<iblist> iblistvect; // vector of lists
  
public:
  struct dboxes {
    ibbox exterior;		// whole region (including boundaries)
    
    ibbox interior;		// interior (without boundaries)
    iblist send_mg_fine;
    iblist send_mg_coarse;
    iblist recv_mg_fine;
    iblist recv_mg_coarse;
    iblistvect send_ref_fine;
    iblistvect send_ref_coarse;
    iblistvect recv_ref_fine;
    iblistvect recv_ref_coarse;
    iblistvect send_sync;	// send while syncing
    iblistvect send_ref_bnd_fine;
    
    ibset boundaries;		// boundaries
    iblistvect recv_sync;	// received while syncing
    iblistvect recv_ref_bnd_coarse; // received from coarser grids
    ibset sync_not;		// not received while syncing (outer boundary of that level)
    ibset recv_not;		// not received while syncing or prolongating (globally outer boundary)
  };
private:
  
  struct dbases {
    ibbox exterior;		// whole region (including boundaries)
    ibbox interior;		// interior (without boundaries)
    ibset boundaries;		// boundaries
  };
  
  typedef vector<dboxes> mboxes; // ... for each multigrid level
  typedef vector<mboxes> cboxes; // ... for each component
  typedef vector<cboxes> rboxes; // ... for each refinement level
  
  typedef vector<dbases> mbases; // ... for each multigrid level
  typedef vector<mbases> rbases; // ... for each refinement level
  
public:				// should be readonly
  
  // Fields
  gh<D> &h;			// hierarchy
  ivect lghosts, ughosts;	// ghost zones
  int prolongation_order;	// order of spatial prolongation operator
  
  rboxes boxes;
  rbases bases;
  
  list<generic_gf<D>*> gfs;
  
public:
  
  // Constructors
  dh (gh<D>& h, const ivect& lghosts, const ivect& ughosts,
      int prolongation_order);
  
  // Destructors
  ~dh ();
  
  // Helpers
  int prolongation_stencil_size () const;
  
  // Modifiers
  void recompose ();
  
  // Grid function management
  void add (generic_gf<D>* f);
  void remove (generic_gf<D>* f);
  
  // Output
  void output (ostream& os) const;
};

template<int D>
inline ostream& operator<< (ostream& os, const dh<D>& d) {
  d.output(os);
  return os;
}



#if defined(TMPL_IMPLICIT)
#  include "dh.cc"
#endif

#endif // DH_HH