This page was generated from mt3d-usgs_example.py. It's also available as a notebook.

MT3D-USGS Example

Demonstrates functionality of the flopy MT3D-USGS module using the ‘Crank-Nicolson’ example distributed with MT3D-USGS.

Problem description:

  • Grid dimensions: 1 Layer, 3 Rows, 650 Columns

  • Stress periods: 3

  • Units are in seconds and meters

  • Flow package: UPW

  • Stress packages: SFR, GHB

  • Solvers: NWT, GCG

[1]:
import os
import string
[2]:
import sys
from pprint import pformat
from tempfile import TemporaryDirectory

import matplotlib as mpl
import matplotlib.pyplot as plt
import numpy as np

import flopy

print(sys.version)
print(f"numpy version: {np.__version__}")
print(f"matplotlib version: {mpl.__version__}")
print(f"flopy version: {flopy.__version__}")
3.12.2 | packaged by conda-forge | (main, Feb 16 2024, 20:50:58) [GCC 12.3.0]
numpy version: 1.26.4
matplotlib version: 3.8.4
flopy version: 3.7.0.dev0
[3]:
temp_dir = TemporaryDirectory()
modelpth = temp_dir.name
modelname = "CrnkNic"
mfexe = "mfnwt"
mtexe = "mt3dusgs"

# Make sure modelpth directory exists
if not os.path.isdir(modelpth):
    os.makedirs(modelpth, exist_ok=True)

# Instantiate MODFLOW object in flopy
mf = flopy.modflow.Modflow(
    modelname=modelname, exe_name=mfexe, model_ws=modelpth, version="mfnwt"
)

Set up model discretization

[4]:
Lx = 650.0
Ly = 15
nrow = 3
ncol = 650
nlay = 1

delr = Lx / ncol
delc = Ly / nrow

xmax = ncol * delr
ymax = nrow * delc

X, Y = np.meshgrid(
    np.linspace(delr / 2, xmax - delr / 2, ncol),
    np.linspace(ymax - delc / 2, 0 + delc / 2, nrow),
)

Instantiate output control (oc) package for MODFLOW-NWT

[5]:
# Output Control: Create a flopy output control object
oc = flopy.modflow.ModflowOc(mf)

Instantiate solver package for MODFLOW-NWT

[6]:
# Newton-Raphson Solver: Create a flopy nwt package object

headtol = 1.0e-4
fluxtol = 5
maxiterout = 5000
thickfact = 1e-06
linmeth = 2
iprnwt = 1
ibotav = 1

nwt = flopy.modflow.ModflowNwt(
    mf,
    headtol=headtol,
    fluxtol=fluxtol,
    maxiterout=maxiterout,
    thickfact=thickfact,
    linmeth=linmeth,
    iprnwt=iprnwt,
    ibotav=ibotav,
    options="SIMPLE",
)

Instantiate discretization (DIS) package for MODFLOW-NWT

[7]:
# The equations for calculating the ground elevation in the 1 Layer CrnkNic model.
# Although Y isn't used, keeping it here for symetry


def topElev(X, Y):
    return 100.0 - (np.ceil(X) - 1) * 0.03


grndElev = topElev(X, Y)
bedRockElev = grndElev - 3.0

Steady = [False, False, False]
nstp = [1, 1, 1]
tsmult = [1.0, 1.0, 1.0]

# Stress periods extend from (12AM-8:29:59AM); (8:30AM-11:30:59AM); (11:31AM-23:59:59PM)
perlen = [30600, 10800, 45000]

# Create the discretization object
# itmuni = 1 (seconds); lenuni = 2 (meters)
dis = flopy.modflow.ModflowDis(
    mf,
    nlay,
    nrow,
    ncol,
    nper=3,
    delr=delr,
    delc=delc,
    top=grndElev,
    botm=bedRockElev,
    laycbd=0,
    itmuni=1,
    lenuni=2,
    steady=Steady,
    nstp=nstp,
    tsmult=tsmult,
    perlen=perlen,
)

Instantiate upstream weighting (UPW) flow package for MODFLOW-NWT

[8]:
# UPW parameters
# UPW must be instantiated after DIS.  Otherwise, during the mf.write_input() procedures,
# flopy will crash.


laytyp = 1
layavg = 2
chani = 1.0
layvka = 1
iphdry = 0
hk = 0.1
hani = 1
vka = 1.0
ss = 0.000001
sy = 0.20
hdry = -888

upw = flopy.modflow.ModflowUpw(
    mf,
    laytyp=laytyp,
    layavg=layavg,
    chani=chani,
    layvka=layvka,
    ipakcb=53,
    hdry=hdry,
    iphdry=iphdry,
    hk=hk,
    hani=hani,
    vka=vka,
    ss=ss,
    sy=sy,
)

Instantiate basic (BAS or BA6) package for MODFLOW-NWT

[9]:
# Create a flopy basic package object
def calc_strtElev(X, Y):
    return 99.5 - (np.ceil(X) - 1) * 0.0001


ibound = np.ones((nlay, nrow, ncol))
ibound[:, 0, :] *= -1
ibound[:, 2, :] *= -1

strtElev = calc_strtElev(X, Y)

bas = flopy.modflow.ModflowBas(mf, ibound=ibound, hnoflo=hdry, strt=strtElev)

Instantiate streamflow routing (SFR2) package for MODFLOW-NWT

[10]:
# Streamflow Routing Package: Try and set up with minimal options in use
# 9 11 IFACE # Data Set 1: ISTCB1  ISTCB2

nstrm = ncol
nss = 6
const = 1.0
dleak = 0.0001
istcb1 = -10
istcb2 = 11
isfropt = 1

segment_data = None
channel_geometry_data = None
channel_flow_data = None
dataset_5 = None
reachinput = True


# The next couple of lines set up the reach_data for the 30x100 hypothetical model.
# Will need to adjust the row based on which grid discretization we're doing.
# Ensure that the stream goes down one of the middle rows of the model.

