"examples/wave/wave-op.py" did not exist on "9cc11ed2ef456b51575f09148b34bf36b7fc2ca0"
Newer
Older
__copyright__ = """
Copyright (C) 2017 Ellis Hoag
Copyright (C) 2017 Andreas Kloeckner
Copyright (C) 2021 University of Illinois Board of Trustees
"""
__license__ = """
Permission is hereby granted, free of charge, to any person obtaining a copy
of this software and associated documentation files (the "Software"), to deal
in the Software without restriction, including without limitation the rights
to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
copies of the Software, and to permit persons to whom the Software is
furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in
all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
THE SOFTWARE.
"""
import pytest
import os
from grudge.array_context import MPIPyOpenCLArrayContext, MPIPytatoArrayContext
logging.basicConfig()
logger.setLevel(logging.INFO)
from pytools.obj_array import flat_obj_array
import grudge.op as op
import grudge.dof_desc as dof_desc
class SimpleTag:
pass
DISTRIBUTED_ACTXS = [MPIPyOpenCLArrayContext, MPIPytatoArrayContext]
def run_test_with_mpi(num_ranks, f, *args):
import pytest
pytest.importorskip("mpi4py")
from pickle import dumps
from base64 import b64encode
invocation_info = b64encode(dumps((f, args))).decode()
from subprocess import check_call
# NOTE: CI uses OpenMPI; -x to pass env vars. MPICH uses -env
check_call([
"mpiexec", "-np", str(num_ranks),
"-x", "RUN_WITHIN_MPI=1",
"-x", f"INVOCATION_INFO={invocation_info}",
sys.executable, __file__])
def run_test_with_mpi_inner():
from pickle import loads
from base64 import b64decode
f, (actx_class, *args) = loads(b64decode(os.environ["INVOCATION_INFO"].encode()))
cl_context = cl.create_some_context()
queue = cl.CommandQueue(cl_context)
if actx_class is MPIPytatoArrayContext:
actx = actx_class(comm, queue, mpi_base_tag=15000)
elif actx_class is MPIPyOpenCLArrayContext:
actx = actx_class(comm, queue, force_device_scalars=True)
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
else:
raise ValueError("unknown actx_class")
f(actx, *args)
# }}}
# {{{ func_comparison
@pytest.mark.parametrize("actx_class", DISTRIBUTED_ACTXS)
@pytest.mark.parametrize("num_ranks", [2])
def test_func_comparison_mpi(actx_class, num_ranks):
run_test_with_mpi(
num_ranks, _test_func_comparison_mpi_communication_entrypoint,
actx_class)
def _test_func_comparison_mpi_communication_entrypoint(actx):
"""Discretize a function, communicate it, check that it matches the
function discretized by the other end.
"""
comm = actx.mpi_communicator
from meshmode.distributed import MPIMeshDistributor, get_partition_by_pymetis
from meshmode.mesh import BTAG_ALL
num_parts = comm.Get_size()
mesh_dist = MPIMeshDistributor(comm)
if mesh_dist.is_mananger_rank():
from meshmode.mesh.generation import generate_regular_rect_mesh
nelements_per_axis=(2,)*2)
part_per_element = get_partition_by_pymetis(mesh, num_parts)
local_mesh = mesh_dist.send_mesh_parts(mesh, part_per_element, num_parts)
else:
local_mesh = mesh_dist.receive_mesh_part()
dcoll = DiscretizationCollection(actx, local_mesh, order=5)
x = actx.thaw(dcoll.nodes())
myfunc = actx.np.sin(np.dot(x, [2, 3]))
from grudge.dof_desc import as_dofdesc
dd_int = as_dofdesc("int_faces")
dd_vol = as_dofdesc("vol")
dd_af = as_dofdesc("all_faces")
all_faces_func = op.project(dcoll, dd_vol, dd_af, myfunc)
int_faces_func = op.project(dcoll, dd_vol, dd_int, myfunc)
bdry_faces_func = op.project(dcoll, BTAG_ALL, dd_af,
op.project(dcoll, dd_vol, BTAG_ALL, myfunc))
Andreas Klöckner
committed
def hopefully_zero():
return (
op.project(
dcoll, "int_faces", "all_faces",
dcoll.opposite_face_connection(
dof_desc.BoundaryDomainTag(
dof_desc.FACE_RESTR_INTERIOR, dof_desc.VTAG_ALL)
)(int_faces_func)
Andreas Klöckner
committed
)
+ sum(op.project(dcoll, tpair.dd, "all_faces", tpair.ext)
for tpair in op.cross_rank_trace_pairs(dcoll, myfunc,
comm_tag=SimpleTag))
) - (all_faces_func - bdry_faces_func)
hopefully_zero_result = actx.compile(hopefully_zero)()
error = actx.to_numpy(flat_norm(hopefully_zero_result, ord=np.inf))
with np.printoptions(threshold=100000000, suppress=True):
logger.debug(hopefully_zero)
logger.info("error: %.5e", error)
# {{{ wave operator
@pytest.mark.parametrize("actx_class", DISTRIBUTED_ACTXS)
@pytest.mark.parametrize("num_ranks", [2])
def test_mpi_wave_op(actx_class, num_ranks):
run_test_with_mpi(num_ranks, _test_mpi_wave_op_entrypoint, actx_class)
def _test_mpi_wave_op_entrypoint(actx, visualize=False):
comm = actx.mpi_communicator
from meshmode.distributed import MPIMeshDistributor, get_partition_by_pymetis
if mesh_dist.is_mananger_rank():
from meshmode.mesh.generation import generate_regular_rect_mesh
mesh = generate_regular_rect_mesh(a=(-0.5,)*dim,
b=(0.5,)*dim,
nelements_per_axis=(16,)*dim)
part_per_element = get_partition_by_pymetis(mesh, num_parts)
local_mesh = mesh_dist.send_mesh_parts(mesh, part_per_element, num_parts)
else:
local_mesh = mesh_dist.receive_mesh_part()
dcoll = DiscretizationCollection(actx, local_mesh, order=order)
def source_f(actx, dcoll, t=0):
source_center = np.array([0.1, 0.22, 0.33])[:dcoll.dim]
source_width = 0.05
source_omega = 3
nodes = actx.thaw(dcoll.nodes())
source_center_dist = flat_obj_array(
[nodes[i] - source_center[i] for i in range(dcoll.dim)]
)
return (
* actx.np.exp(
-np.dot(source_center_dist, source_center_dist)
/ source_width**2
)
)
Andreas Klöckner
committed
from grudge.models.wave import WeakWaveOperator
wave_op = WeakWaveOperator(
dcoll,
0.1,
source_f=source_f,
dirichlet_tag=BTAG_NONE,
neumann_tag=BTAG_NONE,
radiation_tag=BTAG_ALL,
flux_type="upwind",
comm_tag=SimpleTag,
)
fields = flat_obj_array(
dcoll.zeros(actx),
[dcoll.zeros(actx) for i in range(dcoll.dim)]
)
wave_op.estimate_rk4_timestep(actx, dcoll, fields=fields))
from logpyle import LogManager, \
# NOTE: LogManager hangs when using a file on a shared directory.
logmgr = LogManager(log_filename, "w", comm)
add_run_info(logmgr)
add_general_quantities(logmgr)
logger.info("[%04d] dt %.5e nsteps %4d", i_local_rank, dt, nsteps)
if visualize:
from grudge.shortcuts import make_visualizer
vis = make_visualizer(dcoll)
logmgr.tick_before()
for step in range(nsteps):
t = step*dt
fields = rk4_step(fields, t=t, h=dt, f=compiled_rhs)
fields = actx.thaw(actx.freeze(fields))
norm = actx.to_numpy(op.norm(dcoll, fields, 2))
logger.info("[%04d] t = %.5e |u| = %.5e elapsed %.5e",
step, t, norm, time() - t_last_step)
if visualize:
vis.write_parallel_vtk_file(
comm,
f"fld-wave-mpi-{type(actx).__name__}-{{rank:03d}}-{step:04d}.vtu",
[
("u", fields[0]),
("v", fields[1:]),
]
)
assert norm < 1
t_last_step = time()
logmgr.tick_after()
logmgr.tick_before()
logger.info("Rank %d exiting", i_local_rank)
run_test_with_mpi_inner()
elif len(sys.argv) > 1:
exec(sys.argv[1])
from pytest import main
main([__file__])