Skip to content

Commit cd2bd2d

Browse files
switch to direct mesh access
This patch disables all Precice mappings and replaces it by direct mesh access in both directions.
1 parent 54ad7be commit cd2bd2d

File tree

2 files changed

+32
-37
lines changed

2 files changed

+32
-37
lines changed

partitioned-heat-conduction-direct/nutils/heat.py

+21-14
Original file line numberDiff line numberDiff line change
@@ -20,7 +20,7 @@ def main(side='Dirichlet', n=10, degree=1, timestep=.1, alpha=3., beta=1.3):
2020
# define the Nutils mesh
2121
domain, geom = mesh.rectilinear([x_grid, y_grid])
2222
coupling_boundary = domain.boundary['right' if side == 'Dirichlet' else 'left']
23-
coupling_sample = coupling_boundary.sample('gauss', degree=degree * 2)
23+
read_sample = coupling_boundary.sample('gauss', degree=degree * 2)
2424

2525
# Nutils namespace
2626
ns = function.Namespace()
@@ -33,7 +33,7 @@ def main(side='Dirichlet', n=10, degree=1, timestep=.1, alpha=3., beta=1.3):
3333
ns.flux = 'basis_n ?fluxdofs_n' # heat flux
3434
ns.f = 'beta - 2 - 2 alpha' # rhs
3535
ns.uexact = '1 + x_0 x_0 + alpha x_1 x_1 + beta ?t' # analytical solution
36-
ns.readbasis = coupling_sample.basis()
36+
ns.readbasis = read_sample.basis()
3737
ns.readfunc = 'readbasis_n ?readdata_n'
3838

3939
# define the weak form
@@ -45,18 +45,27 @@ def main(side='Dirichlet', n=10, degree=1, timestep=.1, alpha=3., beta=1.3):
4545
else 'top,bottom,right'].integral('(u - uexact)^2 d:x' @ ns, degree=degree * 2)
4646

4747
if side == 'Dirichlet':
48-
sqr += coupling_sample.integral('(u - readfunc)^2 d:x' @ ns)
48+
sqr += read_sample.integral('(u - readfunc)^2 d:x' @ ns)
4949
else:
50-
res += coupling_sample.integral('basis_n readfunc d:x' @ ns)
50+
res += read_sample.integral('basis_n readfunc d:x' @ ns)
5151

5252
# preCICE setup
5353
interface = precice.Interface(side, "../precice-config.xml", 0, 1)
54-
mesh_id = interface.get_mesh_id(side + "-Mesh")
55-
vertex_ids = interface.set_mesh_vertices(mesh_id, coupling_sample.eval(ns.x))
54+
55+
mesh_id_read = interface.get_mesh_id("Dirichlet-Mesh" if side == "Dirichlet" else "Neumann-Mesh")
56+
mesh_id_write = interface.get_mesh_id("Neumann-Mesh" if side == "Dirichlet" else "Dirichlet-Mesh")
57+
58+
vertex_ids_read = interface.set_mesh_vertices(mesh_id_read, read_sample.eval(ns.x))
59+
interface.set_mesh_access_region(mesh_id_write, [.9, 1.1, -.1, 1.1])
60+
61+
precice_dt = interface.initialize()
62+
63+
vertex_ids_write, coords = interface.get_mesh_vertices_and_ids(mesh_id_write)
64+
write_sample = domain.locate(ns.x, coords, eps=1e-10)
5665
precice_write = functools.partial(interface.write_block_scalar_data,
57-
interface.get_data_id("Temperature" if side == "Neumann" else "Heat-Flux", mesh_id), vertex_ids)
66+
interface.get_data_id("Heat-Flux" if side == "Dirichlet" else "Temperature", mesh_id_write), vertex_ids_write)
5867
precice_read = functools.partial(interface.read_block_scalar_data,
59-
interface.get_data_id("Heat-Flux" if side == "Neumann" else "Temperature", mesh_id), vertex_ids)
68+
interface.get_data_id("Temperature" if side == "Dirichlet" else "Heat-Flux", mesh_id_read), vertex_ids_read)
6069

6170
# helper functions to project heat flux to coupling boundary
6271
if side == 'Dirichlet':
@@ -76,13 +85,11 @@ def main(side='Dirichlet', n=10, degree=1, timestep=.1, alpha=3., beta=1.3):
7685
fluxsqr = domain.boundary['right'].boundary['top,bottom'].integral('(flux - uexact_,0)^2 d:x' @ ns, degree=degree*2)
7786
fluxcons = solver.optimize('fluxdofs', fluxsqr, droptol=1e-10, constrain=np.choose(np.isnan(rightcons), [np.nan, 0.]))
7887
# fluxcons is NaN in dofs that are supported on ONLY the right boundary
79-
fluxres = coupling_sample.integral('basis_n flux d:x' @ ns) - res
80-
81-
precice_dt = interface.initialize()
88+
fluxres = read_sample.integral('basis_n flux d:x' @ ns) - res
8289