strmBed_Elev = 98.75 - (np.ceil(X[1, :]) - 1) * 0.0001

s1 = "k,i,j,iseg,ireach,rchlen,strtop,slope,strthick,strhc1\n"
iseg = 0
irch = 0

for y in range(ncol):
    if y <= 37:
        if iseg == 0:
            irch = 1
        else:
            irch += 1
        iseg = 1
        strhc1 = 1.0e-10
    elif y <= 104:
        if iseg == 1:
            irch = 1
        else:
            irch += 1
        iseg = 2
        strhc1 = 1.0e-10
    elif y <= 280:
        if iseg == 2:
            irch = 1
        else:
            irch += 1
        iseg = 3
        strhc1 = 2.946219199e-6
    elif y <= 432:
        if iseg == 3:
            irch = 1
        else:
            irch += 1
        iseg = 4
        strhc1 = 1.375079882e-6
    elif y <= 618:
        if iseg == 4:
            irch = 1
        else:
            irch += 1
        iseg = 5
        strhc1 = 1.764700062e-6
    else:
        if iseg == 5:
            irch = 1
        else:
            irch += 1
        iseg = 6
        strhc1 = 1e-10

    # remember that lay, row, col need to be zero-based and are adjusted accordingly by flopy
    #    layer +  row      +      col     +       iseg      +      irch      +     rchlen      +            strtop          +       slope       +     strthick    +     strmbed K
    s1 += f"0,{1}"
    s1 += f",{y}"
    s1 += f",{iseg}"
    s1 += f",{irch}"
    s1 += f",{delr}"
    s1 += f",{strmBed_Elev[y]}"
    s1 += f",{0.0001}"
    s1 += f",{0.50}"
    s1 += f",{strhc1}\n"


fpth = os.path.join(modelpth, "s1.csv")
f = open(fpth, "w")
f.write(s1)
f.close()

dtype = [
    ("k", "<i4"),
    ("i", "<i4"),
    ("j", "<i4"),
    ("iseg", "<i4"),
    ("ireach", "<f8"),
    ("rchlen", "<f8"),
    ("strtop", "<f8"),
    ("slope", "<f8"),
    ("strthick", "<f8"),
    ("strhc1", "<f8"),
]

f = open(fpth, "rb")

reach_data = np.genfromtxt(f, delimiter=",", names=True, dtype=dtype)
f.close()

s2 = "nseg,icalc,outseg,iupseg,nstrpts,   flow,runoff,etsw,pptsw,        roughch,        roughbk,cdpth,fdpth,awdth,bwdth,width1,width2\n \
      1,    1,     2,     0,      0, 0.0125,   0.0, 0.0,  0.0, 0.082078856000, 0.082078856000,  0.0,  0.0,  0.0,  0.0,   1.5,   1.5\n \
      2,    1,     3,     0,      0,    0.0,   0.0, 0.0,  0.0, 0.143806300000, 0.143806300000,  0.0,  0.0,  0.0,  0.0,   1.5,   1.5\n \
      3,    1,     4,     0,      0,    0.0,   0.0, 0.0,  0.0, 0.104569661821, 0.104569661821,  0.0,  0.0,  0.0,  0.0,   1.5,   1.5\n \
      4,    1,     5,     0,      0,    0.0,   0.0, 0.0,  0.0, 0.126990045841, 0.126990045841,  0.0,  0.0,  0.0,  0.0,   1.5,   1.5\n \
      5,    1,     6,     0,      0,    0.0,   0.0, 0.0,  0.0, 0.183322283828, 0.183322283828,  0.0,  0.0,  0.0,  0.0,   1.5,   1.5\n \
      6,    1,     0,     0,      0,    0.0,   0.0, 0.0,  0.0, 0.183322283828, 0.183322283828,  0.0,  0.0,  0.0,  0.0,   1.5,   1.5"

fpth = os.path.join(modelpth, "s2.csv")
f = open(fpth, "w")
f.write(s2)
f.close()

f = open(fpth, "rb")

segment_data = np.genfromtxt(f, delimiter=",", names=True)
f.close()

# Be sure to convert segment_data to a dictionary keyed on stress period.
segment_data = np.atleast_1d(segment_data)
segment_data = {0: segment_data, 1: segment_data, 2: segment_data}

# There are 3 stress periods
dataset_5 = {0: [nss, 0, 0], 1: [nss, 0, 0], 2: [nss, 0, 0]}

sfr = flopy.modflow.ModflowSfr2(
    mf,
    nstrm=nstrm,
    nss=nss,
    const=const,
    dleak=dleak,
    isfropt=isfropt,
    istcb2=0,
    reachinput=True,
    reach_data=reach_data,
    dataset_5=dataset_5,
    segment_data=segment_data,
    channel_geometry_data=channel_geometry_data,
)

Instantiate gage package for use with MODFLOW-NWT package

[11]:
gages = [
    [1, 38, 61, 1],
    [2, 67, 62, 1],
    [3, 176, 63, 1],
    [4, 152, 64, 1],
    [5, 186, 65, 1],
    [6, 31, 66, 1],
]
files = [
    "CrnkNic.gage",
    "CrnkNic.gag1",
    "CrnkNic.gag2",
    "CrnkNic.gag3",
    "CrnkNic.gag4",
    "CrnkNic.gag5",
    "CrnkNic.gag6",
]
gage = flopy.modflow.ModflowGage(
    mf, numgage=6, gage_data=gages, filenames=files
)

Instantiate linkage with mass transport routing (LMT) package for MODFLOW-NWT (generates linker file)

[12]:
lmt = flopy.modflow.ModflowLmt(
    mf,
    output_file_name="CrnkNic.ftl",
    output_file_header="extended",
    output_file_format="formatted",
    package_flows=["sfr"],
)

Write the MODFLOW input files

[13]:
mf.write_input()

