gel.scripts.get_u_from_gel_mesh

Interfaces with fmtrack GPR-interpolation to obtain target displacements.

An FM-TRACK conda environment must be available on the system with the name "fmtrack". Will use a GPR model saved to a directory with a specific format readable by a subprocess that invokes get_displacements_from_gpr.py.

  1"""Interfaces with fmtrack GPR-interpolation to obtain target displacements.
  2
  3An FM-TRACK conda environment must be available on the system with the
  4name "fmtrack". Will use a GPR model saved to a directory with a
  5specific format readable by a subprocess that invokes
  6`get_displacements_from_gpr.py`.
  7"""
  8from dolfin import *
  9import argparse
 10import os
 11import subprocess
 12import io
 13import numpy as np
 14
 15
 16def get_exp_u_xdmf(
 17        cell_data_dir,
 18        gpr_dir,
 19        outfile
 20    ):
 21    """Creates subprocess to perform GPR-interpolation on the mesh.
 22
 23    * `cell_data_dir`: str path to directory with information needed
 24    to create `gel.geometry.Geometry`
 25    * `gpr_dir`: str path to directory with GPR model files
 26    * `outfile`: str path to .xdmf file that will contain full-shape
 27    experimental, target displacements "u". If None, will write to
 28    `cell_data_dir` as `u_experiment.xdmf`
 29
 30    Side-effects: writes new file `outfile`
 31    """
 32    # Gel Volume Mesh
 33    mesh = Mesh()
 34    with XDMFFile(
 35        os.path.join(cell_data_dir, "reference_domain.xdmf")
 36    ) as infile:
 37        infile.read(mesh)
 38
 39    # Get vertices
 40    gel_vertices = mesh.coordinates()
 41
 42    #
 43    # Pass vertices to GPR subprocess
 44    #
 45    # Write to IO object
 46    s = io.BytesIO()
 47    np.savetxt(s, gel_vertices)
 48
 49    # Spawn process to use GPR model
 50    script_path = os.path.join(
 51        os.path.dirname(os.path.realpath(__file__)),
 52        "get_displacement_from_gpr.py"
 53    )
 54    conda_env = os.path.join(
 55        os.path.dirname(os.environ["CONDA_PREFIX"]),
 56        "fmtrack"
 57    )
 58    script_cmdline = (
 59        f"conda run --no-capture-output -p {conda_env} {script_path}"
 60        f" -v /dev/stdin -d /dev/stdout -g "+'"'+f"{gpr_dir}"+'"'
 61    )
 62    completed_proc = subprocess.run(
 63        script_cmdline,
 64        shell=True,
 65        input=s.getvalue(),
 66        capture_output=True
 67    )
 68
 69    # Read in from IO object
 70    gel_displacements = np.loadtxt(io.BytesIO(completed_proc.stdout))
 71    print(completed_proc.stderr)
 72
 73    # Convert to FE function, then output
 74    V = VectorFunctionSpace(mesh, "P", 1, dim=3)
 75    v2d = vertex_to_dof_map(V)
 76    u = Function(V)
 77
 78    gel_displacements_flat = gel_displacements.flatten()
 79    u.vector()[v2d] = gel_displacements_flat
 80
 81    V = VectorFunctionSpace(mesh, "Lagrange", 1)
 82    u = project(u, V, solver_type='cg', preconditioner_type='amg')
 83    u.rename("u", "displacement")
 84
 85    if outfile is None:
 86        outfile = os.path.join(cell_data_dir, "u_experiment.xdmf")
 87    u_outfile = XDMFFile(outfile)
 88    u_outfile.write_checkpoint(u, "u", 0) # Not appending
 89
 90
 91def get_u_main():
 92    """The function invoked by the command. Parses arguments and passes
 93    to `get_exp_u_xdmf`.
 94    """
 95    parser = argparse.ArgumentParser(
 96        description="Interface for cross-conda-environment "
 97        "communication between newestfenics and fmtrack in order to "
 98        "create a full-shape .xdmf 1st order Lagrange representation "
 99        "of experimental displacents by means of interpolating the GPR "
100        "model. This component is to be run in the newestfenics "
101        "environment."
102    )
103    parser.add_argument(
104        "-c",
105        type=str,
106        metavar="CELL_DATA_DIR",
107        help="directory containing gel geometry"
108    )
109    parser.add_argument(
110        "-g",
111        type=str,
112        metavar="GPR_DIR",
113        help="directory with GPR model files, i.e. gp_U_cleaned.sav etc."
114    )
115    parser.add_argument(
116        "-o",
117        type=str,
118        metavar="OUT_FILE",
119        default=None,
120        help="output full-shape 1st order Lagrange .xdmf file with "
121        "displacement 'u'"
122    )
123    args = parser.parse_args()
124
125    get_exp_u_xdmf(args.c, args.g, args.o)
126
127
128if __name__=="__main__":
129    get_u_main()
def get_exp_u_xdmf(cell_data_dir, gpr_dir, outfile):
17def get_exp_u_xdmf(
18        cell_data_dir,
19        gpr_dir,
20        outfile
21    ):
22    """Creates subprocess to perform GPR-interpolation on the mesh.
23
24    * `cell_data_dir`: str path to directory with information needed
25    to create `gel.geometry.Geometry`
26    * `gpr_dir`: str path to directory with GPR model files
27    * `outfile`: str path to .xdmf file that will contain full-shape
28    experimental, target displacements "u". If None, will write to
29    `cell_data_dir` as `u_experiment.xdmf`
30
31    Side-effects: writes new file `outfile`
32    """
33    # Gel Volume Mesh
34    mesh = Mesh()
35    with XDMFFile(
36        os.path.join(cell_data_dir, "reference_domain.xdmf")
37    ) as infile:
38        infile.read(mesh)
39
40    # Get vertices
41    gel_vertices = mesh.coordinates()
42
43    #
44    # Pass vertices to GPR subprocess
45    #
46    # Write to IO object
47    s = io.BytesIO()
48    np.savetxt(s, gel_vertices)
49
50    # Spawn process to use GPR model
51    script_path = os.path.join(
52        os.path.dirname(os.path.realpath(__file__)),
53        "get_displacement_from_gpr.py"
54    )
55    conda_env = os.path.join(
56        os.path.dirname(os.environ["CONDA_PREFIX"]),
57        "fmtrack"
58    )
59    script_cmdline = (
60        f"conda run --no-capture-output -p {conda_env} {script_path}"
61        f" -v /dev/stdin -d /dev/stdout -g "+'"'+f"{gpr_dir}"+'"'
62    )
63    completed_proc = subprocess.run(
64        script_cmdline,
65        shell=True,
66        input=s.getvalue(),
67        capture_output=True
68    )
69
70    # Read in from IO object
71    gel_displacements = np.loadtxt(io.BytesIO(completed_proc.stdout))
72    print(completed_proc.stderr)
73
74    # Convert to FE function, then output
75    V = VectorFunctionSpace(mesh, "P", 1, dim=3)
76    v2d = vertex_to_dof_map(V)
77    u = Function(V)
78
79    gel_displacements_flat = gel_displacements.flatten()
80    u.vector()[v2d] = gel_displacements_flat
81
82    V = VectorFunctionSpace(mesh, "Lagrange", 1)
83    u = project(u, V, solver_type='cg', preconditioner_type='amg')
84    u.rename("u", "displacement")
85
86    if outfile is None:
87        outfile = os.path.join(cell_data_dir, "u_experiment.xdmf")
88    u_outfile = XDMFFile(outfile)
89    u_outfile.write_checkpoint(u, "u", 0) # Not appending

