aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/gdata.cc
blob: e438a314729884916afc02cdc5586d93bb79a8dd (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
// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/gdata.cc,v 1.23 2003/10/14 16:39:16 schnetter Exp $

#include <assert.h>

#include <iostream>

#include "cctk.h"

#include "bbox.hh"
#include "defs.hh"
#include "dist.hh"
#include "vect.hh"

#include "gdata.hh"

using namespace std;



// Constructors
template<int D>
gdata<D>::gdata (const int varindex_)
  : varindex(varindex_),
    _has_storage(false)
{ }

// Destructors
template<int D>
gdata<D>::~gdata () { }



// Data manipulators
template<int D>
void gdata<D>::copy_from (const gdata* src, const ibbox& box)
{
  assert (has_storage() && src->has_storage());
  assert (all(box.lower()>=extent().lower()
	      && box.lower()>=src->extent().lower()));
  assert (all(box.upper()<=extent().upper()
	      && box.upper()<=src->extent().upper()));
  assert (all(box.stride()==extent().stride()
	      && box.stride()==src->extent().stride()));
  assert (all((box.lower()-extent().lower())%box.stride() == 0
	      && (box.lower()-src->extent().lower())%box.stride() == 0));
  
  if (box.empty()) return;
  
  if (proc() == src->proc()) {
    // copy on same processor
    
    int rank;
    MPI_Comm_rank (dist::comm, &rank);
    if (rank == proc()) {
      copy_from_innerloop (src, box);
    }
    
  } else {
    
    // copy to different processor
    gdata* const tmp = make_typed(varindex);
    tmp->allocate (box, src->proc());
    tmp->copy_from (src, box);
    tmp->change_processor (proc());
    copy_from (tmp, box);
    delete tmp;
    
  }
}



template<int D>
void gdata<D>
::interpolate_from (const vector<const gdata*> srcs,
		    const vector<CCTK_REAL> times,
		    const ibbox& box, const CCTK_REAL time,
		    const int order_space,
		    const int order_time)
{
  assert (has_storage());
  assert (all(box.lower()>=extent().lower()));
  assert (all(box.upper()<=extent().upper()));
  assert (all(box.stride()==extent().stride()));
  assert (all((box.lower()-extent().lower())%box.stride() == 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()));
    assert (all(box.upper()<=srcs[t]->extent().upper()));
  }
  
  assert (! box.empty());
  if (box.empty()) return;
  
  if (proc() == srcs[0]->proc()) {
    // interpolate on same processor
    
    int rank;
    MPI_Comm_rank (dist::comm, &rank);
    if (rank == proc()) {
      interpolate_from_innerloop
	(srcs, times, box, time, order_space, order_time);
    }
    
  } else {
    // interpolate from other processor
    
    gdata* const tmp = make_typed(varindex);
    tmp->allocate (box, srcs[0]->proc());
    tmp->interpolate_from (srcs, times, box, time, order_space, order_time);
    tmp->change_processor (proc());
    copy_from (tmp, box);
    delete tmp;
    
  }
}



template class gdata<3>;