# run the model
success, buff = mf.run_model(silent=True, report=True)
assert success, pformat(buff)

Now draft up MT3D-USGS input files.

[14]:
# Instantiate MT3D-USGS object in flopy
mt = flopy.mt3d.Mt3dms(
    modflowmodel=mf,
    modelname=modelname,
    model_ws=modelpth,
    version="mt3d-usgs",
    namefile_ext="mtnam",
    exe_name=mtexe,
    ftlfilename="CrnkNic.ftl",
    ftlfree=True,
)

Instantiate basic transport (BTN) package for MT3D-USGS

[15]:
btn = flopy.mt3d.Mt3dBtn(
    mt,
    sconc=3.7,
    ncomp=1,
    prsity=0.2,
    cinact=-1.0,
    thkmin=0.001,
    nprs=-1,
    nprobs=10,
    chkmas=True,
    nprmas=10,
    dt0=180,
    mxstrn=2500,
)

Instantiate advection (ADV) package for MT3D-USGS

[16]:
adv = flopy.mt3d.Mt3dAdv(mt, mixelm=0, percel=1.00, mxpart=5000, nadvfd=1)

Instatiate generalized conjugate gradient solver (GCG) package for MT3D-USGS

[17]:
rct = flopy.mt3d.Mt3dRct(mt, isothm=0, ireact=100, igetsc=0, rc1=0.01)
[18]:
gcg = flopy.mt3d.Mt3dGcg(
    mt, mxiter=10, iter1=50, isolve=3, ncrs=0, accl=1, cclose=1e-6, iprgcg=1
)

Instantiate source-sink mixing (SSM) package for MT3D-USGS

[19]:
# For SSM, need to set the constant head boundary conditions to the ambient concentration
# for all 1,300 constant head boundary cells.

itype = flopy.mt3d.Mt3dSsm.itype_dict()
ssm_data = {}
ssm_data[0] = [(0, 0, 0, 3.7, itype["CHD"])]
ssm_data[0].append((0, 2, 0, 3.7, itype["CHD"]))
for i in [0, 2]:
    for j in range(1, ncol):
        ssm_data[0].append((0, i, j, 3.7, itype["CHD"]))

ssm = flopy.mt3d.Mt3dSsm(mt, stress_period_data=ssm_data)

Instantiate streamflow transport (SFT) package for MT3D-USGS

[20]:
dispsf = []
for y in range(ncol):
    if y <= 37:
        dispsf.append(0.12)
    elif y <= 104:
        dispsf.append(0.15)
    elif y <= 280:
        dispsf.append(0.24)
    elif y <= 432:
        dispsf.append(0.31)
    elif y <= 618:
        dispsf.append(0.40)
    else:
        dispsf.append(0.40)

# Enter a list of the observation points
# Each observation is taken as the last reach within the first 5 segments

seg_len = np.unique(reach_data["iseg"], return_counts=True)
obs_sf = np.cumsum(seg_len[1])
obs_sf = obs_sf.tolist()

# The last reach is not an observation point, therefore drop
obs_sf.pop(-1)

# In the first and last stress periods, concentration at the headwater is 3.7
sf_stress_period_data = {0: [0, 0, 3.7], 1: [0, 0, 11.4], 2: [0, 0, 3.7]}

gage_output = [None, None, "CrnkNic.sftobs"]

sft = flopy.mt3d.Mt3dSft(
    mt,
    nsfinit=650,
    mxsfbc=650,
    icbcsf=81,
    ioutobs=82,
    isfsolv=1,
    cclosesf=1.0e-6,
    mxitersf=10,
    crntsf=1.0,
    iprtxmd=0,
    coldsf=3.7,
    dispsf=dispsf,
    nobssf=5,
    obs_sf=obs_sf,
    sf_stress_period_data=sf_stress_period_data,
    filenames=gage_output,
)

sft.dispsf[0].format.fortran = "(10E15.6)"

Write the MT3D-USGS input files and run the model.

