aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/dggf.hh
blob: 76727afee2919654e4b4260d30d2e6485795453a (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
/***************************************************************************
                          dggf.hh  -  Dimension 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/Attic/dggf.hh,v 1.2 2001/07/02 13:22:13 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 DGGF_HH
#define DGGF_HH

#include <assert.h>

#include <iostream>
#include <string>

#include "defs.hh"
#include "dgdata.hh"
#include "th.hh"

using namespace std;



// Forward declaration
class dimgeneric_gf;

// Output
ostream& operator<< (ostream& os, const dimgeneric_gf& f);



// A generic grid function without type information
class dimgeneric_gf {

public:				// should be readonly

  // Fields
  string name;

  th &t;			// time hierarchy
  int tmin, tmax;		// timelevels

public:

  // Constructors
  dimgeneric_gf (const string name, th& t, 
		 const int tmin, const int tmax);
  
  // Destructors
  virtual ~dimgeneric_gf ();

  // Comparison
  bool operator== (const dimgeneric_gf& f) const;
  
  
  
  // Modifiers
  virtual void recompose () = 0;
  
  // Cycle the time levels by rotating the data sets
  virtual void cycle (int rl, int c, int ml) = 0;
  
  
  
  // TODO:
  
  // are these necessary in dimgeneric_gf, or is it sufficient to have
  // them in generic_gf<D>?
  
  // do i really want all dims, or do i just want to disregard the z
  // dim in carpet?
  
  // i likely also have to make a dimgeneric_data class.
  
  // 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
  virtual void copy (int tl, int rl, int c, int ml) = 0;

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

  // Prolongate the boundaries of a component
  virtual void ref_bnd_prolongate (int tl, int rl, int c, int ml,
				   int order_space=1, int order_time=1) = 0;
  
  // Restrict a multigrid level
  virtual void mg_restrict (int tl, int rl, int c, int ml, int order_space=1)
    = 0;

  // Prolongate a multigrid level
  virtual void mg_prolongate (int tl, int rl, int c, int ml, int order_space=1)
    = 0;

  // Restrict a refinement level
  virtual void ref_restrict (int tl, int rl, int c, int ml, int order_space=1)
    = 0;

  // Prolongate a refinement level
  virtual void ref_prolongate (int tl, int rl, int c, int ml,
			       int order_space=1) = 0;
  
  
  
  // Access to the data
  virtual const dimgeneric_data* operator() (int tl, int rl, int c, int ml)
    const = 0;
  
  virtual dimgeneric_data* operator() (int tl, int rl, int c, int ml) = 0;
  
  
  
  // Output
  virtual ostream& output (ostream& os) const = 0;
};



inline ostream& operator<< (ostream& os, const dimgeneric_gf& f) {
  return f.output(os);
}



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

#endif // DGGF_HH