ImpactX
Loading...
Searching...
No Matches
InitMeshRefinement.H
Go to the documentation of this file.
1/* Copyright 2022-2023 The Regents of the University of California, through Lawrence
2 * Berkeley National Laboratory (subject to receipt of any required
3 * approvals from the U.S. Dept. of Energy). All rights reserved.
4 *
5 * This file is part of ImpactX.
6 *
7 * Authors: Axel Huebl
8 * License: BSD-3-Clause-LBNL
9 */
10#pragma once
11
13
15
16#include <AMReX.H>
17#include <AMReX_BLProfiler.H>
18#include <AMReX_ParmParse.H>
19#include <AMReX_REAL.H>
20#include <AMReX_Utility.H>
21
22#include <limits>
23#include <stdexcept>
24#include <string>
25#include <vector>
26
27
29{
30 amrex::Vector<amrex::Real>
32 {
33 amrex::ParmParse pp_algo("algo");
34 amrex::ParmParse pp_amr("amr");
35 amrex::ParmParse pp_geometry("geometry");
36
37 auto space_charge = get_space_charge_algo();
38
39 std::string poisson_solver = "fft";
40 pp_algo.queryAdd("poisson_solver", poisson_solver);
41
42 int max_level = 0;
43 pp_amr.queryWithParser("max_level", max_level);
44
45 if (max_level > 1 && space_charge != SpaceChargeAlgo::False)
46 throw std::runtime_error(
47 "Mesh-refinement (amr.max_level>=0) is only supported with "
48 "space charge modeling (algo.space_charge=1).");
49
50 // The box is expanded beyond the min and max of the particle beam.
51 amrex::Vector<amrex::Real> prob_relative(max_level + 1, 1.0);
52 prob_relative[0] = 3.0; // top/bottom pad the beam on the lowest level by default by its width
53 pp_geometry.queryarr("prob_relative", prob_relative);
54
55 if (prob_relative[0] < 3.0 && space_charge != SpaceChargeAlgo::False && poisson_solver == "multigrid")
57 "ImpactX::read_mr_prob_relative",
58 "Dynamic resizing of the mesh uses a geometry.prob_relative "
59 "padding of less than 3 for level 0. This might result in boundary "
60 "artifacts for space charge calculation. "
61 "There is no minimum good value for this parameter, consider "
62 "doing a convergence test.",
63 ablastr::warn_manager::WarnPriority::high
64 );
65
66 if (prob_relative[0] < 1.0)
67 throw std::runtime_error("geometry.prob_relative must be >= 1.0 (the beam size) on the coarsest level");
68
69 // check that prob_relative[0] > prob_relative[1] > prob_relative[2] ...
70 amrex::Real last_lev_rel = std::numeric_limits<amrex::Real>::max();
71 for (int lev = 0; lev <= max_level; ++lev) {
72 amrex::Real const prob_relative_lvl = prob_relative[lev];
73 if (prob_relative_lvl <= 0.0)
74 throw std::runtime_error("geometry.prob_relative must be strictly positive for all levels");
75 if (prob_relative_lvl > last_lev_rel)
76 throw std::runtime_error("geometry.prob_relative must be descending over refinement levels");
77
78 last_lev_rel = prob_relative_lvl;
79 }
80
81 return prob_relative;
82 }
83} // namespace impactx::initialization
int queryWithParser(const char *name, bool &ref) const
int queryarr(const char *name, std::vector< int > &ref, int start_ix=FIRST, int num_val=ALL) const
int queryAdd(const char *name, T &ref)
void WMRecordWarning(const std::string &topic, const std::string &text, const WarnPriority &priority=WarnPriority::medium)
Definition AmrCoreData.cpp:18
amrex::Vector< amrex::Real > read_mr_prob_relative()
Definition InitMeshRefinement.H:31
SpaceChargeAlgo get_space_charge_algo()
Definition Algorithms.cpp:24