[21]:
mt.write_input()
mt.run_model(silent=True, report=True)
[21]:
(False,
 ['',
  ' MT3D-USGS - Modular 3D Multi-Species Transport Model [Ver 1.1.0]   ',
  ' and based on MT3DMS. MT3D-USGS developed in cooperation by ',
  ' S.S. Papadopulos & Associates and the U.S. Geological Survey',
  '',
  ' Using NAME File: CrnkNic.mtnam',
  '',
  ' STRESS PERIOD NO.    1',
  '',
  ' TIME STEP NO.    1',
  ' FROM TIME =   0.0000     TO    30600.    ',
  '',
  ' Transport Step:    1   Step Size:   180.0     Total Elapsed Time:   180.00    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =  0.4865      [K,I,J]    1    2  289',
  ' Outer Iter.  1  Inner Iter.  2:  Max. DC =  0.5960E-07  [K,I,J]    1    1   36',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    2   Step Size:   180.0     Total Elapsed Time:   360.00    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =  0.9385      [K,I,J]    1    2  394',
  ' Outer Iter.  1  Inner Iter.  2:  Max. DC =  0.1788E-06  [K,I,J]    1    2  294',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =  0.3934E-05  [K,I,J]    1    2  294',
  ' Outer Iter.  2  Inner Iter.  2:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  3  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    3   Step Size:   180.0     Total Elapsed Time:   540.00    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =  0.7602      [K,I,J]    1    1    1',
  ' Outer Iter.  1  Inner Iter.  2:  Max. DC =  0.1788E-06  [K,I,J]    1    2    2',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =  0.8884E-06  [K,I,J]    1    1  642',
  ' Transport Step:    4   Step Size:   180.0     Total Elapsed Time:   720.00    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =  0.3366E-03  [K,I,J]    1    2    1',
  ' Outer Iter.  1  Inner Iter.  2:  Max. DC =  0.4657E-09  [K,I,J]    1    2    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    5   Step Size:   180.0     Total Elapsed Time:   900.00    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    6   Step Size:   180.0     Total Elapsed Time:   1080.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    7   Step Size:   180.0     Total Elapsed Time:   1260.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    8   Step Size:   180.0     Total Elapsed Time:   1440.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    9   Step Size:   180.0     Total Elapsed Time:   1620.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   10   Step Size:   180.0     Total Elapsed Time:   1800.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   11   Step Size:   180.0     Total Elapsed Time:   1980.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   12   Step Size:   180.0     Total Elapsed Time:   2160.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   13   Step Size:   180.0     Total Elapsed Time:   2340.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   14   Step Size:   180.0     Total Elapsed Time:   2520.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   15   Step Size:   180.0     Total Elapsed Time:   2700.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   16   Step Size:   180.0     Total Elapsed Time:   2880.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   17   Step Size:   180.0     Total Elapsed Time:   3060.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   18   Step Size:   180.0     Total Elapsed Time:   3240.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   19   Step Size:   180.0     Total Elapsed Time:   3420.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   20   Step Size:   180.0     Total Elapsed Time:   3600.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   21   Step Size:   180.0     Total Elapsed Time:   3780.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   22   Step Size:   180.0     Total Elapsed Time:   3960.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   23   Step Size:   180.0     Total Elapsed Time:   4140.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   24   Step Size:   180.0     Total Elapsed Time:   4320.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   25   Step Size:   180.0     Total Elapsed Time:   4500.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   26   Step Size:   180.0     Total Elapsed Time:   4680.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   27   Step Size:   180.0     Total Elapsed Time:   4860.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   28   Step Size:   180.0     Total Elapsed Time:   5040.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   29   Step Size:   180.0     Total Elapsed Time:   5220.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   30   Step Size:   180.0     Total Elapsed Time:   5400.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   31   Step Size:   180.0     Total Elapsed Time:   5580.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   32   Step Size:   180.0     Total Elapsed Time:   5760.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   33   Step Size:   180.0     Total Elapsed Time:   5940.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   34   Step Size:   180.0     Total Elapsed Time:   6120.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   35   Step Size:   180.0     Total Elapsed Time:   6300.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   36   Step Size:   180.0     Total Elapsed Time:   6480.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   37   Step Size:   180.0     Total Elapsed Time:   6660.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   38   Step Size:   180.0     Total Elapsed Time:   6840.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   39   Step Size:   180.0     Total Elapsed Time:   7020.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   40   Step Size:   180.0     Total Elapsed Time:   7200.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   41   Step Size:   180.0     Total Elapsed Time:   7380.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   42   Step Size:   180.0     Total Elapsed Time:   7560.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   43   Step Size:   180.0     Total Elapsed Time:   7740.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   44   Step Size:   180.0     Total Elapsed Time:   7920.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   45   Step Size:   180.0     Total Elapsed Time:   8100.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   46   Step Size:   180.0     Total Elapsed Time:   8280.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   47   Step Size:   180.0     Total Elapsed Time:   8460.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   48   Step Size:   180.0     Total Elapsed Time:   8640.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   49   Step Size:   180.0     Total Elapsed Time:   8820.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   50   Step Size:   180.0     Total Elapsed Time:   9000.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   51   Step Size:   180.0     Total Elapsed Time:   9180.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   52   Step Size:   180.0     Total Elapsed Time:   9360.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   53   Step Size:   180.0     Total Elapsed Time:   9540.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   54   Step Size:   180.0     Total Elapsed Time:   9720.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   55   Step Size:   180.0     Total Elapsed Time:   9900.0    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   56   Step Size:   180.0     Total Elapsed Time:   10080.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   57   Step Size:   180.0     Total Elapsed Time:   10260.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   58   Step Size:   180.0     Total Elapsed Time:   10440.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   59   Step Size:   180.0     Total Elapsed Time:   10620.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   60   Step Size:   180.0     Total Elapsed Time:   10800.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   61   Step Size:   180.0     Total Elapsed Time:   10980.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   62   Step Size:   180.0     Total Elapsed Time:   11160.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   63   Step Size:   180.0     Total Elapsed Time:   11340.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   64   Step Size:   180.0     Total Elapsed Time:   11520.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   65   Step Size:   180.0     Total Elapsed Time:   11700.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   66   Step Size:   180.0     Total Elapsed Time:   11880.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   67   Step Size:   180.0     Total Elapsed Time:   12060.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   68   Step Size:   180.0     Total Elapsed Time:   12240.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   69   Step Size:   180.0     Total Elapsed Time:   12420.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   70   Step Size:   180.0     Total Elapsed Time:   12600.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   71   Step Size:   180.0     Total Elapsed Time:   12780.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   72   Step Size:   180.0     Total Elapsed Time:   12960.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   73   Step Size:   180.0     Total Elapsed Time:   13140.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   74   Step Size:   180.0     Total Elapsed Time:   13320.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   75   Step Size:   180.0     Total Elapsed Time:   13500.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   76   Step Size:   180.0     Total Elapsed Time:   13680.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   77   Step Size:   180.0     Total Elapsed Time:   13860.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   78   Step Size:   180.0     Total Elapsed Time:   14040.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   79   Step Size:   180.0     Total Elapsed Time:   14220.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   80   Step Size:   180.0     Total Elapsed Time:   14400.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   81   Step Size:   180.0     Total Elapsed Time:   14580.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   82   Step Size:   180.0     Total Elapsed Time:   14760.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   83   Step Size:   180.0     Total Elapsed Time:   14940.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   84   Step Size:   180.0     Total Elapsed Time:   15120.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   85   Step Size:   180.0     Total Elapsed Time:   15300.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   86   Step Size:   180.0     Total Elapsed Time:   15480.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   87   Step Size:   180.0     Total Elapsed Time:   15660.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   88   Step Size:   180.0     Total Elapsed Time:   15840.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   89   Step Size:   180.0     Total Elapsed Time:   16020.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   90   Step Size:   180.0     Total Elapsed Time:   16200.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   91   Step Size:   180.0     Total Elapsed Time:   16380.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   92   Step Size:   180.0     Total Elapsed Time:   16560.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   93   Step Size:   180.0     Total Elapsed Time:   16740.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   94   Step Size:   180.0     Total Elapsed Time:   16920.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   95   Step Size:   180.0     Total Elapsed Time:   17100.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   96   Step Size:   180.0     Total Elapsed Time:   17280.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   97   Step Size:   180.0     Total Elapsed Time:   17460.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   98   Step Size:   180.0     Total Elapsed Time:   17640.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   99   Step Size:   180.0     Total Elapsed Time:   17820.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  100   Step Size:   180.0     Total Elapsed Time:   18000.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  101   Step Size:   180.0     Total Elapsed Time:   18180.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  102   Step Size:   180.0     Total Elapsed Time:   18360.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  103   Step Size:   180.0     Total Elapsed Time:   18540.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  104   Step Size:   180.0     Total Elapsed Time:   18720.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  105   Step Size:   180.0     Total Elapsed Time:   18900.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  106   Step Size:   180.0     Total Elapsed Time:   19080.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  107   Step Size:   180.0     Total Elapsed Time:   19260.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  108   Step Size:   180.0     Total Elapsed Time:   19440.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  109   Step Size:   180.0     Total Elapsed Time:   19620.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  110   Step Size:   180.0     Total Elapsed Time:   19800.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  111   Step Size:   180.0     Total Elapsed Time:   19980.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  112   Step Size:   180.0     Total Elapsed Time:   20160.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  113   Step Size:   180.0     Total Elapsed Time:   20340.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  114   Step Size:   180.0     Total Elapsed Time:   20520.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  115   Step Size:   180.0     Total Elapsed Time:   20700.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  116   Step Size:   180.0     Total Elapsed Time:   20880.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  117   Step Size:   180.0     Total Elapsed Time:   21060.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  118   Step Size:   180.0     Total Elapsed Time:   21240.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  119   Step Size:   180.0     Total Elapsed Time:   21420.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  120   Step Size:   180.0     Total Elapsed Time:   21600.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  121   Step Size:   180.0     Total Elapsed Time:   21780.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  122   Step Size:   180.0     Total Elapsed Time:   21960.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  123   Step Size:   180.0     Total Elapsed Time:   22140.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  124   Step Size:   180.0     Total Elapsed Time:   22320.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  125   Step Size:   180.0     Total Elapsed Time:   22500.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  126   Step Size:   180.0     Total Elapsed Time:   22680.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  127   Step Size:   180.0     Total Elapsed Time:   22860.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  128   Step Size:   180.0     Total Elapsed Time:   23040.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  129   Step Size:   180.0     Total Elapsed Time:   23220.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  130   Step Size:   180.0     Total Elapsed Time:   23400.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  131   Step Size:   180.0     Total Elapsed Time:   23580.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  132   Step Size:   180.0     Total Elapsed Time:   23760.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  133   Step Size:   180.0     Total Elapsed Time:   23940.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  134   Step Size:   180.0     Total Elapsed Time:   24120.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  135   Step Size:   180.0     Total Elapsed Time:   24300.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  136   Step Size:   180.0     Total Elapsed Time:   24480.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  137   Step Size:   180.0     Total Elapsed Time:   24660.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  138   Step Size:   180.0     Total Elapsed Time:   24840.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  139   Step Size:   180.0     Total Elapsed Time:   25020.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  140   Step Size:   180.0     Total Elapsed Time:   25200.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  141   Step Size:   180.0     Total Elapsed Time:   25380.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  142   Step Size:   180.0     Total Elapsed Time:   25560.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  143   Step Size:   180.0     Total Elapsed Time:   25740.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  144   Step Size:   180.0     Total Elapsed Time:   25920.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  145   Step Size:   180.0     Total Elapsed Time:   26100.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  146   Step Size:   180.0     Total Elapsed Time:   26280.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  147   Step Size:   180.0     Total Elapsed Time:   26460.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  148   Step Size:   180.0     Total Elapsed Time:   26640.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  149   Step Size:   180.0     Total Elapsed Time:   26820.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  150   Step Size:   180.0     Total Elapsed Time:   27000.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  151   Step Size:   180.0     Total Elapsed Time:   27180.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  152   Step Size:   180.0     Total Elapsed Time:   27360.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  153   Step Size:   180.0     Total Elapsed Time:   27540.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  154   Step Size:   180.0     Total Elapsed Time:   27720.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  155   Step Size:   180.0     Total Elapsed Time:   27900.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  156   Step Size:   180.0     Total Elapsed Time:   28080.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  157   Step Size:   180.0     Total Elapsed Time:   28260.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  158   Step Size:   180.0     Total Elapsed Time:   28440.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  159   Step Size:   180.0     Total Elapsed Time:   28620.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  160   Step Size:   180.0     Total Elapsed Time:   28800.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  161   Step Size:   180.0     Total Elapsed Time:   28980.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  162   Step Size:   180.0     Total Elapsed Time:   29160.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  163   Step Size:   180.0     Total Elapsed Time:   29340.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  164   Step Size:   180.0     Total Elapsed Time:   29520.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  165   Step Size:   180.0     Total Elapsed Time:   29700.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  166   Step Size:   180.0     Total Elapsed Time:   29880.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  167   Step Size:   180.0     Total Elapsed Time:   30060.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  168   Step Size:   180.0     Total Elapsed Time:   30240.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  169   Step Size:   180.0     Total Elapsed Time:   30420.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:  170   Step Size:   180.0     Total Elapsed Time:   30600.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  '',
  ' STRESS PERIOD NO.    2',
  '',
  ' TIME STEP NO.    1',
  ' FROM TIME =   30600.     TO    41400.    ',
  '',
  ' Transport Step:    1   Step Size:   180.0     Total Elapsed Time:   30780.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =  0.2062E-04  [K,I,J]    1    1    1',
  ' Outer Iter.  1  Inner Iter.  2:  Max. DC =  0.6548E-10  [K,I,J]    1    2    3',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    2   Step Size:   180.0     Total Elapsed Time:   30960.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    3   Step Size:   180.0     Total Elapsed Time:   31140.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    4   Step Size:   180.0     Total Elapsed Time:   31320.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    5   Step Size:   180.0     Total Elapsed Time:   31500.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    6   Step Size:   180.0     Total Elapsed Time:   31680.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    7   Step Size:   180.0     Total Elapsed Time:   31860.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    8   Step Size:   180.0     Total Elapsed Time:   32040.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    9   Step Size:   180.0     Total Elapsed Time:   32220.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   10   Step Size:   180.0     Total Elapsed Time:   32400.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   11   Step Size:   180.0     Total Elapsed Time:   32580.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   12   Step Size:   180.0     Total Elapsed Time:   32760.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   13   Step Size:   180.0     Total Elapsed Time:   32940.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   14   Step Size:   180.0     Total Elapsed Time:   33120.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   15   Step Size:   180.0     Total Elapsed Time:   33300.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   16   Step Size:   180.0     Total Elapsed Time:   33480.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   17   Step Size:   180.0     Total Elapsed Time:   33660.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   18   Step Size:   180.0     Total Elapsed Time:   33840.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   19   Step Size:   180.0     Total Elapsed Time:   34020.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   20   Step Size:   180.0     Total Elapsed Time:   34200.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   21   Step Size:   180.0     Total Elapsed Time:   34380.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   22   Step Size:   180.0     Total Elapsed Time:   34560.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   23   Step Size:   180.0     Total Elapsed Time:   34740.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   24   Step Size:   180.0     Total Elapsed Time:   34920.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   25   Step Size:   180.0     Total Elapsed Time:   35100.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   26   Step Size:   180.0     Total Elapsed Time:   35280.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   27   Step Size:   180.0     Total Elapsed Time:   35460.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   28   Step Size:   180.0     Total Elapsed Time:   35640.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   29   Step Size:   180.0     Total Elapsed Time:   35820.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   30   Step Size:   180.0     Total Elapsed Time:   36000.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   31   Step Size:   180.0     Total Elapsed Time:   36180.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   32   Step Size:   180.0     Total Elapsed Time:   36360.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   33   Step Size:   180.0     Total Elapsed Time:   36540.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   34   Step Size:   180.0     Total Elapsed Time:   36720.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   35   Step Size:   180.0     Total Elapsed Time:   36900.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   36   Step Size:   180.0     Total Elapsed Time:   37080.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   37   Step Size:   180.0     Total Elapsed Time:   37260.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   38   Step Size:   180.0     Total Elapsed Time:   37440.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   39   Step Size:   180.0     Total Elapsed Time:   37620.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   40   Step Size:   180.0     Total Elapsed Time:   37800.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   41   Step Size:   180.0     Total Elapsed Time:   37980.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   42   Step Size:   180.0     Total Elapsed Time:   38160.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   43   Step Size:   180.0     Total Elapsed Time:   38340.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   44   Step Size:   180.0     Total Elapsed Time:   38520.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   45   Step Size:   180.0     Total Elapsed Time:   38700.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   46   Step Size:   180.0     Total Elapsed Time:   38880.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   47   Step Size:   180.0     Total Elapsed Time:   39060.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   48   Step Size:   180.0     Total Elapsed Time:   39240.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   49   Step Size:   180.0     Total Elapsed Time:   39420.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   50   Step Size:   180.0     Total Elapsed Time:   39600.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   51   Step Size:   180.0     Total Elapsed Time:   39780.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   52   Step Size:   180.0     Total Elapsed Time:   39960.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   53   Step Size:   180.0     Total Elapsed Time:   40140.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   54   Step Size:   180.0     Total Elapsed Time:   40320.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   55   Step Size:   180.0     Total Elapsed Time:   40500.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   56   Step Size:   180.0     Total Elapsed Time:   40680.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   57   Step Size:   180.0     Total Elapsed Time:   40860.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   58   Step Size:   180.0     Total Elapsed Time:   41040.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   59   Step Size:   180.0     Total Elapsed Time:   41220.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   60   Step Size:   180.0     Total Elapsed Time:   41400.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  '',
  ' STRESS PERIOD NO.    3',
  '',
  ' TIME STEP NO.    1',
  ' FROM TIME =   41400.     TO    86400.    ',
  '',
  ' Transport Step:    1   Step Size:   180.0     Total Elapsed Time:   41580.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    2   Step Size:   180.0     Total Elapsed Time:   41760.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    3   Step Size:   180.0     Total Elapsed Time:   41940.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    4   Step Size:   180.0     Total Elapsed Time:   42120.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    5   Step Size:   180.0     Total Elapsed Time:   42300.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    6   Step Size:   180.0     Total Elapsed Time:   42480.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    7   Step Size:   180.0     Total Elapsed Time:   42660.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    8   Step Size:   180.0     Total Elapsed Time:   42840.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:    9   Step Size:   180.0     Total Elapsed Time:   43020.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   10   Step Size:   180.0     Total Elapsed Time:   43200.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   11   Step Size:   180.0     Total Elapsed Time:   43380.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   12   Step Size:   180.0     Total Elapsed Time:   43560.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   13   Step Size:   180.0     Total Elapsed Time:   43740.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   14   Step Size:   180.0     Total Elapsed Time:   43920.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   15   Step Size:   180.0     Total Elapsed Time:   44100.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   16   Step Size:   180.0     Total Elapsed Time:   44280.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   17   Step Size:   180.0     Total Elapsed Time:   44460.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   18   Step Size:   180.0     Total Elapsed Time:   44640.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   19   Step Size:   180.0     Total Elapsed Time:   44820.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   20   Step Size:   180.0     Total Elapsed Time:   45000.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   21   Step Size:   180.0     Total Elapsed Time:   45180.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   22   Step Size:   180.0     Total Elapsed Time:   45360.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   23   Step Size:   180.0     Total Elapsed Time:   45540.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   24   Step Size:   180.0     Total Elapsed Time:   45720.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   25   Step Size:   180.0     Total Elapsed Time:   45900.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   26   Step Size:   180.0     Total Elapsed Time:   46080.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   27   Step Size:   180.0     Total Elapsed Time:   46260.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   28   Step Size:   180.0     Total Elapsed Time:   46440.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   29   Step Size:   180.0     Total Elapsed Time:   46620.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   30   Step Size:   180.0     Total Elapsed Time:   46800.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   31   Step Size:   180.0     Total Elapsed Time:   46980.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   32   Step Size:   180.0     Total Elapsed Time:   47160.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   33   Step Size:   180.0     Total Elapsed Time:   47340.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   34   Step Size:   180.0     Total Elapsed Time:   47520.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   35   Step Size:   180.0     Total Elapsed Time:   47700.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   36   Step Size:   180.0     Total Elapsed Time:   47880.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   37   Step Size:   180.0     Total Elapsed Time:   48060.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   38   Step Size:   180.0     Total Elapsed Time:   48240.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   39   Step Size:   180.0     Total Elapsed Time:   48420.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   40   Step Size:   180.0     Total Elapsed Time:   48600.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   41   Step Size:   180.0     Total Elapsed Time:   48780.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   42   Step Size:   180.0     Total Elapsed Time:   48960.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   43   Step Size:   180.0     Total Elapsed Time:   49140.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   44   Step Size:   180.0     Total Elapsed Time:   49320.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   45   Step Size:   180.0     Total Elapsed Time:   49500.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   46   Step Size:   180.0     Total Elapsed Time:   49680.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   47   Step Size:   180.0     Total Elapsed Time:   49860.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   48   Step Size:   180.0     Total Elapsed Time:   50040.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   49   Step Size:   180.0     Total Elapsed Time:   50220.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   50   Step Size:   180.0     Total Elapsed Time:   50400.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   51   Step Size:   180.0     Total Elapsed Time:   50580.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   52   Step Size:   180.0     Total Elapsed Time:   50760.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   53   Step Size:   180.0     Total Elapsed Time:   50940.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   54   Step Size:   180.0     Total Elapsed Time:   51120.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   55   Step Size:   180.0     Total Elapsed Time:   51300.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   56   Step Size:   180.0     Total Elapsed Time:   51480.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   57   Step Size:   180.0     Total Elapsed Time:   51660.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   58   Step Size:   180.0     Total Elapsed Time:   51840.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   59   Step Size:   180.0     Total Elapsed Time:   52020.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   60   Step Size:   180.0     Total Elapsed Time:   52200.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   61   Step Size:   180.0     Total Elapsed Time:   52380.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   62   Step Size:   180.0     Total Elapsed Time:   52560.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   63   Step Size:   180.0     Total Elapsed Time:   52740.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   64   Step Size:   180.0     Total Elapsed Time:   52920.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   65   Step Size:   180.0     Total Elapsed Time:   53100.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   66   Step Size:   180.0     Total Elapsed Time:   53280.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   67   Step Size:   180.0     Total Elapsed Time:   53460.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   68   Step Size:   180.0     Total Elapsed Time:   53640.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   69   Step Size:   180.0     Total Elapsed Time:   53820.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   70   Step Size:   180.0     Total Elapsed Time:   54000.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   71   Step Size:   180.0     Total Elapsed Time:   54180.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   72   Step Size:   180.0     Total Elapsed Time:   54360.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   73   Step Size:   180.0     Total Elapsed Time:   54540.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   74   Step Size:   180.0     Total Elapsed Time:   54720.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   75   Step Size:   180.0     Total Elapsed Time:   54900.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   76   Step Size:   180.0     Total Elapsed Time:   55080.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   77   Step Size:   180.0     Total Elapsed Time:   55260.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   78   Step Size:   180.0     Total Elapsed Time:   55440.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   79   Step Size:   180.0     Total Elapsed Time:   55620.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   80   Step Size:   180.0     Total Elapsed Time:   55800.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   81   Step Size:   180.0     Total Elapsed Time:   55980.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   82   Step Size:   180.0     Total Elapsed Time:   56160.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   83   Step Size:   180.0     Total Elapsed Time:   56340.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   84   Step Size:   180.0     Total Elapsed Time:   56520.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   85   Step Size:   180.0     Total Elapsed Time:   56700.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   86   Step Size:   180.0     Total Elapsed Time:   56880.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   87   Step Size:   180.0     Total Elapsed Time:   57060.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   88   Step Size:   180.0     Total Elapsed Time:   57240.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   89   Step Size:   180.0     Total Elapsed Time:   57420.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   90   Step Size:   180.0     Total Elapsed Time:   57600.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   91   Step Size:   180.0     Total Elapsed Time:   57780.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   92   Step Size:   180.0     Total Elapsed Time:   57960.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Transport Step:   93   Step Size:   180.0     Total Elapsed Time:   58140.    ',
  ' Outer Iter.  1  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ' Outer Iter.  2  Inner Iter.  1:  Max. DC =   0.000      [K,I,J]    1    1    1',
  ...])

