interpol.h 8.3 KB
Newer Older
Phillip Lino Rall's avatar
Phillip Lino Rall committed
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
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
/**********************************************************************************
 * Copyright 2010 Christoph Pflaum 
 * 		Department Informatik Lehrstuhl 10 - Systemsimulation
 *		Friedrich-Alexander Universität Erlangen-Nürnberg
 * 
 * Licensed under the Apache License, Version 2.0 (the "License");
 * you may not use this file except in compliance with the License.
 * You may obtain a copy of the License at
 *
 * http://www.apache.org/licenses/LICENSE-2.0
 *
 * Unless required by applicable law or agreed to in writing, software
 * distributed under the License is distributed on an "AS IS" BASIS,
 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 * See the License for the specific language governing permissions and
 * limitations under the License.
 **********************************************************************************/

// ------------------------------------------------------------
//
// interpolate.h
//
// ------------------------------------------------------------

#ifndef INTERPOL_H
#define INTERPOL_H

/////////////////////////////////////////////////////////////
// 1. Interpolate from  blockgrid to rectangular blockgrid
// 2. Interpolate from  blockgrid  to  blockgrid
/////////////////////////////////////////////////////////////


/////////////////////////////////////////////////////////////
// 1. Interpolate from  blockgrid to rectangular blockgrid
/////////////////////////////////////////////////////////////


/*
data structure on structured grid
ny-1  *  *  *  *  *  *
      *  +  +  +  +  *
3     *  +  +  +  +  *
2     *  +  +  +  +  *
1     *  +  +  +  +  *
0     *  *  *  *  *  *
      0  1  2  3     nx-1
*/



/**
 * interpolates data from blockgrid_ to rectangular block [pWSD, pENT]
**/
class Interpolate_on_structured_grid {
 public:
/**
 * preparation for interpolation
 @param  nx_ number of structured grid points x-direction
 @param  ny_ number of structured grid points y-direction
 @param  nz_ number of structured grid points z-direction
 @param  pWSD Corner WSD of structured grid
 @param  pENT Corner WSD of structured grid
 @param  blockgrid_ of unstructured grid
**/
  Interpolate_on_structured_grid(int nx_, int ny_, int nz_,
				 D3vector pWSD, D3vector pENT,
				 Blockgrid& blockgrid_);
/**
 * interpolates from Variable u(of unstructured blockgrid) to structured data
 @param  u:  Variable on unstructured blockgrid
 @param  data: pointer to structured grid data i+nx*(j+ny*k)
**/
  template <class DTyp>
    void interpolate(Variable<DTyp>& u, DTyp* data, DTyp defaultInterpolation);

/**
 * interpolates from Variable u(of unstructured blockgrid) to structured data
 @param  u:  Variable on unstructured blockgrid
 @param  data: pointer to structured grid data i+nx*(j+ny*k)
**/
  template <class DTyp>
    void interpolate(Cell_variable<DTyp>& u, DTyp* data, DTyp defaultInterpolation);    
    
    
  ~Interpolate_on_structured_grid();
private:
  int* ids_hex;
  int *ids_i, *ids_j, *ids_k;
  int *typ_tet;

  D3vector *lambda;

  int nx, ny, nz;
  double hx, hy, hz;

  Blockgrid*        blockgrid;
  Unstructured_grid*       ug;
};

template <class DTyp>
class Give_corner_data_of_cube {
 public:
  Give_corner_data_of_cube(Variable<DTyp>& u_, int hex_, int i_, int j_, int k_);
  
  DTyp WSD() { return u->template Give_data<hexahedronEl>(id, i,   j,   k,   Nx, Ny); }
  DTyp ESD() { return u->template Give_data<hexahedronEl>(id, i+1, j,   k,   Nx, Ny); }
  DTyp WND() { return u->template Give_data<hexahedronEl>(id, i,   j+1, k,   Nx, Ny); }
  DTyp END() { return u->template Give_data<hexahedronEl>(id, i+1, j+1, k,   Nx, Ny); }
  DTyp WST() { return u->template Give_data<hexahedronEl>(id, i,   j,   k+1, Nx, Ny); }
  DTyp EST() { return u->template Give_data<hexahedronEl>(id, i+1, j,   k+1, Nx, Ny); }
  DTyp WNT() { return u->template Give_data<hexahedronEl>(id, i,   j+1, k+1, Nx, Ny); }
  DTyp ENT() { return u->template Give_data<hexahedronEl>(id, i+1, j+1, k+1, Nx, Ny); }

 private:
  int id, i, j, k;
  int Nx, Ny;
  Variable<DTyp>* u;
};

template <class DTyp>
Give_corner_data_of_cube<DTyp>::Give_corner_data_of_cube(Variable<DTyp>& u_,
							 int hex_, int i_, int j_, int k_) {
  Blockgrid*        blockgrid;

  u = &u_;
  blockgrid = u->Give_blockgrid();
  id = hex_; i = i_;  j = j_;  k = k_;
  
  Nx = blockgrid->Give_Nx_hexahedron(id);
  Ny = blockgrid->Give_Ny_hexahedron(id);
}


