p3m3dTwoStreamParallel.cpp 39.4 KB
Newer Older
frey_m's avatar
frey_m committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23
//
// Application p3m3dTwoStreamParallel
//   ./p3m3dTwoStreamParallel 16 16 16 rc alpha dt alpha_amp eps Nsteps
//
//   using the "point" distribution will only place one particle
//
// Copyright (c) 2016, Benjamin Ulmer, ETH Zürich
// All rights reserved
//
// Implemented as part of the Master thesis
// "The P3M Model on Emerging Computer Architectures With Application to Microbunching"
// (http://amas.web.psi.ch/people/aadelmann/ETH-Accel-Lecture-1/projectscompleted/cse/thesisBUlmer.pdf)
//
// This file is part of OPAL.
//
// OPAL is free software: you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// You should have received a copy of the GNU General Public License
// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
//
frey_m's avatar
frey_m committed
24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68
#include "Ippl.h"
#include <string>
#include <vector>

#include <iostream>
#include <cfloat>
#include <fstream>
#include <iomanip>
#include "Particle/BoxParticleCachingPolicy.h"
#include "Particle/PairBuilder/HashPairBuilderPeriodicParallel.h"
#include "Particle/PairBuilder/PairConditions.h"
#include "math.h"

#include <random>

//#include "VTKFieldWriter.hpp"
#include "VTKFieldWriterParallel.hpp"
#include "ChargedParticleFactory.hpp"

#define GS

//#define LANDAU

#ifndef LANDAU
#define SPHERE
//#define TWOSTREAM
#endif

#define CALC_ERRORS
////TEMP debug Variables
double RhoSum=0;
//////////////////////////
// dimension of our positions
const unsigned Dim = 3;

// some typedefs
typedef UniformCartesian<Dim, double>                                 Mesh_t;
typedef BoxParticleCachingPolicy<double, Dim, Mesh_t>                 CachingPolicy_t;
typedef ParticleSpatialLayout<double, Dim, Mesh_t, CachingPolicy_t>   playout_t;
typedef playout_t::SingleParticlePos_t                                Vector_t;
typedef Cell                                                          Center_t;
typedef CenteredFieldLayout<Dim, Mesh_t, Center_t>                    FieldLayout_t;
typedef Field<double, Dim, Mesh_t, Center_t>                          Field_t;
typedef Field<int, Dim, Mesh_t, Center_t>                             IField_t;
typedef Field<Vector_t, Dim, Mesh_t, Center_t>                        VField_t;
gsell's avatar
gsell committed
69
typedef Field<std::complex<double>, Dim, Mesh_t, Center_t>            CxField_t;
frey_m's avatar
frey_m committed
70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89
//typedef FFT<RCTransform, Dim, double>                               FFT_t;
typedef FFT<CCTransform, Dim, double>                                 FFT_t;

typedef IntCIC                                                        IntrplCIC_t;
typedef IntNGP                                                        IntrplNGP_t;
typedef IntTSC                                                        IntrplTSC_t;

typedef UniformCartesian<2, double>                                   Mesh2d_t;
typedef CenteredFieldLayout<2, Mesh2d_t, Center_t>                    FieldLayout2d_t;
typedef Field<double, 2, Mesh2d_t, Center_t>                          Field2d_t;

template<class T>
struct ApplyField;

//This is the periodic Greens function with regularization parameter epsilon.
template<unsigned int Dim>
struct SpecializedGreensFunction { };

template<>
struct SpecializedGreensFunction<3> {
90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109
        template<class T, class FT, class FT2>
                static void calculate(Vektor<T, 3> &hrsq, FT &grn, FT2 *grnI, double alpha,double eps) {
                        double r;
                        grn = grnI[0] * hrsq[0] + grnI[1] * hrsq[1] + grnI[2] * hrsq[2];
                        NDIndex<3> lDomain_m = grn.getLayout().getLocalNDIndex();
                        NDIndex<3> elem;
                        for (int i=lDomain_m[0].min(); i<=lDomain_m[0].max(); ++i) {
                                elem[0]=Index(i,i);
                                for (int j=lDomain_m[1].min(); j<=lDomain_m[1].max(); ++j) {
                                        elem[1]=Index(j,j);
                                        for (int k=lDomain_m[2].min(); k<=lDomain_m[2].max(); ++k) {
                                                elem[2]=Index(k,k);

                                                r = real(sqrt(grn.localElement(elem)));
                                                grn.localElement(elem) = 1./(4.*M_PI)*std::complex<double>(erf(alpha*r)/(r+eps));
                                        }
                                }
                        }
                        //grn[0][0][0] = grn[0][0][1];
                }
frey_m's avatar
frey_m committed
110 111 112 113 114

};


template<class PL>
115
class ChargedParticles : public IpplParticleBase<PL> {
116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380
        public:
                ParticleAttrib<double>          Q;
                ParticleAttrib<double>          m;
                ParticleAttrib<double>          SpaceQ;
                ParticleAttrib<double>          Phi; //electrostatic potential
                ParticleAttrib<Vector_t>        EF;
                ParticleAttrib<Vector_t>        v; //velocity of the particles
                ParticleAttrib<int>     ID; //velocity of the particles
                ParticleAttrib<Vektor<double,2> > Rphase; //velocity of the particles