Compare mt3d-usgs results to an analytical solution

[22]:
# Define a function to read SFT output file
def load_ts_from_SFT_output(fname, nd=1):
    f = open(fname)
    iline = 0
    lst = []
    for line in f:
        if line.strip().split()[0].replace(".", "", 1).isdigit():
            l = line.strip().split()
            t = float(l[0])
            loc = int(l[1])
            conc = float(l[2])
            if loc == nd:
                lst.append([t, conc])

    ts = np.array(lst)
    f.close()
    return ts


# Also define a function to read OTIS output file
def load_ts_from_otis(fname, iobs=1):
    f = open(fname)
    iline = 0
    lst = []
    for line in f:
        l = line.strip().split()
        t = float(l[0])
        val = float(l[iobs])
        lst.append([t, val])

    ts = np.array(lst)
    f.close()
    return ts

Load output from SFT as well as from the OTIS solution

[23]:
# Model output
fname_SFTout = os.path.join(modelpth, "CrnkNic.sftcobs.out")

# Loading MT3D-USGS output
ts1_mt3d = load_ts_from_SFT_output(fname_SFTout, nd=38)
ts2_mt3d = load_ts_from_SFT_output(fname_SFTout, nd=105)
ts3_mt3d = load_ts_from_SFT_output(fname_SFTout, nd=281)
ts4_mt3d = load_ts_from_SFT_output(fname_SFTout, nd=433)
ts5_mt3d = load_ts_from_SFT_output(fname_SFTout, nd=619)

