Skip to content

API of modul tough3

outfile2vtu(path, mesh, outfile='OUTFILE', do_sort=True)

Reads an OUTFILE from a TOUGH3 simulation and append simulation results to the vtu mesh for each timestamp in the OUTFILE. Also writes a pvd file to relate files and timestamps in Paraview.

Parameters

path : string path to TOUGH3 OUTFILE mesh : string name of vtu mesh in path outfile : string name of outfile. Default = "OUTFILE". do_sort : boolean parallel simulations may change the order of elements. For this reason sorting may be a good idea. Default = True.

Returns

None.

Source code in voromesh/tough3/__init__.py
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
def outfile2vtu(path, mesh, outfile="OUTFILE", do_sort=True):
    """
    Reads an OUTFILE from a TOUGH3 simulation and append simulation results
    to the vtu mesh for each timestamp in the OUTFILE. Also writes a pvd file
    to relate files and timestamps in Paraview.

    Parameters
    ----------
    path : string
        path to TOUGH3 OUTFILE
    mesh : string
        name of vtu mesh in path
    outfile : string
        name of outfile. Default = "OUTFILE".
    do_sort : boolean
        parallel simulations may change the order of elements. For this
        reason sorting may be a good idea. Default = True.

    Returns
    -------
    None.

    """

    try:
        path_outfile = os.path.join(path, outfile)
        outputs = toughio.read_output(path_outfile)
    except OSError:
        print("Error reading OUTFILE: " + path_outfile)
        sys.exit()

    try:
        path_meshfile = os.path.join(path, mesh)
        mesh = pv.read(path_meshfile)
    except OSError:
        print("Error reading mesh file: " + path_meshfile)
        sys.exit

    pvdfile = os.path.join(path, "outfile.pvd")
    pvd = pvdwriter(pvdfile)

    # Parallel simulations may change the order of elements.
    if do_sort:
        nn = outputs[-1][3]
        ind = np.argsort(nn)

    i = 0

    for output in outputs:
        for key in output[4]:

            if do_sort:
                mesh[key] = output[4][key][ind]
            else:
                mesh[key] = output[4][key]


        timestamp = output[2]/(365.25*24*60*60)
        filename = "mesh_" + str(i) + ".vtu"
        print(filename + " | "+ str(timestamp))

        pvd.append(filename, timestamp)

        mesh.save(os.path.join(path, filename))
        i = i + 1

    pvd.close()
    print("complete!")

update_gener(path2infile, mesh, materials, wells)

Update the well information in the TOUGH3 INFILE

Parameters

path2infile : string path to TOUGH3 INFILE. mesh : pyvista.UnstructuredGrid with cell data see write_mesh function. materials : dict see write_mesh function. wells : dict Define inflow and outflows. Key of dict is the material. The flow is splitted according to the volume fraction, when more than one cell has the same material. Value of dict is a dict with flow components. See TOUGH3 manual and TOUGH3 EOS manuals for details. Example: wells = {"WELL": {"COM1": 28.53, "COM3": 1e-6}}

Returns

None.

Source code in voromesh/tough3/__init__.py
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
def update_gener(path2infile, mesh, materials, wells):
    """
    Update the well information in the TOUGH3 INFILE

    Parameters
    ----------
    path2infile : string
        path to TOUGH3 INFILE.
    mesh : pyvista.UnstructuredGrid with cell data
        see write_mesh function.
    materials : dict
        see write_mesh function.
    wells : dict
        Define inflow and outflows.
        Key of dict is the material. The flow is splitted according to the
        volume fraction, when more than one cell has the same material.
        Value of dict is a dict with flow components. See TOUGH3 manual and
        TOUGH3 EOS manuals for details.
        Example: wells = {"WELL": {"COM1": 28.53, "COM3": 1e-6}}


    Returns
    -------
    None.

    """

    material = [materials[i] for i in mesh["material"]]
    volume = mesh["Volume"]

    _update_gener(material, volume, wells, path2infile,
                  names=False, verbose=False)

write_mesh(path, mesh, materials)

Writes a mesh in TOUGH2/TOUGH3 format.

Parameters

path : string Path for writing the TOUGH3 mesh. mesh : pyvista.UnstructuredGrid with cell data (see below) mesh["Volume"] : cell volumes mesh["material"] : cell material as integer mesh["initial_condition"] : cell initial conditions as list find more information in TOUGH EOS manuals materials : dict Relate numbers in mesh["material"] to material names from the TOUGH INFILE. Example: materials = {1: "WELL", 2: "SAND"}

Returns

None.

Source code in voromesh/tough3/__init__.py
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
def write_mesh(path, mesh, materials):
    """
    Writes a mesh in TOUGH2/TOUGH3 format.

    Parameters
    ----------
    path : string
        Path for writing the TOUGH3 mesh.
    mesh : pyvista.UnstructuredGrid with cell data (see below)
        mesh["Volume"] :  cell volumes
        mesh["material"] : cell material as integer
        mesh["initial_condition"] : cell initial conditions as list
                                    find more information in TOUGH EOS manuals
    materials : dict
        Relate numbers in mesh["material"] to material names from the TOUGH
        INFILE. Example: materials = {1: "WELL", 2: "SAND"}

    Returns
    -------
    None.

    """
    dist, areai, areao, betax_ = calc_conne(mesh)

    volume = mesh["Volume"]                             # 1
    material = [materials[i] for i in mesh["material"]]  # 2

    area = np.zeros(np.size(volume))
    for key in areao:
        area[key] = areao[key]                          # 3

    centers = mesh.cell_centers().points                # 4

    ne = list(areai.keys())                             # 5

    d1 = 0.5 * np.fromiter(dist.values(), dtype=float)  # 7
    d2 = d1                                             # 8

    areax = np.fromiter(areai.values(), dtype=float)    # 9

    betax = np.fromiter(betax_.values(), dtype=float)   # 10

    incon = mesh["initial_condition"]                   # 11

    # Nur Unterschiede in h (1,2 -> 1) und v (3)
    isot = np.ones(np.size(areax))                      # 6
    ind = np.abs(betax) > 0.5
    isot[ind] = 3

    _write_mesh(path, volume, material, area, centers,
                ne, isot, d1, d2, areax, betax, incon)