                ChargedParticles(PL* pl, Vektor<double,3> nr, e_dim_tag /*decomp*/[Dim],Vektor<double,3> extend_l_, Vektor<double,3> extend_r_, Vektor<int,3> Nx_, Vektor<int,3> Nv_, Vektor<double,3> Vmax_) :
                        IpplParticleBase<PL>(pl),
                        nr_m(nr),
                        extend_l(extend_l_),
                        extend_r(extend_r_),
                        Nx(Nx_), Nv(Nv_), Vmax(Vmax_)
        {
                this->addAttribute(Q);
                this->addAttribute(m);
                this->addAttribute(SpaceQ);
                this->addAttribute(Phi);
                this->addAttribute(EF);
                this->addAttribute(v);
                this->addAttribute(ID);
                this->addAttribute(Rphase);

                for (unsigned int i = 0; i < 2 * Dim; ++i) {
                        //use periodic boundary conditions for the particles
                        this->getBConds()[i] = ParticlePeriodicBCond;
                        //boundary conditions used for interpolation kernels allow writes to ghost cells
                        bc_m[i] = new ParallelInterpolationFace<double, Dim, Mesh_t, Center_t>(i);
                        //std periodic boundary conditions for gradient computations etc.
                        vbc_m[i] = new ParallelPeriodicFace<Vector_t, Dim, Mesh_t, Center_t>(i);
                        bcp_m[i] = new ParallelPeriodicFace<double, Dim, Mesh_t, Center_t>(i);
                }

                for (unsigned int d = 0;d<Dim;++d) {
                        rmax_m[d] = extend_r[d];
                        rmin_m[d] = extend_l[d];
                }
                //Initialize the meshes and layouts for 2D phase space interpolation:
                double spacings[2] = {(extend_r[2]-extend_l[2])/(Nx[2]),2.*Vmax[2]/(Nv[2])};
                Vektor<double,2> origin;

                origin(0) = extend_l[2]; origin(1) = -Vmax[2];
                Index I(Nx[2]+1); Index J(Nv[2]+1);
                domain2d_m[0]=I; domain2d_m[1]=J;

                mesh2d_m=Mesh2d_t(domain2d_m, spacings, origin);
                layout2d_m = new FieldLayout2d_t(mesh2d_m);

                BConds<double,2,UniformCartesian<2,double>,Cell> BC;
                BC[0] = new ParallelInterpolationFace<double,2,Mesh2d_t,Cell>(0);
                BC[1] = new ParallelInterpolationFace<double,2,Mesh2d_t,Cell>(1);
                BC[2] = new ParallelInterpolationFace<double,2,Mesh2d_t,Cell>(2);
                BC[3] = new ParallelInterpolationFace<double,2,Mesh2d_t,Cell>(3);

                //set origin and spacing is needed for correct results, even if mesh was created with these paremeters ?!
                mesh2d_m.set_meshSpacing(&(spacings[0]));
                mesh2d_m.set_origin(origin);

                domain2d_m = layout2d_m->getDomain();
                //f_m is used for twostream instability as 2D phase space mesh
                f_m.initialize(mesh2d_m, *layout2d_m, GuardCellSizes<2>(1),BC);

                domain_m = this->getFieldLayout().getDomain();
                lDomain_m = this->getFieldLayout().getLocalNDIndex(); // local domain

                //initialize the FFT
                bool compressTemps = true;
                fft_m = new FFT_t(domain_m,compressTemps);

                fft_m->setDirectionName(+1, "forward");
                fft_m->setDirectionName(-1, "inverse");
                INFOMSG("INIT FFT DONE"<<endl);
        }

                inline const Mesh_t& getMesh() const { return this->getLayout().getLayout().getMesh(); }

                inline Mesh_t& getMesh() { return this->getLayout().getLayout().getMesh(); }

                inline const FieldLayout_t& getFieldLayout() const {
                        return dynamic_cast<FieldLayout_t&>( this->getLayout().getLayout().getFieldLayout());
                }

                inline FieldLayout_t& getFieldLayout() {
                        return dynamic_cast<FieldLayout_t&>(this->getLayout().getLayout().getFieldLayout());
                }

                void update()
                {
                        //should only be needed if meshspacing changes -----------
                        for (unsigned int d = 0;d<Dim;++d) {
                                hr_m[d] = (extend_r[d] - extend_l[d]) / (nr_m[d]);
                        }
                        this->getMesh().set_meshSpacing(&(hr_m[0]));
                        this->getMesh().set_origin(extend_l);
                        //--------------------------------------------------------

                        //init resets the meshes to 0 ?!
                        rhocmpl_m.initialize(getMesh(), getFieldLayout(), GuardCellSizes<Dim>(1));
                        grncmpl_m.initialize(getMesh(), getFieldLayout(), GuardCellSizes<Dim>(1));
                        rho_m.initialize(getMesh(), getFieldLayout(), GuardCellSizes<Dim>(1),bc_m);
                        phi_m.initialize(getMesh(), getFieldLayout(), GuardCellSizes<Dim>(1),bcp_m);
                        eg_m.initialize(getMesh(), getFieldLayout(), GuardCellSizes<Dim>(1), vbc_m);

                        domain_m = this->getFieldLayout().getDomain();
                        lDomain_m = this->getFieldLayout().getLocalNDIndex();
                        /*
                        std::cout << "******************* local domain is " << lDomain_m << std::endl;
                        std::cout << "******************* global domain is " << domain_m << std::endl;
                        */

                        IpplParticleBase<PL>::update();
                }

