aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib
diff options
context:
space:
mode:
authorschnetter <>2002-09-25 13:49:00 +0000
committerschnetter <>2002-09-25 13:49:00 +0000
commit2fa5508433b6112d92820d72b60eac704eb829fd (patch)
tree97ad3e4d47a332592047641058e1b86f8a921391 /Carpet/CarpetLib
parentb0fafc198e4a2ca84d29d2f8817726f73a4c6779 (diff)
Changed the variable types for the time from int to real. This allows
Changed the variable types for the time from int to real. This allows more flexibility for time interpolation. The time interpolators now accept an explicit argument for the current time, instead of using the current time from the time hierarchy. This allows arbitrary current times, which are necessary for the intermediate time steps for time interpolators. The times that are interpolated _from_ still come from the time hierarchy. darcs-hash:20020925134915-07bb3-26927a425ca5c0b52a2ca02bc54673c4b9e1781e.gz
Diffstat (limited to 'Carpet/CarpetLib')
-rw-r--r--Carpet/CarpetLib/src/data.cc72
-rw-r--r--Carpet/CarpetLib/src/data.hh8
-rw-r--r--Carpet/CarpetLib/src/dggf.cc10
-rw-r--r--Carpet/CarpetLib/src/dggf.hh19
-rw-r--r--Carpet/CarpetLib/src/gdata.cc15
-rw-r--r--Carpet/CarpetLib/src/gdata.hh12
-rw-r--r--Carpet/CarpetLib/src/ggf.cc91
-rw-r--r--Carpet/CarpetLib/src/ggf.hh23
-rw-r--r--Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F7715
-rw-r--r--Carpet/CarpetLib/src/prolongate_3d_real8_2tl_o3.F7715
-rw-r--r--Carpet/CarpetLib/src/prolongate_3d_real8_3tl.F7719
-rw-r--r--Carpet/CarpetLib/src/prolongate_3d_real8_3tl_o3.F7719
-rw-r--r--Carpet/CarpetLib/src/th.cc12
-rw-r--r--Carpet/CarpetLib/src/th.hh22
14 files changed, 202 insertions, 150 deletions
diff --git a/Carpet/CarpetLib/src/data.cc b/Carpet/CarpetLib/src/data.cc
index ee6a08c8e..edfb202a5 100644
--- a/Carpet/CarpetLib/src/data.cc
+++ b/Carpet/CarpetLib/src/data.cc
@@ -5,7 +5,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/data.cc,v 1.19 2002/08/30 16:03:19 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/data.cc,v 1.20 2002/09/25 15:49:15 schnetter Exp $
***************************************************************************/
@@ -197,8 +197,8 @@ void data<T,D>
template<class T, int D>
void data<T,D>
::interpolate_from_innerloop (const vector<const generic_data<D>*> gsrcs,
- const vector<int> tls,
- const ibbox& box, const int tl,
+ const vector<CCTK_REAL> times,
+ const ibbox& box, const CCTK_REAL time,
const int order_space,
const int order_time)
{
@@ -209,7 +209,7 @@ void data<T,D>
assert (all((box.lower()-extent().lower())%box.stride() == 0));
vector<const data*> srcs(gsrcs.size());
for (int t=0; t<(int)srcs.size(); ++t) srcs[t] = (const data*)gsrcs[t];
- assert (srcs.size() == tls.size() && srcs.size()>0);
+ assert (srcs.size() == times.size() && srcs.size()>0);
for (int t=0; t<(int)srcs.size(); ++t) {
assert (srcs[t]->has_storage());
assert (all(box.lower()>=srcs[t]->extent().lower()));
@@ -229,7 +229,7 @@ void data<T,D>
src_fac[t] = 1;
for (int tt=0; tt<(int)src_fac.size(); ++tt) {
if (tt!=t) {
- src_fac[t] *= (T)(tl - tls[tt]) / (T)(tls[t] - tls[tt]);
+ src_fac[t] *= (T)(time - times[tt]) / (T)(times[t] - times[tt]);
}
}
}
@@ -366,40 +366,40 @@ extern "C" {
const int regbbox[3][3]);
void CCTK_FCALL CCTK_FNAME(prolongate_3d_real8_2tl)
- (const CCTK_REAL8* src1, const int& t1,
- const CCTK_REAL8* src2, const int& t2,
+ (const CCTK_REAL8* src1, const CCTK_REAL8& t1,
+ const CCTK_REAL8* src2, const CCTK_REAL8& t2,
const int& srciext, const int& srcjext, const int& srckext,
- CCTK_REAL8* dst, const int& t,
+ CCTK_REAL8* dst, const CCTK_REAL8& t,
const int& dstiext, const int& dstjext, const int& dstkext,
const int srcbbox[3][3],
const int dstbbox[3][3],
const int regbbox[3][3]);
void CCTK_FCALL CCTK_FNAME(prolongate_3d_real8_2tl_o3)
- (const CCTK_REAL8* src1, const int& t1,
- const CCTK_REAL8* src2, const int& t2,
+ (const CCTK_REAL8* src1, const CCTK_REAL8& t1,
+ const CCTK_REAL8* src2, const CCTK_REAL8& t2,
const int& srciext, const int& srcjext, const int& srckext,
- CCTK_REAL8* dst, const int& t,
+ CCTK_REAL8* dst, const CCTK_REAL8& t,
const int& dstiext, const int& dstjext, const int& dstkext,
const int srcbbox[3][3],
const int dstbbox[3][3],
const int regbbox[3][3]);
void CCTK_FCALL CCTK_FNAME(prolongate_3d_real8_3tl)
- (const CCTK_REAL8* src1, const int& t1,
- const CCTK_REAL8* src2, const int& t2,
- const CCTK_REAL8* src3, const int& t3,
+ (const CCTK_REAL8* src1, const CCTK_REAL8& t1,
+ const CCTK_REAL8* src2, const CCTK_REAL8& t2,
+ const CCTK_REAL8* src3, const CCTK_REAL8& t3,
const int& srciext, const int& srcjext, const int& srckext,
- CCTK_REAL8* dst, const int& t,
+ CCTK_REAL8* dst, const CCTK_REAL8& t,
const int& dstiext, const int& dstjext, const int& dstkext,
const int srcbbox[3][3],
const int dstbbox[3][3],
const int regbbox[3][3]);
void CCTK_FCALL CCTK_FNAME(prolongate_3d_real8_3tl_o3)
- (const CCTK_REAL8* src1, const int& t1,
- const CCTK_REAL8* src2, const int& t2,
- const CCTK_REAL8* src3, const int& t3,
+ (const CCTK_REAL8* src1, const CCTK_REAL8& t1,
+ const CCTK_REAL8* src2, const CCTK_REAL8& t2,
+ const CCTK_REAL8* src3, const CCTK_REAL8& t3,
const int& srciext, const int& srcjext, const int& srckext,
- CCTK_REAL8* dst, const int& t,
+ CCTK_REAL8* dst, const CCTK_REAL8& t,
const int& dstiext, const int& dstjext, const int& dstkext,
const int srcbbox[3][3],
const int dstbbox[3][3],
@@ -410,8 +410,8 @@ extern "C" {
template<>
void data<CCTK_REAL8,3>
::interpolate_from_innerloop (const vector<const generic_data<3>*> gsrcs,
- const vector<int> tls,
- const ibbox& box, const int tl,
+ const vector<CCTK_REAL> times,
+ const ibbox& box, const CCTK_REAL time,
const int order_space,
const int order_time)
{
@@ -422,7 +422,7 @@ void data<CCTK_REAL8,3>
assert (all((box.lower()-extent().lower())%box.stride() == 0));
vector<const data*> srcs(gsrcs.size());
for (int t=0; t<(int)srcs.size(); ++t) srcs[t] = (const data*)gsrcs[t];
- assert (srcs.size() == tls.size() && srcs.size()>0);
+ assert (srcs.size() == times.size() && srcs.size()>0);
for (int t=0; t<(int)srcs.size(); ++t) {
assert (srcs[t]->has_storage());
assert (all(box.lower()>=srcs[t]->extent().lower()));
@@ -505,20 +505,20 @@ void data<CCTK_REAL8,3>
case 0:
case 1:
CCTK_FNAME(prolongate_3d_real8_2tl)
- ((const CCTK_REAL8*)srcs[0]->storage(), tls[0],
- (const CCTK_REAL8*)srcs[1]->storage(), tls[1],
+ ((const CCTK_REAL8*)srcs[0]->storage(), times[0],
+ (const CCTK_REAL8*)srcs[1]->storage(), times[1],
srcshp[0], srcshp[1], srcshp[2],
- (CCTK_REAL8*)storage(), tl,
+ (CCTK_REAL8*)storage(), time,
dstshp[0], dstshp[1], dstshp[2],
srcbbox, dstbbox, regbbox);
break;
case 2:
case 3:
CCTK_FNAME(prolongate_3d_real8_2tl_o3)
- ((const CCTK_REAL8*)srcs[0]->storage(), tls[0],
- (const CCTK_REAL8*)srcs[1]->storage(), tls[1],
+ ((const CCTK_REAL8*)srcs[0]->storage(), times[0],
+ (const CCTK_REAL8*)srcs[1]->storage(), times[1],
srcshp[0], srcshp[1], srcshp[2],
- (CCTK_REAL8*)storage(), tl,
+ (CCTK_REAL8*)storage(), time,
dstshp[0], dstshp[1], dstshp[2],
srcbbox, dstbbox, regbbox);
break;
@@ -533,22 +533,22 @@ void data<CCTK_REAL8,3>
case 0:
case 1:
CCTK_FNAME(prolongate_3d_real8_3tl)
- ((const CCTK_REAL8*)srcs[0]->storage(), tls[0],
- (const CCTK_REAL8*)srcs[1]->storage(), tls[1],
- (const CCTK_REAL8*)srcs[2]->storage(), tls[2],
+ ((const CCTK_REAL8*)srcs[0]->storage(), times[0],
+ (const CCTK_REAL8*)srcs[1]->storage(), times[1],
+ (const CCTK_REAL8*)srcs[2]->storage(), times[2],
srcshp[0], srcshp[1], srcshp[2],
- (CCTK_REAL8*)storage(), tl,
+ (CCTK_REAL8*)storage(), time,
dstshp[0], dstshp[1], dstshp[2],
srcbbox, dstbbox, regbbox);
break;
case 2:
case 3:
CCTK_FNAME(prolongate_3d_real8_3tl_o3)
- ((const CCTK_REAL8*)srcs[0]->storage(), tls[0],
- (const CCTK_REAL8*)srcs[1]->storage(), tls[1],
- (const CCTK_REAL8*)srcs[2]->storage(), tls[2],
+ ((const CCTK_REAL8*)srcs[0]->storage(), times[0],
+ (const CCTK_REAL8*)srcs[1]->storage(), times[1],
+ (const CCTK_REAL8*)srcs[2]->storage(), times[2],
srcshp[0], srcshp[1], srcshp[2],
- (CCTK_REAL8*)storage(), tl,
+ (CCTK_REAL8*)storage(), time,
dstshp[0], dstshp[1], dstshp[2],
srcbbox, dstbbox, regbbox);
break;
diff --git a/Carpet/CarpetLib/src/data.hh b/Carpet/CarpetLib/src/data.hh
index 983f7f141..468b7cced 100644
--- a/Carpet/CarpetLib/src/data.hh
+++ b/Carpet/CarpetLib/src/data.hh
@@ -5,7 +5,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/data.hh,v 1.10 2002/05/05 22:16:59 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/data.hh,v 1.11 2002/09/25 15:49:15 schnetter Exp $
***************************************************************************/
@@ -26,6 +26,8 @@
#include <iostream>
#include <string>
+#include "cctk.h"
+
#include "defs.hh"
#include "dist.hh"
#include "bbox.hh"
@@ -94,8 +96,8 @@ public:
void copy_from_innerloop (const generic_data<D>* gsrc,
const ibbox& box);
void interpolate_from_innerloop (const vector<const generic_data<D>*> gsrcs,
- const vector<int> tls,
- const ibbox& box, const int tl,
+ const vector<CCTK_REAL> times,
+ const ibbox& box, const CCTK_REAL time,
const int order_space,
const int order_time);
diff --git a/Carpet/CarpetLib/src/dggf.cc b/Carpet/CarpetLib/src/dggf.cc
index 0091fe800..c41e6a61a 100644
--- a/Carpet/CarpetLib/src/dggf.cc
+++ b/Carpet/CarpetLib/src/dggf.cc
@@ -6,7 +6,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/Attic/dggf.cc,v 1.3 2002/05/05 22:17:00 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/Attic/dggf.cc,v 1.4 2002/09/25 15:49:15 schnetter Exp $
***************************************************************************/
@@ -22,12 +22,8 @@
#include <assert.h>
#include <stdlib.h>
-// #include <iostream>
-// #include <string>
-
-// #include "defs.hh"
-// #include "dh.hh"
-// #include "th.hh"
+#include "defs.hh"
+#include "th.hh"
#include "dggf.hh"
diff --git a/Carpet/CarpetLib/src/dggf.hh b/Carpet/CarpetLib/src/dggf.hh
index 926d76f00..193903f70 100644
--- a/Carpet/CarpetLib/src/dggf.hh
+++ b/Carpet/CarpetLib/src/dggf.hh
@@ -6,7 +6,7 @@
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.4 2002/05/05 22:17:01 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/Attic/dggf.hh,v 1.5 2002/09/25 15:49:15 schnetter Exp $
***************************************************************************/
@@ -27,6 +27,8 @@
#include <iostream>
#include <string>
+#include "cctk.h"
+
#include "defs.hh"
#include "dgdata.hh"
#include "th.hh"
@@ -101,19 +103,24 @@ public:
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) = 0;
+ virtual void ref_bnd_prolongate (int tl, int rl, int c, int ml,
+ CCTK_REAL time) = 0;
// Restrict a multigrid level
- virtual void mg_restrict (int tl, int rl, int c, int ml) = 0;
+ virtual void mg_restrict (int tl, int rl, int c, int ml,
+ CCTK_REAL time) = 0;
// Prolongate a multigrid level
- virtual void mg_prolongate (int tl, int rl, int c, int ml) = 0;
+ virtual void mg_prolongate (int tl, int rl, int c, int ml,
+ CCTK_REAL time) = 0;
// Restrict a refinement level
- virtual void ref_restrict (int tl, int rl, int c, int ml) = 0;
+ virtual void ref_restrict (int tl, int rl, int c, int ml,
+ CCTK_REAL time) = 0;
// Prolongate a refinement level
- virtual void ref_prolongate (int tl, int rl, int c, int ml) = 0;
+ virtual void ref_prolongate (int tl, int rl, int c, int ml,
+ CCTK_REAL time) = 0;
diff --git a/Carpet/CarpetLib/src/gdata.cc b/Carpet/CarpetLib/src/gdata.cc
index 6e7df5472..76184b3d1 100644
--- a/Carpet/CarpetLib/src/gdata.cc
+++ b/Carpet/CarpetLib/src/gdata.cc
@@ -5,7 +5,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/gdata.cc,v 1.19 2002/05/05 22:17:01 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/gdata.cc,v 1.20 2002/09/25 15:49:16 schnetter Exp $
***************************************************************************/
@@ -22,6 +22,8 @@
#include <iostream>
+#include "cctk.h"
+
#include "bbox.hh"
#include "defs.hh"
#include "dist.hh"
@@ -85,8 +87,8 @@ void generic_data<D>::copy_from (const generic_data* src, const ibbox& box)
template<int D>
void generic_data<D>
::interpolate_from (const vector<const generic_data*> srcs,
- const vector<int> tls,
- const ibbox& box, const int tl,
+ const vector<CCTK_REAL> times,
+ const ibbox& box, const CCTK_REAL time,
const int order_space,
const int order_time)
{
@@ -95,7 +97,7 @@ void generic_data<D>
assert (all(box.upper()<=extent().upper()));
assert (all(box.stride()==extent().stride()));
assert (all((box.lower()-extent().lower())%box.stride() == 0));
- assert (srcs.size() == tls.size() && srcs.size()>0);
+ assert (srcs.size() == times.size() && srcs.size()>0);
for (int t=0; t<(int)srcs.size(); ++t) {
assert (srcs[t]->has_storage());
assert (all(box.lower()>=srcs[t]->extent().lower()));
@@ -108,7 +110,8 @@ void generic_data<D>
int rank;
MPI_Comm_rank (dist::comm, &rank);
if (rank == proc()) {
- interpolate_from_innerloop (srcs, tls, box, tl, order_space, order_time);
+ interpolate_from_innerloop
+ (srcs, times, box, time, order_space, order_time);
}
} else {
@@ -116,7 +119,7 @@ void generic_data<D>
generic_data* const tmp = make_typed();
tmp->allocate (box, srcs[0]->proc());
- tmp->interpolate_from (srcs, tls, box, tl, order_space, order_time);
+ tmp->interpolate_from (srcs, times, box, time, order_space, order_time);
tmp->change_processor (proc());
copy_from (tmp, box);
delete tmp;
diff --git a/Carpet/CarpetLib/src/gdata.hh b/Carpet/CarpetLib/src/gdata.hh
index 78ca1c935..66eeb00dd 100644
--- a/Carpet/CarpetLib/src/gdata.hh
+++ b/Carpet/CarpetLib/src/gdata.hh
@@ -5,7 +5,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/gdata.hh,v 1.14 2002/05/05 22:17:02 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/gdata.hh,v 1.15 2002/09/25 15:49:16 schnetter Exp $
***************************************************************************/
@@ -27,6 +27,8 @@
#include <iostream>
#include <string>
+#include "cctk.h"
+
#include "defs.hh"
#include "dgdata.hh"
#include "dist.hh"
@@ -99,8 +101,8 @@ public:
// Data manipulators
void copy_from (const generic_data* src, const ibbox& box);
void interpolate_from (const vector<const generic_data*> srcs,
- const vector<int> tls,
- const ibbox& box, const int tl,
+ const vector<CCTK_REAL> times,
+ const ibbox& box, const CCTK_REAL time,
const int order_space,
const int order_time);
protected:
@@ -108,8 +110,8 @@ protected:
copy_from_innerloop (const generic_data* src, const ibbox& box) = 0;
virtual void
interpolate_from_innerloop (const vector<const generic_data*> srcs,
- const vector<int> tls,
- const ibbox& box, const int tl,
+ const vector<CCTK_REAL> times,
+ const ibbox& box, const CCTK_REAL time,
const int order_space,
const int order_time) = 0;
diff --git a/Carpet/CarpetLib/src/ggf.cc b/Carpet/CarpetLib/src/ggf.cc
index fa24943c5..ae733aa35 100644
--- a/Carpet/CarpetLib/src/ggf.cc
+++ b/Carpet/CarpetLib/src/ggf.cc
@@ -6,7 +6,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/ggf.cc,v 1.18 2002/08/30 16:03:20 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/ggf.cc,v 1.19 2002/09/25 15:49:16 schnetter Exp $
***************************************************************************/
@@ -25,6 +25,8 @@
#include <iostream>
#include <string>
+#include "cctk.h"
+
#include "defs.hh"
#include "dh.hh"
#include "th.hh"
@@ -102,7 +104,8 @@ void generic_gf<D>::recompose () {
// Initialise from coarser level, if possible
// TODO: init only un-copied regions
if (rl>0) {
- ref_prolongate (tl,rl,c,ml);
+ const CCTK_REAL time = t.time(tl,rl,ml);
+ ref_prolongate (tl,rl,c,ml,time);
} // if rl
// Copy from old storage, if possible
@@ -144,7 +147,8 @@ void generic_gf<D>::recompose () {
sync (tl,rl,c,ml);
// TODO: assert that reflevel 0 boundaries are copied
if (rl>0) {
- ref_bnd_prolongate (tl,rl,c,ml);
+ const CCTK_REAL time = t.time(tl,rl,ml);
+ ref_bnd_prolongate (tl,rl,c,ml,time);
} // if rl
} // for ml
} // for c
@@ -189,9 +193,8 @@ void generic_gf<D>::copytoprevs (int rl, int c, int ml) {
assert (c>=0 && c<h.components(rl));
assert (ml>=0 && ml<h.mglevels(rl,c));
for (int tl=tmin; tl<=tmax-1; ++tl) {
- storage[tl-tmin][rl][c][ml]->copy_from
-(storage[tmax-tmin][rl][c][ml], storage[tmax-tmin][rl][c][ml]->extent
-());
+ storage[tl-tmin][rl][c][ml]->copy_from
+ (storage[tmax-tmin][rl][c][ml], storage[tmax-tmin][rl][c][ml]->extent());
}
}
@@ -288,7 +291,8 @@ template<int D>
void generic_gf<D>::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)
+ const ibbox dh<D>::dboxes::* send_list,
+ CCTK_REAL time)
{
assert (tl1>=tmin && tl1<=tmax);
assert (rl1>=0 && rl1<h.reflevels());
@@ -302,15 +306,15 @@ void generic_gf<D>::intercat (int tl1, int rl1, int c1, int ml1,
assert (ml2>=0 && ml2<h.mglevels(rl2,c2));
vector<const generic_data<D>*> gsrcs(tl2s.size());
- vector<int> tls(tl2s.size());
+ vector<CCTK_REAL> times(tl2s.size());
for (int i=0; i<(int)gsrcs.size(); ++i) {
assert (rl2<(int)storage[tl2s[i]-tmin].size());
assert (c2<(int)storage[tl2s[i]-tmin][rl2].size());
assert (ml2<(int)storage[tl2s[i]-tmin][rl2][ml2].size());
gsrcs[i] = storage[tl2s[i]-tmin][rl2][c2][ml2];
- tls[i] = t.time(tl2s[i],rl2,ml2);
+ times[i] = t.time(tl2s[i],rl2,ml2);
}
- const int tl = t.time(tl1,rl1,ml1);
+// const CCTK_REAL time = t.time(tl1,rl1,ml1);
const ibbox recv = d.boxes[rl1][c1][ml1].*recv_list;
const ibbox send = d.boxes[rl2][c2][ml2].*send_list;
@@ -318,7 +322,7 @@ void generic_gf<D>::intercat (int tl1, int rl1, int c1, int ml1,
// interpolate the content
assert (recv==send);
storage[tl1-tmin][rl1][c1][ml1]->interpolate_from
- (gsrcs, tls, recv, tl,
+ (gsrcs, times, recv, time,
d.prolongation_order_space, prolongation_order_time);
}
@@ -327,7 +331,8 @@ template<int D>
void generic_gf<D>::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)
+ const iblist dh<D>::dboxes::* send_list,
+ const CCTK_REAL time)
{
assert (tl1>=tmin && tl1<=tmax);
assert (rl1>=0 && rl1<h.reflevels());
@@ -341,15 +346,15 @@ void generic_gf<D>::intercat (int tl1, int rl1, int c1, int ml1,
assert (ml2>=0 && ml2<h.mglevels(rl2,c2));
vector<const generic_data<D>*> gsrcs(tl2s.size());
- vector<int> tls(tl2s.size());
+ vector<CCTK_REAL> times(tl2s.size());
for (int i=0; i<(int)gsrcs.size(); ++i) {
assert (rl2<(int)storage[tl2s[i]-tmin].size());
assert (c2<(int)storage[tl2s[i]-tmin][rl2].size());
assert (ml2<(int)storage[tl2s[i]-tmin][rl2][ml2].size());
gsrcs[i] = storage[tl2s[i]-tmin][rl2][c2][ml2];
- tls[i] = t.time(tl2s[i],rl2,ml2);
+ times[i] = t.time(tl2s[i],rl2,ml2);
}
- const int tl = t.time(tl1,rl1,ml1);
+// const CCTK_REAL time = t.time(tl1,rl1,ml1);
const iblist recv = d.boxes[rl1][c1][ml1].*recv_list;
const iblist send = d.boxes[rl2][c2][ml2].*send_list;
@@ -360,7 +365,7 @@ void generic_gf<D>::intercat (int tl1, int rl1, int c1, int ml1,
// (use the send boxes for communication)
// interpolate the content
storage[tl1-tmin][rl1][c1][ml1]->interpolate_from
- (gsrcs, tls, *r, tl,
+ (gsrcs, times, *r, time,
d.prolongation_order_space, prolongation_order_time);
}
}
@@ -370,7 +375,8 @@ template<int D>
void generic_gf<D>::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)
+ const iblistvect dh<D>::dboxes::* send_listvect,
+ const CCTK_REAL time)
{
assert (tl1>=tmin && tl1<=tmax);
assert (rl1>=0 && rl1<h.reflevels());
@@ -385,15 +391,15 @@ void generic_gf<D>::intercat (int tl1, int rl1, int c1, int ml1,
assert (ml2>=0 && ml2<h.mglevels(rl2,c2));
vector<const generic_data<D>*> gsrcs(tl2s.size());
- vector<int> tls(tl2s.size());
+ vector<CCTK_REAL> times(tl2s.size());
for (int i=0; i<(int)gsrcs.size(); ++i) {
assert (rl2<(int)storage[tl2s[i]-tmin].size());
assert (c2<(int)storage[tl2s[i]-tmin][rl2].size());
assert (ml2<(int)storage[tl2s[i]-tmin][rl2][ml2].size());
gsrcs[i] = storage[tl2s[i]-tmin][rl2][c2][ml2];
- tls[i] = t.time(tl2s[i],rl2,ml2);
+ times[i] = t.time(tl2s[i],rl2,ml2);
}
- const int tl = t.time(tl1,rl1,ml1);
+// const CCTK_REAL time = t.time(tl1,rl1,ml1);
const iblist recv = (d.boxes[rl1][c1][ml1].*recv_listvect)[c2];
const iblist send = (d.boxes[rl2][c2][ml2].*send_listvect)[c1];
@@ -404,7 +410,7 @@ void generic_gf<D>::intercat (int tl1, int rl1, int c1, int ml1,
// (use the send boxes for communication)
// interpolate the content
storage[tl1-tmin][rl1][c1][ml1]->interpolate_from
- (gsrcs, tls, *r, tl,
+ (gsrcs, times, *r, time,
d.prolongation_order_space, prolongation_order_time);
}
}
@@ -414,7 +420,8 @@ void generic_gf<D>::intercat (int tl1, int rl1, int c1, int ml1,
// Copy a component from the next time level
template<int D>
-void generic_gf<D>::copy (int tl, int rl, int c, int ml) {
+void generic_gf<D>::copy (int tl, int rl, int c, int ml)
+{
// Copy
copycat (tl ,rl,c,ml, &dh<D>::dboxes::exterior,
tl+1,rl, ml, &dh<D>::dboxes::exterior);
@@ -422,7 +429,8 @@ void generic_gf<D>::copy (int tl, int rl, int c, int ml) {
// Synchronise the boundaries a component
template<int D>
-void generic_gf<D>::sync (int tl, int rl, int c, int ml) {
+void generic_gf<D>::sync (int tl, int rl, int c, int ml)
+{
// Copy
copycat (tl,rl,c,ml, &dh<D>::dboxes::recv_sync,
tl,rl, ml, &dh<D>::dboxes::send_sync);
@@ -430,7 +438,9 @@ void generic_gf<D>::sync (int tl, int rl, int c, int ml) {
// Prolongate the boundaries of a component
template<int D>
-void generic_gf<D>::ref_bnd_prolongate (int tl, int rl, int c, int ml) {
+void generic_gf<D>::ref_bnd_prolongate (int tl, int rl, int c, int ml,
+ CCTK_REAL time)
+{
// Interpolate
assert (rl>=1);
vector<int> tl2s;
@@ -439,49 +449,57 @@ void generic_gf<D>::ref_bnd_prolongate (int tl, int rl, int c, int ml) {
tl2s.resize(prolongation_order_time+1);
for (int i=0; i<=prolongation_order_time; ++i) tl2s[i] = tmax - i;
intercat (tl ,rl ,c,ml, &dh<D>::dboxes::recv_ref_bnd_coarse,
- tl2s,rl-1, ml, &dh<D>::dboxes::send_ref_bnd_fine);
+ tl2s,rl-1, ml, &dh<D>::dboxes::send_ref_bnd_fine,
+ time);
}
// Restrict a multigrid level
template<int D>
-void generic_gf<D>::mg_restrict (int tl, int rl, int c, int ml)
+void generic_gf<D>::mg_restrict (int tl, int rl, int c, int ml,
+ CCTK_REAL time)
{
// Require same times
assert (t.get_time(rl,ml) == t.get_time(rl,ml-1));
const vector<int> tl2s(1,tl);
intercat (tl ,rl,c,ml, &dh<D>::dboxes::recv_mg_coarse,
- tl2s,rl, ml-1, &dh<D>::dboxes::send_mg_fine);
+ tl2s,rl, ml-1, &dh<D>::dboxes::send_mg_fine,
+ time);
}
// Prolongate a multigrid level
template<int D>
-void generic_gf<D>::mg_prolongate (int tl, int rl, int c, int ml)
+void generic_gf<D>::mg_prolongate (int tl, int rl, int c, int ml,
+ CCTK_REAL time)
{
// Require same times
assert (t.get_time(rl,ml) == t.get_time(rl,ml+1));
const vector<int> tl2s(1,tl);
intercat (tl ,rl,c,ml, &dh<D>::dboxes::recv_mg_coarse,
- tl2s,rl, ml+1, &dh<D>::dboxes::send_mg_fine);
+ tl2s,rl, ml+1, &dh<D>::dboxes::send_mg_fine,
+ time);
}
// Restrict a refinement level
template<int D>
-void generic_gf<D>::ref_restrict (int tl, int rl, int c, int ml)
+void generic_gf<D>::ref_restrict (int tl, int rl, int c, int ml,
+ CCTK_REAL time)
{
// Require same times
// SHH: removed assert and added warning
- if (t.get_time(rl,ml) != t.get_time(rl+1,ml)) {
- printf("WARNING: Time on rl %d is %d, time on rl %d is %d\n",rl,t.get_time(rl,ml),rl+1,t.get_time(rl+1,ml));
- }
+// if (t.get_time(rl,ml) != t.get_time(rl+1,ml)) {
+// printf("WARNING: Time on rl %d is %d, time on rl %d is %d\n",rl,t.get_time(rl,ml),rl+1,t.get_time(rl+1,ml));
+// }
//assert (t.get_time(rl,ml) == t.get_time(rl+1,ml));
const vector<int> tl2s(1,tl);
intercat (tl ,rl ,c,ml, &dh<D>::dboxes::recv_ref_fine,
- tl2s,rl+1, ml, &dh<D>::dboxes::send_ref_coarse);
+ tl2s,rl+1, ml, &dh<D>::dboxes::send_ref_coarse,
+ time);
}
// Prolongate a refinement level
template<int D>
-void generic_gf<D>::ref_prolongate (int tl, int rl, int c, int ml)
+void generic_gf<D>::ref_prolongate (int tl, int rl, int c, int ml,
+ CCTK_REAL time)
{
assert (rl>=1);
vector<int> tl2s;
@@ -490,7 +508,8 @@ void generic_gf<D>::ref_prolongate (int tl, int rl, int c, int ml)
tl2s.resize(prolongation_order_time+1);
for (int i=0; i<=prolongation_order_time; ++i) tl2s[i] = tmax - i;
intercat (tl ,rl ,c,ml, &dh<D>::dboxes::recv_ref_coarse,
- tl2s,rl-1, ml, &dh<D>::dboxes::send_ref_fine);
+ tl2s,rl-1, ml, &dh<D>::dboxes::send_ref_fine,
+ time);
}
diff --git a/Carpet/CarpetLib/src/ggf.hh b/Carpet/CarpetLib/src/ggf.hh
index a49bb1bad..66012e469 100644
--- a/Carpet/CarpetLib/src/ggf.hh
+++ b/Carpet/CarpetLib/src/ggf.hh
@@ -6,7 +6,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/ggf.hh,v 1.11 2002/07/18 14:33:50 shawley Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/ggf.hh,v 1.12 2002/09/25 15:49:16 schnetter Exp $
***************************************************************************/
@@ -28,6 +28,8 @@
#include <string>
#include <vector>
+#include "cctk.h"
+
#include "defs.hh"
#include "dggf.hh"
#include "dh.hh"
@@ -138,19 +140,22 @@ protected:
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);
+ 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);
+ 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);
+ const iblistvect dh<D>::dboxes::* send_listvect,
+ CCTK_REAL time);
@@ -169,19 +174,19 @@ public:
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);
+ 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);
+ 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);
+ 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);
+ 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);
+ void ref_prolongate (int tl, int rl, int c, int ml, CCTK_REAL time);
diff --git a/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77 b/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
index 6cb18abbb..46979c219 100644
--- a/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
+++ b/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
@@ -1,5 +1,5 @@
c -*-Fortran-*-
-c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77,v 1.5 2002/01/09 17:45:41 schnetter Exp $
+c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77,v 1.6 2002/09/25 15:49:16 schnetter Exp $
#include "cctk.h"
@@ -15,12 +15,12 @@ c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d
integer srciext, srcjext, srckext
CCTK_REAL8 src1(srciext,srcjext,srckext)
- integer t1
+ CCTK_REAL8 t1
CCTK_REAL8 src2(srciext,srcjext,srckext)
- integer t2
+ CCTK_REAL8 t2
integer dstiext, dstjext, dstkext
CCTK_REAL8 dst(dstiext,dstjext,dstkext)
- integer t
+ CCTK_REAL8 t
c bbox(:,1) is lower boundary (inclusive)
c bbox(:,2) is upper boundary (inclusive)
c bbox(:,3) is stride
@@ -105,9 +105,12 @@ c Linear (first order) interpolation
if (t1.eq.t2) then
call CCTK_WARN (0, "Internal error: arrays have same time")
end if
+ if (t.lt.min(t1,t2) .or. t.gt.max(t1,t2)) then
+ call CCTK_WARN (0, "Internal error: extrapolation")
+ end if
- s1fac = (t - t2) * one / (t1 - t2)
- s2fac = (t - t1) * one / (t2 - t1)
+ s1fac = (t - t2) / (t1 - t2)
+ s2fac = (t - t1) / (t2 - t1)
diff --git a/Carpet/CarpetLib/src/prolongate_3d_real8_2tl_o3.F77 b/Carpet/CarpetLib/src/prolongate_3d_real8_2tl_o3.F77
index 00a286952..145b7aeb9 100644
--- a/Carpet/CarpetLib/src/prolongate_3d_real8_2tl_o3.F77
+++ b/Carpet/CarpetLib/src/prolongate_3d_real8_2tl_o3.F77
@@ -1,5 +1,5 @@
c -*-Fortran-*-
-c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_2tl_o3.F77,v 1.7 2002/01/09 17:45:41 schnetter Exp $
+c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_2tl_o3.F77,v 1.8 2002/09/25 15:49:16 schnetter Exp $
#include "cctk.h"
@@ -15,12 +15,12 @@ c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d
integer srciext, srcjext, srckext
CCTK_REAL8 src1(srciext,srcjext,srckext)
- integer t1
+ CCTK_REAL8 t1
CCTK_REAL8 src2(srciext,srcjext,srckext)
- integer t2
+ CCTK_REAL8 t2
integer dstiext, dstjext, dstkext
CCTK_REAL8 dst(dstiext,dstjext,dstkext)
- integer t
+ CCTK_REAL8 t
c bbox(:,1) is lower boundary (inclusive)
c bbox(:,2) is upper boundary (inclusive)
c bbox(:,3) is stride
@@ -114,9 +114,12 @@ c Linear (first order) interpolation
if (t1.eq.t2) then
call CCTK_WARN (0, "Internal error: arrays have same time")
end if
+ if (t.lt.min(t1,t2) .or. t.gt.max(t1,t2)) then
+ call CCTK_WARN (0, "Internal error: extrapolation")
+ end if
- s1fac = (t - t2) * one / (t1 - t2)
- s2fac = (t - t1) * one / (t2 - t1)
+ s1fac = (t - t2) / (t1 - t2)
+ s2fac = (t - t1) / (t2 - t1)
diff --git a/Carpet/CarpetLib/src/prolongate_3d_real8_3tl.F77 b/Carpet/CarpetLib/src/prolongate_3d_real8_3tl.F77
index f9c451070..6a41b64c4 100644
--- a/Carpet/CarpetLib/src/prolongate_3d_real8_3tl.F77
+++ b/Carpet/CarpetLib/src/prolongate_3d_real8_3tl.F77
@@ -1,5 +1,5 @@
c -*-Fortran-*-
-c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_3tl.F77,v 1.3 2002/01/09 17:45:41 schnetter Exp $
+c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_3tl.F77,v 1.4 2002/09/25 15:49:16 schnetter Exp $
#include "cctk.h"
@@ -15,14 +15,14 @@ c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d
integer srciext, srcjext, srckext
CCTK_REAL8 src1(srciext,srcjext,srckext)
- integer t1
+ CCTK_REAL8 t1
CCTK_REAL8 src2(srciext,srcjext,srckext)
- integer t2
+ CCTK_REAL8 t2
CCTK_REAL8 src3(srciext,srcjext,srckext)
- integer t3
+ CCTK_REAL8 t3
integer dstiext, dstjext, dstkext
CCTK_REAL8 dst(dstiext,dstjext,dstkext)
- integer t
+ CCTK_REAL8 t
c bbox(:,1) is lower boundary (inclusive)
c bbox(:,2) is upper boundary (inclusive)
c bbox(:,3) is stride
@@ -107,10 +107,13 @@ c Quadratic (second order) interpolation
if (t1.eq.t2 .or. t1.eq.t3 .or. t2.eq.t3) then
call CCTK_WARN (0, "Internal error: arrays have same time")
end if
+ if (t.lt.min(t1,t2,t3) .or. t.gt.max(t1,t2,t3)) then
+ call CCTK_WARN (0, "Internal error: extrapolation")
+ end if
- s1fac = (t - t2) * (t - t3) * one / ((t1 - t2) * (t1 - t3))
- s2fac = (t - t1) * (t - t3) * one / ((t2 - t1) * (t2 - t3))
- s3fac = (t - t1) * (t - t2) * one / ((t3 - t1) * (t3 - t2))
+ s1fac = (t - t2) * (t - t3) / ((t1 - t2) * (t1 - t3))
+ s2fac = (t - t1) * (t - t3) / ((t2 - t1) * (t2 - t3))
+ s3fac = (t - t1) * (t - t2) / ((t3 - t1) * (t3 - t2))
diff --git a/Carpet/CarpetLib/src/prolongate_3d_real8_3tl_o3.F77 b/Carpet/CarpetLib/src/prolongate_3d_real8_3tl_o3.F77
index 2ce036fea..b36f51439 100644
--- a/Carpet/CarpetLib/src/prolongate_3d_real8_3tl_o3.F77
+++ b/Carpet/CarpetLib/src/prolongate_3d_real8_3tl_o3.F77
@@ -1,5 +1,5 @@
c -*-Fortran-*-
-c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_3tl_o3.F77,v 1.7 2002/01/09 17:45:41 schnetter Exp $
+c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d_real8_3tl_o3.F77,v 1.8 2002/09/25 15:49:16 schnetter Exp $
#include "cctk.h"
@@ -15,14 +15,14 @@ c $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/prolongate_3d
integer srciext, srcjext, srckext
CCTK_REAL8 src1(srciext,srcjext,srckext)
- integer t1
+ CCTK_REAL8 t1
CCTK_REAL8 src2(srciext,srcjext,srckext)
- integer t2
+ CCTK_REAL8 t2
CCTK_REAL8 src3(srciext,srcjext,srckext)
- integer t3
+ CCTK_REAL8 t3
integer dstiext, dstjext, dstkext
CCTK_REAL8 dst(dstiext,dstjext,dstkext)
- integer t
+ CCTK_REAL8 t
c bbox(:,1) is lower boundary (inclusive)
c bbox(:,2) is upper boundary (inclusive)
c bbox(:,3) is stride
@@ -116,10 +116,13 @@ c Quadratic (second order) interpolation
if (t1.eq.t2 .or. t1.eq.t3 .or. t2.eq.t3) then
call CCTK_WARN (0, "Internal error: arrays have same time")
end if
+ if (t.lt.min(t1,t2,t3) .or. t.gt.max(t1,t2,t3)) then
+ call CCTK_WARN (0, "Internal error: extrapolation")
+ end if
- s1fac = (t - t2) * (t - t3) * one / ((t1 - t2) * (t1 - t3))
- s2fac = (t - t1) * (t - t3) * one / ((t2 - t1) * (t2 - t3))
- s3fac = (t - t1) * (t - t2) * one / ((t3 - t1) * (t3 - t2))
+ s1fac = (t - t2) * (t - t3) / ((t1 - t2) * (t1 - t3))
+ s2fac = (t - t1) * (t - t3) / ((t2 - t1) * (t2 - t3))
+ s3fac = (t - t1) * (t - t2) / ((t3 - t1) * (t3 - t2))
diff --git a/Carpet/CarpetLib/src/th.cc b/Carpet/CarpetLib/src/th.cc
index 7d22c8829..7dc292585 100644
--- a/Carpet/CarpetLib/src/th.cc
+++ b/Carpet/CarpetLib/src/th.cc
@@ -6,7 +6,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/th.cc,v 1.7 2002/05/05 22:17:03 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/th.cc,v 1.8 2002/09/25 15:49:17 schnetter Exp $
***************************************************************************/
@@ -20,9 +20,12 @@
***************************************************************************/
#include <assert.h>
+#include <math.h>
#include <iostream>
+#include "cctk.h"
+
#include "defs.hh"
#include "dggh.hh"
@@ -33,7 +36,7 @@ using namespace std;
// Constructors
-th::th (dimgeneric_gh* h, const int basedelta)
+th::th (dimgeneric_gh* h, const CCTK_REAL basedelta)
: h(h), delta(basedelta) {
h->add(this);
}
@@ -49,7 +52,7 @@ void th::recompose () {
deltas.resize(h->reflevels());
for (int rl=0; rl<h->reflevels(); ++rl) {
const int old_mglevels = times[rl].size();
- int mgtime;
+ CCTK_REAL mgtime;
// Select default time
if (old_mglevels==0 && rl==0) {
mgtime = 0;
@@ -64,7 +67,8 @@ void th::recompose () {
if (rl==0 && ml==0) {
deltas[rl][ml] = delta;
} else if (ml==0) {
- assert (deltas[rl-1][ml] % h->reffact == 0);
+// assert (deltas[rl-1][ml] % h->reffact == 0);
+ assert (fabs(fmod(deltas[rl-1][ml], h->reffact)) < 1e-10);
deltas[rl][ml] = deltas[rl-1][ml] / h->reffact;
} else {
deltas[rl][ml] = deltas[rl][ml-1] * h->mgfact;
diff --git a/Carpet/CarpetLib/src/th.hh b/Carpet/CarpetLib/src/th.hh
index 0affc5e5f..61d489d6a 100644
--- a/Carpet/CarpetLib/src/th.hh
+++ b/Carpet/CarpetLib/src/th.hh
@@ -6,7 +6,7 @@
copyright : (C) 2000 by Erik Schnetter
email : schnetter@astro.psu.edu
- $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/th.hh,v 1.7 2002/06/07 20:20:43 schnetter Exp $
+ $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/th.hh,v 1.8 2002/09/25 15:49:17 schnetter Exp $
***************************************************************************/
@@ -27,6 +27,8 @@
#include <iostream>
#include <vector>
+#include "cctk.h"
+
#include "defs.hh"
#include "dggh.hh"
@@ -52,14 +54,14 @@ public: // should be readonly
private:
- int delta; // time step
- vector<vector<int> > times; // current times
- vector<vector<int> > deltas; // time steps
+ CCTK_REAL delta; // time step
+ vector<vector<CCTK_REAL> > times; // current times
+ vector<vector<CCTK_REAL> > deltas; // time steps
public:
// Constructors
- th (dimgeneric_gh* h, const int basedelta);
+ th (dimgeneric_gh* h, const CCTK_REAL basedelta);
// Destructors
~th ();
@@ -68,13 +70,13 @@ public:
void recompose ();
// Time management
- int get_time (const int rl, const int ml) const {
+ CCTK_REAL get_time (const int rl, const int ml) const {
assert (rl>=0 && rl<h->reflevels());
assert (ml>=0 && ml<h->mglevels(rl,0));
return times[rl][ml];
}
- void set_time (const int rl, const int ml, const int t) {
+ void set_time (const int rl, const int ml, const CCTK_REAL t) {
assert (rl>=0 && rl<h->reflevels());
assert (ml>=0 && ml<h->mglevels(rl,0));
times[rl][ml] = t;
@@ -84,19 +86,19 @@ public:
set_time(rl,ml, get_time(rl,ml) + get_delta(rl,ml));
}
- int get_delta (const int rl, const int ml) const {
+ CCTK_REAL get_delta (const int rl, const int ml) const {
assert (rl>=0 && rl<h->reflevels());
assert (ml>=0 && ml<h->mglevels(rl,0));
return deltas[rl][ml];
}
- void set_delta (const int rl, const int ml, const int dt) {
+ void set_delta (const int rl, const int ml, const CCTK_REAL dt) {
assert (rl>=0 && rl<h->reflevels());
assert (ml>=0 && ml<h->mglevels(rl,0));
deltas[rl][ml] = dt;
}
- int time (const int tl, const int rl, const int ml) const {
+ CCTK_REAL time (const int tl, const int rl, const int ml) const {
assert (rl>=0 && rl<h->reflevels());
assert (ml>=0 && ml<h->mglevels(rl,0));
return get_time(rl, ml) + tl * get_delta(rl, ml);