pydfnWorks
python wrapper for dfnWorks
dfn2graph.py
Go to the documentation of this file.
1 import networkx as nx
2 import numpy as np
3 import json
4 
8 from networkx.readwrite import json_graph
9 
10 import matplotlib
11 matplotlib.use('Agg')
12 
13 import matplotlib.pylab as plt
14 from itertools import islice
15 
16 
17 def create_graph(self, graph_type, inflow, outflow):
18  """Header function to create a graph based on a DFN
19 
20  Parameters
21  ----------
22  self : object
23  DFN Class object
24  graph_type : string
25  Option for what graph representation of the DFN is requested. Currently supported are fracture, intersection, and bipartitie
26  inflow : string
27  Name of inflow boundary (connect to source)
28  outflow : string
29  Name of outflow boundary (connect to target)
30 
31  Returns
32  -------
33  G : NetworkX Graph
34  Graph based on DFN
35 
36  Notes
37  -----
38 
39 """
40 
41  if graph_type == "fracture":
42  G = create_fracture_graph(inflow, outflow)
43  elif graph_type == "intersection":
44  G = create_intersection_graph(inflow, outflow)
45  elif graph_type == "bipartite":
46  G = create_bipartite_graph(inflow, outflow)
47  else:
48  print("ERROR! Unknown graph type")
49  return []
50  return G
51 
52 
54  outflow,
55  topology_file="connectivity.dat",
56  fracture_info="fracture_info.dat"):
57  """ Create a graph based on topology of network. Fractures
58  are represented as nodes and if two fractures intersect
59  there is an edge between them in the graph.
60 
61  Source and Target node are added to the graph.
62 
63  Parameters
64  ----------
65  inflow : string
66  Name of inflow boundary (connect to source)
67  outflow : string
68  Name of outflow boundary (connect to target)
69  topology_file : string
70  Name of adjacency matrix file for a DFN default=connectivity.dat
71  fracture_infor : str
72  filename for fracture information
73 
74  Returns
75  -------
76  G : NetworkX Graph
77  NetworkX Graph where vertices in the graph correspond to fractures and edges indicated two fractures intersect
78 
79  Notes
80  -----
81  """
82  print("--> Loading Graph based on topology in " + topology_file)
83  G = nx.Graph(representation="fracture")
84  with open(topology_file, "r") as infile:
85  for i, line in enumerate(infile):
86  conn = [int(n) for n in line.split()]
87  for j in conn:
88  G.add_edge(i + 1, j)
89 
90  inflow_filename = inflow + ".dat"
91  outflow_filename = outflow + ".dat"
92  inflow = np.genfromtxt(inflow_filename).astype(int)
93  outflow = np.genfromtxt(outflow_filename).astype(int)
94 
95  try:
96  if len(inflow) > 1:
97  inflow = list(inflow)
98  except:
99  inflow = [inflow.tolist()]
100 
101  try:
102  if len(outflow) > 1:
103  outflow = list(outflow)
104  except:
105  outflow = [outflow.tolist()]
106 
107  G.add_node('s')
108  G.add_node('t')
109  G.add_edges_from(zip(['s'] * (len(inflow)), inflow))
110  G.add_edges_from(zip(outflow, ['t'] * (len(outflow))))
111  add_perm(G, fracture_info)
112  print("--> Graph loaded")
113  return G
114 
115 
116 def boundary_index(bc_name):
117  """Determines boundary index in intersections_list.dat from name
118 
119  Parameters
120  ----------
121  bc_name : string
122  Boundary condition name
123 
124  Returns
125  -------
126  bc_index : int
127  integer indexing of cube faces
128 
129  Notes
130  -----
131  top = 1
132  bottom = 2
133  left = 3
134  front = 4
135  right = 5
136  back = 6
137  """
138  bc_dict = {
139  "top": -1,
140  "bottom": -2,
141  "left": -3,
142  "front": -4,
143  "right": -5,
144  "back": -6
145  }
146  try:
147  return bc_dict[bc_name]
148  except:
149  error = "Unknown boundary condition: %s\nExiting\n" % bc
150  sys.stderr.write(error)
151  sys.exit(1)
152 
153 
155  outflow,
156  intersection_file="intersection_list.dat",
157  fracture_info="fracture_info.dat"):
158  """ Create a graph based on topology of network.
159  Edges are represented as nodes and if two intersections
160  are on the same fracture, there is an edge between them in the graph.
161 
162  Source and Target node are added to the graph.
163 
164  Parameters
165  ----------
166  inflow : string
167  Name of inflow boundary
168  outflow : string
169  Name of outflow boundary
170  intersection_file : string
171  File containing intersection information
172  File Format:
173  fracture 1, fracture 2, x center, y center, z center, intersection length
174 
175  fracture_infor : str
176  filename for fracture information
177  Returns
178  -------
179  G : NetworkX Graph
180  Vertices have attributes x,y,z location and length. Edges has attribute length
181 
182  Notes
183  -----
184  Aperture and Perm on edges can be added using add_app and add_perm functions
185  """
186 
187  print("Creating Graph Based on DFN")
188  print("Intersections being mapped to nodes and fractures to edges")
189  inflow_index = boundary_index(inflow)
190  outflow_index = boundary_index(outflow)
191 
192  f = open(intersection_file)
193  f.readline()
194  frac_edges = []
195  for line in f:
196  frac_edges.append(line.rstrip().split())
197  f.close()
198 
199  # Tag mapping
200  G = nx.Graph(representation="intersection")
201  remove_list = []
202 
203  # each edge in the DFN is a node in the graph
204  for i in range(len(frac_edges)):
205  f1 = int(frac_edges[i][0])
206  keep = True
207  if frac_edges[i][1] == 's' or frac_edges[i][1] == 't':
208  f2 = frac_edges[i][1]
209  elif int(frac_edges[i][1]) > 0:
210  f2 = int(frac_edges[i][1])
211  elif int(frac_edges[i][1]) == inflow_index:
212  f2 = 's'
213  elif int(frac_edges[i][1]) == outflow_index:
214  f2 = 't'
215  elif int(frac_edges[i][1]) < 0:
216  keep = False
217 
218  if keep:
219  # note fractures of the intersection
220  G.add_node(i, frac=(f1, f2))
221  # keep intersection location and length
222  G.nodes[i]['x'] = float(frac_edges[i][2])
223  G.nodes[i]['y'] = float(frac_edges[i][3])
224  G.nodes[i]['z'] = float(frac_edges[i][4])
225  G.nodes[i]['length'] = float(frac_edges[i][5])
226 
227  nodes = list(nx.nodes(G))
228  f1 = nx.get_node_attributes(G, 'frac')
229  # identify which edges are on whcih fractures
230  for i in nodes:
231  e = set(f1[i])
232  for j in nodes:
233  if i != j:
234  tmp = set(f1[j])
235  x = e.intersection(tmp)
236  if len(x) > 0:
237  x = list(x)[0]
238  # Check for Boundary Intersections
239  # This stops boundary fractures from being incorrectly
240  # connected
241  # If not, add edge between
242  if x != 's' and x != 't':
243  xi = G.nodes[i]['x']
244  yi = G.nodes[i]['y']
245  zi = G.nodes[i]['z']
246 
247  xj = G.nodes[j]['x']
248  yj = G.nodes[j]['y']
249  zj = G.nodes[j]['z']
250 
251  distance = np.sqrt((xi - xj)**2 + (yi - yj)**2 +
252  (zi - zj)**2)
253  G.add_edge(i, j, frac=x, length=distance)
254 
255  # Add Sink and Source nodes
256  G.add_node('s')
257  G.add_node('t')
258 
259  for i in nodes:
260  e = set(f1[i])
261  if len(e.intersection(set('s'))) > 0 or len(e.intersection(set(
262  [-1]))) > 0:
263  G.add_edge(i, 's', frac='s', length=0.0)
264  if len(e.intersection(set('t'))) > 0 or len(e.intersection(set(
265  [-2]))) > 0:
266  G.add_edge(i, 't', frac='t', length=0.0)
267  add_perm(G, fracture_info)
268  print("Graph Construction Complete")
269  return G
270 
271 
273  outflow,
274  intersection_list='intersection_list.dat',
275  fracture_info='fracture_info.dat'):
276  """Creates a bipartite graph of the DFN.
277  Nodes are in two sets, fractures and intersections, with edges connecting them.
278 
279 
280  Parameters
281  ----------
282  inflow : str
283  name of inflow boundary
284  outflow : str
285  name of outflow boundary
286  intersection_list: str
287  filename of intersections generated from DFN
288  fracture_infor : str
289  filename for fracture information
290 
291  Returns
292  -------
293  B : NetworkX Graph
294 
295  Notes
296  -----
297  See Hyman et al. 2018 "Identifying Backbones in Three-Dimensional Discrete Fracture Networks: A Bipartite Graph-Based Approach" SIAM Multiscale Modeling and Simulation for more details
298 """
299 
300  print("--> Creating Bipartite Graph")
301 
302  # generate sequential letter sequence as ids for fractures
303  # e..g aaaaa aaaaab aaaaac
304  from itertools import product
305 
306  def intersection_id_generator(length=5):
307  chars = 'abcdefghijklmnopqrstuvwxyz'
308  for p in product(chars, repeat=length):
309  yield (''.join(p))
310 
311  B = nx.Graph(representation="bipartite")
312  # keep track of the sets of fractures and intersections
313  B.fractures = set()
314  B.intersections = set()
315  intersection_id = intersection_id_generator()
316 
317  inflow_index = boundary_index(inflow)
318  outflow_index = boundary_index(outflow)
319 
320  with open(intersection_list) as f:
321  header = f.readline()
322  data = f.read().strip()
323  for line in data.split('\n'):
324  fracture1, fracture2, x, y, z, length = line.split(' ')
325  fracture1 = int(fracture1)
326  fracture2 = int(fracture2)
327  if fracture2 < 0:
328  if fracture2 == inflow_index:
329  fracture2 = 's'
330  elif fracture2 == outflow_index:
331  fracture2 = 't'
332  intersection = next(intersection_id)
333  # add intersection node explicitly to include intersection properties
334  B.add_node(intersection,
335  x=float(x),
336  y=float(y),
337  z=float(z),
338  length=float(length))
339  B.intersections.add(intersection)
340 
341  B.add_edge(intersection, fracture1, frac=fracture1)
342  B.fractures.add(fracture1)
343  if fracture2 > 0 or fracture2 == 's' or fracture2 == 't':
344  B.add_edge(intersection, fracture2, frac=fracture2)
345  B.fractures.add(fracture2)
346 
347  # add source and sink for intersections so they will appear in intersection projection
348  B.add_edge('intersection_s', 's')
349  B.add_edge('intersection_t', 't')
350 
351  # add fracture info
352  with open(fracture_info) as f:
353  header = f.readline()
354  data = f.read().strip()
355  for fracture, line in enumerate(data.split('\n'), 1):
356  c, perm, aperture = line.split(' ')
357  B.nodes[fracture]['perm'] = float(perm)
358  B.nodes[fracture]['aperture'] = float(aperture)
359 
360  print("--> Complete")
361 
362  return B
363 
364 
365 def add_fracture_source(self, G, source):
366  """Returns the k shortest paths in a graph
367 
368  Parameters
369  ----------
370  G : NetworkX Graph
371  NetworkX Graph based on a DFN
372  source_list : list
373  list of integers corresponding to fracture numbers
374  remove_old_source: bool
375  remove old source from the graph
376 
377  Returns
378  -------
379  G : NetworkX Graph
380 
381  Notes
382  -----
383  bipartite graph not supported
384 
385  """
386 
387  if not type(source) == list:
388  source = [source]
389 
390  print("--> Adding new source connections")
391  print("--> Warning old source will be removed!!!")
392 
393  if G.graph['representation'] == "fracture":
394  # removing old source term and all connections
395  G.remove_node('s')
396  # add new source node
397  G.add_node('s')
398 
399  G.nodes['s']['perm'] = 1.0
400  G.nodes['s']['iperm'] = 1.0
401 
402  for u in source:
403  G.add_edge(u, 's')
404 
405  elif G.graph['representation'] == "intersection":
406  # removing old source term and all connections
407  nodes_to_remove = ['s']
408  for u, d in G.nodes(data=True):
409  if u != 's' and u != 't':
410  f1, f2 = d["frac"]
411  #print("node {0}: f1 {1}, f2 {2}".format(u,f1,f2))
412  if f2 == 's':
413  nodes_to_remove.append(u)
414 
415  print("--> Removing nodes: ", nodes_to_remove)
416  G.remove_nodes_from(nodes_to_remove)
417 
418  # add new source node
419  G.add_node('s')
420  for u, d in G.nodes(data=True):
421  if u != 's' and u != 't':
422  f1 = d["frac"][0]
423  f2 = d["frac"][1]
424  if f1 in source:
425  print(
426  "--> Adding edge between {0} and new source / fracture {1}"
427  .format(u, f1))
428  G.add_edge(u, 's', frac=f1, length=0., perm=1., iperm=1.)
429  elif f2 in source:
430  print(
431  "--> Adding edge between {0} and new source / fracture {1}"
432  .format(u, f2))
433  G.add_edge(u, 's', frac=f2, length=0., perm=1., iperm=1.)
434 
435  elif G.graph['representation'] == "bipartite":
436  print("--> Not supported for bipartite graph")
437  print("--> Returning unchanged graph")
438  return G
439 
440 
441 def add_fracture_target(self, G, target):
442  """Returns the k shortest paths in a graph
443 
444  Parameters
445  ----------
446  G : NetworkX Graph
447  NetworkX Graph based on a DFN
448  target : list
449  list of integers corresponding to fracture numbers
450  Returns
451  -------
452  G : NetworkX Graph
453 
454  Notes
455  -----
456  bipartite graph not supported
457 
458  """
459 
460  if not type(target) == list:
461  source = [target]
462 
463  print("--> Adding new target connections")
464  print("--> Warning old target will be removed!!!")
465 
466  if G.graph['representation'] == "fracture":
467  # removing old target term and all connections
468  G.remove_node('t')
469  # add new target node
470  G.add_node('t')
471 
472  G.nodes['t']['perm'] = 1.0
473  G.nodes['t']['iperm'] = 1.0
474 
475  for u in target:
476  G.add_edge(u, 't')
477 
478  elif G.graph['representation'] == "intersection":
479  # removing old target term and all connections
480  nodes_to_remove = ['t']
481  for u, d in G.nodes(data=True):
482  if u != 's' and u != 't':
483  f1, f2 = d["frac"]
484  #print("node {0}: f1 {1}, f2 {2}".format(u,f1,f2))
485  if f2 == 't':
486  nodes_to_remove.append(u)
487 
488  print("--> Removing nodes: ", nodes_to_remove)
489  G.remove_nodes_from(nodes_to_remove)
490 
491  # add new target node
492  G.add_node('t')
493  for u, d in G.nodes(data=True):
494  if u != 's' and u != 't':
495  f1 = d["frac"][0]
496  f2 = d["frac"][1]
497  if f1 in target:
498  print(
499  "--> Adding edge between {0} and new target / fracture {1}"
500  .format(u, f1))
501  G.add_edge(u, 't', frac=f1, length=0., perm=1., iperm=1.)
502  elif f2 in target:
503  print(
504  "--> Adding edge between {0} and new target / fracture {1}"
505  .format(u, f2))
506  G.add_edge(u, 't', frac=f2, length=0., perm=1., iperm=1.)
507 
508  elif G.graph['representation'] == "bipartite":
509  print("--> Not supported for bipartite graph")
510  print("--> Returning unchanged graph")
511  return G
512 
513 
514 def k_shortest_paths(G, k, source, target, weight):
515  """Returns the k shortest paths in a graph
516 
517  Parameters
518  ----------
519  G : NetworkX Graph
520  NetworkX Graph based on a DFN
521  k : int
522  Number of requested paths
523  source : node
524  Starting node
525  target : node
526  Ending node
527  weight : string
528  Edge weight used for finding the shortest path
529 
530  Returns
531  -------
532  paths : sets of nodes
533  a list of lists of nodes in the k shortest paths
534 
535  Notes
536  -----
537  Edge weights must be numerical and non-negative
538 """
539  return list(
540  islice(nx.shortest_simple_paths(G, source, target, weight=weight), k))
541 
542 
543 def k_shortest_paths_backbone(self, G, k, source='s', target='t', weight=None):
544  """Returns the subgraph made up of the k shortest paths in a graph
545 
546  Parameters
547  ----------
548  G : NetworkX Graph
549  NetworkX Graph based on a DFN
550  k : int
551  Number of requested paths
552  source : node
553  Starting node
554  target : node
555  Ending node
556  weight : string
557  Edge weight used for finding the shortest path
558 
559  Returns
560  -------
561  H : NetworkX Graph
562  Subgraph of G made up of the k shortest paths
563 
564  Notes
565  -----
566  See Hyman et al. 2017 "Predictions of first passage times in sparse discrete fracture networks using graph-based reductions" Physical Review E for more details
567 """
568 
569  print("\n--> Determining %d shortest paths in the network" % k)
570  H = G.copy()
571  k_shortest = set([])
572  for path in k_shortest_paths(G, k, source, target, weight):
573  k_shortest |= set(path)
574  k_shortest.remove('s')
575  k_shortest.remove('t')
576  path_nodes = sorted(list(k_shortest))
577  path_nodes.append('s')
578  path_nodes.append('t')
579  nodes = list(G.nodes())
580  secondary = list(set(nodes) - set(path_nodes))
581  for n in secondary:
582  H.remove_node(n)
583  return H
584  print("--> Complete\n")
585 
586 
587 def pull_source_and_target(nodes, source='s', target='t'):
588  """Removes source and target from list of nodes, useful for dumping subnetworks to file for remeshing
589 
590  Parameters
591  ----------
592  nodes :list
593  List of nodes in the graph
594  source : node
595  Starting node
596  target : node
597  Ending node
598  Returns
599  -------
600  nodes : list
601  List of nodes with source and target nodes removed
602 
603  Notes
604  -----
605 
606 """
607  for node in [source, target]:
608  try:
609  nodes.remove(node)
610  except:
611  pass
612  return nodes
613 
614 
615 def dump_fractures(self, G, filename):
616  """Write fracture numbers assocaited with the graph G out into an ASCII file inputs
617 
618  Parameters
619  ----------
620  self : object
621  DFN Class
622  G : NetworkX graph
623  NetworkX Graph based on the DFN
624  filename : string
625  Output filename
626 
627  Returns
628  -------
629 
630  Notes
631  -----
632  """
633 
634  if G.graph['representation'] == "fracture":
635  nodes = list(G.nodes())
636  elif G.graph['representation'] == "intersection":
637  nodes = []
638  for u, v, d in G.edges(data=True):
639  nodes.append(G[u][v]['frac'])
640  nodes = list(set(nodes))
641  elif G.graph['representation'] == "bipartite":
642  nodes = []
643  for u, v, d in G.edges(data=True):
644  nodes.append(G[u][v]['frac'])
645  nodes = list(set(nodes))
646 
647  nodes = pull_source_and_target(nodes)
648  fractures = [int(i) for i in nodes]
649  fractures = sorted(fractures)
650  print("--> Dumping %s" % filename)
651  np.savetxt(filename, fractures, fmt="%d")
652 
653 
654 def greedy_edge_disjoint(self, G, source='s', target='t', weight='None', k=''):
655  """
656  Greedy Algorithm to find edge disjoint subgraph from s to t.
657  See Hyman et al. 2018 SIAM MMS
658 
659  Parameters
660  ----------
661  self : object
662  DFN Class Object
663  G : NetworkX graph
664  NetworkX Graph based on the DFN
665  source : node
666  Starting node
667  target : node
668  Ending node
669  weight : string
670  Edge weight used for finding the shortest path
671  k : int
672  Number of edge disjoint paths requested
673 
674  Returns
675  -------
676  H : NetworkX Graph
677  Subgraph of G made up of the k shortest of all edge-disjoint paths from source to target
678 
679  Notes
680  -----
681  1. Edge weights must be numerical and non-negative.
682  2. See Hyman et al. 2018 "Identifying Backbones in Three-Dimensional Discrete Fracture Networks: A Bipartite Graph-Based Approach" SIAM Multiscale Modeling and Simulation for more details
683 
684  """
685  print("--> Identifying edge disjoint paths")
686  if G.graph['representation'] != "intersection":
687  print(
688  "--> ERROR!!! Wrong type of DFN graph representation\nRepresentation must be intersection\nReturning Empty Graph\n"
689  )
690  return nx.Graph()
691  Gprime = G.copy()
692  Hprime = nx.Graph()
693  Hprime.graph['representation'] = G.graph['representation']
694  cnt = 0
695 
696  # if a number of paths in not provided k will equal the min cut between s and t
697  min_cut = len(nx.minimum_edge_cut(G, 's', 't'))
698  if k == '' or k > min_cut:
699  k = min_cut
700 
701  while nx.has_path(Gprime, source, target):
702  path = nx.shortest_path(Gprime, source, target, weight=weight)
703  H = Gprime.subgraph(path)
704  Hprime.add_edges_from(H.edges(data=True))
705  Gprime.remove_edges_from(list(H.edges()))
706 
707  cnt += 1
708  if cnt > k:
709  break
710  print("--> Complete")
711  return Hprime
712 
713 
714 def plot_graph(self, G, source='s', target='t', output_name="dfn_graph"):
715  """ Create a png of a graph with source nodes colored blue, target red, and all over nodes black
716 
717  Parameters
718  ----------
719  G : NetworkX graph
720  NetworkX Graph based on the DFN
721  source : node
722  Starting node
723  target : node
724  Ending node
725  output_name : string
726  Name of output file (no .png)
727 
728  Returns
729  -------
730 
731  Notes
732  -----
733  Image is written to output_name.png
734 
735  """
736  print("\n--> Plotting Graph")
737  print("--> Output file: %s.png" % output_name)
738  # get positions for all nodes
739  pos = nx.spring_layout(G)
740  nodes = list(G.nodes)
741  # draw nodes
742  nx.draw_networkx_nodes(G,
743  pos,
744  nodelist=nodes,
745  node_color='k',
746  node_size=10,
747  alpha=1.0)
748  nx.draw_networkx_nodes(G,
749  pos,
750  nodelist=[source],
751  node_color='b',
752  node_size=50,
753  alpha=1.0)
754  nx.draw_networkx_nodes(G,
755  pos,
756  nodelist=[target],
757  node_color='r',
758  node_size=50,
759  alpha=1.0)
760 
761  nx.draw_networkx_edges(G, pos, width=1.0, alpha=0.5)
762  plt.axis('off')
763  plt.tight_layout()
764  plt.savefig(output_name + ".png")
765  plt.clf()
766  print("--> Plotting Graph Complete\n")
767 
768 
769 def dump_json_graph(self, G, name):
770  """Write graph out in json format
771 
772  Parameters
773  ----------
774  self : object
775  DFN Class
776  G :networkX graph
777  NetworkX Graph based on the DFN
778  name : string
779  Name of output file (no .json)
780 
781  Returns
782  -------
783 
784  Notes
785  -----
786 
787 """
788  print("--> Dumping Graph into file: " + name + ".json")
789  jsondata = json_graph.node_link_data(G)
790  with open(name + '.json', 'w') as fp:
791  json.dump(jsondata, fp)
792  print("--> Complete")
793 
794 
795 def load_json_graph(self, name):
796  """ Read in graph from json format
797 
798  Parameters
799  ----------
800  self : object
801  DFN Class
802  name : string
803  Name of input file (no .json)
804 
805  Returns
806  -------
807  G :networkX graph
808  NetworkX Graph based on the DFN
809 """
810 
811  print("Loading Graph in file: " + name + ".json")
812  fp = open(name + '.json')
813  G = json_graph.node_link_graph(json.load(fp))
814  print("Complete")
815  return G
816 
817 
818 def add_perm(G, fracture_info="fracture_info.dat"):
819  """ Add fracture permeability to Graph. If Graph representation is
820  fracture, then permeability is a node attribute. If graph representation
821  is intersection, then permeability is an edge attribute
822 
823 
824  Parameters
825  ----------
826  G :networkX graph
827  NetworkX Graph based on the DFN
828 
829  fracture_infor : str
830  filename for fracture information
831  Returns
832  -------
833 
834  Notes
835  -----
836 
837 """
838 
839  perm = np.genfromtxt(fracture_info, skip_header=1)[:, 1]
840  if G.graph['representation'] == "fracture":
841  nodes = list(nx.nodes(G))
842  for n in nodes:
843  if n != 's' and n != 't':
844  G.nodes[n]['perm'] = perm[n - 1]
845  G.nodes[n]['iperm'] = 1.0 / perm[n - 1]
846  else:
847  G.nodes[n]['perm'] = 1.0
848  G.nodes[n]['iperm'] = 1.0
849 
850  elif G.graph['representation'] == "intersection":
851  edges = list(nx.edges(G))
852  for u, v in edges:
853  x = G[u][v]['frac']
854  if x != 's' and x != 't':
855  G[u][v]['perm'] = perm[x - 1]
856  G[u][v]['iperm'] = 1.0 / perm[x - 1]
857  else:
858  G[u][v]['perm'] = 1.0
859  G[u][v]['iperm'] = 1.0
860  elif G.graph['representation'] == "bipartite":
861  # add fracture info
862  with open(fracture_info) as f:
863  header = f.readline()
864  data = f.read().strip()
865  for fracture, line in enumerate(data.split('\n'), 1):
866  c, perm, aperture = line.split(' ')
867  G.nodes[fracture]['perm'] = float(perm)
868  G.nodes[fracture]['iperm'] = 1.0 / float(perm)
869  G.nodes[fracture]['aperture'] = float(aperture)
870 
871 
872 def add_area(G, fracture_info="fracture_info.dat"):
873  ''' Read Fracture aperture from fracture_info.dat and
874  load on the edges in the graph. Graph must be intersection to node
875  representation
876 
877  Parameters
878  ----------
879  G : NetworkX Graph
880  networkX graph
881  fracture_info : str
882  filename for fracture information
883 
884  Returns
885  -------
886  None
887 '''
888 
889  aperture = np.genfromtxt(fracture_info, skip_header=1)[:, 2]
890  edges = list(nx.edges(G))
891  for u, v in edges:
892  x = G.edges[u, v]['frac']
893  if x != 's' and x != 't':
894  G.edges[u,
895  v]['area'] = aperture[x - 1] * (G.nodes[u]['length'] +
896  G.nodes[v]['length']) / 2.0
897  else:
898  G.edges[u, v]['area'] = 1.0
899  return
900 
901 
902 def add_weight(G):
903  '''Compute weight w = K*A/L associated with each edge
904  Parameters
905  ----------
906  G : NetworkX Graph
907  networkX graph
908 
909  Returns
910  -------
911  None
912 '''
913  edges = list(nx.edges(G))
914  for u, v in edges:
915  if G.edges[u, v]['length'] > 0:
916  G.edges[u, v]['weight'] = G.edges[u, v]['perm'] * G.edges[
917  u, v]['area'] / G.edges[u, v]['length']
918  return
def add_fracture_source(self, G, source)
Definition: dfn2graph.py:365
def k_shortest_paths_backbone(self, G, k, source='s', target='t', weight=None)
Definition: dfn2graph.py:543
def create_graph(self, graph_type, inflow, outflow)
Definition: dfn2graph.py:17
def dump_json_graph(self, G, name)
Definition: dfn2graph.py:769
def k_shortest_paths(G, k, source, target, weight)
Definition: dfn2graph.py:514
def dump_fractures(self, G, filename)
Definition: dfn2graph.py:615
def pull_source_and_target(nodes, source='s', target='t')
Definition: dfn2graph.py:587
def add_perm(G, fracture_info="fracture_info.dat")
Definition: dfn2graph.py:818
def create_intersection_graph(inflow, outflow, intersection_file="intersection_list.dat", fracture_info="fracture_info.dat")
Definition: dfn2graph.py:157
def add_area(G, fracture_info="fracture_info.dat")
Definition: dfn2graph.py:872
def create_fracture_graph(inflow, outflow, topology_file="connectivity.dat", fracture_info="fracture_info.dat")
Definition: dfn2graph.py:56
def load_json_graph(self, name)
Definition: dfn2graph.py:795
def plot_graph(self, G, source='s', target='t', output_name="dfn_graph")
Definition: dfn2graph.py:714
def create_bipartite_graph(inflow, outflow, intersection_list='intersection_list.dat', fracture_info='fracture_info.dat')
Definition: dfn2graph.py:275
def greedy_edge_disjoint(self, G, source='s', target='t', weight='None', k='')
Definition: dfn2graph.py:654
def add_fracture_target(self, G, target)
Definition: dfn2graph.py:441