                void interpolate_distribution(Vektor<double,3> dx, Vektor<double,3> dv){
                        f_m=0;
                        for (unsigned i=0; i<this->getLocalNum(); ++i) {
                                SpaceQ[i]=Q[i]/(dx[0]*dx[1]*dv[0]*dv[1]);
                                Rphase[i]=Vektor<double,2>(this->R[i][2],v[i][2]);
                        }
                        //this->SpaceQ.scatter(this->f_m, this->Rphase, IntrplCIC_t());
                        this->Q.scatter(this->f_m, this->Rphase, IntrplCIC_t());
                        std::cout << "scatter done" << std::endl;
                }

                void calc_kinetic_energy() {
                        double loc_kinetic_energy=0;
                        double v2;
                        for (unsigned i=0; i<this->getLocalNum(); ++i) {
                                v2 = (v[i][0]*v[i][0]+v[i][1]*v[i][1]+v[i][2]*v[i][2]);
                                loc_kinetic_energy+=0.5*v2;
                        }
                        kinetic_energy=0;
                        Ippl::Comm->barrier();
                        reduce(loc_kinetic_energy, kinetic_energy, OpAddAssign());
                }

                void  calc_field_energy() {
                        NDIndex<3> elem;
                        double cell_volume = hr_m[0]*hr_m[1]*hr_m[2];
                        field_energy=0;
                        field_energy=0.5*cell_volume*sum(dot(eg_m,eg_m));

                        rhomax=max(rho_m)/(hr_m[0]*hr_m[1]*hr_m[2]);
                        //rhomax=max(rho_m);
                        integral_phi_m=0.5*sum(rho_m*phi_m);
                }

                void  calc_potential_energy() {
                        double loc_potential_energy=0;
                        for (unsigned i=0; i<this->getLocalNum(); ++i) {
                                loc_potential_energy+=0.5*(Q[i])*Phi[i];
                        }
                        potential_energy=0;
                        Ippl::Comm->barrier();
                        reduce(loc_potential_energy, potential_energy, OpAddAssign());

                }

                void calc_Amplitude_E(){
                        //computes the maximum amplitude in the electric field
                        AmplitudeEfield=max(sqrt(dot(eg_m,eg_m)));
                        eg_m=eg_m*Vektor<double,3>(0,0,1);
                        AmplitudeEFz=max(sqrt(dot(eg_m,eg_m)));
                }


                void calculatePairForces(double interaction_radius, double eps, double alpha)
                {
                        if (interaction_radius>0){
                                HashPairBuilderPeriodicParallel< ChargedParticles<playout_t> > HPB(*this);
                                HPB.for_each(RadiusCondition<double, Dim>(interaction_radius), ApplyField<double>(-1,interaction_radius,eps,alpha),extend_l, extend_r);
                        }
                }

