aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/ggf.hh
blob: 66012e469118ab6aeaadd9d6567d40295e381e72 (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
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
/***************************************************************************
                          ggf.hh  -  Generic Grid Function
                          grid function without type information
                             -------------------
    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/ggf.hh,v 1.12 2002/09/25 15:49:16 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 GGF_HH
#define GGF_HH

#include <assert.h>

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

#include "cctk.h"

#include "defs.hh"
#include "dggf.hh"
#include "dh.hh"
#include "gdata.hh"
#include "gh.hh"
#include "th.hh"

using namespace std;



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

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



// A generic grid function without type information
template<int D>
class generic_gf: public dimgeneric_gf {

  // Types

  typedef vect<int,D>    ivect;
  typedef bbox<int,D>    ibbox;
  typedef bboxset<int,D> ibset;
  typedef list<ibbox>    iblist;
  typedef vector<iblist> iblistvect;

  typedef generic_data<D>* tdata; // data ...
  typedef vector<tdata>    mdata; // ... for each multigrid level
  typedef vector<mdata>    cdata; // ... for each component
  typedef vector<cdata>    rdata; // ... for each refinement level
  typedef vector<rdata>    fdata; // ... for each time level

public:				// should be readonly

  // Fields

  gh<D> &h;			// grid hierarchy
  dh<D> &d;			// data hierarchy

protected:
  fdata storage;		// storage

public:

  // Constructors
  generic_gf (const string name, th& t, dh<D>& d,
              const int tmin, const int tmax,
	      const int prolongation_order_time);

  // Destructors
  virtual ~generic_gf ();

  // Comparison
  bool operator== (const generic_gf<D>& f) const;



  // Modifiers
  void recompose ();

  // Cycle the time levels by rotating the data sets
  void cycle (int rl, int c, int ml);
  
  // Flip the time levels by exchanging the data sets
  void flip (int rl, int c, int ml);
  
  // Copy data from current time level to all previous time levels
  void copytoprevs (int rl, int c, int ml);
  
  
  // Helpers
  
protected:
  
  virtual generic_data<D>* typed_data() = 0;
  
  
  
  // Operations
  
protected:
  
  // Copy a region
  void copycat (int tl1, int rl1, int c1, int ml1,
		const ibbox dh<D>::dboxes::* recv_list,
		int tl2, int rl2, int ml2,
		const ibbox dh<D>::dboxes::* send_list);

  // Copy regions
  void copycat (int tl1, int rl1, int c1, int ml1,
		const iblist dh<D>::dboxes::* recv_list,
		int tl2, int rl2, int ml2,
		const iblist dh<D>::dboxes::* send_list);

  // Copy regions
  void copycat (int tl1, int rl1, int c1, int ml1,
		const iblistvect dh<D>::dboxes::* recv_listvect,
		int tl2, int rl2, int ml2,
		const iblistvect dh<D>::dboxes::* send_listvect);
  
  // Interpolate a region
  void intercat (int tl1, int rl1, int c1, int ml1,
		 const ibbox dh<D>::dboxes::* recv_list,
		 const vector<int> tl2s, int rl2, int ml2,
		 const ibbox dh<D>::dboxes::* send_list,
		 CCTK_REAL time);

  // Interpolate regions
  void intercat (int tl1, int rl1, int c1, int ml1,
		 const iblist dh<D>::dboxes::* recv_list,
		 const vector<int> tl2s, int rl2, int ml2,
		 const iblist dh<D>::dboxes::* send_list,
		 CCTK_REAL time);

  // Interpolate regions
  void intercat (int tl1, int rl1, int c1, int ml1,
		 const iblistvect dh<D>::dboxes::* recv_listvect,
		 const vector<int> tl2s, int rl2, int ml2,
		 const iblistvect dh<D>::dboxes::* send_listvect,
		 CCTK_REAL time);



public:

  // The grid boundaries have to be updated after calling mg_restrict,
  // mg_prolongate, ref_restrict, or ref_prolongate.

  // "Updating" means here that the boundaries have to be
  // synchronised.  They don't need to be prolongated.

  // Copy a component from the next time level
  void copy (int tl, int rl, int c, int ml);

  // Synchronise the boundaries of a component
  void sync (int tl, int rl, int c, int ml);

  // Prolongate the boundaries of a component
  void ref_bnd_prolongate (int tl, int rl, int c, int ml, CCTK_REAL time);

  // Restrict a multigrid level
  void mg_restrict (int tl, int rl, int c, int ml, CCTK_REAL time);

  // Prolongate a multigrid level
  void mg_prolongate (int tl, int rl, int c, int ml, CCTK_REAL time);

  // Restrict a refinement level
  void ref_restrict (int tl, int rl, int c, int ml, CCTK_REAL time);

  // Prolongate a refinement level
  void ref_prolongate (int tl, int rl, int c, int ml, CCTK_REAL time);
  
  
  
  // Access to the data
  virtual const generic_data<D>* operator() (int tl, int rl, int c, int ml)
    const = 0;
  
  virtual generic_data<D>* operator() (int tl, int rl, int c, int ml) = 0;
  
  
  
  // Output
  virtual ostream& output (ostream& os) const = 0;
};



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



#endif // GGF_HH