template <class DTyp>
void Interpolate_on_structured_grid::interpolate(Variable<DTyp>& u, DTyp* data, DTyp defaultInterpolation) {
  int i,j,k, id_hex, typ;
  int ind_global;


  for(int id=0;id<ug->Give_number_hexahedra();++id)
    u.template Update<hexahedronEl>(id);

  for(int ks = 0; ks < nz;++ks)
    for(int js = 0; js < ny;++js)
      for(int is = 0; is < nx;++is) {
	ind_global = is+nx*(js+ny*ks);
	i = ids_i[ind_global];
	j = ids_j[ind_global];
	k = ids_k[ind_global];
	id_hex = ids_hex[ind_global];

	/*
	// test : das geht nicht::
	if(id_hex < 0) {
	  ind_global = is+nx*((ny-1-js)+ny*ks);
	  i = ids_i[ind_global];
	  j = ids_j[ind_global];
	  k = ids_k[ind_global];
	  id_hex = ids_hex[ind_global];
	}
*/		
	if(id_hex < 0) data[ind_global] = defaultInterpolation;
	else {
	  Give_corner_data_of_cube<DTyp> du(u, id_hex, i,j,k);

	  typ = typ_tet[ind_global];

	
	  /*
	  if( ind_global == 3)
	    cout << " \n ind_global: " << is + nx * (js + ny * ks)
	         << " typ: " << typ
	         << " id_hex " << id_hex
	         << " i: " << i
	         << " j: " << j
	         << " k: " << k
	         << endl;
  	  */

	  if(typ==0) data[ind_global] = interpolate_in_tet(lambda[ind_global],
							   du.WND(),du.WNT(),du.WST(),du.EST());
	  if(typ==1) data[ind_global] = interpolate_in_tet(lambda[ind_global],
							   du.EST(),du.WND(),du.WST(),du.ESD());
	  if(typ==2) data[ind_global] = interpolate_in_tet(lambda[ind_global],
							   du.WND(),du.WSD(),du.WST(),du.ESD());
	  if(typ==3) data[ind_global] = interpolate_in_tet(lambda[ind_global],
							   du.EST(),du.WND(),du.ESD(),du.END());
	  if(typ==4) data[ind_global] = interpolate_in_tet(lambda[ind_global],
							   du.ENT(),du.WNT(),du.EST(),du.END());
	  if(typ==5) data[ind_global] = interpolate_in_tet(lambda[ind_global],
							   du.WNT(),du.WND(),du.EST(),du.END());
	}
      }
}



template <class DTyp>
void Interpolate_on_structured_grid::interpolate(Cell_variable<DTyp>& u, DTyp* data, DTyp defaultInterpolation) {
  int i,j,k, id_hex, typ;
  int ind_global;


  for(int id=0;id<ug->Give_number_hexahedra();++id)
    u.template Update<hexahedronEl>(id);

  for(int ks = 0; ks < nz;++ks)
    for(int js = 0; js < ny;++js)
      for(int is = 0; is < nx;++is) {
	ind_global = is+nx*(js+ny*ks);
	i = ids_i[ind_global];
	j = ids_j[ind_global];
	k = ids_k[ind_global];
	id_hex = ids_hex[ind_global];
		
	if(id_hex < 0) data[ind_global] = defaultInterpolation;
	else {
	   Blockgrid* blockgrid = u.getBlockgrid();
           int  Nx = blockgrid->Give_Nx_hexahedron(id_hex);
           int  Ny = blockgrid->Give_Ny_hexahedron(id_hex);
           int  Nz = blockgrid->Give_Nz_hexahedron(id_hex);
	    
	   data[ind_global] = u.Give_cell_hexahedra(id_hex,Ind_loc_matrix_hexahedra(i,j,k),i,j,k,Nx,Ny);
	}
      }
}

/////////////////////////////////////////////////////////////
// 2. Interpolate from  blockgrid  to  blockgrid
/////////////////////////////////////////////////////////////

/**
 * interpolates data from blockgrid_ to rectangular block [pWSD, pENT]
**/
class Interpolate_on_block_grid {
 public:
/**
 * preparation for interpolation
 @param  nx_ number of structured grid points x-direction
 @param  ny_ number of structured grid points y-direction
 @param  nz_ number of structured grid points z-direction
 @param  blockgrid_ of unstructured grid
**/
  Interpolate_on_block_grid(int nx_, int ny_, int nz_,
		            Blockgrid* blockgrid_from, Blockgrid* blockgrid_to_);
  ~Interpolate_on_block_grid();
  
  /**
   * interpoliert Daten. Falls an einem Punkt nicht interpoliert werden kann
   * da U_from keine Zelle hat, dann wird
   *     defaultInterpolation
   * verwendet.
   ***/ 
  void interpolate(Variable<double>* U_from, Variable<double>* U_to, double defaultInterpolation = 0.0);
    		
  double evaluate(double coord_x, double coord_y, double coord_z);
		  
private:
  int nx, ny, nz;
  double hx, hy, hz;
  Interpolate_on_structured_grid* interpolatorStructured;
  double* data;
  Blockgrid* blockgrid_to;
  
  D3vector pWSD;
  D3vector pENT;
  
};


#endif // INTERPOL_H