                //compute Greens function without making use of index computations
                void calcGrealSpace(double alpha, double eps) {
                        NDIndex<3> elem, elem_mirr,elem_all;
                        Vector_t hrsq(hr_m * hr_m);
                        int N =  domain_m[0].max();
                        INFOMSG("domainlength = " << domain_m[0].max() << " "<<domain_m[1].max()<< " " << domain_m[2].max() << endl);
                        INFOMSG("hr_m = " << hr_m << endl);
                        //Loop over first octant of domain
                        for (int i=lDomain_m[0].min(); i<=(lDomain_m[0].max())/2.+1; ++i) {
                                elem[0]=Index(i,i);
                                elem_mirr[0]=Index(N-i+1,N-i+1);
                                for (int j=lDomain_m[1].min(); j<=(lDomain_m[1].max())/2.+1; ++j) {
                                        elem[1]=Index(j,j);
                                        elem_mirr[1]=Index(N-j+1,N-j+1);
                                        for (int k=lDomain_m[2].min(); k<=(lDomain_m[2].max())/2.+1; ++k) {
                                                elem[2]=Index(k,k);
                                                elem_mirr[2]=Index(N-k+1,N-k+1);

                                                double val = 0;
                                                double r = 0;
                                                r = sqrt((i)*(i)*hr_m[0]*hr_m[0]+(j)*(j)*hr_m[1]*hr_m[1]+(k)*(k)*hr_m[2]*hr_m[2]);
                                                val = 1./(4.*M_PI)*erf(alpha*r)/(r+eps);

                                                elem_all=elem;
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);
                                                //mirror on x axis
                                                elem_all[0]=elem_mirr[0];
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);
                                                //mirror on y axis
                                                elem_all[0]=elem[0];
                                                elem_all[1]=elem_mirr[1];
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);
                                                //mirror on z axis
                                                elem_all[1]=elem[1];
                                                elem_all[2]=elem_mirr[2];
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);
                                                //mirror on x,y axis
                                                elem_all[0]=elem_mirr[0];
                                                elem_all[1]=elem_mirr[1];
                                                elem_all[2]=elem[2];
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);
                                                //mirror on y,z axis
                                                elem_all[0]=elem[0];
                                                elem_all[1]=elem_mirr[1];
                                                elem_all[2]=elem_mirr[2];
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);
                                                //mirror on z,x axis
                                                elem_all[0]=elem_mirr[0];
                                                elem_all[1]=elem[1];
                                                elem_all[2]=elem_mirr[2];
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);
                                                //mirror on x,y,z axis
                                                elem_all[0]=elem_mirr[0];
                                                elem_all[1]=elem_mirr[1];
                                                elem_all[2]=elem_mirr[2];
                                                grncmpl_m.localElement(elem_all) = std::complex<double>(val);
                                                //rhocmpl_m.localElement(elem_all) *= std::complex<double>(val);

                                        }
                                }
                        }

                }



                void calculateGridForces(double /*interaction_radius*/, double alpha, double eps, bool GcalcKSpace, int it=0)
                {
                        //this->Q.scatter(this->rho_m, this->R, IntrplTSC_t());
                        rho_m[lDomain_m]=0; //!!!!!! there has to be a better way than setting rho to 0 every time
                        this->Q.scatter(this->rho_m, this->R, IntrplCIC_t());
                        //this->Q.scatter(this->rho_m, this->R, IntrplNGP_t());
                        //dumpVTKScalar(rho_m,this,it,"RhoInterpol");

                        //rhocmpl_m[domain_m] = rho_m[domain_m];
                        rhocmpl_m[lDomain_m] = rho_m[lDomain_m]/(hr_m[0]*hr_m[1]*hr_m[2]);
                        RhoSum=sum(real(rhocmpl_m));

                        std::cout << "total charge in densitty field before ion subtraction is" << sum(real(rhocmpl_m))<< std::endl;
                        std::cout << "max total charge in densitty field before ion subtraction is" << max(real(rhocmpl_m)) << std::endl;
                        //subtract the background charge of the ions
frey_m's avatar
frey_m committed
381
#ifndef SPHERE
382 383 384
                        //rhocmpl_m[domain_m]=rho_m[domain_m]+hr_m[0]*hr_m[1]*hr_m[2];
                        //rhocmpl_m[domain_m]=1./(nr_m[0]*nr_m[1]*nr_m[2])+rho_m[domain_m];
                        rhocmpl_m[lDomain_m]=1.+rhocmpl_m[lDomain_m];
frey_m's avatar
frey_m committed
385 386 387
#endif

#ifdef SPHERE
388 389
                        //rhocmpl_m[domain_m]=1./(nr_m[0]*nr_m[1]*nr_m[2])+rho_m[domain_m];
                        rhocmpl_m[lDomain_m]=1./(hr_m[0]*hr_m[1]*hr_m[2]*(nr_m[0]*nr_m[1]*nr_m[2]))+rhocmpl_m[lDomain_m];
frey_m's avatar
frey_m committed
390
#endif
391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508
                        std::cout << "total charge in densitty field after ion subtraction is" << sum(real(rhocmpl_m)) << std::endl;

                        dumpVTKScalar(rhocmpl_m,this,it,"RhoInterpolCompl");

                        //compute rhoHat and store in rhocmpl_m
                        fft_m->transform("inverse", rhocmpl_m);

                        if (GcalcKSpace)
                                std::cout << "NOT IMPLEMENTED" << std::endl;
                        else {
                                //compute G in real space and store in grncmpl_m
                                //calcGrealSpace(alpha,eps);

                                ////////compute G with Index Magic///////////////////
                                // Fields used to eliminate excess calculation in greensFunction()
                                IField_t grnIField_m[3];

                                // mesh and layout objects for rho_m
                                Mesh_t *mesh_m = &(getMesh());
                                FieldLayout_t *layout_m = &(getFieldLayout());

                                //This loop stores in grnIField_m[i] the index of the ith dimension mirrored at the central axis. e.g. grnIField_m[0]=[(0 1 2 3 ... 3 2 1) ; (0 1 2 3 ... 3 2 1; ...)]
                                for (int i = 0; i < 3; ++i) {
                                        grnIField_m[i].initialize(*mesh_m, *layout_m);
                                        grnIField_m[i][domain_m] = where(lt(domain_m[i], nr_m[i]/2),
                                                        domain_m[i] * domain_m[i],
                                                        (nr_m[i]-domain_m[i]) *
                                                        (nr_m[i]-domain_m[i]));
                                }

                                Vector_t hrsq(hr_m * hr_m);
                                SpecializedGreensFunction<3>::calculate(hrsq, grncmpl_m, grnIField_m, alpha,eps);
                                /////////////////////////////////////////////////
                                //dumpVTKScalar(grncmpl_m,this,it,"GRealSpace");

                                //transform G -> Ghat and store in grncmpl_m
                                fft_m->transform("inverse", grncmpl_m);
                                //grncmpl_m[0][0][0]=std::complex<double>(0);
                                //multiply in fourier space and obtain PhiHat in rhocmpl_m
                                rhocmpl_m *= grncmpl_m;
                        }

                        //compute electrostatic potential Phi in real space by FFT PhiHat -> Phi and store it in rhocmpl_m
                        fft_m->transform("forward", rhocmpl_m);

                        //take only the real part and store in phi_m (has periodic bc instead of interpolation bc)
                        phi_m[lDomain_m] = real(rhocmpl_m[lDomain_m])*hr_m[0]*hr_m[1]*hr_m[2];
                        //dumpVTKScalar(phi_m,this,it,"PhiRealSpace");
                        //dumpVTKScalar( phi_m, this,it, "Phi_m") ;

                        //compute Electric field on the grid by -Grad(Phi) store in eg_m
                        eg_m = -Grad1Ord(phi_m, eg_m);
                        //std::cout << "Sum of eg_m =" << sum(dot(eg_m,eg_m)) << std::endl;
                        //dumpVTKVector(eg_m, this,it+1,"EField");

                        //interpolate the electric field to the particle positions
                        //eg_m=eg_m*(hr_m[0]*hr_m[1]*hr_m[2]);
                        EF.gather(eg_m, this->R,  IntrplCIC_t());
                        //interpolate electrostatic potenital to the particle positions
                        Phi.gather(phi_m, this->R, IntrplCIC_t());
                }


                Vector_t getRmin() {
                        return this->rmin_m;
                }
                Vector_t getRmax() {
                        return this->rmax_m;
                }

                Vector_t get_hr() { return hr_m;}

                //private:
                BConds<double, Dim, Mesh_t, Center_t> bc_m;
                BConds<double, Dim, Mesh_t, Center_t> bcp_m;
                BConds<Vector_t, Dim, Mesh_t, Center_t> vbc_m;

                CxField_t rhocmpl_m;
                CxField_t grncmpl_m;

                Field_t rho_m;
                Field_t phi_m;

                VField_t eg_m;

                Vektor<int,Dim> nr_m;
                Vector_t hr_m;
                Vector_t rmax_m;
                Vector_t rmin_m;
                Vektor<double,3> extend_l;
                Vektor<double,3> extend_r;
                Mesh_t *mesh_m;
                FieldLayout_t *layout_m;
                NDIndex<Dim> domain_m;
                NDIndex<Dim> lDomain_m;

                double kinetic_energy;
                double field_energy;
                double field_energy_gather;
                double integral_phi_m;
                double potential_energy;
                double AmplitudeEfield;
                double AmplitudeEFz;
                double total_charge;
                double rhomax;

                FFT_t *fft_m;

                e_dim_tag decomp_m[Dim];

                Vektor<int,3> Nx;
                Vektor<int,3> Nv;
                Vektor<double,3> Vmax;
                //Fields for tracking distribution function
                Field2d_t f_m;
                Mesh2d_t mesh2d_m;
                NDIndex<2> domain2d_m;
                FieldLayout2d_t *layout2d_m;