# OTIS results located here
fname_OTIS = "../../examples/data/mt3d_test/mfnwt_mt3dusgs/sft_crnkNic/OTIS_solution.out"

# Loading OTIS output
ts1_Otis = load_ts_from_otis(fname_OTIS, 1)
ts2_Otis = load_ts_from_otis(fname_OTIS, 2)
ts3_Otis = load_ts_from_otis(fname_OTIS, 3)
ts4_Otis = load_ts_from_otis(fname_OTIS, 4)
ts5_Otis = load_ts_from_otis(fname_OTIS, 5)

Set up some plotting functions

[24]:
def set_plot_params():
    import matplotlib as mpl
    from matplotlib.font_manager import FontProperties

    mpl.rcParams["font.sans-serif"] = "Arial"
    mpl.rcParams["font.serif"] = "Times"
    mpl.rcParams["font.cursive"] = "Zapf Chancery"
    mpl.rcParams["font.fantasy"] = "Comic Sans MS"
    mpl.rcParams["font.monospace"] = "Courier New"
    mpl.rcParams["pdf.compression"] = 0
    mpl.rcParams["pdf.fonttype"] = 42

    ticksize = 10
    mpl.rcParams["legend.fontsize"] = 7
    mpl.rcParams["axes.labelsize"] = 12
    mpl.rcParams["xtick.labelsize"] = ticksize
    mpl.rcParams["ytick.labelsize"] = ticksize
    return


