// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/ggf.hh,v 1.14 2003/05/02 14:23:12 schnetter Exp $ #ifndef GGF_HH #define GGF_HH #include #include #include #include #include "cctk.h" #include "defs.hh" #include "dh.hh" #include "gdata.hh" #include "gh.hh" #include "th.hh" using namespace std; // Forward declaration template class ggf; // Output template ostream& operator<< (ostream& os, const ggf& f); // A generic grid function without type information template class ggf { // Types typedef vect ivect; typedef bbox ibbox; typedef bboxset ibset; typedef list iblist; typedef vector iblistvect; typedef gdata* tdata; // data ... typedef vector mdata; // ... for each multigrid level typedef vector cdata; // ... for each component typedef vector rdata; // ... for each refinement level typedef vector fdata; // ... for each time level public: // should be readonly // Fields string name; th &t; // time hierarchy int tmin, tmax; // timelevels int prolongation_order_time; // order of temporal prolongation operator gh &h; // grid hierarchy dh &d; // data hierarchy protected: fdata storage; // storage public: // Constructors ggf (const string name, th& t, dh& d, const int tmin, const int tmax, const int prolongation_order_time); // Destructors virtual ~ggf (); // Comparison bool operator== (const ggf& f) const; // Modifiers void recompose (const int initialise_upto = -1); // 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 gdata* typed_data() = 0; // Operations protected: // Copy a region void copycat (int tl1, int rl1, int c1, int ml1, const ibbox dh::dboxes::* recv_list, int tl2, int rl2, int ml2, const ibbox dh::dboxes::* send_list); // Copy regions void copycat (int tl1, int rl1, int c1, int ml1, const iblist dh::dboxes::* recv_list, int tl2, int rl2, int ml2, const iblist dh::dboxes::* send_list); // Copy regions void copycat (int tl1, int rl1, int c1, int ml1, const iblistvect dh::dboxes::* recv_listvect, int tl2, int rl2, int ml2, const iblistvect dh::dboxes::* send_listvect); // Interpolate a region void intercat (int tl1, int rl1, int c1, int ml1, const ibbox dh::dboxes::* recv_list, const vector tl2s, int rl2, int ml2, const ibbox dh::dboxes::* send_list, CCTK_REAL time); // Interpolate regions void intercat (int tl1, int rl1, int c1, int ml1, const iblist dh::dboxes::* recv_list, const vector tl2s, int rl2, int ml2, const iblist dh::dboxes::* send_list, CCTK_REAL time); // Interpolate regions void intercat (int tl1, int rl1, int c1, int ml1, const iblistvect dh::dboxes::* recv_listvect, const vector tl2s, int rl2, int ml2, const iblistvect dh::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 gdata* operator() (int tl, int rl, int c, int ml) const = 0; virtual gdata* operator() (int tl, int rl, int c, int ml) = 0; // Output virtual ostream& output (ostream& os) const = 0; }; template inline ostream& operator<< (ostream& os, const ggf& f) { return f.output(os); } #endif // GGF_HH