frey_m's avatar
frey_m committed
509 510 511 512
};

template<class T>
struct ApplyField {
513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547
        ApplyField(T c, double r, double epsilon, double alpha) : C(c), R(r), eps(epsilon), a(alpha) {}
        void operator()(std::size_t i, std::size_t j, ChargedParticles<playout_t> &P,Vektor<double,3> &shift) const
        {
                Vector_t diff = P.R[i] - (P.R[j]+shift);
                double sqr = 0;

                for (unsigned d = 0; d<Dim; ++d)
                        sqr += diff[d]*diff[d];

                //compute r with softening parameter, unsoftened r is obtained by sqrt(sqr)
                if(sqr!=0) {
                        double r = std::sqrt(sqr+eps*eps);

                        //for order two transition
                        if (P.Q[i]!=0 && P.Q[j]!=0) {
                                //compute potential energy
                                double phi =1./(4.*M_PI)*(1.-erf(a*sqrt(sqr)))/r;

                                //compute force
                                Vector_t Fij = 1./(4.*M_PI)*C*(diff/sqrt(sqr))*((2.*a*exp(-a*a*sqr))/(sqrt(M_PI)*r)+(1.-erf(a*sqrt(sqr)))/(r*r));

                                //Actual Force is F_ij multiplied by Qi*Qj
                                //The electrical field on particle i is E=F/q_i and hence:
                                P.EF[i] -= P.Q[j]*Fij;
                                P.EF[j] += P.Q[i]*Fij;
                                //update potential per particle
                                P.Phi[i] += P.Q[j]*phi;
                                P.Phi[j] +=     P.Q[i]*phi;
                        }
                }
        }
        T C;
        double R;
        double eps;
        double a;
frey_m's avatar
frey_m committed
548 549 550
};