8390
# write initial data
8491
if interface.is_action_required(precice.action_write_initial_data()):
85-
precice_write(coupling_sample.eval(0.))
92+
precice_write(write_sample.eval(0.))
8693
interface.mark_action_fulfilled(precice.action_write_initial_data())
8794

8895
interface.initialize_data()
@@ -130,9 +137,9 @@ def main(side='Dirichlet', n=10, degree=1, timestep=.1, alpha=3., beta=1.3):
130137
if interface.is_write_data_required(dt):
131138
if side == 'Dirichlet':
132139
fluxdofs = solver.solve_linear('fluxdofs', fluxres, arguments=dict(lhs0=lhs0, lhs=lhs, dt=dt, t=t), constrain=fluxcons)
133-
write_data = coupling_sample.eval('flux' @ ns, fluxdofs=fluxdofs)
140+
write_data = write_sample.eval('flux' @ ns, fluxdofs=fluxdofs)
134141
else:
135-
write_data = coupling_sample.eval('u' @ ns, lhs=lhs)
142+
write_data = write_sample.eval('u' @ ns, lhs=lhs)
136143
precice_write(write_data)
137144

138145
# do the coupling

partitioned-heat-conduction-direct/precice-config.xml

+11-23
Original file line numberDiff line numberDiff line change
@@ -7,8 +7,8 @@
77
enabled="true" />
88
</log>
99

10-
<solver-interface dimensions="2">
11-
<data:scalar name="Temperature" />
10+
<solver-interface dimensions="2" experimental="yes">
11+
<data:scalar name="Temperature" />
1212
<data:scalar name="Heat-Flux" />
1313

1414
<mesh name="Dirichlet-Mesh">
@@ -23,28 +23,16 @@
2323

2424
<participant name="Dirichlet">
2525
<use-mesh name="Dirichlet-Mesh" provide="yes" />
26-
<use-mesh name="Neumann-Mesh" from="Neumann" />
27-
<write-data name="Heat-Flux" mesh="Dirichlet-Mesh" />
26+
<use-mesh name="Neumann-Mesh" from="Neumann" direct-access="true" />
27+
<write-data name="Heat-Flux" mesh="Neumann-Mesh" />
2828
<read-data name="Temperature" mesh="Dirichlet-Mesh" />
29-
<mapping:rbf-thin-plate-splines
30-
direction="read"
31-
from="Neumann-Mesh"
32-
to="Dirichlet-Mesh"
33-
constraint="consistent"
34-
x-dead="true" />
3529
</participant>
3630

3731
<participant name="Neumann">
3832
<use-mesh name="Neumann-Mesh" provide="yes" />
39-
<use-mesh name="Dirichlet-Mesh" from="Dirichlet" />
40-
<write-data name="Temperature" mesh="Neumann-Mesh" />
33+
<use-mesh name="Dirichlet-Mesh" from="Dirichlet" direct-access="true" />
34+
<write-data name="Temperature" mesh="Dirichlet-Mesh" />
4135
<read-data name="Heat-Flux" mesh="Neumann-Mesh" />
42-
<mapping:rbf-thin-plate-splines
43-
direction="read"
44-
from="Dirichlet-Mesh"
45-
to="Neumann-Mesh"
46-
constraint="consistent"
47-
x-dead="true" />
4836
</participant>
4937

5038
<m2n:sockets from="Dirichlet" to="Neumann" exchange-directory=".." />
@@ -54,17 +42,17 @@
5442
<max-time value="1.0" />
5543
<time-window-size value="0.1" />
5644
<max-iterations value="100" />
57-
<exchange data="Heat-Flux" mesh="Dirichlet-Mesh" from="Dirichlet" to="Neumann" />
45+
<exchange data="Heat-Flux" mesh="Neumann-Mesh" from="Dirichlet" to="Neumann" />
5846
<exchange
5947
data="Temperature"
60-
mesh="Neumann-Mesh"
48+
mesh="Dirichlet-Mesh"
6149
from="Neumann"
6250
to="Dirichlet"
6351
initialize="true" />
64-
<relative-convergence-measure data="Heat-Flux" mesh="Dirichlet-Mesh" limit="1e-5" />
65-
<relative-convergence-measure data="Temperature" mesh="Neumann-Mesh" limit="1e-5" />
52+
<relative-convergence-measure data="Heat-Flux" mesh="Neumann-Mesh" limit="1e-5" suffices="yes" />
53+
<relative-convergence-measure data="Temperature" mesh="Dirichlet-Mesh" limit="1e-5" suffices="yes" />
6654
<acceleration:IQN-ILS>
67-
<data name="Temperature" mesh="Neumann-Mesh" />
55+
<data name="Temperature" mesh="Dirichlet-Mesh" />
6856
<initial-relaxation value="0.1" />
6957
<max-used-iterations value="10" />
7058
<time-windows-reused value="5" />

0 commit comments

Comments
 (0)