OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
FMGPoissonSolver.h
Go to the documentation of this file.
1//
2// Class FMGPoissonSolver
3// Interface to the Fortran based AMR Poisson multigrid solver of BoxLib.
4// The Fortran solver is part of AMReX till version 18.07
5// (Link: https://github.com/AMReX-Codes/amrex/archive/18.07.tar.gz)
6//
7// Copyright (c) 2016 - 2020, Matthias Frey, Paul Scherrer Institute, Villigen PSI, Switzerland
8// All rights reserved.
9//
10// Implemented as part of the PhD thesis
11// "Precise Simulations of Multibunches in High Intensity Cyclotrons"
12//
13// This file is part of OPAL.
14//
15// OPAL is free software: you can redistribute it and/or modify
16// it under the terms of the GNU General Public License as published by
17// the Free Software Foundation, either version 3 of the License, or
18// (at your option) any later version.
19//
20// You should have received a copy of the GNU General Public License
21// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
22//
23#ifndef FMG_POISSON_SOLVER_H_
24#define FMG_POISSON_SOLVER_H_
25
27#include "Amr/AmrBoxLib.h"
28
29#include <AMReX_MultiFabUtil.H>
30#include <AMReX_BLFort.H>
31#include <AMReX_FMultiGrid.H>
32
33class FMGPoissonSolver : public AmrPoissonSolver< AmrBoxLib > {
34
35private:
37 typedef amrex::Vector<AmrBoxLib::AmrField_t*> AmrFieldContainer_pt;
43
44public:
45
53 FMGPoissonSolver(AmrBoxLib* itsAmrObject_p);
54
69 unsigned short baseLevel,
70 unsigned short finestLevel,
71 bool prevAsGuess = true);
72
73 double getXRangeMin(unsigned short level = 0);
74 double getXRangeMax(unsigned short level = 0);
75 double getYRangeMin(unsigned short level = 0);
76 double getYRangeMax(unsigned short level = 0);
77 double getZRangeMin(unsigned short level = 0);
78 double getZRangeMax(unsigned short level = 0);
79
80
85 Inform &print(Inform &os) const;
86
87
88private:
89
95 void initParameters_m();
96
109 double solveWithF90_m(const AmrFieldContainer_pt& rho,
110 const AmrFieldContainer_pt& phi,
111 const amrex::Vector< AmrFieldContainer_pt >& grad_phi_edge,
112 const GeomContainer_t& geom,
113 int baseLevel,
114 int finestLevel);
115
121 /*
122 void interpolate_m(AmrScalarFieldContainer_t& phi, const GeomContainer_t& geom,
123 double l0norm, int finestLevel);
124 */
125private:
126 int bc_m[2*AMREX_SPACEDIM];
127 double reltol_m;
128 double abstol_m;
129};
130
131
133 return fs.print(os);
134}
135
136#endif
Inform & operator<<(Inform &os, const FMGPoissonSolver &fs)
FRONT * fs
Definition: hypervolume.cpp:59
amr::AmrGeometry_t AmrGeometry_t
Definition: AmrBoxLib.h:52
amr::AmrProcMap_t AmrProcMap_t
Definition: AmrBoxLib.h:51
amr::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
Definition: AmrBoxLib.h:43
amr::AmrGrid_t AmrGrid_t
Definition: AmrBoxLib.h:50
amr::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
Definition: AmrBoxLib.h:42
amr::AmrGeomContainer_t AmrGeomContainer_t
Definition: AmrBoxLib.h:44
amrex::Vector< AmrBoxLib::AmrField_t * > AmrFieldContainer_pt
double abstol_m
Absolute tolerance for solver.
Inform & print(Inform &os) const
double getZRangeMax(unsigned short level=0)
double reltol_m
Relative tolearance for solver.
double getYRangeMin(unsigned short level=0)
AmrBoxLib::AmrGeomContainer_t GeomContainer_t
void solve(AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, unsigned short baseLevel, unsigned short finestLevel, bool prevAsGuess=true)
AmrBoxLib::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
FMGPoissonSolver(AmrBoxLib *itsAmrObject_p)
AmrBoxLib::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
AmrBoxLib::AmrProcMap_t AmrProcMap_t
double getXRangeMax(unsigned short level=0)
double solveWithF90_m(const AmrFieldContainer_pt &rho, const AmrFieldContainer_pt &phi, const amrex::Vector< AmrFieldContainer_pt > &grad_phi_edge, const GeomContainer_t &geom, int baseLevel, int finestLevel)
AmrBoxLib::AmrGeometry_t AmrGeometry_t
double getZRangeMin(unsigned short level=0)
double getXRangeMin(unsigned short level=0)
double getYRangeMax(unsigned short level=0)
int bc_m[2 *AMREX_SPACEDIM]
Boundary conditions.
AmrBoxLib::AmrGrid_t AmrGrid_t
Definition: Inform.h:42