OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
EllipticDomain.cpp
Go to the documentation of this file.
1 //
2 // Class EllipticDomain
3 // This class provides an elliptic beam pipe. The mesh adapts to the bunch size
4 // in the longitudinal direction. At the intersection of the mesh with the
5 // beam pipe, three stencil interpolation methods are available.
6 //
7 // Copyright (c) 2008, Yves Ineichen, ETH Zürich,
8 // 2013 - 2015, Tülin Kaman, Paul Scherrer Institut, Villigen PSI, Switzerland
9 // 2017 - 2020, Paul Scherrer Institut, Villigen PSI, Switzerland
10 // All rights reserved
11 //
12 // Implemented as part of the master thesis
13 // "A Parallel Multigrid Solver for Beam Dynamics"
14 // and the paper
15 // "A fast parallel Poisson solver on irregular domains applied to beam dynamics simulations"
16 // (https://doi.org/10.1016/j.jcp.2010.02.022)
17 //
18 // This file is part of OPAL.
19 //
20 // OPAL is free software: you can redistribute it and/or modify
21 // it under the terms of the GNU General Public License as published by
22 // the Free Software Foundation, either version 3 of the License, or
23 // (at your option) any later version.
24 //
25 // You should have received a copy of the GNU General Public License
26 // along with OPAL. If not, see <https://www.gnu.org/licenses/>.
27 //
28 #include "Solvers/EllipticDomain.h"
29 
30 #include <map>
31 #include <cmath>
32 #include <iostream>
33 
34 //FIXME: ORDER HOW TO TRAVERSE NODES IS FIXED, THIS SHOULD BE MORE GENERIC! (PLACES MARKED)
35 
36 
38  std::string interpl)
39  : RegularDomain(nr, hr, interpl)
40 {
41  Vector_t min(-bgeom->getA(), -bgeom->getB(), bgeom->getS());
42  Vector_t max( bgeom->getA(), bgeom->getB(), bgeom->getS() + bgeom->getLength());
45  setMinMaxZ(min[2], max[2]);
46 }
47 
49  //nothing so far
50 }
51 
52 // for this geometry we only have to calculate the intersection on
53 // one x-y-plane
54 // for the moment we center the ellipse around the center of the grid
56  // there is nothing to be done if the mesh spacings in x and y have not changed
57  if (hr[0] == getHr()[0] &&
58  hr[1] == getHr()[1])
59  {
60  hasGeometryChanged_m = false;
61  return;
62  }
63 
64  setHr(hr);
65  hasGeometryChanged_m = true;
66  //reset number of points inside domain
67 
68  // clear previous coordinate maps
69  idxMap_m.clear();
70  coordMap_m.clear();
71  //clear previous intersection points
72  intersectYDir_m.clear();
73  intersectXDir_m.clear();
74 
75  // build a index and coordinate map
76  int idx = 0;
77  int x, y;
78 
79  int nxy = 0;
80 
81  /* we need to scan the full x-y-plane on all cores
82  * in order to figure out the number of valid
83  * grid points per plane --> otherwise we might
84  * get not unique global indices in the Tpetra::CrsMatrix
85  */
86  for (x = 0; x < nr_m[0]; ++x) {
87  for (y = 0; y < nr_m[1]; ++y) {
88  if (isInside(x, y, 1)) {
89  idxMap_m[toCoordIdx(x, y)] = idx;
90  coordMap_m[idx++] = toCoordIdx(x, y);
91  nxy++;
92  }
93  }
94  }
95 
96  setNumXY(nxy);
97 
98  switch (interpolationMethod_m) {
99  case CONSTANT:
100  break;
101  case LINEAR:
102  case QUADRATIC:
103 
104  double smajsq = getXRangeMax() * getXRangeMax();
105  double sminsq = getYRangeMax() * getYRangeMax();
106  double yd = 0.0;
107  double xd = 0.0;
108  double pos = 0.0;
109 
110  // calculate intersection with the ellipse
111  for (x = localId[0].first(); x <= localId[0].last(); x++) {
112  pos = getXRangeMin() + hr_m[0] * (x + 0.5);
113  if (std::abs(pos) >= getXRangeMax())
114  {
115  intersectYDir_m.insert(std::pair<int, double>(x, 0));
116  intersectYDir_m.insert(std::pair<int, double>(x, 0));
117  } else {
118  yd = std::abs(std::sqrt(sminsq - sminsq * pos * pos / smajsq));
119  intersectYDir_m.insert(std::pair<int, double>(x, yd));
120  intersectYDir_m.insert(std::pair<int, double>(x, -yd));
121  }
122 
123  }
124 
125  for (y = localId[0].first(); y < localId[1].last(); y++) {
126  pos = getYRangeMin() + hr_m[1] * (y + 0.5);
127  if (std::abs(pos) >= getYRangeMax())
128  {
129  intersectXDir_m.insert(std::pair<int, double>(y, 0));
130  intersectXDir_m.insert(std::pair<int, double>(y, 0));
131  } else {
132  xd = std::abs(std::sqrt(smajsq - smajsq * pos * pos / sminsq));
133  intersectXDir_m.insert(std::pair<int, double>(y, xd));
134  intersectXDir_m.insert(std::pair<int, double>(y, -xd));
135  }
136  }
137  }
138 }
139 
140 void EllipticDomain::linearInterpolation(int x, int y, int z, StencilValue_t& value,
141  double &scaleFactor) const
142 {
143  scaleFactor = 1.0;
144 
145  double cx = getXRangeMin() + hr_m[0] * (x + 0.5);
146  double cy = getYRangeMin() + hr_m[1] * (y + 0.5);
147 
148  double dx = 0.0;
149  std::multimap<int, double>::const_iterator it = intersectXDir_m.find(y);
150 
151  if (cx < 0)
152  ++it;
153  dx = it->second;
154 
155  double dy = 0.0;
156  it = intersectYDir_m.find(x);
157  if (cy < 0)
158  ++it;
159  dy = it->second;
160 
161 
162  double dw = hr_m[0];
163  double de = hr_m[0];
164  double dn = hr_m[1];
165  double ds = hr_m[1];
166  value.center = 0.0;
167 
168  // we are a right boundary point
169  if (!isInside(x + 1, y, z)) {
170  value.center += 1.0 / ((dx - cx) * de);
171  value.east = 0.0;
172  } else {
173  value.center += 1.0 / (de * de);
174  value.east = -1.0 / (de * de);
175  }
176 
177  // we are a left boundary point
178  if (!isInside(x - 1, y, z)) {
179  value.center += 1.0 / ((std::abs(dx) - std::abs(cx)) * dw);
180  value.west = 0.0;
181  } else {
182  value.center += 1.0 / (dw * dw);
183  value.west = -1.0 / (dw * dw);
184  }
185 
186  // we are a upper boundary point
187  if (!isInside(x, y + 1, z)) {
188  value.center += 1.0 / ((dy - cy) * dn);
189  value.north = 0.0;
190  } else {
191  value.center += 1.0 / (dn * dn);
192  value.north = -1.0 / (dn * dn);
193  }
194 
195  // we are a lower boundary point
196  if (!isInside(x, y - 1, z)) {
197  value.center += 1.0 / ((std::abs(dy) - std::abs(cy)) * ds);
198  value.south = 0.0;
199  } else {
200  value.center += 1.0 / (ds * ds);
201  value.south = -1.0 / (ds * ds);
202  }
203 
204  // handle boundary condition in z direction
205  value.front = -1.0 / (hr_m[2] * hr_m[2]);
206  value.back = -1.0 / (hr_m[2] * hr_m[2]);
207  value.center += 2.0 / (hr_m[2] * hr_m[2]);
208  robinBoundaryStencil(z, value.front, value.back, value.center);
209 }
210 
211 void EllipticDomain::quadraticInterpolation(int x, int y, int z,
212  StencilValue_t& value,
213  double &scaleFactor) const
214 {
215  scaleFactor = 1.0;
216 
217  double cx = (x - (nr_m[0] - 1) / 2.0) * hr_m[0];
218  double cy = (y - (nr_m[1] - 1) / 2.0) * hr_m[1];
219 
220  // since every vector for elliptic domains has ALWAYS size 2 we
221  // can catch all cases manually
222  double dx = 0.0;
223  std::multimap<int, double>::const_iterator it = intersectXDir_m.find(y);
224  if (cx < 0)
225  ++it;
226  dx = it->second;
227 
228  double dy = 0.0;
229  it = intersectYDir_m.find(x);
230  if (cy < 0)
231  ++it;
232  dy = it->second;
233 
234  double dw = hr_m[0];
235  double de = hr_m[0];
236  double dn = hr_m[1];
237  double ds = hr_m[1];
238 
239  value.west = 0.0;
240  value.east = 0.0;
241  value.south = 0.0;
242  value.north = 0.0;
243  value.center = 0.0;
244 
245  // we are a right boundary point
246  if (!isInside(x + 1, y, z)) {
247  double s = dx - cx;
248  value.center -= 2.0 * (s - 1.0) / (s * de);
249  value.east = 0.0;
250  value.west += (s - 1.0) / ((s + 1.0) * de);
251  } else {
252  value.center += 1.0 / (de * de);
253  value.east = -1.0 / (de * de);
254  }
255 
256  // we are a left boundary point
257  if (!isInside(x - 1, y, z)) {
258  double s = std::abs(dx) - std::abs(cx);
259  value.center -= 2.0 * (s - 1.0) / (s * de);
260  value.west = 0.0;
261  value.east += (s - 1.0) / ((s + 1.0) * de);
262  } else {
263  value.center += 1.0 / (dw * dw);
264  value.west = -1.0 / (dw * dw);
265  }
266 
267  // we are a upper boundary point
268  if (!isInside(x, y + 1, z)) {
269  double s = dy - cy;
270  value.center -= 2.0 * (s - 1.0) / (s * dn);
271  value.north = 0.0;
272  value.south += (s - 1.0) / ((s + 1.0) * dn);
273  } else {
274  value.center += 1.0 / (dn * dn);
275  value.north = -1.0 / (dn * dn);
276  }
277 
278  // we are a lower boundary point
279  if (!isInside(x, y - 1, z)) {
280  double s = std::abs(dy) - std::abs(cy);
281  value.center -= 2.0 * (s - 1.0) / (s * dn);
282  value.south = 0.0;
283  value.north += (s - 1.0) / ((s + 1.0) * dn);
284  } else {
285  value.center += 1.0 / (ds * ds);
286  value.south = -1.0 / (ds * ds);
287  }
288 
289  // handle boundary condition in z direction
290  value.front = -1.0 / (hr_m[2] * hr_m[2]);
291  value.back = -1.0 / (hr_m[2] * hr_m[2]);
292  value.center += 2.0 / (hr_m[2] * hr_m[2]);
293  robinBoundaryStencil(z, value.front, value.back, value.center);
294 }
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
const int nr
Definition: ClassicRandom.h:24
@ QUADRATIC
@ CONSTANT
@ LINEAR
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:84
T::PETE_Expr_t::PETE_Return_t min(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:76
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
void compute(Vector_t hr, NDIndex< 3 > localId)
calculates intersection
int toCoordIdx(int x, int y) const
conversion from (x,y) to index in xy plane
void linearInterpolation(int x, int y, int z, StencilValue_t &value, double &scaleFactor) const override
different interpolation methods for boundary points
void quadraticInterpolation(int x, int y, int z, StencilValue_t &value, double &scaleFactor) const override
EllipticPointList_t intersectYDir_m
all intersection points with grid lines in Y direction
bool isInside(int x, int y, int z) const
queries if a given (x,y,z) coordinate lies inside the domain
EllipticPointList_t intersectXDir_m
all intersection points with grid lines in X direction
EllipticDomain(BoundaryGeometry *bgeom, IntVector_t nr, Vector_t hr, std::string interpl)
IntVector_t nr_m
number of mesh points in each direction
std::map< int, int > coordMap_m
mapping idx -> (x,y,z)
double getXRangeMax() const
void setRangeMin(const Vector_t &min)
double getYRangeMax() const
int interpolationMethod_m
interpolation type
bool hasGeometryChanged_m
flag indicating if geometry has changed for the current time-step
void setMinMaxZ(double minz, double maxz)
double getXRangeMin() const
Vector_t getHr() const
std::map< int, int > idxMap_m
mapping (x,y,z) -> idx
double getYRangeMin() const
Vector_t hr_m
mesh-spacings in each direction
void setRangeMax(const Vector_t &max)
void setHr(Vector_t hr)
void setNumXY(int nxy)
Definition: RegularDomain.h:36
void robinBoundaryStencil(int z, double &F, double &B, double &C) const
function to handle the open boundary condition in longitudinal direction