Navigation

Main Page
Download
Support
Installation
Tutorial
Examples
Reference Manual
   Version 1.0.5new_solver
      Class Index
      File List
   Version 1.0.5dev
   Version 1.0.5b
   Version 1.0.4dev
   Version 1.0.4
Publications


Hosted by Get Ion Beam Simulator at SourceForge.net. Fast, secure and Free Open Source software downloads
particleiterator.hpp
Go to the documentation of this file.
1 
5 /* Copyright (c) 2005-2013 Taneli Kalvas, Tobin Jones. All rights reserved.
6  *
7  * You can redistribute this software and/or modify it under the terms
8  * of the GNU General Public License as published by the Free Software
9  * Foundation; either version 2 of the License, or (at your option)
10  * any later version.
11  *
12  * This library is distributed in the hope that it will be useful, but
13  * WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
15  * General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this library (file "COPYING" included in the package);
19  * if not, write to the Free Software Foundation, Inc., 51 Franklin
20  * Street, Fifth Floor, Boston, MA 02110-1301 USA
21  *
22  * If you have questions about your rights to use or distribute this
23  * software, please contact Berkeley Lab's Technology Transfer
24  * Department at TTD@lbl.gov. Other questions, comments and bug
25  * reports should be sent directly to the author via email at
26  * taneli.kalvas@jyu.fi.
27  *
28  * NOTICE. This software was developed under partial funding from the
29  * U.S. Department of Energy. As such, the U.S. Government has been
30  * granted for itself and others acting on its behalf a paid-up,
31  * nonexclusive, irrevocable, worldwide license in the Software to
32  * reproduce, prepare derivative works, and perform publicly and
33  * display publicly. Beginning five (5) years after the date
34  * permission to assert copyright is obtained from the U.S. Department
35  * of Energy, and subject to any subsequent five (5) year renewals,
36  * the U.S. Government is granted for itself and others acting on its
37  * behalf a paid-up, nonexclusive, irrevocable, worldwide license in
38  * the Software to reproduce, prepare derivative works, distribute
39  * copies to the public, perform publicly and display publicly, and to
40  * permit others to do so.
41  */
42 
43 #ifndef PARTICLEITERATOR_HPP
44 #define PARTICLEITERATOR_HPP 1
45 
46 
47 #include <vector>
48 #include <iostream>
49 #include <algorithm>
50 #include <iomanip>
51 #include <gsl/gsl_odeiv.h>
52 #include <gsl/gsl_poly.h>
53 #include "geometry.hpp"
54 #include "mat3d.hpp"
55 #include "compmath.hpp"
56 #include "trajectory.hpp"
57 #include "particles.hpp"
58 #include "vectorfield.hpp"
59 #include "meshscalarfield.hpp"
60 #include "scharge.hpp"
61 #include "scheduler.hpp"
62 #include "polysolver.hpp"
63 #include "particledatabase.hpp"
64 #include "cfifo.hpp"
65 #include "ibsimu.hpp"
66 
67 
68 //#define DEBUG_PARTICLE_ITERATOR 1
69 
70 
71 #ifdef DEBUG_PARTICLE_ITERATOR
72 #define DEBUG_MESSAGE(x) ibsimu.message(MSG_DEBUG_GENERAL,1) << x
73 #define DEBUG_INC_INDENT() ibsimu.inc_indent()
74 #define DEBUG_DEC_INDENT() ibsimu.dec_indent()
75 #else
76 #define DEBUG_MESSAGE(x) do {} while(0)
77 #define DEBUG_INC_INDENT() do {} while(0)
78 #define DEBUG_DEC_INDENT() do {} while(0)
79 #endif
80 
81 
82 #define COLLISION_EPS 1.0e-6
83 
84 
88  PARTICLE_ITERATOR_ADAPTIVE = 0,
89  PARTICLE_ITERATOR_FIXED_STEP_LEN
90 };
91 
92 
93 
102 template <class PP> class ColData {
103 public:
104  PP _x;
105  int _dir;
110  ColData() : _dir(0) {}
111 
114  ColData( PP x, int dir ) : _x(x), _dir(dir) {}
115 
120  bool operator<( const ColData &cd ) const {
121  return( _x[0] < cd._x[0] );
122  }
123 
132  static void build_coldata_linear( std::vector<ColData> &coldata, const Mesh &mesh,
133  const PP &x1, const PP &x2 ) {
134 
135  DEBUG_MESSAGE( "Building coldata using linear interpolation\n" );
136  DEBUG_INC_INDENT();
137 
138  coldata.clear();
139 
140  for( size_t a = 0; a < PP::dim(); a++ ) {
141 
142  int a1 = (int)floor( (x1[2*a+1]-mesh.origo(a))/mesh.h() );
143  int a2 = (int)floor( (x2[2*a+1]-mesh.origo(a))/mesh.h() );
144  if( a1 > a2 ) {
145  int a = a2;
146  a2 = a1;
147  a1 = a;
148  }
149 
150  for( int b = a1+1; b <= a2; b++ ) {
151 
152  // Save intersection coordinates
153  double K = (b*mesh.h() + mesh.origo(a) - x1[2*a+1]) /
154  (x2[2*a+1] - x1[2*a+1]);
155  if( K < 0.0 ) K = 0.0;
156  else if( K > 1.0 ) K = 1.0;
157 
158  DEBUG_MESSAGE( "Adding point " << x1 + (x2-x1)*K << "\n" );
159 
160  if( x2[2*a+1] > x1[2*a+1] )
161  coldata.push_back( ColData( x1 + (x2-x1)*K, a+1 ) );
162  else
163  coldata.push_back( ColData( x1 + (x2-x1)*K, -a-1 ) );
164  }
165  }
166 
167  // Sort intersections in increasing time order
168  sort( coldata.begin(), coldata.end() );
169 
170  DEBUG_DEC_INDENT();
171  DEBUG_MESSAGE( "Coldata built\n" );
172  }
173 
182  static void build_coldata_poly( std::vector<ColData> &coldata, const Mesh &mesh,
183  const PP &x1, const PP &x2 ) {
184 
185  DEBUG_MESSAGE( "Building coldata using polynomial interpolation\n" );
186  DEBUG_INC_INDENT();
187 
188  coldata.clear();
189 
190  // Construct trajectory representation
191  TrajectoryRep1D *traj = new TrajectoryRep1D[PP::dim()];
192  for( size_t a = 0; a < PP::dim(); a++ ) {
193  traj[a].construct( x2[0]-x1[0],
194  x1[2*a+1], x1[2*a+2],
195  x2[2*a+1], x2[2*a+2] );
196  DEBUG_MESSAGE( "Trajectory polynomial " << a << " order = "
197  << traj[a].get_representation_order() << "\n" );
198  }
199 
200  // Solve trajectory intersections
201  for( size_t a = 0; a < PP::dim(); a++ ) {
202 
203  // Mesh number of x1 (start point)
204  int i = (int)floor( (x1[2*a+1]-mesh.origo(a))/mesh.h() );
205 
206  // Search to negative (dj = -1) and positive (dj = +1) mesh directions
207  for( int dj = -1; dj <= 1; dj += 2 ) {
208  int j = i;
209  if( dj == +1 )
210  j = i+1;
211  int Kcount; // Solution counter
212  double K[3]; // Solution array
213  while( 1 ) {
214 
215  // Intersection point
216  double val = mesh.origo(a) + mesh.h() * j;
217  if( val < mesh.origo(a)-mesh.h() )
218  break;
219  else if( val > mesh.max(a)+mesh.h() )
220  break;
221 
222  DEBUG_MESSAGE( "Searching intersections at coord(" << a << ") = " << val << "\n" );
223 
224  Kcount = traj[a].solve( K, val );
225  if( Kcount == 0 )
226  break; // No valid roots
227 
228  // Save roots to coldata
229  DEBUG_MESSAGE( "Found " << Kcount << " valid roots\n" );
230  DEBUG_INC_INDENT();
231  for( int b = 0; b < Kcount; b++ ) {
232  PP xcol;
233  double x, v;
234  xcol(0) = x1[0] + K[b]*(x2[0]-x1[0]);
235  for( size_t c = 0; c < PP::dim(); c++ ) {
236  traj[c].coord( x, v, K[b] );
237  if( a == c )
238  xcol[2*c+1] = val; // limit numerical inaccuracy
239  else
240  xcol[2*c+1] = x;
241  xcol[2*c+2] = v;
242  }
243  if( mesh.geom_mode() == MODE_CYL )
244  xcol[5] = x1[5] + K[b]*(x2[5]-x1[5]);
245 
246  DEBUG_MESSAGE( "K = " << K[b] << "\n" );
247  DEBUG_MESSAGE( "Adding point " << xcol << "\n" );
248 
249  if( xcol[2*a+2] >= 0.0 )
250  coldata.push_back( ColData( xcol, a+1 ) );
251  else
252  coldata.push_back( ColData( xcol, -a-1 ) );
253  }
254  DEBUG_DEC_INDENT();
255 
256  j += dj;
257  }
258  }
259  }
260 
261  // Sort intersections in increasing time order
262  sort( coldata.begin(), coldata.end() );
263 
264  delete[] traj;
265 
266  DEBUG_DEC_INDENT();
267  DEBUG_MESSAGE( "Coldata built\n" );
268  }
269 
270 };
271 
272 
280 template <class PP> class ParticleIterator {
281 
282  gsl_odeiv_system _system;
283  gsl_odeiv_step *_step;
284  gsl_odeiv_control *_control;
285  gsl_odeiv_evolve *_evolve;
290  scharge_deposition_e _scharge_dep;
291  double _epsabs;
292  double _epsrel;
293  uint32_t _maxsteps;
294  double _maxt;
295  bool _save_points;
296  uint32_t _trajdiv;
298  bool _mirror[6];
299  bool _surface_collision;
300 
301  ParticleIteratorData _pidata;
302  TrajectoryHandlerCallback *_thand_cb;
303  TrajectoryEndCallback *_tend_cb;
305  const TrajectoryEndCallback *_bsup_cb;
306  ParticleDataBase *_pdb;
307  pthread_mutex_t *_scharge_mutex;
309  PP _xi;
311  std::vector<PP> _traj;
312  std::vector<ColData<PP> > _coldata;
313  CFiFo<PP,4> _cdpast;
315  ParticleStatistics _stat;
322  void save_trajectory_point( PP x ) {
323 
324  try {
325  _traj.push_back( x );
326  } catch( std::bad_alloc ) {
327  throw( ErrorNoMem( ERROR_LOCATION, "Out of memory saving trajectory" ) );
328  }
329  }
330 
333  uint32_t get_solid( int i, int j, int k ) {
334  uint32_t node;
335  if( PP::dim() == 2 ) {
336  node = _pidata._geom->mesh( i, j );
337  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
338  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
339  return( node & SMESH_BOUNDARY_NUMBER_MASK );
340  node = _pidata._geom->mesh( i+1, j );
341  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
342  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
343  return( node & SMESH_BOUNDARY_NUMBER_MASK );
344  node = _pidata._geom->mesh( i, j+1 );
345  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
346  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
347  return( node & SMESH_BOUNDARY_NUMBER_MASK );
348  node = _pidata._geom->mesh( i+1, j+1 );
349  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
350  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
351  return( node & SMESH_BOUNDARY_NUMBER_MASK );
352  } else {
353  node = _pidata._geom->mesh( i, j, k );
354  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
355  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
356  return( node & SMESH_BOUNDARY_NUMBER_MASK );
357  node = _pidata._geom->mesh( i+1, j, k );
358  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
359  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
360  return( node & SMESH_BOUNDARY_NUMBER_MASK );
361  node = _pidata._geom->mesh( i, j+1, k );
362  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
363  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
364  return( node & SMESH_BOUNDARY_NUMBER_MASK );
365  node = _pidata._geom->mesh( i+1, j+1, k );
366  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
367  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
368  return( node & SMESH_BOUNDARY_NUMBER_MASK );
369  node = _pidata._geom->mesh( i, j, k+1 );
370  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
371  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
372  return( node & SMESH_BOUNDARY_NUMBER_MASK );
373  node = _pidata._geom->mesh( i+1, j, k+1 );
374  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
375  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
376  return( node & SMESH_BOUNDARY_NUMBER_MASK );
377  node = _pidata._geom->mesh( i, j+1, k+1 );
378  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
379  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
380  return( node & SMESH_BOUNDARY_NUMBER_MASK );
381  node = _pidata._geom->mesh( i+1, j+1, k+1 );
382  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
383  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
384  return( node & SMESH_BOUNDARY_NUMBER_MASK );
385  }
386  return( 0 );
387  }
388 
394  bool check_collision_surface( Particle<PP> &particle, const PP &x1, const PP &x2,
395  PP &status_x, const int32_t i[3] ) {
396 
397  DEBUG_MESSAGE( "Checking collisions with solids - surface check\n" );
398  DEBUG_INC_INDENT();
399 
400  if( i[0] < 0 )
401  return( true );
402  else if( i[0] >= (int32_t)_pidata._geom->size(0)-1 )
403  return( true );
404  else if( i[1] < 0 )
405  return( true );
406  else if( i[1] >= (int32_t)_pidata._geom->size(1)-1 )
407  return( true );
408  else if( i[2] < 0 )
409  return( true );
410  else if( i[2] >= (int32_t)_pidata._geom->size(2)-1 )
411  return( true );
412 
413  Vec3D v1 = x1.location();
414  Vec3D v2 = x2.location();
415 
416  DEBUG_MESSAGE( "v1 = " << v1 << "\n" );
417  DEBUG_MESSAGE( "v2 = " << v2 << "\n" );
418 
419  int32_t tric = _pidata._geom->surface_trianglec( i[0], i[1], i[2] );
420  int32_t ptr = _pidata._geom->surface_triangle_ptr( i[0], i[1], i[2] );
421 
422  DEBUG_MESSAGE( "tric = " << tric << "\n" );
423 
424  // Go through surface triangles at mesh cube i
425  for( int32_t a = 0; a < tric; a++ ) {
426  const VTriangle &tri = _pidata._geom->surface_triangle( ptr+a );
427  const Vec3D &va = _pidata._geom->surface_vertex( tri[0] );
428  const Vec3D &vb = _pidata._geom->surface_vertex( tri[1] );
429  const Vec3D &vc = _pidata._geom->surface_vertex( tri[2] );
430 
431  DEBUG_MESSAGE( "a = " << a << "\n" );
432  DEBUG_MESSAGE( "tri[" << a << "][0] = " << va << "\n" );
433  DEBUG_MESSAGE( "tri[" << a << "][1] = " << vb << "\n" );
434  DEBUG_MESSAGE( "tri[" << a << "][2] = " << vc << "\n" );
435 
436  // Solve for intersection between trajectory segment and surface triangle
437  // r1 + (r2-r1)*K[0] = ra + (rb-ra)*K[1] + (rc-ra)*K[2]
438  Mat3D m( v2[0]-v1[0], -vb[0]+va[0], -vc[0]+va[0],
439  v2[1]-v1[1], -vb[1]+va[1], -vc[1]+va[1],
440  v2[2]-v1[2], -vb[2]+va[2], -vc[2]+va[2] );
441  double mdet = m.determinant();
442  if( mdet == 0.0 )
443  continue;
444  Mat3D minv = m.inverse( mdet );
445  Vec3D off( -v1[0]+va[0], -v1[1]+va[1], -v1[2]+va[2] );
446  Vec3D K = minv*off;
447  double K3 = K[1]+K[2];
448  DEBUG_MESSAGE( "K = " << K << "\n" );
449 
450  // Check for intersection at valid ranges
451  // Allow COLLISION_EPS amount of overlap, double inclusion is
452  // not an issue here, missing an intersection is a problem.
453  if( K[0] > -COLLISION_EPS && K[0] < 1.0+COLLISION_EPS &&
454  K[1] > -COLLISION_EPS && K[1] < 1.0+COLLISION_EPS &&
455  K[2] > -COLLISION_EPS && K[2] < 1.0+COLLISION_EPS &&
456  K3 > -COLLISION_EPS && K3 < 1.0+COLLISION_EPS ) {
457 
458  DEBUG_MESSAGE( "Intersection found\n" );
459 
460  // Found intersection, set collision coordinates
461  for( uint32_t b = 0; b < PP::size(); b++ )
462  status_x[b] = x1[b] + K[0]*(x2[b]-x1[b]);
463 
464  // Remove all points from trajectory after time status_x[0].
465  // Does this ever happen???
466  for( int32_t b = _traj.size()-1; b > 0; b-- ) {
467  if( _traj[b][0] > status_x[0] )
468  _traj.pop_back();
469  else
470  break;
471  }
472 
473  // Update status and statistics
474  particle.set_status( PARTICLE_COLL );
475  uint32_t solid = get_solid( i[0], i[1], i[2] );
476  _stat.add_bound_collision( solid, particle.IQ() );
477 
478  if( _tsur_cb )
479  (*_tsur_cb)( &particle, &status_x, ptr+a, K[1], K[2] );
480 
481  DEBUG_MESSAGE( "Solid collision detected\n" );
482  DEBUG_DEC_INDENT();
483 
484  return( false );
485  }
486  }
487 
488  DEBUG_MESSAGE( "No collisions\n" );
489  DEBUG_DEC_INDENT();
490 
491  return( true );
492  }
493 
502  bool check_collision_solid( Particle<PP> &particle, const PP &x1, const PP &x2,
503  PP &status_x ) {
504 
505  DEBUG_MESSAGE( "Checking collisions with solids - inside check\n" );
506  DEBUG_INC_INDENT();
507 
508  // If inside solid, bracket for collision point
509  Vec3D v2 = x2.location();
510  int32_t bound = _pidata._geom->inside( v2 );
511  if( bound < 7 ) {
512  DEBUG_MESSAGE( "No collisions\n" );
513  DEBUG_DEC_INDENT();
514  return( true );
515  }
516  Vec3D vc;
517  Vec3D v1 = x1.location();
518  double K = _pidata._geom->bracket_surface( bound, v2, v1, vc );
519 
520  // Calculate new PP
521  for( size_t a = 0; a < PP::size(); a++ )
522  status_x[a] = x2[a] + K*(x1[a]-x2[a]);
523 
524  // Remove all points from trajectory after time status_x[0].
525  for( size_t a = _traj.size()-1; a > 0; a-- ) {
526  if( _traj[a][0] > status_x[0] )
527  _traj.pop_back();
528  else
529  break;
530  }
531 
532  // Save last trajectory point and update status
533  //save_trajectory_point( status_x );
534  particle.set_status( PARTICLE_COLL );
535 
536  // Update collision statistics for boundary
537  _stat.add_bound_collision( bound, particle.IQ() );
538 
539  DEBUG_MESSAGE( "Solid collision detected\n" );
540  DEBUG_DEC_INDENT();
541 
542  return( false );
543  }
544 
549  bool check_collision( Particle<PP> &particle, const PP &x1, const PP &x2,
550  PP &status_x, const int i[3] ) {
551 
552  if( _surface_collision )
553  return( check_collision_surface( particle, x1, x2, status_x, i ) );
554  return( check_collision_solid( particle, x1, x2, status_x ) );
555  }
556 
557 
565  void handle_mirror( size_t c, int i[3], size_t a, int border, PP &x2 ) {
566 
567  DEBUG_MESSAGE( "Mirror trajectory\n" );
568  DEBUG_INC_INDENT();
569 
570  double xmirror;
571  if( border < 0 ) {
572  xmirror = _pidata._geom->origo(a);
573  i[a] = -i[a]-1;
574  } else {
575  xmirror = _pidata._geom->max(a);
576  i[a] = 2*_pidata._geom->size(a)-i[a]-3;
577  }
578 
579  DEBUG_MESSAGE( "xmirror = " << xmirror << "\n" );
580  DEBUG_MESSAGE( "i = (" << i[0] << ", " << i[1] << ", " << i[2] << ")\n" );
581  DEBUG_MESSAGE( "xi = " << _xi << "\n" );
582 
583  // Check if found edge at first encounter
584  bool caught_at_boundary = false;
585  if( _coldata[c]._dir == border*((int)a+1) &&
586  ( i[a] == 0 || i[a] == (int)_pidata._geom->size(a)-2 ) ) {
587  caught_at_boundary = true;
588  DEBUG_MESSAGE( "caught_at_boundary\n" );
589  }
590 
591  // Mirror traj back to _xi
592  if( caught_at_boundary ) {
593  save_trajectory_point( _coldata[c]._x );
594  } else {
595  for( int b = _traj.size()-1; b > 0; b-- ) {
596  if( _traj[b][0] >= _xi[0] ) {
597 
598  DEBUG_MESSAGE( "mirroring traj[" << b << "] = " << _traj[b] << "\n" );
599  _traj[b][2*a+1] = 2.0*xmirror - _traj[b][2*a+1];
600  _traj[b][2*a+2] *= -1.0;
601  } else
602  break;
603  }
604  }
605 
606  // Mirror rest of the coldata
607  for( size_t b = c; b < _coldata.size(); b++ ) {
608  if( (size_t)abs(_coldata[b]._dir) == a+1 )
609  _coldata[b]._dir *= -1;
610  _coldata[b]._x[2*a+1] = 2.0*xmirror - _coldata[b]._x[2*a+1];
611  _coldata[b]._x[2*a+2] *= -1.0;
612  }
613 
614  if( caught_at_boundary )
615  save_trajectory_point( _coldata[c]._x );
616 
617  // Mirror calculation point
618  x2[2*a+1] = 2.0*xmirror - x2[2*a+1];
619  x2[2*a+2] *= -1.0;
620 
621  // Coordinates changed, reset integrator
622  gsl_odeiv_step_reset( _step );
623  gsl_odeiv_evolve_reset( _evolve );
624 
625  DEBUG_DEC_INDENT();
626  }
627 
628 
629  void handle_collision( Particle<PP> &particle, uint32_t bound, size_t c, PP &status_x ) {
630 
631  DEBUG_MESSAGE( "Handle collision\n" );
632 
633  //save_trajectory_point( _coldata[c]._x );
634  status_x = _coldata[c]._x;
635  particle.set_status( PARTICLE_OUT );
636  _stat.add_bound_collision( bound, particle.IQ() );
637  }
638 
639 
642  bool is_solid( int i, int j ) {
643  uint32_t node = _pidata._geom->mesh_check( i, j );
644  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
645  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
646  return( true );
647  return( false );
648  }
649 
652  bool is_solid( int i, int j, int k ) {
653  uint32_t node = _pidata._geom->mesh_check( i, j, k );
654  if( (node & SMESH_NODE_ID_MASK) == SMESH_NODE_ID_DIRICHLET &&
655  (node & SMESH_BOUNDARY_NUMBER_MASK) >= 7 )
656  return( true );
657  return( false );
658  }
659 
668  bool handle_trajectory_advance( Particle<PP> &particle, size_t c, int i[3], PP &x2 ) {
669 
670  bool surface_collision = false;
671  DEBUG_MESSAGE( "Handle trajectory advance\n" );
672  DEBUG_INC_INDENT();
673 
674  // Check for collisions with solids and advance coordinates i.
675  if( PP::dim() == 2 ) {
676  if( _coldata[c]._dir == -1 ) {
677  if( (is_solid(i[0],i[1]) || is_solid(i[0], i[1]+1)) &&
678  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
679  surface_collision = true;
680  i[0]--;
681  } else if( _coldata[c]._dir == +1 ) {
682  if( (is_solid(i[0]+1,i[1]) || is_solid(i[0]+1,i[1]+1)) &&
683  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
684  surface_collision = true;
685  i[0]++;
686  } else if( _coldata[c]._dir == -2 ) {
687  if( (is_solid(i[0],i[1]) || is_solid(i[0]+1,i[1])) &&
688  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
689  surface_collision = true;
690  i[1]--;
691  } else {
692  if( (is_solid(i[0], i[1]+1) || is_solid(i[0]+1,i[1]+1)) &&
693  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
694  surface_collision = true;
695  i[1]++;
696  }
697  } else if( PP::dim() == 3 ) {
698  if( _coldata[c]._dir == -1 ) {
699  if( (is_solid(i[0], i[1], i[2] ) ||
700  is_solid(i[0], i[1]+1,i[2] ) ||
701  is_solid(i[0], i[1], i[2]+1) ||
702  is_solid(i[0], i[1]+1,i[2]+1)) &&
703  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
704  surface_collision = true;
705  i[0]--;
706  } else if( _coldata[c]._dir == +1 ) {
707  if( (is_solid(i[0]+1,i[1], i[2] ) ||
708  is_solid(i[0]+1,i[1]+1,i[2] ) ||
709  is_solid(i[0]+1,i[1], i[2]+1) ||
710  is_solid(i[0]+1,i[1]+1,i[2]+1)) &&
711  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
712  surface_collision = true;
713  i[0]++;
714  } else if( _coldata[c]._dir == -2 ) {
715  if( (is_solid(i[0], i[1],i[2] ) ||
716  is_solid(i[0]+1,i[1],i[2] ) ||
717  is_solid(i[0], i[1],i[2]+1) ||
718  is_solid(i[0]+1,i[1],i[2]+1)) &&
719  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
720  surface_collision = true;
721  i[1]--;
722  } else if( _coldata[c]._dir == +2 ) {
723  if( (is_solid(i[0], i[1]+1,i[2] ) ||
724  is_solid(i[0]+1,i[1]+1,i[2] ) ||
725  is_solid(i[0], i[1]+1,i[2]+1) ||
726  is_solid(i[0]+1,i[1]+1,i[2]+1)) &&
727  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
728  surface_collision = true;
729  i[1]++;
730  } else if( _coldata[c]._dir == -3 ) {
731  if( (is_solid(i[0], i[1], i[2]) ||
732  is_solid(i[0]+1,i[1], i[2]) ||
733  is_solid(i[0], i[1]+1,i[2]) ||
734  is_solid(i[0]+1,i[1]+1,i[2])) &&
735  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
736  surface_collision = true;
737  i[2]--;
738  } else {
739  if( (is_solid(i[0], i[1], i[2]+1) ||
740  is_solid(i[0]+1,i[1], i[2]+1) ||
741  is_solid(i[0], i[1]+1,i[2]+1) ||
742  is_solid(i[0]+1,i[1]+1,i[2]+1)) &&
743  !check_collision( particle, _xi, _coldata[c]._x, x2, i ) )
744  surface_collision = true;
745  i[2]++;
746  }
747  } else {
748  throw( Error( ERROR_LOCATION, "unsupported dimension number" ) );
749  }
750 
751  if( surface_collision ) {
752  DEBUG_MESSAGE( "Surface collision!\n" );
753  DEBUG_DEC_INDENT();
754  return( false );
755  }
756 
757  // Check for collisions/mirroring with simulation boundary. Here
758  // coordinates i are already advanced to next mesh.
759  for( size_t a = 0; a < PP::dim(); a++ ) {
760 
761  if( i[a] < 0 ) {
762  DEBUG_MESSAGE( "Boundary collision at boundary " << 2*a+0 << "\n" );
763  if( _mirror[2*a] )
764  handle_mirror( c, i, a, -1, x2 );
765  else {
766  handle_collision( particle, 1+2*a, c, x2 );
767  DEBUG_DEC_INDENT();
768  return( false );
769  }
770  } else if( i[a] >= (int32_t)(_pidata._geom->size(a)-1) ) {
771  DEBUG_MESSAGE( "Boundary collision at boundary " << 2*a+1 << "\n" );
772  if( _mirror[2*a+1] )
773  handle_mirror( c, i, a, +1, x2 );
774  else {
775  handle_collision( particle, 2+2*a, c, x2 );
776  DEBUG_DEC_INDENT();
777  return( false );
778  }
779  }
780  }
781 
782  DEBUG_MESSAGE( "No collision.\n" );
783  DEBUG_DEC_INDENT();
784  return( true );
785  }
786 
792  bool limit_trajectory_advance( const PP &x1, PP &x2 ) {
793 
794  bool touched = false;
795  for( size_t a = 0; a < PP::dim(); a++ ) {
796 
797  double lim1 = _pidata._geom->origo(a) -
798  (_pidata._geom->size(a)-1)*_pidata._geom->h();
799  double lim2 = _pidata._geom->origo(a) +
800  2*(_pidata._geom->size(a)-1)*_pidata._geom->h();
801 
802  if( x2[2*a+1] < lim1 ) {
803  double K = (lim1 - x1[2*a+1]) / (x2[2*a+1] - x1[2*a+1]);
804  x2 = x1 + K*(x2-x1);
805  touched = true;
806  DEBUG_MESSAGE( "Limiting step to " << x2 << "\n" );
807  } else if(x2[2*a+1] > lim2 ) {
808  double K = (lim2 - x1[2*a+1]) / (x2[2*a+1] - x1[2*a+1]);
809  x2 = x1 + K*(x2-x1);
810  touched = true;
811  DEBUG_MESSAGE( "Limiting step to " << x2 << "\n" );
812  }
813  }
814 
815  return( touched );
816  }
817 
822  void build_coldata( bool force_linear, const PP &x1, const PP &x2 ) {
823 
824  try {
825  if( _intrp == TRAJECTORY_INTERPOLATION_POLYNOMIAL && !force_linear )
826  ColData<PP>::build_coldata_poly( _coldata, *_pidata._geom, x1, x2 );
827  else
828  ColData<PP>::build_coldata_linear( _coldata, *_pidata._geom, x1, x2 );
829  } catch( std::bad_alloc ) {
830  throw( ErrorNoMem( ERROR_LOCATION, "out of memory building collision data" ) );
831  }
832  }
833 
834 
854  bool handle_trajectory( Particle<PP> &particle, const PP &x1, PP &x2,
855  bool force_linear, bool first_step ) {
856 
857  DEBUG_MESSAGE( "Handle trajectory from x1 to x2:\n" <<
858  " x1 = " << x1 << "\n" <<
859  " x2 = " << x2 << "\n" );
860  DEBUG_INC_INDENT();
861 
862  // Limit trajectory advance to double the simulation box
863  // If limitation done, force to linear interpolation
864  if( limit_trajectory_advance( x1, x2 ) )
865  force_linear = true;
866 
867  build_coldata( force_linear, x1, x2 );
868 
869  // TODO
870  // Remove entrance to geometry if coming from outside or make
871  // code to skip the collision detection for these particles
872 
873  // No intersections, nothing to do
874  if( _coldata.size() == 0 ) {
875  DEBUG_MESSAGE( "No coldata\n" );
876  DEBUG_DEC_INDENT();
877  return( true );
878  }
879 
880  // Starting mesh index
881  int i[3] = {0, 0, 0};
882  for( size_t cdir = 0; cdir < PP::dim(); cdir++ )
883  i[cdir] = (int)floor( (x1[2*cdir+1]-_pidata._geom->origo(cdir))*_pidata._geom->div_h() );
884 
885  // Process intersection points
886  DEBUG_MESSAGE( "Process coldata points:\n" );
887  DEBUG_INC_INDENT();
888  for( size_t a = 0; a < _coldata.size(); a++ ) {
889 
890  if( _save_points )
891  save_trajectory_point( _coldata[a]._x );
892 
893  DEBUG_MESSAGE( "Coldata " << a << "\n" <<
894  " x = " << _coldata[a]._x << "\n" <<
895  " i = (" << std::setw(3) << i[0] << ", "
896  << std::setw(3) << i[1] << ", "
897  << std::setw(3) << i[2] << ")\n" <<
898  " dir = " << _coldata[a]._dir << "\n" );
899 
900  // Update space charge for mesh volume i
901  if( _scharge_dep == SCHARGE_DEPOSITION_LINEAR && _pidata._scharge ) {
902  _cdpast.push( _coldata[a]._x );
903  scharge_add_from_trajectory_linear( *_pidata._scharge, _scharge_mutex, particle.IQ(),
904  _coldata[a]._dir, _cdpast, i );
905  }
906 
907  // Advance particle in mesh, update i, check for possible
908  // collisions and do mirroring for trajectory and coldata.
909  handle_trajectory_advance( particle, a, i, x2 );
910 
911  if( _scharge_dep == SCHARGE_DEPOSITION_PIC && _pidata._scharge ) {
912  scharge_add_from_trajectory_pic( *_pidata._scharge, _scharge_mutex, particle.IQ(),
913  _xi, _coldata[a]._x );
914  }
915 
916  // Call trajectory handler callback
917  if( _thand_cb )
918  (*_thand_cb)( &particle, &_coldata[a]._x, &x2 );
919 
920  // Clear coldata and exit if particle collided.
921  if( particle.get_status() != PARTICLE_OK ) {
922  save_trajectory_point( x2 );
923  _coldata.clear();
924  DEBUG_DEC_INDENT();
925  DEBUG_MESSAGE( "Interrupting coldata processing\n" );
926  DEBUG_DEC_INDENT();
927  return( false );
928  }
929 
930  // Update last accepted intersection point xi.
931  _xi = _coldata[a]._x;
932  }
933 
934  _coldata.clear();
935  DEBUG_DEC_INDENT();
936  DEBUG_MESSAGE( "Coldata done\n" );
937  DEBUG_DEC_INDENT();
938  return( true );
939  }
940 
941 
944  bool axis_mirror_required( const PP &x2 ) {
945  return( _pidata._geom->geom_mode() == MODE_CYL &&
946  x2[4] < 0.0 &&
947  x2[3] <= 0.01*_pidata._geom->h() &&
948  x2[3]*fabs(x2[5]) <= 1.0e-9*fabs(x2[4]) );
949 
950  }
951 
952 
958  bool handle_axis_mirror_step( Particle<PP> &particle, const PP &x1, PP &x2 ) {
959 
960  // Get acceleration at x2
961  double dxdt[5];
962  PP::get_derivatives( x2[0], &x2[1], dxdt, (void *)&_pidata );
963 
964  // Calculate crossover point assuming zero acceleration in
965  // r-direction and constant acceleration in x-direction
966  double dt = -x2[3]/x2[4];
967  PP xc;
968  xc.clear();
969  xc[0] = x2[0]+dt;
970  xc[1] = x2[1]+(x2[2]+0.5*dxdt[1]*dt)*dt;
971  xc[2] = x2[2];
972  xc[3] = x2[3]+x2[4]*dt;
973  xc[4] = x2[4];
974  xc[5] = x2[5];
975 
976  // Mirror x2 to x3
977  PP x3 = 2*xc - x2;
978  x3[3] *= -1.0;
979  x3[4] *= -1.0;
980  x3[5] *= -1.0;
981 
982  DEBUG_MESSAGE( "Particle mirror:\n" <<
983  " x1: " << x1 << "\n" <<
984  " x2: " << x2 << "\n" <<
985  " xc: " << xc << "\n" <<
986  " x3: " << x3 << "\n" );
987 
988  // Handle step with linear interpolation to avoid going to r<=0
989  if( !handle_trajectory( particle, x2, x3, true, false ) )
990  return( false ); // Particle done
991 
992  // Save trajectory calculation points
993  save_trajectory_point( x2 );
994  save_trajectory_point( xc );
995  xc[4] *= -1.0;
996  xc[5] *= -1.0;
997  save_trajectory_point( xc );
998 
999  // Next step not a continuation of previous one, reset
1000  // integrator
1001  gsl_odeiv_step_reset( _step );
1002  gsl_odeiv_evolve_reset( _evolve );
1003 
1004  // Continue iteration at mirrored point
1005  x2 = x3;
1006  return( true );
1007  }
1008 
1022  bool check_particle_definition( PP &x ) {
1023 
1024  DEBUG_MESSAGE( "Particle defined at x = " << x << "\n" );
1025  DEBUG_INC_INDENT();
1026 
1027  // Check for NaN
1028  for( size_t a = 0; a < PP::size(); a++ ) {
1029  if( comp_isnan( x[a] ) ) {
1030  DEBUG_MESSAGE( "Particle coordinate NaN\n" );
1031  DEBUG_DEC_INDENT();
1032  return( false );
1033  }
1034  }
1035 
1036  // Check if inside solids or outside simulation geometry box.
1037  if( _surface_collision ) {
1038  if( _pidata._geom->surface_inside( x.location() ) ) {
1039  DEBUG_MESSAGE( "Particle inside solid or outside simulation\n" );
1040  DEBUG_DEC_INDENT();
1041  return( false );
1042  }
1043  } else {
1044  if( _pidata._geom->inside( x.location() ) ) {
1045  DEBUG_MESSAGE( "Particle inside solid or outside simulation\n" );
1046  DEBUG_DEC_INDENT();
1047  return( false );
1048  }
1049  }
1050 
1051  // Check if particle on simulation geometry border and directed outwards
1052  /*
1053  for( size_t a = 0; a < PP::dim(); a++ ) {
1054  if( x[2*a+1] == _pidata._geom->origo(a) && x[2*a+2] < 0.0 ) {
1055  if( _mirror[2*a] ) {
1056  x[2*a+2] *= -1.0;
1057 #ifdef DEBUG_PARTICLE_ITERATOR
1058  ibsimu.message(MSG_DEBUG_GENERAL,1) << "Mirroring to:\n";
1059  ibsimu.message(MSG_DEBUG_GENERAL,1) << " x = " << x << "\n";
1060 #endif
1061  } else {
1062  return( false );
1063  }
1064 
1065  } else if( x[2*a+1] == _pidata._geom->max(a) & x[2*a+2] > 0.0 ) {
1066  if( _mirror[2*a+1] ) {
1067  x[2*a+2] *= -1.0;
1068 #ifdef DEBUG_PARTICLE_ITERATOR
1069  ibsimu.message(MSG_DEBUG_GENERAL,1) << "Mirroring to:\n";
1070  ibsimu.message(MSG_DEBUG_GENERAL,1) << " x = " << x << "\n";
1071 #endif
1072  } else {
1073  return( false );
1074  }
1075  }
1076  }
1077 
1078 #ifdef DEBUG_PARTICLE_ITERATOR
1079  ibsimu.message(MSG_DEBUG_GENERAL,1) << "Definition ok\n";
1080 #endif
1081 
1082  */
1083  DEBUG_MESSAGE( "Ok\n" );
1084  DEBUG_DEC_INDENT();
1085  return( true );
1086  }
1087 
1088  double calculate_dt( const PP &x, const double *dxdt ) {
1089 
1090  double spd = 0.0, acc = 0.0;
1091 
1092  for( size_t a = 0; a < (PP::size()-1)/2; a++ ) {
1093  //ibsimu.message(MSG_DEBUG_GENERAL,1) << "spd += " << dxdt[2*a]*dxdt[2*a] << "\n";
1094  spd += dxdt[2*a]*dxdt[2*a];
1095  //ibsimu.message(MSG_DEBUG_GENERAL,1) << "acc += " << dxdt[2*a+1]*dxdt[2*a+1] << "\n";
1096  acc += dxdt[2*a+1]*dxdt[2*a+1];
1097  }
1098  if( _pidata._geom->geom_mode() == MODE_CYL ) {
1099  //ibsimu.message(MSG_DEBUG_GENERAL,1) << "MODE_CYL\n";
1100  //ibsimu.message(MSG_DEBUG_GENERAL,1) << "spd += " << x[3]*x[3]*x[5]*x[5] << "\n";
1101  spd += x[3]*x[3]*x[5]*x[5];
1102  //ibsimu.message(MSG_DEBUG_GENERAL,1) << "acc += " << x[3]*x[3]*dxdt[4]*dxdt[4] << "\n";
1103  acc += x[3]*x[3]*dxdt[4]*dxdt[4];
1104  }
1105  //ibsimu.message(MSG_DEBUG_GENERAL,1) << "spd = " << sqrt(spd) << "\n";
1106  //ibsimu.message(MSG_DEBUG_GENERAL,1) << "acc = " << sqrt(acc) << "\n";
1107  spd = _pidata._geom->h() / sqrt(spd);
1108  acc = sqrt( 2.0*_pidata._geom->h() / sqrt(acc) );
1109 
1110  return( spd < acc ? spd : acc );
1111  }
1112 
1113 public:
1114 
1142  ParticleIterator( particle_iterator_type_e type, double epsabs, double epsrel,
1144  uint32_t maxsteps, double maxt, bool save_points,
1145  uint32_t trajdiv, bool mirror[6], MeshScalarField *scharge,
1146  pthread_mutex_t *scharge_mutex,
1147  const VectorField *efield, const VectorField *bfield,
1148  const Geometry *geom )
1149  : _type(type), _intrp(intrp), _scharge_dep(scharge_dep), _epsabs(epsabs), _epsrel(epsrel),
1150  _maxsteps(maxsteps), _maxt(maxt), _save_points(save_points), _trajdiv(trajdiv),
1151  _surface_collision(false), _pidata(scharge,efield,bfield,geom),
1152  _thand_cb(0), _tend_cb(0), _tsur_cb(0), _bsup_cb(0), _pdb(0), _scharge_mutex(scharge_mutex),
1153  _stat(geom->number_of_boundaries()) {
1154 
1155  // Initialize mirroring
1156  _mirror[0] = mirror[0];
1157  _mirror[1] = mirror[1];
1158  _mirror[2] = mirror[2];
1159  _mirror[3] = mirror[3];
1160  _mirror[4] = mirror[4];
1161  _mirror[5] = mirror[5];
1162 
1163  // Initialize system of ordinary differential equations (ODE)
1164  _system.jacobian = NULL;
1165  _system.params = (void *)&_pidata;
1166  _system.function = PP::get_derivatives;
1167  _system.dimension = PP::size()-1; // Time is not part of differential equation dimensions
1168 
1169  // Make scale
1170  // 2D: x vx y vy
1171  // Cyl: x vx r vr omega
1172  // 3D: x vx y vy z vz
1173  double scale_abs[PP::size()-1];
1174  for( uint32_t a = 0; a < (uint32_t)PP::size()-2; a+=2 ) {
1175  scale_abs[a+0] = 1.0;
1176  scale_abs[a+1] = 1.0e6;
1177  }
1178  if( _pidata._geom->geom_mode() == MODE_CYL )
1179  scale_abs[4] = 1.0;
1180 
1181  // Initialize ODE solver
1182  _step = gsl_odeiv_step_alloc( gsl_odeiv_step_rkck, _system.dimension );
1183  //_control = gsl_odeiv_control_standard_new( _epsabs, _epsrel, 1.0, 1.0 );
1184  _control = gsl_odeiv_control_scaled_new( _epsabs, _epsrel, 1.0, 1.0, scale_abs, PP::size()-1 );
1185  _evolve = gsl_odeiv_evolve_alloc( _system.dimension );
1186  }
1187 
1188 
1192  gsl_odeiv_evolve_free( _evolve );
1193  gsl_odeiv_control_free( _control );
1194  gsl_odeiv_step_free( _step );
1195  }
1196 
1197 
1200  void set_surface_collision( bool surface_collision ) {
1201  if( surface_collision && _pidata._geom->geom_mode() == MODE_2D )
1202  throw( Error( ERROR_LOCATION, "2D surface collision not supported" ) );
1203  if( surface_collision && !_pidata._geom->surface_built() )
1204  throw( Error( ERROR_LOCATION, "surface model not built" ) );
1205  _surface_collision = surface_collision;
1206  }
1207 
1208 
1212  _thand_cb = thand_cb;
1213  }
1214 
1215 
1219  _tend_cb = tend_cb;
1220  _pdb = pdb;
1221  }
1222 
1223 
1227  _tsur_cb = tsur_cb;
1228  }
1229 
1230 
1234  _pidata.set_bfield_suppression_callback( bsup_cb );
1235  }
1236 
1239  void set_relativistic( bool enable ) {
1240  _pidata.set_relativistic( enable );
1241  }
1242 
1245  const ParticleStatistics &get_statistics( void ) const {
1246  return( _stat );
1247  }
1248 
1249 
1258  void operator()( Particle<PP> *particle, uint32_t pi ) {
1259 
1260  DEBUG_MESSAGE( "Calculating particle " << pi << "\n" );
1261  DEBUG_INC_INDENT();
1262 
1263  // Check particle status
1264  if( particle->get_status() != PARTICLE_OK ) {
1265  DEBUG_DEC_INDENT();
1266  return;
1267  }
1268 
1269  // Copy starting point to x and
1270  PP x = particle->x();
1271 
1272  // Check particle definition
1273  if( !check_particle_definition( x ) ) {
1274  particle->set_status( PARTICLE_BADDEF );
1275  _stat.inc_end_baddef();
1276  DEBUG_DEC_INDENT();
1277  return;
1278  }
1279  particle->x() = x;
1280 
1281  // Reset trajectory and save first trajectory point.
1282  _traj.clear();
1283  save_trajectory_point( x );
1284  _pidata._qm = particle->qm();
1285  _xi = x;
1286  if( _scharge_dep == SCHARGE_DEPOSITION_LINEAR ) {
1287  _cdpast.clear();
1288  _cdpast.push( x );
1289  }
1290 
1291  // Reset integrator
1292  gsl_odeiv_step_reset( _step );
1293  gsl_odeiv_evolve_reset( _evolve );
1294 
1295  // Make initial guess for step size
1296  double dxdt[PP::size()-1];
1297  PP::get_derivatives( 0.0, &x[1], dxdt, (void *)&_pidata );
1298  double dt = calculate_dt( x, dxdt );
1299 
1300  DEBUG_MESSAGE( "Starting values for iteration:\n" <<
1301  " dxdt = " );
1302  for( size_t a = 0; a < PP::size()-1; a++ )
1303  DEBUG_MESSAGE( dxdt[a] << " " );
1304  DEBUG_MESSAGE( "\n" <<
1305  " dt = " << dt << "\n" );
1306  DEBUG_MESSAGE( "Starting iteration\n" );
1307 
1308  // Iterate ODEs until maximum steps are done, time is used
1309  // or particle collides.
1310  PP x2;
1311  size_t nstp = 0; // Steps taken
1312  while( nstp < _maxsteps && x[0] < _maxt ) {
1313 
1314  // Take a step.
1315  x2 = x;
1316 
1317  DEBUG_MESSAGE( "\nStep *** *** *** ****\n" <<
1318  " x = " << x2 << "\n" <<
1319  " dt = " << dt << " (proposed)\n" );
1320 
1321  while( true ) {
1322  int retval = gsl_odeiv_evolve_apply( _evolve, _control, _step, &_system,
1323  &x2[0], _maxt, &dt, &x2[1] );
1324  if( retval == IBSIMU_DERIV_ERROR ) {
1325  DEBUG_MESSAGE( "Step rejected\n" <<
1326  " x2 = " << x2 << "\n" <<
1327  " dt = " << dt << "\n" );
1328  x2[0] = x[0]; // Reset time (this shouldn't be necessary - there
1329  // is a bug in GSL-1.12, report has been sent)
1330  dt *= 0.5;
1331  if( dt == 0.0 )
1332  throw( Error( ERROR_LOCATION, "too small step size" ) );
1333  //nstp++;
1334  continue;
1335  } else if( retval == GSL_SUCCESS ) {
1336  break;
1337  } else {
1338  throw( Error( ERROR_LOCATION, "gsl_odeiv_evolve_apply failed" ) );
1339  }
1340  }
1341 
1342  // Check step count number and step size validity
1343  if( nstp >= _maxsteps )
1344  break;
1345  if( x2[0] == x[0] )
1346  throw( Error( ERROR_LOCATION, "too small step size" ) );
1347 
1348  // Increase step count.
1349  nstp++;
1350 
1351  DEBUG_MESSAGE( "Step accepted from x1 to x2:\n" <<
1352  " x1 = " << x << "\n" <<
1353  " x2 = " << x2 << "\n" );
1354 
1355  // Handle collisions and space charge of step.
1356  if( !handle_trajectory( *particle, x, x2, false, nstp == 0 ) ) {
1357  x = x2;
1358  break; // Particle done
1359  }
1360 
1361  // Check if particle mirroring is required to avoid
1362  // singularity at symmetry axis.
1363  if( axis_mirror_required( x2 ) ) {
1364  if( !handle_axis_mirror_step( *particle, x, x2 ) )
1365  break; // Particle done
1366  }
1367 
1368  // Propagate coordinates
1369  x = x2;
1370 
1371  // Save trajectory point
1372  save_trajectory_point( x2 );
1373  }
1374 
1375  DEBUG_MESSAGE( "Done iterating\n" );
1376 
1377  // Check if step count or time limited
1378  if( nstp == _maxsteps ) {
1379  particle->set_status( PARTICLE_NSTP );
1380  _stat.inc_end_step();
1381  } else if( x[0] >= _maxt ) {
1382  particle->set_status( PARTICLE_TIME );
1383  _stat.inc_end_time();
1384  }
1385 
1386  // Save step count
1387  _stat.inc_sum_steps( nstp );
1388 
1389  // Save trajectory of current particle
1390  if( _trajdiv != 0 && pi % _trajdiv == 0 )
1391  particle->copy_trajectory( _traj );
1392 
1393  // Save last particle location
1394  particle->x() = x;
1395 
1396  // Call trajectory end callback
1397  if( _tend_cb )
1398  (*_tend_cb)( particle, _pdb );
1399 
1400  DEBUG_DEC_INDENT();
1401  }
1402 
1403 };
1404 
1405 
1406 #endif
1407 
uint32_t inside(const Vec3D &x) const
Return if point is inside solids.
Definition: geometry.cpp:332
Geometry definition
double div_h(void) const
Returns reciprocal of mesh cell size (1/h).
Definition: mesh.hpp:150
static void build_coldata_poly(std::vector< ColData > &coldata, const Mesh &mesh, const PP &x1, const PP &x2)
Find mesh intersections of polynomially interpolated particle trajectory segment. ...
Definition: particleiterator.hpp:182
Basic error class.
Definition: error.hpp:142
void set_relativistic(bool enable)
Set relativistic particle iteration.
Definition: particleiterator.hpp:1239
Trajectory interpolation solver.
void set_trajectory_surface_collision_callback(TrajectorySurfaceCollisionCallback *tsur_cb)
Set trajectory surface collision callback.
Definition: particleiterator.hpp:1226
uint32_t mesh_check(int32_t i) const
Returns number from solid mesh array.
Definition: geometry.cpp:367
Vector field.
Definition: vectorfield.hpp:56
Vertex-based triangle representation.
Definition: vtriangle.hpp:56
bool surface_built(void) const
Is the solid surface representation built?
Definition: geometry.hpp:503
Mat3D inverse(void) const
Return inverse matrix.
Definition: mat3d.cpp:83
MeshScalarField * _scharge
Space charge field or NULL.
Definition: particles.hpp:915
void set_trajectory_handler_callback(TrajectoryHandlerCallback *thand_cb)
Set trajectory handler callback.
Definition: particleiterator.hpp:1211
double bracket_surface(uint32_t n, const Vec3D &xin, const Vec3D &xout, Vec3D &xsurf) const
Find solid n surface location by bracketing.
Definition: geometry.cpp:412
double determinant(void) const
Return determinant of matrix.
Definition: mat3d.cpp:74
Particle and particle point objects
Temporary data bundle for particle iterators.
Definition: particles.hpp:914
Mesh geometry definion.
Definition: mesh.hpp:67
Polynomial interpolation.
Definition: types.hpp:154
Vec3D origo(void) const
Returns origo vector of geometry.
Definition: mesh.hpp:128
void set_bfield_suppression_callback(const CallbackFunctorD_V *bsup_cb)
Set B-field potential dependent suppression callback.
Definition: particles.hpp:930
Trajectory representation between two calculated points in 1d.
Definition: trajectory.hpp:65
Deposition to nodes as a linear function of distance to closet trajectory segment.
Definition: types.hpp:163
2D geometry
Definition: types.hpp:61
First-in first-out container with cyclic memory.
Trajectory end callback.
Definition: particledatabase.hpp:71
Particle databases
Polynomial solver.
const VTriangle & surface_triangle(int32_t a) const
Return reference to surface triangle a.
Definition: geometry.hpp:540
void set_relativistic(bool enable)
Set relativistic particle iteration.
Definition: particles.hpp:936
Particle iterator class for continuous Vlasov-type iteration.
Definition: particleiterator.hpp:280
Scalar field class.
Definition: meshscalarfield.hpp:70
Compatibility math.
Trajectory handler callback.
Definition: particledatabase.hpp:58
static void build_coldata_linear(std::vector< ColData > &coldata, const Mesh &mesh, const PP &x1, const PP &x2)
Find mesh intersections of linearly interpolated particle trajectory segment.
Definition: particleiterator.hpp:132
Geometry defining class.
Definition: geometry.hpp:179
Particle-in-cell type deposition to neighbouring nodes in each cell.
Definition: types.hpp:162
const uint32_t & mesh(int32_t i) const
Returns a const reference to solid mesh array.
Definition: geometry.hpp:399
PP _x
Mesh intersection coordinates.
Definition: particleiterator.hpp:104
const Vec3D & surface_vertex(int32_t a) const
Return reference to surface vertex a.
Definition: geometry.hpp:534
#define ERROR_LOCATION
Macro for setting error location when throwing errors.
Definition: error.hpp:72
Job scheduler for parallel processing.
ParticleIterator(particle_iterator_type_e type, double epsabs, double epsrel, trajectory_interpolation_e intrp, scharge_deposition_e scharge_dep, uint32_t maxsteps, double maxt, bool save_points, uint32_t trajdiv, bool mirror[6], MeshScalarField *scharge, pthread_mutex_t *scharge_mutex, const VectorField *efield, const VectorField *bfield, const Geometry *geom)
Constructor for new particle iterator.
Definition: particleiterator.hpp:1142
double h(void) const
Returns mesh cell size.
Definition: mesh.hpp:146
PP & x()
Return reference to coordinate data.
Definition: particles.hpp:798
Space charge deposition functions.
particle_status_e get_status()
Return particle status.
Definition: particles.hpp:695
Definition: callback.hpp:61
~ParticleIterator()
Destructor.
Definition: particleiterator.hpp:1191
uint32_t surface_trianglec(void) const
Return total surface triangle count.
Definition: geometry.hpp:528
double qm() const
Return charge per mass ratio (q/m) [C/kg].
Definition: particles.hpp:718
void set_bfield_suppression_callback(const CallbackFunctorD_V *bsup_cb)
Set B-field potential dependent suppression callback.
Definition: particleiterator.hpp:1233
void operator()(Particle< PP > *particle, uint32_t pi)
Iterate a particle from start to end.
Definition: particleiterator.hpp:1258
Particle class in some geometry.
Definition: particles.hpp:739
ColData(PP x, int dir)
Constructor for collision at x into direction dir.
Definition: particleiterator.hpp:114
void set_status(particle_status_e status)
Set particle status.
Definition: particles.hpp:699
scharge_deposition_e
Space charge depostition type.
Definition: types.hpp:161
Particle iteration statistics.
Definition: particlestatistics.hpp:57
Vec3D max(void) const
Returns vector pointing to the last mesh point opposite of origo.
Definition: mesh.hpp:137
Mesh intersection (collision) coordinate data
Definition: particleiterator.hpp:102
Error class for memory allocation errors.
Definition: error.hpp:185
uint32_t surface_triangle_ptr(int32_t i, int32_t j, int32_t k) const
Return index of surface triangle at mesh cube (i,j,k).
Definition: geometry.hpp:546
const ParticleStatistics & get_statistics(void) const
Get particle iterator statistics.
Definition: particleiterator.hpp:1245
particle_iterator_type_e
Particle iterator type.
Definition: particleiterator.hpp:87
void coord(double &x, double &v, double K)
Calculate location x and velocity v at parametric time K.
Definition: trajectory.cpp:152
Three-by-three matrices.
double _qm
Precalculated q/m.
Definition: particles.hpp:919
bool operator<(const ColData &cd) const
Compare coldata entry times.
Definition: particleiterator.hpp:120
double IQ() const
Return current or charge carried by trajectory or particle cloud [A/C].
Definition: particles.hpp:706
int solve(double K[3], double x, int extrapolate=0)
Solves for trajectory intersection with location.
Definition: trajectory.cpp:186
void construct(double dt, double x1, double v1, double x2, double v2, trajectory_rep_e force=TRAJ_EMPTY)
Construct representation of trajectory from (x1,v1) to (x2,v2) in time dt.
Definition: trajectory.cpp:53
Int3D size(void) const
Returns size array of geometry.
Definition: mesh.hpp:116
const Geometry * _geom
Geometry.
Definition: particles.hpp:918
void scharge_add_from_trajectory_pic(MeshScalarField &scharge, pthread_mutex_t *mutex, double I, const ParticleP2D &x1, const ParticleP2D &x2)
Function for adding charge to space charge density map from particle trajectory in 2d simulation...
Definition: scharge.cpp:171
void set_surface_collision(bool surface_collision)
Enable/disable surface collision model.
Definition: particleiterator.hpp:1200
Three-by-three matrix.
Definition: mat3d.hpp:58
Cylindrically symmetric geometry.
Definition: types.hpp:62
Three dimensional vector.
Definition: vec3d.hpp:58
void copy_trajectory(const std::vector< PP > &traj)
Define trajectory by copying.
Definition: particles.hpp:822
void set_trajectory_end_callback(TrajectoryEndCallback *tend_cb, ParticleDataBase *pdb)
Set trajectory end callback.
Definition: particleiterator.hpp:1218
ColData()
Default constructor.
Definition: particleiterator.hpp:110
Trajectory surface collision callback.
Definition: particledatabase.hpp:86
int _dir
Direction of particle at intersection. i: -1/+1, j: -2/+2, k: -3:/+3.
Definition: particleiterator.hpp:105
geom_mode_e geom_mode(void) const
Returns geometry mode.
Definition: mesh.hpp:108
Particle database base class.
Definition: particledatabase.hpp:191
trajectory_interpolation_e
Trajectory interpolation type.
Definition: types.hpp:153
uint32_t surface_inside(const Vec3D &x) const
Finds if point is inside surface triangulation.
Definition: geometry.cpp:2224
Mesh based scalar fields.
Ion Beam Simulator global settings.
Vector field base.


Reference manual for Ion Beam Simulator 1.0.5new_solver
Generated by Doxygen 1.8.5 on Tue May 19 2015 09:15:42.