File:/home/sbrandt/cactus/Cactus/arrangements/Carpet/CarpetLib/src/ggf.hh
1:#ifndef GGF_HH
2:#define GGF_HH
3:
4:#include <cctk.h>
5:
6:#include <cassert>
7:#include <iostream>
8:#include <string>
9:#include <vector>
10:
11:#include "defs.hh"
12:#include "dh.hh"
13:#include "gdata.hh"
14:#include "gh.hh"
15:#include "th.hh"
16:
17:using namespace std;
18:
19:
20:
21:// Forward declaration
22:class ggf;
23:
24:// Output
25:ostream& operator<< (ostream& os, const ggf& f);
26:
27:
28:
29:// A generic grid function without type information
30:class ggf {
31:  
32:  static list<ggf*> allggf;
33:  list<ggf*>::iterator allggfi;
34:
35:  // Types
36:  typedef vector <pseudoregion_t> pvect;
37:  typedef vector <sendrecv_pseudoregion_t> srpvect;
38:  
39:  typedef gdata*        tdata;  // data ...
40:  typedef vector<tdata> fdata;  // ... for each time level
41:  typedef vector<fdata> cdata;  // ... for each local component
42:  typedef vector<cdata> rdata;  // ... for each refinement level
43:  typedef vector<rdata> mdata;  // ... for each multigrid level
44:  
45:public:    // should be readonly
46:  
47:  // Fields
48:  const int varindex;           // Cactus variable index
49:  const operator_type transport_operator;
50:  
51:  const th &t;                  // time hierarchy
52:  const int prolongation_order_time; // order of temporal prolongation operator
53:  
54:  const gh &h;                  // grid hierarchy
55:  dh &d;   // data hierarchy
56:  dh::ggf_handle dh_handle;
57:
58:protected:
59:  vector<vector<int> > timelevels_; // time levels [ml][rl]
60:
61:  mdata storage;  // storage
62:  
63:public:
64:  const int vectorlength;       // vector length
65:  const int vectorindex;        // index of *this
66:  ggf* const vectorleader;      // first vector element
67:  
68:private:
69:  mdata oldstorage;             // temporary storage
70:  
71:public:
72:
73:  // Constructors
74:  ggf (const int varindex, const operator_type transport_operator,
75:       th& t, dh& d,
76:       const int prolongation_order_time,
77:       const int vectorlength, const int vectorindex,
78:       ggf* const vectorleader);
79:
80:  // Destructors
81:  virtual ~ggf ();
82:
83:  // Comparison
84:  bool operator== (const ggf& f) const CCTK_MEMBER_ATTRIBUTE_PURE;
85:  
86:  // Querying
87:  int timelevels (int const ml, int const rl) const
88:  {
89:    return timelevels_.AT(ml).AT(rl);
90:  }
91:  
92:  
93:  
94:  // Modifiers
95:  void set_timelevels (int ml, int rl, int new_timelevels);
96:
97:  void recompose_crop ();
98:  void recompose_allocate (int rl);
99:  void recompose_fill (comm_state& state, int rl, bool do_prolongate);
100:  void recompose_free_old (int rl);
101:  void recompose_free (int rl);
102:
103:  // Cycle the time levels by rotating the data sets
104:  void cycle_all (int rl, int ml);
105:
106:  // Uncycle the time levels by rotating the data sets
107:  void uncycle_all (int rl, int ml);
108:  
109:  // Flip the time levels by exchanging the data sets
110:  void flip_all (int rl, int ml);
111:
112:  // Fill all time levels from the current time level
113:  void fill_all (int rl, int ml);
114:
115:  // The grid boundaries have to be updated after calling mg_restrict,
116:  // mg_prolongate, ref_restrict, or ref_prolongate.
117:
118:  // "Updating" means here that the boundaries have to be
119:  // synchronised. They don't need to be prolongated.
120:
121:  // Synchronise the boundaries of a component
122:  void sync_all (comm_state& state, int tl, int rl, int ml);
123:
124:  // Prolongate the boundaries of a component
125:  void ref_bnd_prolongate_all (comm_state& state,
126:                               int tl, int rl, int ml, CCTK_REAL time);
127:
128:  // Restrict a multigrid level
129:  void mg_restrict_all (comm_state& state,
130:                        int tl, int rl, int ml, CCTK_REAL time);
131:
132:  // Prolongate a multigrid level
133:  void mg_prolongate_all (comm_state& state,
134:                          int tl, int rl, int ml, CCTK_REAL time);
135:
136:  // Restrict a refinement level
137:  void ref_restrict_all (comm_state& state,
138:                         int tl, int rl, int ml);
139:
140:  // Prolongate a refinement level
141:  void ref_prolongate_all (comm_state& state,
142:                           int tl, int rl, int ml, CCTK_REAL time);
143:
144:  // Reflux a refinement level
145:  void ref_reflux_all (comm_state& state,
146:                       int tl, int rl, int ml, int dir, int face);
147:
148:  // Reflux a refinement level
149:  void ref_reflux_prolongate_all (comm_state& state,
150:                                  int tl, int rl, int ml, int dir, int face);
151:  
152:  
153:  
154:  // Helpers
155:  
156:  virtual gdata* typed_data (int tl, int rl, int lc, int ml) const = 0;
157:  virtual gdata* new_typed_data () const = 0;
158:  
159:  
160:  
161:protected:
162:  
163:  // Transfer regions
164:  void
165:  transfer_from_all (comm_state & state,
166:                     int tl1, int rl1, int ml1,
167:                     srpvect const dh::fast_dboxes::* sendrecvs,
168:                     vector<int> const & tl2s, int rl2, int ml2,
169:                     CCTK_REAL const & time,
170:                     bool use_old_storage = false,
171:                     bool flip_send_recv = false,
172:                     islab const *restrict slabinfo = NULL);
173:  
174:  void
175:  transfer_from_all (comm_state & state,
176:                     int tl1, int rl1, int ml1,
177:                     srpvect const dh::fast_dboxes::* sendrecvs,
178:                     int tl2, int rl2, int ml2,
179:                     bool use_old_storage = false,
180:                     bool flip_send_recv = false,
181:                     islab const *restrict slabinfo = NULL)
182:  {
183:    vector <int> tl2s(1);
184:    tl2s.AT(0) = tl2;
185:    CCTK_REAL const time = t.get_time (ml2,rl2,tl2);
186:    transfer_from_all (state,
187:                       tl1, rl1, ml1,
188:                       sendrecvs,
189:                       tl2s, rl2, ml2,
190:                       time,
191:                       use_old_storage,
192:                       flip_send_recv,
193:                       slabinfo);
194:  }
195:
196:
197:
198:public:
199:  
200:  // Access to the data
201:  const gdata* data_pointer (int const tl, int const rl, int const lc, int const ml) const
202:  {
203:    assert (rl>=0 and rl<h.reflevels());
204:    assert (lc>=0 and lc<h.local_components(rl));
205:    assert (ml>=0 and ml<h.mglevels());
206:    assert (tl>=0 and tl<timelevels(ml, rl));
207:    return storage.AT(ml).AT(rl).AT(lc).AT(tl);
208:  }
209:  gdata* data_pointer (int const tl, int const rl, int const lc, int const ml)
210:  {
211:    assert (rl>=0 and rl<h.reflevels());
212:    assert (lc>=0 and lc<h.local_components(rl));
213:    assert (ml>=0 and ml<h.mglevels());
214:    assert (tl>=0 and tl<timelevels(ml, rl));
215:    return storage.AT(ml).AT(rl).AT(lc).AT(tl);
216:  }
217:  
218:#if 0
219:  virtual const gdata* operator() (int tl, int rl, int lc, int ml) const
220:    CCTK_MEMBER_ATTRIBUTE_PURE = 0;
221:  virtual gdata* operator() (int tl, int rl, int lc, int ml)
222:    CCTK_MEMBER_ATTRIBUTE_PURE = 0;
223:#endif
224:  
225:  
226:  
227:  // Output
228:  virtual size_t memory () const CCTK_MEMBER_ATTRIBUTE_PURE = 0;
229:  static size_t allmemory () CCTK_MEMBER_ATTRIBUTE_PURE;
230:  virtual ostream& output (ostream& os) const = 0;
231:
232:private:
233:  ggf ();                       // canonical default construtor
234:  ggf (const ggf &);            // canonical copy construtor
235:  ggf & operator= (const ggf &); // canonical copy
236:
237:};
238:
239:
240:
241:inline size_t memoryof (ggf const & f)
242:{
243:  return f.memory ();
244:}
245:
246:inline ostream& operator<< (ostream& os, const ggf& f)
247:{
248:  return f.output(os);
249:}
250:
251:
252:
253:#endif // GGF_HH