aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/data.cc
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@aei.mpg.de>2005-06-07 16:08:00 +0000
committerErik Schnetter <schnetter@aei.mpg.de>2005-06-07 16:08:00 +0000
commitf05030a3f6a1470372ce29dba62cb2a167d1fd7a (patch)
treec0036aa9df6ed5927a80a6cbe4139ece9ae6ac67 /Carpet/CarpetLib/src/data.cc
parentbd0a846069dc0e91461bc20996358b3ff4ad0339 (diff)
CarpetLib: Add new transport operator type op_copy
Add a new transport operator type op_copy. This "prolongation" operator does not interpolate in time, but rather copies the from the newest time level instead. Such grid functions need only one time level. This is intended for prolongating or restricing grid functions like the ADM constraints; if done properly; they will have the same values on the coarse and fine grids. (However, this does not work for the ADM constraints, because such grid functions still need to be set in EVOL, not in ANALYSIS.) darcs-hash:20050607160833-891bb-cfd1c7630f8996606328d7c7e9fe326561106aba.gz
Diffstat (limited to 'Carpet/CarpetLib/src/data.cc')
-rw-r--r--Carpet/CarpetLib/src/data.cc77
1 files changed, 71 insertions, 6 deletions
diff --git a/Carpet/CarpetLib/src/data.cc b/Carpet/CarpetLib/src/data.cc
index c8dca0d70..1d7de42e0 100644
--- a/Carpet/CarpetLib/src/data.cc
+++ b/Carpet/CarpetLib/src/data.cc
@@ -992,6 +992,7 @@ void data<T>
switch (transport_operator) {
+ case op_copy:
case op_Lagrange:
case op_TVD:
case op_ENO:
@@ -1036,6 +1037,60 @@ void data<T>
box, sext, dext);
switch (transport_operator) {
+ case op_copy:
+ assert (times.size() == 1);
+ assert (srcs.size()>=1);
+ switch (order_space) {
+ case 0:
+ case 1:
+ if (all (sext.stride() == dext.stride() * 2)) {
+ CCTK_FNAME(prolongate_3d_real8_rf2)
+ ((const CCTK_REAL8*)srcs[0]->storage(),
+ srcshp[0], srcshp[1], srcshp[2],
+ (CCTK_REAL8*)storage(),
+ dstshp[0], dstshp[1], dstshp[2],
+ srcbbox, dstbbox, regbbox);
+ } else {
+ CCTK_FNAME(prolongate_3d_real8)
+ ((const CCTK_REAL8*)srcs[0]->storage(),
+ srcshp[0], srcshp[1], srcshp[2],
+ (CCTK_REAL8*)storage(),
+ dstshp[0], dstshp[1], dstshp[2],
+ srcbbox, dstbbox, regbbox);
+ }
+ break;
+ case 2:
+ case 3:
+ if (all (sext.stride() == dext.stride() * 2)) {
+ CCTK_FNAME(prolongate_3d_real8_o3_rf2)
+ ((const CCTK_REAL8*)srcs[0]->storage(),
+ srcshp[0], srcshp[1], srcshp[2],
+ (CCTK_REAL8*)storage(),
+ dstshp[0], dstshp[1], dstshp[2],
+ srcbbox, dstbbox, regbbox);
+ } else {
+ CCTK_FNAME(prolongate_3d_real8_o3)
+ ((const CCTK_REAL8*)srcs[0]->storage(),
+ srcshp[0], srcshp[1], srcshp[2],
+ (CCTK_REAL8*)storage(),
+ dstshp[0], dstshp[1], dstshp[2],
+ srcbbox, dstbbox, regbbox);
+ }
+ break;
+ case 4:
+ case 5:
+ CCTK_FNAME(prolongate_3d_real8_o5)
+ ((const CCTK_REAL8*)srcs[0]->storage(),
+ srcshp[0], srcshp[1], srcshp[2],
+ (CCTK_REAL8*)storage(),
+ dstshp[0], dstshp[1], dstshp[2],
+ srcbbox, dstbbox, regbbox);
+ break;
+ default:
+ assert (0);
+ }
+ break;
+
case op_Lagrange:
switch (order_time) {
@@ -1396,12 +1451,22 @@ void data<CCTK_REAL8>
min_time = min(min_time, times[tl]);
max_time = max(max_time, times[tl]);
}
- if (time < min_time - eps || time > max_time + eps) {
- ostringstream buf;
- buf << "Internal error: extrapolation in time."
- << " time=" << time
- << " times=" << times;
- CCTK_WARN (0, buf.str().c_str());
+ if (transport_operator != op_copy) {
+ if (time < min_time - eps || time > max_time + eps) {
+ ostringstream buf;
+ buf << "Internal error: extrapolation in time."
+ << " time=" << time
+ << " times=" << times;
+ CCTK_WARN (0, buf.str().c_str());
+ }
+ } else {
+ if (time > max_time + eps) {
+ ostringstream buf;
+ buf << "Internal error: extrapolation into the future."
+ << " time=" << time
+ << " times=" << times;
+ CCTK_WARN (0, buf.str().c_str());
+ }
}
}