Creates subprocess to perform GPR-interpolation on the mesh.

  • cell_data_dir: str path to directory with information needed to create gel.geometry.Geometry
  • gpr_dir: str path to directory with GPR model files
  • outfile: str path to .xdmf file that will contain full-shape experimental, target displacements "u". If None, will write to cell_data_dir as u_experiment.xdmf

Side-effects: writes new file outfile

def get_u_main():
 92def get_u_main():
 93    """The function invoked by the command. Parses arguments and passes
 94    to `get_exp_u_xdmf`.
 95    """
 96    parser = argparse.ArgumentParser(
 97        description="Interface for cross-conda-environment "
 98        "communication between newestfenics and fmtrack in order to "
 99        "create a full-shape .xdmf 1st order Lagrange representation "
100        "of experimental displacents by means of interpolating the GPR "
101        "model. This component is to be run in the newestfenics "
102        "environment."
103    )
104    parser.add_argument(
105        "-c",
106        type=str,
107        metavar="CELL_DATA_DIR",
108        help="directory containing gel geometry"
109    )
110    parser.add_argument(
111        "-g",
112        type=str,
113        metavar="GPR_DIR",
114        help="directory with GPR model files, i.e. gp_U_cleaned.sav etc."
115    )
116    parser.add_argument(
117        "-o",
118        type=str,
119        metavar="OUT_FILE",
120        default=None,
121        help="output full-shape 1st order Lagrange .xdmf file with "
122        "displacement 'u'"
123    )
124    args = parser.parse_args()
125
126    get_exp_u_xdmf(args.c, args.g, args.o)

The function invoked by the command. Parses arguments and passes to get_exp_u_xdmf.