def set_sizexaxis(a, fmt, sz):
    success = 0
    x = a.get_xticks()
    #   print x
    xc = np.chararray(len(x), itemsize=16)
    for i in range(0, len(x)):
        text = fmt % (x[i])
        xc[i] = string.strip(string.ljust(text, 16))
    #   print xc
    a.set_xticklabels(xc, size=sz)
    success = 1
    return success


def set_sizeyaxis(a, fmt, sz):
    success = 0
    y = a.get_yticks()
    #   print y
    yc = np.chararray(len(y), itemsize=16)
    for i in range(0, len(y)):
        text = fmt % (y[i])
        yc[i] = string.strip(string.ljust(text, 16))
    #   print yc
    a.set_yticklabels(yc, size=sz)
    success = 1
    return success

Compare output:

[25]:
# set up figure
try:
    plt.close("all")
except:
    pass

set_plot_params()

fig = plt.figure(figsize=(6, 4), facecolor="w")
ax = fig.add_subplot(1, 1, 1)

ax.plot(ts1_Otis[:, 0], ts1_Otis[:, 1], "k-", linewidth=1.0)
ax.plot(ts2_Otis[:, 0], ts2_Otis[:, 1], "b-", linewidth=1.0)
ax.plot(ts3_Otis[:, 0], ts3_Otis[:, 1], "r-", linewidth=1.0)
ax.plot(ts4_Otis[:, 0], ts4_Otis[:, 1], "g-", linewidth=1.0)
ax.plot(ts5_Otis[:, 0], ts5_Otis[:, 1], "c-", linewidth=1.0)

