Extension: Geometries
#include "dg/geometries/geometries.h"
fieldaligned_tmp.h
Go to the documentation of this file.
1#pragma once
2#include <cmath>
3#include <array>
4#include <cusp/csr_matrix.h>
5
6#include "dg/algorithm.h"
7#include "magnetic_field.h"
8#include "fluxfunctions.h"
9#include "curvilinear.h"
10
11namespace dg{
12namespace geo{
13
17{
18 einsPlus = 0,
19 einsPlusT,
20 einsMinus,
22 zeroPlus,
23 zeroMinus,
24 zeroPlusT,
27};
28
31typedef ONE FullLimiter;
32
35typedef ZERO NoLimiter;
37namespace detail{
38
39
40struct DSFieldCylindrical3
41{
42 DSFieldCylindrical3( const dg::geo::CylindricalVectorLvl0& v): m_v(v){}
43 void operator()( double t, const std::array<double,3>& y,
44 std::array<double,3>& yp) const {
45 double R = y[0], Z = y[1];
46 double vz = m_v.z()(R, Z);
47 yp[0] = m_v.x()(R, Z)/vz;
48 yp[1] = m_v.y()(R, Z)/vz;
49 yp[2] = 1./vz;
50 }
51 private:
53};
54
55struct DSFieldCylindrical4
56{
57 DSFieldCylindrical4( const dg::geo::CylindricalVectorLvl1& v): m_v(v){}
58 void operator()( double t, const std::array<double,3>& y,
59 std::array<double,3>& yp) const {
60 double R = y[0], Z = y[1];
61 double vx = m_v.x()(R,Z);
62 double vy = m_v.y()(R,Z);
63 double vz = m_v.z()(R,Z);
64 double divvvz = m_v.divvvz()(R,Z);
65 yp[0] = vx/vz;
66 yp[1] = vy/vz;
67 yp[2] = divvvz*y[2];
68 }
69
70 private:
72};
73
74struct DSField
75{
76 DSField() = default;
77 //z component of v may not vanish
78 DSField( const dg::geo::CylindricalVectorLvl1& v,
79 const dg::aGeometry2d& g ):
80 m_g(g)
81 {
82 dg::HVec v_zeta, v_eta;
83 dg::pushForwardPerp( v.x(), v.y(), v_zeta, v_eta, g);
84 dg::HVec vx = dg::pullback( v.x(), g);
85 dg::HVec vy = dg::pullback( v.y(), g);
86 dg::HVec vz = dg::pullback( v.z(), g);
87 dg::HVec divvvz = dg::pullback( v.divvvz(), g);
88 dg::blas1::pointwiseDivide(v_zeta, vz, v_zeta);
89 dg::blas1::pointwiseDivide(v_eta, vz, v_eta);
90 dzetadphi_ = dg::forward_transform( v_zeta, g );
91 detadphi_ = dg::forward_transform( v_eta, g );
92 dvdphi_ = dg::forward_transform( divvvz, g );
93 }
94 //interpolate the vectors given in the constructor on the given point
95 void operator()(double t, const std::array<double,3>& y, std::array<double,3>& yp) const
96 {
97 // shift point into domain
98 yp[0] = interpolate(dg::lspace, dzetadphi_, y[0], y[1], *m_g);
99 yp[1] = interpolate(dg::lspace, detadphi_, y[0], y[1], *m_g);
100 yp[2] = interpolate(dg::lspace, dvdphi_, y[0], y[1], *m_g)*y[2];
101 }
102 private:
103 thrust::host_vector<double> dzetadphi_, detadphi_, dvdphi_;
105};
106
107//used in constructor of Fieldaligned
108template<class real_type>
109void integrate_all_fieldlines2d( const dg::geo::CylindricalVectorLvl1& vec,
110 const dg::aRealGeometry2d<real_type>& grid_field,
111 const dg::aRealTopology2d<real_type>& grid_evaluate,
112 std::array<thrust::host_vector<real_type>,3>& yp,
113 const thrust::host_vector<double>& vol0,
114 thrust::host_vector<real_type>& yp2b,
115 thrust::host_vector<bool>& in_boxp,
116 real_type deltaPhi, real_type eps)
117{
118 //grid_field contains the global geometry for the field and the boundaries
119 //grid_evaluate contains the points to actually integrate
120 std::array<thrust::host_vector<real_type>,3> y{
121 dg::evaluate( dg::cooX2d, grid_evaluate),
122 dg::evaluate( dg::cooY2d, grid_evaluate),
123 vol0
124 };
125 yp.fill(dg::evaluate( dg::zero, grid_evaluate));
126 //construct field on high polynomial grid, then integrate it
127 dg::geo::detail::DSField field;
128 if( !dynamic_cast<const dg::CartesianGrid2d*>( &grid_field))
129 field = dg::geo::detail::DSField( vec, grid_field);
130
131 //field in case of cartesian grid
132 dg::geo::detail::DSFieldCylindrical4 cyl_field(vec);
133 const unsigned size = grid_evaluate.size();
135 "Dormand-Prince-7-4-5", std::array<real_type,3>{0,0,0});
137 if( dynamic_cast<const dg::CartesianGrid2d*>( &grid_field))
139 cyl_field, dg::pid_control, dg::fast_l2norm, eps, 1e-10);
140 else
142 field, dg::pid_control, dg::fast_l2norm, eps, 1e-10);
143
144 for( unsigned i=0; i<size; i++)
145 {
146 std::array<real_type,3> coords{y[0][i],y[1][i],y[2][i]}, coordsP;
147 //x,y,s
148 real_type phi1 = deltaPhi;
149 odeint.set_dt( deltaPhi/2.);
150 odeint.integrate( 0, coords, phi1, coordsP);
151 yp[0][i] = coordsP[0], yp[1][i] = coordsP[1], yp[2][i] = coordsP[2];
152 }
153 yp2b.assign( grid_evaluate.size(), deltaPhi); //allocate memory for output
154 in_boxp.resize( yp2b.size());
155 //Now integrate again but this time find the boundary distance
156 for( unsigned i=0; i<size; i++)
157 {
158 std::array<real_type,3> coords{y[0][i],y[1][i],y[2][i]}, coordsP;
159 in_boxp[i] = grid_field.contains( yp[0][i], yp[1][i]) ? true : false;
160 if( false == in_boxp[i])
161 {
162 //x,y,s
163 real_type phi1 = deltaPhi;
164 odeint.integrate_in_domain( 0., coords, phi1, coordsP, 0., (const
165 dg::aRealTopology2d<real_type>&)grid_field, eps);
166 yp2b[i] = phi1;
167 }
168 }
169}
170
171
172}//namespace detail
174
181struct WallFieldlineDistance : public aCylindricalFunctor<WallFieldlineDistance>
182{
199 const dg::aRealTopology2d<double>& domain,
200 double maxPhi, double eps, std::string type) :
201 m_domain( domain), m_cyl_field(vec),
202 m_deltaPhi( maxPhi), m_eps( eps), m_type(type)
203 {
204 if( m_type != "phi" && m_type != "s")
205 throw std::runtime_error( "Distance type "+m_type+" not recognized!\n");
206 }
217 double do_compute( double R, double Z) const
218 {
219 std::array<double,3> coords{ R, Z, 0}, coordsP(coords);
220 // determine sign
221 m_cyl_field( 0., coords, coordsP);
222 double sign = coordsP[2] > 0 ? +1. : -1.;
223 double phi1 = sign*m_deltaPhi; // we integrate negative ...
224 try{
226 "Dormand-Prince-7-4-5", coords);
228 m_cyl_field, dg::pid_control, dg::fast_l2norm, m_eps,
229 1e-10);
230 odeint.integrate_in_domain( 0., coords, phi1, coordsP, 0.,
231 m_domain, m_eps);
232 //integration
233 }catch (std::exception& e)
234 {
235 // if not possible the distance is large
236 //std::cerr << e.what();
237 phi1 = sign*m_deltaPhi;
238 coordsP[2] = 1e6*phi1;
239 }
240 if( m_type == "phi")
241 return sign*phi1;
242 return coordsP[2];
243 }
244
245 private:
246 const dg::Grid2d m_domain;
247 dg::geo::detail::DSFieldCylindrical3 m_cyl_field;
248 double m_deltaPhi, m_eps;
249 std::string m_type;
250};
251
264struct WallFieldlineCoordinate : public aCylindricalFunctor<WallFieldlineCoordinate>
265{
269 const dg::aRealTopology2d<double>& domain,
270 double maxPhi, double eps, std::string type) :
271 m_domain( domain), m_cyl_field(vec),
272 m_deltaPhi( maxPhi), m_eps( eps), m_type(type)
273 {
274 if( m_type != "phi" && m_type != "s")
275 throw std::runtime_error( "Distance type "+m_type+" not recognized!\n");
276 }
277 double do_compute( double R, double Z) const
278 {
279 double phiP = m_deltaPhi, phiM = -m_deltaPhi;
280 std::array<double,3> coords{ R, Z, 0}, coordsP(coords), coordsM(coords);
281 // determine sign
282 m_cyl_field( 0., coords, coordsP);
283 double sign = coordsP[2] > 0 ? +1. : -1.;
284 try{
286 dg::Adaptive<dg::ERKStep<std::array<double,3>>>(
287 "Dormand-Prince-7-4-5", coords), m_cyl_field,
288 dg::pid_control, dg::fast_l2norm, m_eps, 1e-10);
289 odeint.integrate_in_domain( 0., coords, phiP, coordsP, 0.,
290 m_domain, m_eps);
291 odeint.integrate_in_domain( 0., coords, phiM, coordsM, 0.,
292 m_domain, m_eps);
293 }catch (std::exception& e)
294 {
295 // if not possible the distance is large
296 phiP = m_deltaPhi;
297 coordsP[2] = 1e6*phiP;
298 phiM = -m_deltaPhi;
299 coordsM[2] = 1e6*phiM;
300 }
301 if( m_type == "phi")
302 return sign*(-phiP-phiM)/(phiP-phiM);
303 double sP = coordsP[2], sM = coordsM[2];
304 double value = sign*(-sP-sM)/(sP-sM);
305 if( (phiM <= -m_deltaPhi) && (phiP >= m_deltaPhi))
306 return 0.; //return exactly zero if sheath not reached
307 if( (phiM <= -m_deltaPhi))
308 return value*sign > 0 ? value : 0.; // avoid false negatives
309 if( (phiP >= m_deltaPhi))
310 return value*sign < 0 ? value : 0.; // avoid false positives
311 return value;
312 }
313
314 private:
315 const dg::Grid2d m_domain;
316 dg::geo::detail::DSFieldCylindrical3 m_cyl_field;
317 double m_deltaPhi, m_eps;
318 std::string m_type;
319};
320
321
360
372template<class ProductGeometry, class IMatrix, class container >
373struct Fieldaligned
374{
375
381 template <class Limiter>
383 const ProductGeometry& grid,
386 Limiter limit = FullLimiter(),
387 double eps = 1e-5,
388 unsigned mx=10, unsigned my=10,
389 double deltaPhi = -1,
390 std::string interpolation_method = "dg",
391 bool benchmark=true
392 ):
393 Fieldaligned( dg::geo::createBHat(vec),
394 grid, bcx, bcy, limit, eps, mx, my, deltaPhi, interpolation_method)
395 {
396 }
397
401 template <class Limiter>
403 const ProductGeometry& grid,
406 Limiter limit = FullLimiter(),
407 double eps = 1e-5,
408 unsigned mx=10, unsigned my=10,
409 double deltaPhi = -1,
410 std::string interpolation_method = "dg",
411 bool benchmark=true
412 );
418 template<class ...Params>
419 void construct( Params&& ...ps)
420 {
421 //construct and swap
422 *this = Fieldaligned( std::forward<Params>( ps)...);
423 }
424
425 dg::bc bcx()const{
426 return m_bcx;
427 }
428 dg::bc bcy()const{
429 return m_bcy;
430 }
431
432
442 void set_boundaries( dg::bc bcz, double left, double right)
443 {
444 m_bcz = bcz;
445 const dg::Grid1d g2d( 0, 1, 1, m_perp_size);
446 m_left = dg::evaluate( dg::CONSTANT(left), g2d);
447 m_right = dg::evaluate( dg::CONSTANT(right),g2d);
448 }
449
459 void set_boundaries( dg::bc bcz, const container& left, const container& right)
460 {
461 m_bcz = bcz;
462 m_left = left;
463 m_right = right;
464 }
465
476 void set_boundaries( dg::bc bcz, const container& global, double scal_left, double scal_right)
477 {
478 dg::split( global, m_f, *m_g);
479 dg::blas1::axpby( scal_left, m_f[0], 0, m_left);
480 dg::blas1::axpby( scal_right, m_f[m_Nz-1], 0, m_right);
481 m_bcz = bcz;
482 }
483
492 void operator()(enum whichMatrix which, const container& in, container& out);
493
494 double deltaPhi() const{return m_deltaPhi;}
497 const container& hbm()const {
498 return m_hbm;
499 }
502 const container& hbp()const {
503 return m_hbp;
504 }
507 const container& sqrtG()const {
508 return m_G;
509 }
512 const container& sqrtGm()const {
513 return m_Gm;
514 }
517 const container& sqrtGp()const {
518 return m_Gp;
519 }
522 const container& bphi()const {
523 return m_bphi;
524 }
527 const container& bphiM()const {
528 return m_bphiM;
529 }
532 const container& bphiP()const {
533 return m_bphiP;
534 }
537 const container& bbm()const {
538 return m_bbm;
539 }
542 const container& bbo()const {
543 return m_bbo;
544 }
547 const container& bbp()const {
548 return m_bbp;
549 }
551 const ProductGeometry& grid()const{return *m_g;}
552
566 container interpolate_from_coarse_grid( const ProductGeometry& grid_coarse, const container& coarse);
576 void integrate_between_coarse_grid( const ProductGeometry& grid_coarse, const container& coarse, container& out );
577
578
607 template< class BinaryOp, class UnaryOp>
608 container evaluate( BinaryOp binary, UnaryOp unary,
609 unsigned p0, unsigned rounds) const;
610
611 std::string method() const{return m_interpolation_method;}
612
613 private:
614 void ePlus( enum whichMatrix which, const container& in, container& out);
615 void eMinus(enum whichMatrix which, const container& in, container& out);
616 void zero( enum whichMatrix which, const container& in, container& out);
617 IMatrix m_plus, m_minus, m_plusT, m_minusT; //2d interpolation matrices
618 container m_tmp; // 3d size
619 IMatrix m_back, m_forw;
620 IMatrix m_stencil, m_stencilY;
621 container m_hbm, m_hbp; //3d size
622 container m_G, m_Gm, m_Gp; // 3d size
623 container m_bphi, m_bphiM, m_bphiP; // 3d size
624 container m_bbm, m_bbp, m_bbo; //3d size masks
625
626 container m_left, m_right, m_temp2d; //perp_size
627 container m_limiter; //perp_size
628 container m_ghostM, m_ghostP; //perp_size
629 unsigned m_Nz, m_perp_size;
630 dg::bc m_bcx, m_bcy, m_bcz;
631 std::vector<dg::View<const container>> m_f;
632 std::vector<dg::View< container>> m_temp, m_temp2;
635 double m_deltaPhi;
636 std::string m_interpolation_method;
637
638 bool m_have_adjoint = false;
639 void updateAdjoint( )
640 {
641 m_plusT = dg::transpose( m_plus);
642 m_minusT = dg::transpose( m_minus);
643 m_have_adjoint = true;
644 }
645};
646
648
650template<class Geometry, class IMatrix, class container>
651template <class Limiter>
654 const Geometry& grid,
655 dg::bc bcx, dg::bc bcy, Limiter limit, double eps,
656 unsigned mx, unsigned my, double deltaPhi, std::string interpolation_method, bool benchmark)
657{
658 m_stencil = dg::create::limiter_stencil( dg::coo3d::x, grid, grid.bcx());
659 m_stencilY = dg::create::limiter_stencil( dg::coo3d::y, grid, grid.bcy());
661 if( (grid.bcx() == PER && bcx != PER) || (grid.bcx() != PER && bcx == PER) )
662 throw( dg::Error(dg::Message(_ping_)<<"Fieldaligned: Got conflicting periodicity in x. The grid says "<<bc2str(grid.bcx())<<" while the parameter says "<<bc2str(bcx)));
663 if( (grid.bcy() == PER && bcy != PER) || (grid.bcy() != PER && bcy == PER) )
664 throw( dg::Error(dg::Message(_ping_)<<"Fieldaligned: Got conflicting boundary conditions in y. The grid says "<<bc2str(grid.bcy())<<" while the parameter says "<<bc2str(bcy)));
665 m_Nz=grid.Nz(), m_bcx = bcx, m_bcy = bcy, m_bcz=grid.bcz();
666 m_g.reset(grid);
667 if( deltaPhi <=0) deltaPhi = grid.hz();
668 //if( interpolation_method != "dg")
669 {
670 m_back = dg::create::inv_backproject( grid); //from equidist to dg
671 m_forw = dg::create::backproject( grid); // from dg to equidist
672 }
674 dg::Timer t;
675 if( benchmark) t.tic();
676 dg::ClonePtr<dg::aGeometry2d> grid_original( grid.perp_grid()) ;
677 dg::ClonePtr<dg::aGeometry3d> grid_original3d( grid) ;
678 dg::ClonePtr<dg::aGeometry2d> grid_transform( grid_original) ;
679 dg::ClonePtr<dg::aGeometry3d> grid_transform3d( grid) ;
680 //if( interpolation_method == "dg-const")
681 //{
682 // grid_transform->set( 1, grid.gx().size(), grid.gy().size());
683 // grid_transform3d->set( 1, grid.gx().size(), grid.gy().size(), grid.gz().size());
684 //}
685 dg::ClonePtr<dg::aGeometry2d> grid_magnetic = grid_original;//INTEGRATE HIGH ORDER GRID
686 grid_magnetic->set( grid_original->n() < 3 ? 4 : 7, grid_magnetic->Nx(), grid_magnetic->Ny());
688 if( interpolation_method == "dg-const" || interpolation_method == "dg")
689 {
690 grid_fine = grid_original;
691 if( interpolation_method == "dg-const")
692 {
693 //mx *= grid.n(), my *= grid.n();
694 //grid_original->set( 1, grid.gx().size(), grid.gy().size());
695 //grid_original3d->set( 1, grid.gx().size(), grid.gy().size(), grid.gz().size());
696 //m_inv_linear = dg::create::inv_fem_linear2const2d( *grid_original3d);
697 m_forw = dg::create::smoothing( *grid_transform);
698 }
699 grid_fine->multiplyCellNumbers((double)mx, (double)my);
700 }
701 else
702 {
703 grid_original->set( 1, grid.gx().size(), grid.gy().size());
704 grid_original3d->set( 1, grid.gx().size(), grid.gy().size(), grid.gz().size());
705
706 if( interpolation_method == "linear")
707 {
708 m_inv_linear = dg::create::inv_fem_mass2d( *grid_original3d);
709 dg::FemRefinement fem_refX( mx), fem_refY(my);
710 dg::Grid2d tmp( *grid_original);
711 dg::CartesianRefinedGrid2d ref( fem_refX, fem_refY, tmp.x0(), tmp.x1(),
712 tmp.y0(), tmp.y1(), tmp.n(), tmp.Nx(), tmp.Ny(),
713 tmp.bcx(), tmp.bcy());
714 grid_fine = ref;
715 }
716 if ( interpolation_method == "linear-const")
717 {
718 grid_fine = grid_original;
719 grid_fine->multiplyCellNumbers((double)mx, (double)my);
720 m_inv_linear = dg::create::inv_fem_linear2const2d( *grid_original3d);
721 }
722 }
723 if( benchmark)
724 {
725 t.toc();
726 std::cout << "# DS: High order grid gen took: "<<t.diff()<<"\n";
727 t.tic();
728 }
730 std::array<thrust::host_vector<double>,3> yp_trafo, ym_trafo, yp, ym;
731 thrust::host_vector<bool> in_boxp, in_boxm;
732 thrust::host_vector<double> hbp, hbm;
733 thrust::host_vector<double> vol = dg::tensor::volume(grid_transform3d->metric()), vol2d0;
734 auto vol2d = dg::split( vol, *grid_transform3d);
735 dg::assign( vol2d[0], vol2d0);
736 detail::integrate_all_fieldlines2d( vec, *grid_magnetic, *grid_transform,
737 yp_trafo, vol2d0, hbp, in_boxp, deltaPhi, eps);
738 detail::integrate_all_fieldlines2d( vec, *grid_magnetic, *grid_transform,
739 ym_trafo, vol2d0, hbm, in_boxm, -deltaPhi, eps);
740 dg::HVec Xf = dg::pullback( dg::cooX2d, *grid_fine);
741 dg::HVec Yf = dg::pullback( dg::cooY2d, *grid_fine);
742 {
744 *grid_transform, dg::NEU, dg::NEU, grid_transform->n() < 3 ? "cubic" : "dg"); //INTERPOLATE TO FINE GRID
745 yp.fill(dg::evaluate( dg::zero, *grid_fine));
746 ym = yp;
747 for( int i=0; i<2; i++) //only R and Z get interpolated
748 {
749 dg::blas2::symv( interpolate, yp_trafo[i], yp[i]);
750 dg::blas2::symv( interpolate, ym_trafo[i], ym[i]);
751 }
752 } // release memory for interpolate matrix
753 if( benchmark)
754 {
755 t.toc();
756 std::cout << "# DS: Computing all points took: "<<t.diff()<<"\n";
757 t.tic();
758 }
760 dg::IHMatrix plusFine, minusFine, plus, minus;
761 if( interpolation_method == "dg-const" || interpolation_method == "dg")
762 {
763 plusFine = dg::create::interpolation( yp[0], yp[1],
764 *grid_transform, bcx, bcy, "dg");
765 minusFine = dg::create::interpolation( ym[0], ym[1],
766 *grid_transform, bcx, bcy, "dg");
767 }
768 else
769 {
770 dg::IHMatrix plusFineTmp = dg::create::interpolation( yp[0], yp[1],
771 *grid_original, bcx, bcy, "linear");
772 dg::IHMatrix minusFineTmp = dg::create::interpolation( ym[0], ym[1],
773 *grid_original, bcx, bcy, "linear");
774 dg::IHMatrix forw = dg::create::backproject( *grid_transform); // from dg to equidist
775 cusp::multiply( plusFineTmp, forw, plusFine);
776 cusp::multiply( minusFineTmp, forw, minusFine);
777 }
778 if( mx == 1 && my == 1)
779 {
780 plus = plusFine;
781 minus = minusFine;
782 }
783 else
784 {
786 if( interpolation_method == "dg-const")
787 {
788 //dg::IHMatrix proj = dg::create::transformation( *grid_original, *grid_fine);
789 //auto back = dg::create::inv_backproject( *grid_transform);
790 //cusp::multiply( back, proj, projection);
791 projection = dg::create::projection( *grid_original, *grid_fine);
792 }
793 else if( interpolation_method == "dg")
794 projection = dg::create::projection( *grid_original, *grid_fine);
795 else if( interpolation_method == "linear")
796 {
798 Yf, *grid_original, dg::NEU, dg::NEU, "linear"));
801 dg::IHMatrix tmp;
802 cusp::multiply( projection, Wf, tmp);
803 cusp::multiply( Vf, tmp, Wf);
804 auto back = dg::create::inv_backproject( *grid_transform);
805 cusp::multiply( back, Wf, projection);
806 }
807 else
808 {
809 dg::IHMatrix proj = dg::create::projection( *grid_original, *grid_fine);
810 auto back = dg::create::inv_backproject( *grid_transform);
811 cusp::multiply( back, proj, projection);
812 }
813 cusp::multiply( projection, plusFine, plus);
814 cusp::multiply( projection, minusFine, minus);
815 }
816 if( benchmark)
817 {
818 t.toc();
819 std::cout << "# DS: Multiplication PI took: "<<t.diff()<<"\n";
820 }
821 dg::blas2::transfer( plus, m_plus);
822 dg::blas2::transfer( minus, m_minus);
824 dg::HVec hbphi( yp_trafo[2]), hbphiP(hbphi), hbphiM(hbphi);
825 hbphi = dg::pullback( vec.z(), *grid_transform);
826 //this is a pullback bphi( R(zeta, eta), Z(zeta, eta)):
827 if( dynamic_cast<const dg::CartesianGrid2d*>( grid_transform.get()))
828 {
829 for( unsigned i=0; i<hbphiP.size(); i++)
830 {
831 hbphiP[i] = vec.z()(yp_trafo[0][i], yp_trafo[1][i]);
832 hbphiM[i] = vec.z()(ym_trafo[0][i], ym_trafo[1][i]);
833 }
834 }
835 else
836 {
837 dg::HVec Ihbphi = dg::pullback( vec.z(), *grid_magnetic);
838 dg::HVec Lhbphi = dg::forward_transform( Ihbphi, *grid_magnetic);
839 for( unsigned i=0; i<yp_trafo[0].size(); i++)
840 {
841 hbphiP[i] = dg::interpolate( dg::lspace, Lhbphi, yp_trafo[0][i],
842 yp_trafo[1][i], *grid_magnetic);
843 hbphiM[i] = dg::interpolate( dg::lspace, Lhbphi, ym_trafo[0][i],
844 ym_trafo[1][i], *grid_magnetic);
845 }
846 }
847 dg::assign3dfrom2d( hbphi, m_bphi, *grid_transform3d);
848 dg::assign3dfrom2d( hbphiM, m_bphiM, *grid_transform3d);
849 dg::assign3dfrom2d( hbphiP, m_bphiP, *grid_transform3d);
850
851 dg::assign3dfrom2d( yp_trafo[2], m_Gp, *grid_transform3d);
852 dg::assign3dfrom2d( ym_trafo[2], m_Gm, *grid_transform3d);
853 // The weights don't matter since they fall out in Div and Lap anyway
854 // But they are good for testing
855 m_G = vol;
856 container weights = dg::create::weights( *grid_transform3d);
858 dg::blas1::pointwiseDot( m_Gp, weights, m_Gp);
859 dg::blas1::pointwiseDot( m_Gm, weights, m_Gm);
860
861 dg::assign( dg::evaluate( dg::zero, *grid_transform3d), m_hbm);
862 m_f = dg::split( (const container&)m_hbm, *grid_transform3d);
863 m_temp = dg::split( m_hbm, *grid_transform3d);
864 m_temp2 = dg::split( m_hbm, *grid_transform3d);
865 dg::assign3dfrom2d( hbp, m_hbp, *grid_transform3d);
866 dg::assign3dfrom2d( hbm, m_hbm, *grid_transform3d);
867 dg::blas1::scal( m_hbm, -1.);
868
870 thrust::host_vector<double> bbm( in_boxp.size(),0.), bbo(bbm), bbp(bbm);
871 for( unsigned i=0; i<in_boxp.size(); i++)
872 {
873 if( !in_boxp[i] && !in_boxm[i])
874 bbo[i] = 1.;
875 else if( !in_boxp[i] && in_boxm[i])
876 bbp[i] = 1.;
877 else if( in_boxp[i] && !in_boxm[i])
878 bbm[i] = 1.;
879 // else all are 0
880 }
881 dg::assign3dfrom2d( bbm, m_bbm, *grid_transform3d);
882 dg::assign3dfrom2d( bbo, m_bbo, *grid_transform3d);
883 dg::assign3dfrom2d( bbp, m_bbp, *grid_transform3d);
884
885 m_deltaPhi = deltaPhi; // store for evaluate
886 m_tmp = m_bbm;
887 m_interpolation_method = interpolation_method;
889
890 m_perp_size = grid_transform->size();
891 dg::assign( dg::pullback(limit, *grid_transform), m_limiter);
892 dg::assign( dg::evaluate(dg::zero, *grid_transform), m_left);
893 m_temp2d = m_ghostM = m_ghostP = m_right = m_left;
894}
895
896
897template<class G, class I, class container>
899 const G& grid, const container& in)
900{
901 //I think we need grid as input to split input vector and we need to interpret
902 //the grid nodes as node centered not cell-centered!
903 //idea: apply I+/I- cphi - 1 times in each direction and then apply interpolation formula
904 assert( m_g->Nz() % grid.Nz() == 0);
905 unsigned Nz_coarse = grid.Nz(), Nz = m_g->Nz();
906 unsigned cphi = Nz / Nz_coarse;
907
908 container out = dg::evaluate( dg::zero, *m_g);
909 container helper = dg::evaluate( dg::zero, *m_g);
910 dg::split( helper, m_temp, *m_g);
911 std::vector<dg::View< container>> out_split = dg::split( out, *m_g);
912 std::vector<dg::View< const container>> in_split = dg::split( in, grid);
913 for ( int i=0; i<(int)Nz_coarse; i++)
914 {
915 //1. copy input vector to appropriate place in output
916 dg::blas1::copy( in_split[i], out_split[i*cphi]);
917 dg::blas1::copy( in_split[i], m_temp[i*cphi]);
918 }
919 //Step 1 needs to finish so that m_temp contains values everywhere
920 //2. Now apply plus and minus T to fill in the rest
921 for ( int i=0; i<(int)Nz_coarse; i++)
922 {
923 for( int j=1; j<(int)cphi; j++)
924 {
926 dg::blas2::symv( m_minus, out_split[i*cphi+j-1], out_split[i*cphi+j]);
928 dg::blas2::symv( m_plus, m_temp[(i*cphi+cphi+1-j)%Nz], m_temp[i*cphi+cphi-j]);
929 }
930 }
931 //3. Now add up with appropriate weights
932 for( int i=0; i<(int)Nz_coarse; i++)
933 for( int j=1; j<(int)cphi; j++)
934 {
935 double alpha = (double)(cphi-j)/(double)cphi;
936 double beta = (double)j/(double)cphi;
937 dg::blas1::axpby( alpha, out_split[i*cphi+j], beta, m_temp[i*cphi+j], out_split[i*cphi+j]);
938 }
939 return out;
940}
941template<class G, class I, class container>
942void Fieldaligned<G, I,container>::integrate_between_coarse_grid( const G& grid, const container& in, container& out)
943{
944 assert( m_g->Nz() % grid.Nz() == 0);
945 unsigned Nz_coarse = grid.Nz(), Nz = m_g->Nz();
946 unsigned cphi = Nz / Nz_coarse;
947
948 out = in;
949 container helperP( in), helperM(in), tempP(in), tempM(in);
950
951 //1. Apply plus and minus T and sum up
952 for( int j=1; j<(int)cphi; j++)
953 {
955 dg::blas2::symv( m_minus, helperP, tempP);
956 dg::blas1::axpby( (double)(cphi-j)/(double)cphi, tempP, 1., out );
957 helperP.swap(tempP);
959 dg::blas2::symv( m_plus, helperM, tempM);
960 dg::blas1::axpby( (double)(cphi-j)/(double)cphi, tempM, 1., out );
961 helperM.swap(tempM);
962 }
963 dg::blas1::scal( out, 1./(double)cphi);
964}
965
966template<class G, class I, class container>
967void Fieldaligned<G, I, container >::operator()(enum whichMatrix which, const container& f, container& fe)
968{
969 if( which == einsPlus || which == einsMinusT ) ePlus( which, f, fe);
970 else if(which == einsMinus || which == einsPlusT ) eMinus( which, f, fe);
971 else if(which == zeroMinus || which == zeroPlus ||
972 which == zeroMinusT|| which == zeroPlusT ||
973 which == zeroForw ) zero( which, f, fe);
974 if( m_interpolation_method == "dg-const")
975 {
976 //dg::split( m_tmp, m_temp2, *m_g);
977 //dg::split( fe, m_temp, *m_g);
978 //for( unsigned i0 = 0; i0 <m_Nz; i0++)
979 // dg::blas2::symv( m_forw, m_temp[i0], m_temp2[i0]);
980 //m_tmp.swap(fe);
981 dg::blas2::stencil( dg::CSRSlopeLimiter<double>(10), m_stencil, fe, m_tmp);
982 dg::blas2::stencil( dg::CSRSlopeLimiter<double>(10), m_stencilY, m_tmp, fe);
983 }
984}
985
986template< class G, class I, class container>
987void Fieldaligned<G, I, container>::zero( enum whichMatrix which,
988 const container& f, container& f0)
989{
990 dg::split( f, m_f, *m_g);
991 dg::split( f0, m_temp, *m_g);
992 //1. compute 2d interpolation in every plane and store in m_temp
993 for( unsigned i0=0; i0<m_Nz; i0++)
994 {
995 if(which == zeroPlus)
996 dg::blas2::symv( m_plus, m_f[i0], m_temp[i0]);
997 else if(which == zeroMinus)
998 dg::blas2::symv( m_minus, m_f[i0], m_temp[i0]);
999 else if(which == zeroPlusT)
1000 {
1001 if( ! m_have_adjoint) updateAdjoint( );
1002 dg::blas2::symv( m_plusT, m_f[i0], m_temp[i0]);
1003 }
1004 else if(which == zeroMinusT)
1005 {
1006 if( ! m_have_adjoint) updateAdjoint( );
1007 dg::blas2::symv( m_minusT, m_f[i0], m_temp[i0]);
1008 }
1009 }
1010 if( which == zeroForw)
1011 {
1012 if ( m_interpolation_method == "linear-const" )
1013 {
1014 dg::blas2::symv( m_forw, f, f0);
1015 dg::blas2::symv( m_inv_linear.tri(), f0, m_tmp);
1016 dg::blas2::symv( m_back, m_tmp, f0);
1017 }
1018 else if ( m_interpolation_method == "linear" )
1019 {
1020 dg::blas2::symv( m_forw, f, f0);
1021 dg::blas2::symv( m_inv_linear.tri(), f0, m_tmp);
1022 dg::blas2::symv( m_back, m_tmp, f0);
1023 }
1024 else
1025 dg::blas1::copy( f, f0);
1026 }
1027}
1028template< class G, class I, class container>
1029void Fieldaligned<G, I, container>::ePlus( enum whichMatrix which,
1030 const container& f, container& fpe)
1031{
1032 dg::split( f, m_f, *m_g);
1033 dg::split( fpe, m_temp, *m_g);
1034 //1. compute 2d interpolation in every plane and store in m_temp
1035 for( unsigned i0=0; i0<m_Nz; i0++)
1036 {
1037 unsigned ip = (i0==m_Nz-1) ? 0:i0+1;
1038 if(which == einsPlus)
1039 dg::blas2::symv( m_plus, m_f[ip], m_temp[i0]);
1040 else if(which == einsMinusT)
1041 {
1042 if( ! m_have_adjoint) updateAdjoint( );
1043 dg::blas2::symv( m_minusT, m_f[ip], m_temp[i0]);
1044 }
1045 }
1046 //2. apply right boundary conditions in last plane
1047 unsigned i0=m_Nz-1;
1048 if( m_bcz != dg::PER)
1049 {
1050 if( m_bcz == dg::DIR || m_bcz == dg::NEU_DIR)
1051 dg::blas1::axpby( 2, m_right, -1., m_f[i0], m_ghostP);
1052 if( m_bcz == dg::NEU || m_bcz == dg::DIR_NEU)
1053 dg::blas1::axpby( m_deltaPhi, m_right, 1., m_f[i0], m_ghostP);
1054 //interlay ghostcells with periodic cells: L*g + (1-L)*fpe
1055 dg::blas1::axpby( 1., m_ghostP, -1., m_temp[i0], m_ghostP);
1056 dg::blas1::pointwiseDot( 1., m_limiter, m_ghostP, 1., m_temp[i0]);
1057 }
1058}
1059
1060template< class G, class I, class container>
1061void Fieldaligned<G, I, container>::eMinus( enum whichMatrix which,
1062 const container& f, container& fme)
1063{
1064 dg::split( f, m_f, *m_g);
1065 dg::split( fme, m_temp, *m_g);
1066 //1. compute 2d interpolation in every plane and store in m_temp
1067 for( unsigned i0=0; i0<m_Nz; i0++)
1068 {
1069 unsigned im = (i0==0) ? m_Nz-1:i0-1;
1070 if(which == einsPlusT)
1071 {
1072 if( ! m_have_adjoint) updateAdjoint( );
1073 dg::blas2::symv( m_plusT, m_f[im], m_temp[i0]);
1074 }
1075 else if (which == einsMinus)
1076 dg::blas2::symv( m_minus, m_f[im], m_temp[i0]);
1077 }
1078 //2. apply left boundary conditions in first plane
1079 unsigned i0=0;
1080 if( m_bcz != dg::PER)
1081 {
1082 if( m_bcz == dg::DIR || m_bcz == dg::DIR_NEU)
1083 dg::blas1::axpby( 2., m_left, -1., m_f[i0], m_ghostM);
1084 if( m_bcz == dg::NEU || m_bcz == dg::NEU_DIR)
1085 dg::blas1::axpby( -m_deltaPhi, m_left, 1., m_f[i0], m_ghostM);
1086 //interlay ghostcells with periodic cells: L*g + (1-L)*fme
1087 dg::blas1::axpby( 1., m_ghostM, -1., m_temp[i0], m_ghostM);
1088 dg::blas1::pointwiseDot( 1., m_limiter, m_ghostM, 1., m_temp[i0]);
1089 }
1090}
1091
1092template<class G, class I, class container>
1093template< class BinaryOp, class UnaryOp>
1094container Fieldaligned<G, I,container>::evaluate( BinaryOp binary,
1095 UnaryOp unary, unsigned p0, unsigned rounds) const
1096{
1097 //idea: simply apply I+/I- enough times on the init2d vector to get the result in each plane
1098 //unary function is always such that the p0 plane is at x=0
1099 assert( p0 < m_g->Nz());
1100 const dg::ClonePtr<aGeometry2d> g2d = m_g->perp_grid();
1101 container init2d = dg::pullback( binary, *g2d);
1102 container zero2d = dg::evaluate( dg::zero, *g2d);
1103
1104 container temp(init2d), tempP(init2d), tempM(init2d);
1105 container vec3d = dg::evaluate( dg::zero, *m_g);
1106 std::vector<container> plus2d(m_Nz, zero2d), minus2d(plus2d), result(plus2d);
1107 unsigned turns = rounds;
1108 if( turns ==0) turns++;
1109 //first apply Interpolation many times, scale and store results
1110 for( unsigned r=0; r<turns; r++)
1111 for( unsigned i0=0; i0<m_Nz; i0++)
1112 {
1113 dg::blas1::copy( init2d, tempP);
1114 dg::blas1::copy( init2d, tempM);
1115 unsigned rep = r*m_Nz + i0;
1116 for(unsigned k=0; k<rep; k++)
1117 {
1119 dg::blas2::symv( m_minus, tempP, temp);
1120 temp.swap( tempP);
1122 dg::blas2::symv( m_plus, tempM, temp);
1123 temp.swap( tempM);
1124 }
1125 dg::blas1::scal( tempP, unary( (double)rep*m_deltaPhi ) );
1126 dg::blas1::scal( tempM, unary( -(double)rep*m_deltaPhi ) );
1127 dg::blas1::axpby( 1., tempP, 1., plus2d[i0]);
1128 dg::blas1::axpby( 1., tempM, 1., minus2d[i0]);
1129 }
1130 //now we have the plus and the minus filaments
1131 if( rounds == 0) //there is a limiter
1132 {
1133 for( unsigned i0=0; i0<m_Nz; i0++)
1134 {
1135 int idx = (int)i0 - (int)p0;
1136 if(idx>=0)
1137 result[i0] = plus2d[idx];
1138 else
1139 result[i0] = minus2d[abs(idx)];
1140 thrust::copy( result[i0].begin(), result[i0].end(), vec3d.begin() + i0*m_perp_size);
1141 }
1142 }
1143 else //sum up plus2d and minus2d
1144 {
1145 for( unsigned i0=0; i0<m_Nz; i0++)
1146 {
1147 unsigned revi0 = (m_Nz - i0)%m_Nz; //reverted index
1148 dg::blas1::axpby( 1., plus2d[i0], 0., result[i0]);
1149 dg::blas1::axpby( 1., minus2d[revi0], 1., result[i0]);
1150 }
1151 dg::blas1::axpby( -1., init2d, 1., result[0]);
1152 for(unsigned i0=0; i0<m_Nz; i0++)
1153 {
1154 int idx = ((int)i0 -(int)p0 + m_Nz)%m_Nz; //shift index
1155 thrust::copy( result[idx].begin(), result[idx].end(), vec3d.begin() + i0*m_perp_size);
1156 }
1157 }
1158 return vec3d;
1159}
1160
1161
1163
1195template<class BinaryOp, class UnaryOp>
1196thrust::host_vector<double> fieldaligned_evaluate(
1197 const aProductGeometry3d& grid,
1198 const CylindricalVectorLvl0& vec,
1199 const BinaryOp& binary,
1200 const UnaryOp& unary,
1201 unsigned p0,
1202 unsigned rounds,
1203 double eps = 1e-5)
1204{
1205 unsigned Nz = grid.Nz();
1206 const dg::ClonePtr<aGeometry2d> g2d = grid.perp_grid();
1207 // Construct for field-aligned output
1208 dg::HVec tempP = dg::evaluate( dg::zero, *g2d), tempM( tempP);
1209 std::vector<dg::HVec> plus2d(Nz, tempP), minus2d(plus2d), result(plus2d);
1210 dg::HVec vec3d = dg::evaluate( dg::zero, grid);
1211 dg::HVec init2d = dg::pullback( binary, *g2d);
1212 std::array<dg::HVec,3> yy0{
1213 dg::pullback( dg::cooX2d, *g2d),
1214 dg::pullback( dg::cooY2d, *g2d),
1215 dg::evaluate( dg::zero, *g2d)}, yy1(yy0), xx0( yy0), xx1(yy0); //s
1216 dg::geo::detail::DSFieldCylindrical3 cyl_field(vec);
1217 double deltaPhi = grid.hz();
1218 double phiM0 = 0., phiP0 = 0.;
1219 unsigned turns = rounds;
1220 if( turns == 0) turns++;
1221 for( unsigned r=0; r<turns; r++)
1222 for( unsigned i0=0; i0<Nz; i0++)
1223 {
1224 unsigned rep = r*Nz + i0;
1225 if( rep == 0)
1226 tempM = tempP = init2d;
1227 else
1228 {
1230 "Dormand-Prince-7-4-5", std::array<double,3>{0,0,0});
1232 cyl_field, dg::pid_control, dg::fast_l2norm, eps,
1233 1e-10);
1234 for( unsigned i=0; i<g2d->size(); i++)
1235 {
1236 // minus direction needs positive integration!
1237 double phiM1 = phiM0 + deltaPhi;
1238 std::array<double,3>
1239 coords0{yy0[0][i],yy0[1][i],yy0[2][i]}, coords1;
1240 odeint.integrate_in_domain( phiM0, coords0, phiM1, coords1,
1241 deltaPhi, *g2d, eps);
1242 yy1[0][i] = coords1[0], yy1[1][i] = coords1[1], yy1[2][i] =
1243 coords1[2];
1244 tempM[i] = binary( yy1[0][i], yy1[1][i]);
1245
1246 // plus direction needs negative integration!
1247 double phiP1 = phiP0 - deltaPhi;
1248 coords0 = std::array<double,3>{xx0[0][i],xx0[1][i],xx0[2][i]};
1249 odeint.integrate_in_domain( phiP0, coords0, phiP1, coords1,
1250 -deltaPhi, *g2d, eps);
1251 xx1[0][i] = coords1[0], xx1[1][i] = coords1[1], xx1[2][i] =
1252 coords1[2];
1253 tempP[i] = binary( xx1[0][i], xx1[1][i]);
1254 }
1255 std::swap( yy0, yy1);
1256 std::swap( xx0, xx1);
1257 phiM0 += deltaPhi;
1258 phiP0 -= deltaPhi;
1259 }
1260 dg::blas1::scal( tempM, unary( -(double)rep*deltaPhi ) );
1261 dg::blas1::scal( tempP, unary( (double)rep*deltaPhi ) );
1262 dg::blas1::axpby( 1., tempM, 1., minus2d[i0]);
1263 dg::blas1::axpby( 1., tempP, 1., plus2d[i0]);
1264 }
1265 //now we have the plus and the minus filaments
1266 if( rounds == 0) //there is a limiter
1267 {
1268 for( unsigned i0=0; i0<Nz; i0++)
1269 {
1270 int idx = (int)i0 - (int)p0;
1271 if(idx>=0)
1272 result[i0] = plus2d[idx];
1273 else
1274 result[i0] = minus2d[abs(idx)];
1275 thrust::copy( result[i0].begin(), result[i0].end(), vec3d.begin() +
1276 i0*g2d->size());
1277 }
1278 }
1279 else //sum up plus2d and minus2d
1280 {
1281 for( unsigned i0=0; i0<Nz; i0++)
1282 {
1283 unsigned revi0 = (Nz - i0)%Nz; //reverted index
1284 dg::blas1::axpby( 1., plus2d[i0], 0., result[i0]);
1285 dg::blas1::axpby( 1., minus2d[revi0], 1., result[i0]);
1286 }
1287 dg::blas1::axpby( -1., init2d, 1., result[0]);
1288 for(unsigned i0=0; i0<Nz; i0++)
1289 {
1290 int idx = ((int)i0 -(int)p0 + Nz)%Nz; //shift index
1291 thrust::copy( result[idx].begin(), result[idx].end(), vec3d.begin()
1292 + i0*g2d->size());
1293 }
1294 }
1295 return vec3d;
1296}
1297
1298}//namespace geo
1299}//namespace dg
void assign(const from_ContainerType &from, ContainerType &to, Params &&... ps)
static DG_DEVICE double cooY2d(double x, double y)
static DG_DEVICE double zero(double x)
static DG_DEVICE double cooX2d(double x, double y)
void copy(const ContainerTypeIn &source, ContainerTypeOut &target)
void plus(ContainerType &x, get_value_type< ContainerType > alpha)
void axpby(get_value_type< ContainerType > alpha, const ContainerType1 &x, get_value_type< ContainerType > beta, ContainerType &y)
void scal(ContainerType &x, get_value_type< ContainerType > alpha)
void pointwiseDivide(get_value_type< ContainerType > alpha, const ContainerType1 &x1, const ContainerType2 &x2, get_value_type< ContainerType > beta, ContainerType &y)
void pointwiseDot(get_value_type< ContainerType > alpha, const ContainerType1 &x1, const ContainerType2 &x2, get_value_type< ContainerType > beta, ContainerType &y)
void transfer(const MatrixType &x, AnotherMatrixType &y)
void symv(MatrixType &&M, const ContainerType1 &x, ContainerType2 &y)
void stencil(FunctorType f, MatrixType &&M, const ContainerType1 &x, ContainerType2 &y)
static std::string bc2str(bc bcx)
thrust::host_vector< real_type > evaluate(UnaryOp f, const RealGrid1d< real_type > &g)
thrust::host_vector< real_type > fem_inv_weights(const RealGrid1d< real_type > &g)
dg::InverseKroneckerTriDiagonal2d< dg::HVec_t< real_type > > inv_fem_mass2d(const aRealTopology3d< real_type > &g)
dg::InverseKroneckerTriDiagonal2d< dg::HVec_t< real_type > > inv_fem_linear2const2d(const aRealTopology3d< real_type > &g)
whichMatrix
Enum for the use in Fieldaligned.
Definition: fieldaligned.h:17
ZERO NoLimiter
No Limiter.
Definition: fieldaligned.h:35
thrust::host_vector< double > fieldaligned_evaluate(const aProductGeometry3d &grid, const CylindricalVectorLvl0 &vec, const BinaryOp &binary, const UnaryOp &unary, unsigned p0, unsigned rounds, double eps=1e-5)
Evaluate a 2d functor and transform to all planes along the fieldlines
Definition: fieldaligned.h:1204
ONE FullLimiter
Full Limiter means there is a limiter everywhere.
Definition: fieldaligned.h:31
@ zeroPlusT
transposed plus interpolation in the current plane
Definition: fieldaligned.h:24
@ einsMinus
minus interpolation in previous plane
Definition: fieldaligned.h:20
@ einsMinusT
transposed minus interpolation in next plane
Definition: fieldaligned.h:21
@ zeroMinusT
transposed minus interpolation in the current plane
Definition: fieldaligned.h:25
@ einsPlusT
transposed plus interpolation in previous plane
Definition: fieldaligned.h:19
@ zeroForw
from dg to transformed coordinates
Definition: fieldaligned.h:26
@ zeroMinus
minus interpolation in the current plane
Definition: fieldaligned.h:23
@ zeroPlus
plus interpolation in the current plane
Definition: fieldaligned.h:22
@ einsPlus
plus interpolation in next plane
Definition: fieldaligned.h:18
MPI_Vector< thrust::host_vector< real_type > > weights(const aRealMPITopology2d< real_type > &g)
cusp::coo_matrix< int, real_type, cusp::host_memory > diagonal(const thrust::host_vector< real_type > &diagonal)
cusp::coo_matrix< int, real_type, cusp::host_memory > interpolation(const thrust::host_vector< real_type > &x, const RealGrid1d< real_type > &g, dg::bc bcx=dg::NEU, std::string method="dg")
real_type interpolate(dg::space sp, const thrust::host_vector< real_type > &v, real_type x, const RealGrid1d< real_type > &g, dg::bc bcx=dg::NEU)
dg::MIHMatrix_t< real_type > projection(const aRealMPITopology2d< real_type > &g_new, const aRealMPITopology2d< real_type > &g_old, std::string method="dg")
CylindricalVectorLvl1 createBHat(const TokamakMagneticField &mag)
Contravariant components of the magnetic unit vector field and its Divergence and derivative in cylin...
Definition: magnetic_field.h:931
get_host_vector< Geometry > volume(const Geometry &g)
thrust::host_vector< real_type > forward_transform(const thrust::host_vector< real_type > &in, const aRealTopology2d< real_type > &g)
void transpose(unsigned nx, unsigned ny, const ContainerType &in, ContainerType &out)
thrust::host_vector< real_type > pullback(const Functor &f, const aRealGeometryX2d< real_type > &g)
void pushForwardPerp(const Functor1 &vR, const Functor2 &vZ, container &vx, container &vy, const Geometry &g)
void assign3dfrom2d(const thrust::host_vector< real_type > &in2d, Container &out, const aRealTopology3d< real_type > &grid)
dg::IHMatrix_t< real_type > backproject(const RealGrid1d< real_type > &g)
dg::IHMatrix_t< real_type > inv_backproject(const RealGrid1d< real_type > &g)
void split(SharedContainer &in, std::vector< View< SharedContainer > > &out, const aRealTopology3d< real_type > &grid)
dg::IHMatrix_t< real_type > limiter_stencil(const RealGrid1d< real_type > &g, dg::bc bound)
ContainerType volume(const SparseTensor< ContainerType > &t)
static auto pid_control
static auto fast_l2norm
IHMatrix_t< double > IHMatrix
thrust::host_vector< double > HVec
void set_dt(value_type dt)
void integrate_in_domain(value_type t0, const ContainerType &u0, value_type &t1, ContainerType &u1, value_type dt, Domain &&domain, value_type eps_root)
double diff() const
void toc()
void tic()
aRealGeometry2d< real_type > * perp_grid() const
bool contains(real_type x, real_type y) const
unsigned size() const
real_type hz() const
unsigned Nz() const
void integrate(value_type t0, const ContainerType &u0, value_type t1, ContainerType &u1)
Definition: fluxfunctions.h:412
This struct bundles a vector field and its divergence.
Definition: fluxfunctions.h:440
const CylindricalFunctor & y() const
y-component of the vector
Definition: fluxfunctions.h:468
const CylindricalFunctor & x() const
x-component of the vector
Definition: fluxfunctions.h:466
const CylindricalFunctor & divvvz() const
Definition: fluxfunctions.h:474
const CylindricalFunctor & z() const
z-component of the vector
Definition: fluxfunctions.h:470
Create and manage interpolation matrices from fieldline integration.
Definition: fieldaligned.h:433
dg::bc bcx() const
Definition: fieldaligned.h:484
const container & hbp() const
Distance between the planes .
Definition: fieldaligned_tmp.h:502
const container & bphi() const
bphi
Definition: fieldaligned_tmp.h:522
container interpolate_from_coarse_grid(const ProductGeometry &grid_coarse, const container &coarse)
Interpolate along fieldlines from a coarse to a fine grid in phi.
void integrate_between_coarse_grid(const ProductGeometry &grid_coarse, const container &coarse, container &out)
Integrate a 2d function on the fine grid.
const container & bbp() const
Mask plus, 1 if fieldline intersects wall in plus direction but not in minus direction,...
Definition: fieldaligned_tmp.h:547
container evaluate(BinaryOp binary, UnaryOp unary, unsigned p0, unsigned rounds) const
Evaluate a 2d functor and transform to all planes along the fieldline
const container & sqrtG() const
Volume form (including weights) .
Definition: fieldaligned_tmp.h:507
void operator()(enum whichMatrix which, const container &in, container &out)
Apply the interpolation to three-dimensional vectors.
const container & sqrtGm() const
Volume form on minus plane (including weights) .
Definition: fieldaligned_tmp.h:512
void set_boundaries(dg::bc bcz, double left, double right)
Set boundary conditions in the limiter region.
Definition: fieldaligned_tmp.h:442
std::string method() const
Definition: fieldaligned_tmp.h:611
const container & hbm() const
Distance between the planes and the boundary .
Definition: fieldaligned_tmp.h:497
void set_boundaries(dg::bc bcz, const container &global, double scal_left, double scal_right)
Set boundary conditions in the limiter region.
Definition: fieldaligned_tmp.h:476
const container & bphiM() const
bphi on minus plane
Definition: fieldaligned_tmp.h:527
dg::bc bcy() const
Definition: fieldaligned.h:487
double deltaPhi() const
Definition: fieldaligned.h:553
Fieldaligned(const dg::geo::TokamakMagneticField &vec, const ProductGeometry &grid, dg::bc bcx=dg::NEU, dg::bc bcy=dg::NEU, Limiter limit=FullLimiter(), double eps=1e-5, unsigned mx=10, unsigned my=10, double deltaPhi=-1, std::string interpolation_method="dg", bool benchmark=true)
Construct from a magnetic field and a grid.
Definition: fieldaligned_tmp.h:382
Fieldaligned()
do not allocate memory; no member call except construct is valid
Definition: fieldaligned_tmp.h:377
void construct(Params &&...ps)
Perfect forward parameters to one of the constructors.
Definition: fieldaligned_tmp.h:419
const container & bbo() const
Mask both, 1 if fieldline intersects wall in plus direction and in minus direction,...
Definition: fieldaligned_tmp.h:542
const container & bbm() const
Mask minus, 1 if fieldline intersects wall in minus direction but not in plus direction,...
Definition: fieldaligned_tmp.h:537
const container & sqrtGp() const
Volume form on plus plane (including weights) .
Definition: fieldaligned_tmp.h:517
const container & bphiP() const
bphi on plus plane
Definition: fieldaligned_tmp.h:532
const ProductGeometry & grid() const
Grid used for construction.
Definition: fieldaligned.h:610
Fieldaligned(const dg::geo::CylindricalVectorLvl1 &vec, const ProductGeometry &grid, dg::bc bcx=dg::NEU, dg::bc bcy=dg::NEU, Limiter limit=FullLimiter(), double eps=1e-5, unsigned mx=10, unsigned my=10, double deltaPhi=-1, std::string interpolation_method="dg", bool benchmark=true)
Construct from a vector field and a grid.
void set_boundaries(dg::bc bcz, const container &left, const container &right)
Set boundary conditions in the limiter region.
Definition: fieldaligned_tmp.h:459
A tokamak field as given by R0, Psi and Ipol plus Meta-data like shape and equilibrium.
Definition: magnetic_field.h:162
WallFieldlineCoordinate(const dg::geo::CylindricalVectorLvl0 &vec, const dg::aRealTopology2d< double > &domain, double maxPhi, double eps, std::string type)
Construct with vector field, domain.
Definition: fieldaligned_tmp.h:267
double do_compute(double R, double Z) const
Definition: fieldaligned_tmp.h:277
WallFieldlineDistance(const dg::geo::CylindricalVectorLvl0 &vec, const dg::aRealTopology2d< double > &domain, double maxPhi, double eps, std::string type)
Construct with vector field, domain.
Definition: fieldaligned_tmp.h:197
double do_compute(double R, double Z) const
Integrate fieldline until wall is reached.
Definition: fieldaligned_tmp.h:217