int main(int argc, char *argv[]){
551 552 553
        Ippl ippl(argc, argv);
        Inform msg(argv[0]);
        Inform msg2all(argv[0],INFORM_ALL_NODES);
frey_m's avatar
frey_m committed
554

555 556
        IpplTimings::TimerRef allTimer = IpplTimings::getTimer("AllTimer");
        IpplTimings::startTimer(allTimer);
frey_m's avatar
frey_m committed
557

558
        Vektor<int,Dim> nr;
frey_m's avatar
frey_m committed
559

560
        unsigned param = 1;
frey_m's avatar
frey_m committed
561

562 563 564 565 566 567 568
        if(Dim == 3) {
                nr = Vektor<int,Dim>(atoi(argv[1]),atoi(argv[2]),atoi(argv[3]));
                param = 4;
        } else {
                nr = Vektor<int,Dim>(atoi(argv[1]),atoi(argv[2]));
                param = 3;
        }
frey_m's avatar
frey_m committed
569

570
        double interaction_radius = atof(argv[param++]);
frey_m's avatar
frey_m committed
571
#ifndef SPHERE
572
        double k = 0.5;
frey_m's avatar
frey_m committed
573 574
#endif

575
        /////////// setup the computational domain boundaries /////////
frey_m's avatar
frey_m committed
576
#ifdef TWOSTREAM
577 578 579 580
        //Vektor<double,Dim> extend_r(.25*M_PI,.25*M_PI,2.*M_PI/k);
        Vektor<double,Dim> extend_r(2.*M_PI/k,2.*M_PI/k,2.*M_PI/k);
        Vektor<double,Dim> extend_l(0,0,0);
        //Vektor<double,Dim> extend_r(1,1,1);
frey_m's avatar
frey_m committed
581 582 583
#endif

#ifdef LANDAU
584 585 586
        Vektor<double,Dim> extend_l(0,0,0);
        Vektor<double,Dim> extend_r(2.*M_PI/k,2.*M_PI/k,2.*M_PI/k);
        //Vektor<double,Dim> extend_r(1,1,2.*M_PI/k);
frey_m's avatar
frey_m committed
587 588 589
#endif

#ifdef SPHERE
590 591 592
        double L = 4.;
        Vektor<double,Dim> extend_l(-L,-L,-L);
        Vektor<double,Dim> extend_r(L,L,L);
frey_m's avatar
frey_m committed
593
#endif
594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631
        ///////////////////////////////////////////////////////////////

        //read the remaining sim params
        double alpha =atof(argv[param++]);
        double dt = atof(argv[param++]);
        double eps = atof(argv[param++]);
        int iterations =  atoi(argv[param++]);
        Vektor<double,3> source(0,0,0);
        source[0]=atof(argv[param++]);
        source[1]=atof(argv[param++]);
        source[2]=atof(argv[param++]);



        ///////// setup the initial layout ///////////////////////////////////////
        e_dim_tag decomp[Dim];
        Mesh_t *mesh;
        FieldLayout_t *FL;
        ChargedParticles<playout_t>  *P;

        NDIndex<Dim> domain;
        for (unsigned i=0; i<Dim; i++)
                domain[i] = domain[i] = Index(nr[i]+1);
        std::cout << "Initially created domain = " << domain << std::endl;

        for (unsigned d=0; d < Dim; ++d)
                decomp[d] = PARALLEL;
        // create mesh and layout objects for this problem domain
        mesh          = new Mesh_t(domain);
        FL            = new FieldLayout_t(*mesh, decomp);
        playout_t* PL = new playout_t(*FL, *mesh);

        PL->setAllCacheDimensions(interaction_radius);
        PL->enableCaching();

        /////////////////////////////////////////////////////////////////////////

        /////////// Create the particle distribution /////////////////////////////////////////////////////
frey_m's avatar
frey_m committed
632
#ifdef TWOSTREAM
633 634 635 636 637 638 639 640 641 642 643 644
        Vektor<int,Dim> Nx(4,4,32);
        Vektor<int,Dim> Nv(8,8,128);
        Vektor<double,Dim> Vmax(6,6,6);
        double ampl_alpha = atof(argv[param++]);
        //refinement factor for mesh in 2d phase space
        int refine = 1;
        P = new ChargedParticles<playout_t>(PL, nr, decomp, extend_l, extend_r,refine*Nx,refine*Nv,Vmax);

        createParticleDistributionTwoStream(P,extend_l,extend_r,Nx,Nv,Vmax,ampl_alpha);

        P->interpolate_distribution((extend_r-extend_l)/(Nx),2.*Vmax/(Nv));
        write_f_field(P->f_m,P,0,"f_m");
frey_m's avatar
frey_m committed
645 646 647
#endif

#ifdef LANDAU
648 649 650 651 652 653 654 655 656 657 658 659
        Vektor<int,Dim> Nx(8,8,8);
        Vektor<int,Dim> Nv(32,32,32);
        Vektor<double,Dim> Vmax(6,6,6);
        double ampl_alpha = atof(argv[param++]);
        //refinement factor for mesh in 2d phase space
        int refine = 1;
        P = new ChargedParticles<playout_t>(PL, nr, decomp, extend_l, extend_r,refine*Nx,refine*Nv,Vmax);
        createParticleDistributionLandau(P,extend_l,extend_r,Nx,Nv,Vmax,ampl_alpha);
        //std::cout << "charge per particle pls: " << std::endl;
        //double qi;
        //std::cin >> qi;
        //createParticleDistribution(P, "even",10000, qi, extend_l,extend_r);
frey_m's avatar
frey_m committed
660 661 662 663
#endif


#ifdef SPHERE
664 665 666 667 668 669 670 671 672 673 674 675 676 677 678
        Vektor<int,Dim> Nx(16,16,16);
        Vektor<int,Dim> Nv(16,16,16);
        Vektor<double,Dim> Vmax(6,6,6);
        //refinement factor for mesh in 2d phase space
        int refine = 4;
        P = new ChargedParticles<playout_t>(PL, nr, decomp, extend_l, extend_r,refine*Nx,refine*Nv,Vmax);
        createParticleDistribution(P,"random", 20000, 0.00005, extend_l,extend_r,source,1.,0);
        //createParticleDistribution(P,"random", 10, 0.1, extend_l,extend_r,source,1.,0);
                //createParticleDistribution(P,"manual", 2, 1, extend_l,extend_r,1.)0
        for (unsigned i=0; i<P->getLocalNum(); ++i) {
                P->ID[i]=int(Ippl::myNode());
        }
        P->update();

        dumpParticlesCSV(P,0);
frey_m's avatar
frey_m committed
679 680
#endif

681
        /////////////////////////////////////////////////////////////////////////////////////////////
frey_m's avatar
frey_m committed
682

683 684 685 686 687
        /////// Print mesh informations ////////////////////////////////////////////////////////////
        INFOMSG(P->getMesh() << endl);
        INFOMSG(P->getFieldLayout() << endl);
        msg << endl << endl;
        Ippl::Comm->barrier();
frey_m's avatar
frey_m committed
688

689
        //dumpParticlesCSV(P,0);
frey_m's avatar
frey_m committed
690

691 692 693
        INFOMSG(P->getMesh() << endl);
        INFOMSG(P->getFieldLayout() << endl);
        msg << endl << endl;
frey_m's avatar
frey_m committed
694

695 696 697 698 699
        msg<<"number of particles = " << endl;
        msg<< P->getTotalNum() << endl;
        msg<<"Total charge Q = " << endl;
        msg<< P->total_charge << endl;
        ////////////////////////////////////////////////////////////////////////////////////////////
frey_m's avatar
frey_m committed
700

701 702 703
        //START TIMESTEPPING LOOP
        msg << "Starting iterations ..." << endl;
        bool GcalcKSpace = false;
frey_m's avatar
frey_m committed
704

705 706 707
        // calculate initial grid forces
        P->calculateGridForces(interaction_radius,alpha,eps,GcalcKSpace,999);
        dumpVTKVector(P->eg_m, P,0,"EFieldAfterPMandPP");
frey_m's avatar
frey_m committed
708

709 710 711 712 713
        //compute quantities to check correctness:
        P->calc_field_energy();
        P->calc_potential_energy();
        P->calc_kinetic_energy();
        writeEnergy(P,0);
frey_m's avatar
frey_m committed
714 715

#ifndef TWOSTREAM
716 717
        P->calc_Amplitude_E();
        writeEamplitude(P,0);
frey_m's avatar
frey_m committed
718
#endif
719 720 721 722 723 724 725
        for (int it=0; it<iterations; it++) {
                // advance the particle positions
                // basic leapfrogging timestep scheme.  velocities are offset
                // by half a timestep from the positions.
                assign(P->R, P->R + dt * P->v);
                // update particle distribution across processors
                P->update();
frey_m's avatar
frey_m committed
726

727 728 729 730
                // compute the electric field
                msg << "calculating grid" << endl;
                IpplTimings::TimerRef gridTimer = IpplTimings::getTimer("GridTimer");
                IpplTimings::startTimer(gridTimer);
frey_m's avatar
frey_m committed
731

732 733
                P->calculateGridForces(interaction_radius,alpha,eps,GcalcKSpace,it);
                IpplTimings::stopTimer(gridTimer);
frey_m's avatar
frey_m committed
734

735
                msg << "calculating pairs" << endl;
frey_m's avatar
frey_m committed
736

737 738
                IpplTimings::TimerRef particleTimer = IpplTimings::getTimer("ParticleTimer");
                IpplTimings::startTimer(particleTimer);
frey_m's avatar
frey_m committed
739

740 741
                P->calculatePairForces(interaction_radius,eps,alpha);
                IpplTimings::stopTimer(particleTimer);
frey_m's avatar
frey_m committed
742

743
                P->update();
frey_m's avatar
frey_m committed
744 745 746 747
#ifdef SPHERE
if (it%1==0)
dumpVTKVector(P->eg_m, P,it+1,"EFieldAfterPMandPP");
#endif
748 749
                //second part of leapfrog: advance velocitites
                //assign(P->P, P->P + dt * P->Q/P->m * P->EF);
frey_m's avatar
frey_m committed
750 751


752
                //assign(P->P, P->P + dt * P->EF);
frey_m's avatar
frey_m committed
753 754

#ifdef SPHERE
755 756 757 758 759 760 761 762 763 764 765
                //Print the particle positions
                if (it%1==0){
                        dumpParticlesCSV(P,it+1);
                }
#endif

                //compute quantities
                P->calc_field_energy();
                P->calc_kinetic_energy();
                P->calc_potential_energy();
                writeEnergy(P,it+1);
frey_m's avatar
frey_m committed
766
#ifndef TWOSTREAM
767 768
                P->calc_Amplitude_E();
                writeEamplitude(P,it+1);
frey_m's avatar
frey_m committed
769 770 771
#endif

#ifdef TWOSTREAM
772 773 774
                msg << "start interpolation to phase space " << endl;
                P->interpolate_distribution((extend_r-extend_l)/(Nx),2.*Vmax/(Nv));
                write_f_field(P->f_m,P,it+1,"f_m");
frey_m's avatar
frey_m committed
775
#endif
776 777 778
                msg << "Finished iteration " << it << endl;
        }
        Ippl::Comm->barrier();
frey_m's avatar
frey_m committed
779

780 781
        msg<<"number of particles = " << endl;
        msg<< P->getTotalNum() << endl;
frey_m's avatar
frey_m committed
782

783
        IpplTimings::stopTimer(allTimer);
frey_m's avatar
frey_m committed
784

785
        IpplTimings::print();
frey_m's avatar
frey_m committed
786 787 788


#ifdef CALC_ERRORS
789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875
//      Vektor<double,3> source(3.5,3.5,3.5);
        //Vektor<double,3> source(0,0,0);
        double sphere_radius = 1;
        std::cout << "total E-Field = " << sum(dot(P->EF,P->EF)) << std::endl;
        std::cout << "total Potential = " << sum(P->Phi) << std::endl;
        std::cout << "source = " << source << std::endl;
        std::cout << "RhoSum = " << RhoSum << std::endl;

        double error=0, errorV = 0;
        double l2_exact=0, l2_V_exact = 0;
        double total_E_exact=0, total_V_exact = 0;
        double total_E=0, total_V = 0;

        double total_charge = 1.;
        double k0 = 1./(4.*M_PI);

        for (unsigned i=0; i<P->getLocalNum(); ++i)
        {
                double radius = std::sqrt(dot(source - P->R[i], source - P->R[i]));
                double E = std::sqrt(dot(P->EF[i], P->EF[i]));
                double V = P->Phi[i];

                double diff=0, diffV = 0;
                double exact=0,exactV = 0;


                if(radius <= sphere_radius)
                {
                        exact = k0*(total_charge*radius/(sphere_radius*sphere_radius*sphere_radius));
                        diff = E - exact;

                        exactV=k0*total_charge/(2.*sphere_radius)*(3.-radius*radius/(sphere_radius*sphere_radius));
                        diffV = V-exactV;

                }
                else
                {
                        if(radius>0){

                                exact = k0*(total_charge*1./(radius*radius));
                                diff = E - exact;

                                exactV = k0*total_charge/radius;
                                diff= V-exactV;
                        }
                }

                total_E+=E;
                total_V+=V;
                total_E_exact+=exact;
                total_V_exact+=exactV;

                error += diff*diff;
                errorV += diffV*diffV;

                l2_exact += exact*exact;
                l2_V_exact += exactV*exactV;
        }
        //reduce all relevant quantities:
        double Error, L2_exact, ErrorV,L2_V_exact, Total_E, Total_V, Total_E_exact, Total_V_exact; // A, Potential_energy;
        reduce(error, Error, OpAddAssign());
        reduce(l2_exact, L2_exact, OpAddAssign());
        reduce(errorV, ErrorV, OpAddAssign());
        reduce(l2_V_exact, L2_V_exact, OpAddAssign());

        reduce(total_E, Total_E, OpAddAssign());
        reduce(total_V, Total_V, OpAddAssign());
        reduce(total_E_exact, Total_E_exact, OpAddAssign());
        reduce(total_V_exact, Total_V_exact, OpAddAssign());

        //reduce(P->potential_energy, Potential_energy, OpAddAssign());

        double Relative_error = std::sqrt(Error)/std::sqrt(L2_exact);
        double Relative_V_error = std::sqrt(ErrorV)/std::sqrt(L2_V_exact);
        double U = k0*3./5.*total_charge*total_charge/sphere_radius; //electric energy stored in solid charged sphere for infinite domain

        //      double Ufinite = 0.0394785; //electric energy stored in solid charged sphere for finite domain 8^3
        //double U = Ufinite;

        if(Ippl::myNode()==0) {
                std::cout << "master node prints: Q = " << total_charge << std::endl;
                std::ofstream ofs;
                ofs.open ("data/statistics.txt", std::ofstream::out | std::ofstream::app);
                //mesh size , n particle, r_cut, alpha, smoothing eps, absolut_err, relative error, relative error in total E-field, absolut_V_err, relative V error, relative error in total V, deviation in sum(U) from solid sphere
                ofs << nr[0] << "," << P->getTotalNum() << "," << interaction_radius << "," << eps << "," << alpha << "," << Error << "," << Relative_error << "," << fabs(Total_E-Total_E_exact)/Total_E_exact <<  "," << ErrorV << "," << Relative_V_error << "," << fabs(Total_V-Total_V_exact)/Total_V_exact << "," << std::abs(P->potential_energy-U)/U << std::endl;
                ofs.close();
        }
frey_m's avatar
frey_m committed
876 877 878 879 880

#endif



881 882 883
        delete P;
        delete FL;
        delete mesh;
frey_m's avatar
frey_m committed
884

885
        return 0;
frey_m's avatar
frey_m committed
886
}