ax.plot(
    (ts1_mt3d[:, 0]) / 3600,
    ts1_mt3d[:, 1],
    "kD",
    markersize=2.0,
    mfc="none",
    mec="k",
)
ax.plot(
    (ts2_mt3d[:, 0]) / 3600,
    ts2_mt3d[:, 1],
    "b*",
    markersize=3.0,
    mfc="none",
    mec="b",
)
ax.plot((ts3_mt3d[:, 0]) / 3600, ts3_mt3d[:, 1], "r+", markersize=3.0)
ax.plot(
    (ts4_mt3d[:, 0]) / 3600,
    ts4_mt3d[:, 1],
    "g^",
    markersize=2.0,
    mfc="none",
    mec="g",
)
ax.plot(
    (ts5_mt3d[:, 0]) / 3600,
    ts5_mt3d[:, 1],
    "co",
    markersize=2.0,
    mfc="none",
    mec="c",
)

# customize plot
ax.set_xlabel("Time, hours")
ax.set_ylabel("Concentration, mg L-1")
ax.set_ylim([3.5, 13])
ticksize = 10

# legend
leg = ax.legend(
    (
        "Otis, Site 1",
        "Otis, Site 2",
        "Otis, Site 3",
        "Otis, Site 4",
        "Otis, Site 5",
        "MT3D-USGS, Site 1",
        "MT3D-USGS, Site 2",
        "MT3D-USGS, Site 3",
        "MT3D-USGS, Site 4",
        "MT3D-USGS, Site 5",
    ),
    loc="upper right",
    labelspacing=0.25,
    columnspacing=1,
    handletextpad=0.5,
    handlelength=2.0,
    numpoints=1,
)
leg._drawFrame = False

plt.show()
../_images/Notebooks_mt3d-usgs_example_46_0.png
[26]:
try:
    # ignore PermissionError on Windows
    temp_dir.cleanup()
except:
    pass