pydfnWorks
python wrapper for dfnWorks
|
Functions | |
def | lagrit2pflotran (self, inp_file='', mesh_type='', hex2tet=False) |
def | zone2ex (self, uge_file='', zone_file='', face='', boundary_cell_area=1.e-1) |
def | write_perms_and_correct_volumes_areas (self) |
def | pflotran (self, transient=False, restart=False, restart_file='') |
def | check_pflotran_convergence (pflotran_input_file) |
def | pflotran_cleanup (self, index_start=0, index_finish=1, filename='') |
def | parse_pflotran_vtk_python (self, grid_vtk_file='') |
Variables | |
inp_file | |
mesh_type | |
uge_file | |
h | |
perm_cell_file | |
aper_cell_file | |
vtk_file | |
functions for using pflotran in dfnworks
def pydfnworks.dfnFlow.pflotran.check_pflotran_convergence | ( | pflotran_input_file | ) |
checks pflotran_input_file.out for convergence Parameters ---------- pflotran_input_file : string pflotran_input_file Returns ---------- bool True if solver converged / False if not Notes ----------
Definition at line 501 of file pflotran.py.
Referenced by pydfnworks.dfnFlow.pflotran.pflotran().
def pydfnworks.dfnFlow.pflotran.lagrit2pflotran | ( | self, | |
inp_file = '' , |
|||
mesh_type = '' , |
|||
hex2tet = False |
|||
) |
Takes output from LaGriT and processes it for use in PFLOTRAN. Calls the functuon write_perms_and_correct_volumes_areas() and zone2ex Parameters -------------- self : object DFN Class inp_file : str Name of the inp (AVS) file produced by LaGriT mesh_type : str The type of mesh hex2tet : bool True if hex mesh elements should be converted to tet elements, False otherwise. Returns -------- None Notes -------- None
Definition at line 14 of file pflotran.py.
def pydfnworks.dfnFlow.pflotran.parse_pflotran_vtk_python | ( | self, | |
grid_vtk_file = '' |
|||
) |
Adds CELL_DATA to POINT_DATA in the VTK output from PFLOTRAN. Parameters ---------- self : object DFN Class grid_vtk_file : string Name of vtk file with mesh. Typically local_dfnFlow_file.vtk Returns -------- None Notes -------- If DFN class does not have a vtk file, inp2vtk_python is called
Definition at line 602 of file pflotran.py.
def pydfnworks.dfnFlow.pflotran.pflotran | ( | self, | |
transient = False , |
|||
restart = False , |
|||
restart_file = '' |
|||
) |
Run PFLOTRAN. Copy PFLOTRAN run file into working directory and run with ncpus Parameters ---------- self : object DFN Class transient : bool Boolean if PFLOTRAN is running in transient mode restart : bool Boolean if PFLOTRAN is restarting from checkpoint restart_file : string Filename of restart file Returns ---------- None Notes ---------- Runs PFLOTRAN Executable, see http://www.pflotran.org/ for details on PFLOTRAN input cards
Definition at line 432 of file pflotran.py.
References pydfnworks.dfnFlow.pflotran.check_pflotran_convergence().
def pydfnworks.dfnFlow.pflotran.pflotran_cleanup | ( | self, | |
index_start = 0 , |
|||
index_finish = 1 , |
|||
filename = '' |
|||
) |
pflotran_cleanup Concatenate PFLOTRAN output files and then delete them Parameters ----------- self : object DFN Class index : int If PFLOTRAN has multiple dumps use this to pick which dump is put into cellinfo.dat and darcyvel.dat Returns ---------- None Notes ---------- Can be run in a loop over all pflotran dumps
Definition at line 541 of file pflotran.py.
def pydfnworks.dfnFlow.pflotran.write_perms_and_correct_volumes_areas | ( | self | ) |
Write permeability values to perm_file, write aperture values to aper_file, and correct volume areas in uge_file Parameters ---------- self : object DFN Class Returns --------- None Notes ---------- Calls executable correct_uge
Definition at line 276 of file pflotran.py.
def pydfnworks.dfnFlow.pflotran.zone2ex | ( | self, | |
uge_file = '' , |
|||
zone_file = '' , |
|||
face = '' , |
|||
boundary_cell_area = 1.e-1 |
|||
) |
Convert zone files from LaGriT into ex format for LaGriT Parameters ----------- self : object DFN Class uge_file : string Name of uge file zone_file : string Name of zone file Face : Face of the plane corresponding to the zone file zone_file : string Name of zone file to work on. Can be 'all' processes all directions, top, bottom, left, right, front, back boundary_cell_area : double should be a large value relative to the mesh size to force pressure boundary conditions. Returns ---------- None Notes ---------- the boundary_cell_area should be a function of h, the mesh resolution
Definition at line 96 of file pflotran.py.
pydfnworks.dfnFlow.pflotran.aper_cell_file |
Definition at line 319 of file pflotran.py.
pydfnworks.dfnFlow.pflotran.h |
Definition at line 220 of file pflotran.py.
pydfnworks.dfnFlow.pflotran.inp_file |
Definition at line 47 of file pflotran.py.
pydfnworks.dfnFlow.pflotran.mesh_type |
Definition at line 58 of file pflotran.py.
pydfnworks.dfnFlow.pflotran.perm_cell_file |
Definition at line 312 of file pflotran.py.
pydfnworks.dfnFlow.pflotran.uge_file |
Definition at line 72 of file pflotran.py.
pydfnworks.dfnFlow.pflotran.vtk_file |
Definition at line 627 of file pflotran.py.