source: trunk/anuga_core/source/anuga/utilities/sww_merge.py @ 9019

Last change on this file since 9019 was 9019, checked in by steve, 11 years ago

Added storing of centroids to the sww file

File size: 25.2 KB
Line 
1"""
2    Merge a list of .sww files together into a single file.
3"""
4
5import numpy as num
6from anuga.utilities.numerical_tools import ensure_numeric
7
8from anuga.file.netcdf import NetCDFFile
9from anuga.config import netcdf_mode_r, netcdf_mode_w, netcdf_mode_a
10from anuga.config import netcdf_float, netcdf_float32, netcdf_int
11from anuga.file.sww import SWW_file, Write_sww
12
13def sww_merge(domain_global_name, np, verbose=False):
14
15    output = domain_global_name+".sww"
16    swwfiles = [ domain_global_name+"_P"+str(np)+"_"+str(v)+".sww" for v in range(np)]
17
18    _sww_merge(swwfiles, output, verbose)
19
20
21def sww_merge_parallel(domain_global_name, np, verbose=False, delete_old=False):
22
23    output = domain_global_name+".sww"
24    swwfiles = [ domain_global_name+"_P"+str(np)+"_"+str(v)+".sww" for v in range(np)]
25
26    fid = NetCDFFile(swwfiles[0], netcdf_mode_r)
27
28    try: # works with netcdf4
29        number_of_volumes = len(fid.dimensions['number_of_volumes'])
30        number_of_points = len(fid.dimensions['number_of_points'])
31    except: # works with scientific.io.netcdf
32        number_of_volumes = int(fid.dimensions['number_of_volumes'])
33        number_of_points = int(fid.dimensions['number_of_points'])
34
35    fid.close()
36
37    if 3*number_of_volumes == number_of_points:
38        _sww_merge_parallel_non_smooth(swwfiles, output, verbose, delete_old)
39    else:
40        _sww_merge_parallel_smooth(swwfiles, output, verbose, delete_old)
41       
42
43def _sww_merge(swwfiles, output, verbose=False):
44    """
45        Merge a list of sww files into a single file.
46       
47        May be useful for parallel runs. Note that colinear points and
48        edges are not merged: there will essentially be multiple meshes within
49        the one sww file.
50       
51        The sww files to be merged must have exactly the same timesteps. Note
52        that some advanced information and custom quantities may not be
53        exported.
54       
55        swwfiles is a list of .sww files to merge.
56        output is the output filename, including .sww extension.
57        verbose True to log output information
58    """
59
60    if verbose:
61        print "MERGING SWW Files"
62       
63    static_quantities = ['elevation']
64    dynamic_quantities = ['stage', 'xmomentum', 'ymomentum']
65   
66    first_file = True
67    tri_offset = 0
68    for filename in swwfiles:
69        if verbose:
70            print 'Reading file ', filename, ':'   
71   
72        fid = NetCDFFile(filename, netcdf_mode_r)
73
74
75
76       
77       
78        tris = fid.variables['volumes'][:]       
79         
80        if first_file:
81            times = fid.variables['time'][:]
82            x = []
83            y = []
84            out_tris = list(tris) 
85            out_s_quantities = {}
86            out_d_quantities = {}
87
88
89            xllcorner = fid.xllcorner
90            yllcorner = fid.yllcorner
91
92            order      = fid.order
93            xllcorner  = fid.xllcorner;
94            yllcorner  = fid.yllcorner ;
95            zone       = fid.zone;
96            false_easting  = fid.false_easting;
97            false_northing = fid.false_northing;
98            datum      = fid.datum;
99            projection = fid.projection;
100
101           
102            for quantity in static_quantities:
103                out_s_quantities[quantity] = []
104
105            # Quantities are stored as a 2D array of timesteps x data.
106            for quantity in dynamic_quantities:
107                out_d_quantities[quantity] = [ [] for _ in range(len(times))]
108                 
109            description = 'merged:' + getattr(fid, 'description')         
110            first_file = False
111        else:
112            for tri in tris:
113                # Advance new tri indices to point at newly appended points.
114                verts = [vertex+tri_offset for vertex in tri]
115                out_tris.append(verts)
116
117
118
119        try: # works with netcdf4
120            num_pts = len(fid.dimensions['number_of_points'])
121        except: # works with scientific.io.netcdf
122            num_pts = int(fid.dimensions['number_of_points'])
123
124        tri_offset += num_pts
125       
126        if verbose:
127            print '  new triangle index offset is ', tri_offset
128           
129        x.extend(list(fid.variables['x'][:]))
130        y.extend(list(fid.variables['y'][:]))
131       
132        # Grow the list of static quantities associated with the x,y points
133        for quantity in static_quantities:
134            out_s_quantities[quantity].extend(fid.variables[quantity][:])
135           
136        #Collate all dynamic quantities according to their timestep
137        for quantity in dynamic_quantities:
138            time_chunks = fid.variables[quantity][:]
139            for i, time_chunk in enumerate(time_chunks):
140                out_d_quantities[quantity][i].extend(time_chunk)           
141   
142    # Mash all points into a single big list   
143    points = [[xx, yy] for xx, yy in zip(x, y)]
144
145    points = num.asarray(points).astype(netcdf_float32)
146
147    fid.close()
148
149    #---------------------------
150    # Write out the SWW file
151    #---------------------------
152
153    if verbose:
154        print 'Writing file ', output, ':'
155    fido = NetCDFFile(output, netcdf_mode_w)
156    sww = Write_sww(static_quantities, dynamic_quantities)
157    sww.store_header(fido, times,
158                             len(out_tris),
159                             len(points),
160                             description=description,
161                             sww_precision=netcdf_float32)
162
163
164
165    from anuga.coordinate_transforms.geo_reference import Geo_reference
166    geo_reference = Geo_reference()
167   
168    sww.store_triangulation(fido, points, out_tris, points_georeference=geo_reference)
169
170    fido.order      = order
171    fido.xllcorner  = xllcorner;
172    fido.yllcorner  = yllcorner ;
173    fido.zone       = zone;
174    fido.false_easting  = false_easting;
175    fido.false_northing = false_northing;
176    fido.datum      = datum;
177    fido.projection = projection;
178       
179    sww.store_static_quantities(fido, verbose=verbose, **out_s_quantities)
180
181    # Write out all the dynamic quantities for each timestep
182    for q in dynamic_quantities:
183        q_values = out_d_quantities[q]
184        for i, time_slice in enumerate(q_values):
185            fido.variables[q][i] = num.array(time_slice, netcdf_float32)
186       
187        # This updates the _range values
188        q_range = fido.variables[q + Write_sww.RANGE][:]
189        q_values_min = num.min(q_values)
190        if q_values_min < q_range[0]:
191            fido.variables[q + Write_sww.RANGE][0] = q_values_min
192        q_values_max = num.max(q_values)
193        if q_values_max > q_range[1]:
194            fido.variables[q + Write_sww.RANGE][1] = q_values_max       
195
196                                       
197    fido.close()
198
199
200def _sww_merge_parallel_smooth(swwfiles, output,  verbose=False, delete_old=False):
201    """
202        Merge a list of sww files into a single file.
203       
204        Use to merge files created by parallel runs.
205
206        The sww files to be merged must have exactly the same timesteps.
207
208        It is assumed that the separate sww files have been stored in non_smooth
209        format.
210
211        Note that some advanced information and custom quantities may not be
212        exported.
213       
214        swwfiles is a list of .sww files to merge.
215        output is the output filename, including .sww extension.
216        verbose True to log output information
217    """
218
219    if verbose:
220        print "MERGING SWW Files"
221       
222   
223    first_file = True
224    tri_offset = 0
225    for filename in swwfiles:
226        if verbose:
227            print 'Reading file ', filename, ':'   
228   
229        fid = NetCDFFile(filename, netcdf_mode_r)
230         
231        if first_file:
232
233            times    = fid.variables['time'][:]
234            n_steps = len(times)
235            #number_of_timesteps = fid.dimensions['number_of_timesteps']
236            #print n_steps, number_of_timesteps
237            starttime = int(fid.starttime)
238           
239            out_s_quantities = {}
240            out_d_quantities = {}
241
242            out_s_c_quantities = {}
243            out_d_c_quantities = {}
244
245
246            xllcorner = fid.xllcorner
247            yllcorner = fid.yllcorner
248
249            number_of_global_triangles = int(fid.number_of_global_triangles)
250            number_of_global_nodes     = int(fid.number_of_global_nodes)
251
252            order      = fid.order
253            xllcorner  = fid.xllcorner;
254            yllcorner  = fid.yllcorner ;
255            zone       = fid.zone;
256            false_easting  = fid.false_easting;
257            false_northing = fid.false_northing;
258            datum      = fid.datum;
259            projection = fid.projection;
260
261            g_volumes = num.zeros((number_of_global_triangles,3),num.int)
262            g_x = num.zeros((number_of_global_nodes,),num.float32)
263            g_y = num.zeros((number_of_global_nodes,),num.float32)
264
265            g_points = num.zeros((number_of_global_nodes,2),num.float32)
266
267            #=====================================
268            # Deal with the vertex based variables
269            #=====================================
270            quantities = set(['elevation', 'friction', 'stage', 'xmomentum',
271                              'ymomentum', 'xvelocity', 'yvelocity', 'height'])
272            variables = set(fid.variables.keys())
273
274            quantities = list(quantities & variables)
275           
276            static_quantities = []
277            dynamic_quantities = []
278
279            for quantity in quantities:
280                # Test if quantity is static
281                if n_steps == fid.variables[quantity].shape[0]:
282                    dynamic_quantities.append(quantity)
283                else:
284                    static_quantities.append(quantity)
285               
286            for quantity in static_quantities:
287                out_s_quantities[quantity] = num.zeros((number_of_global_nodes,),num.float32)
288
289            # Quantities are stored as a 2D array of timesteps x data.
290            for quantity in dynamic_quantities:
291                out_d_quantities[quantity] = \
292                      num.zeros((n_steps,number_of_global_nodes),num.float32)
293
294            #=======================================
295            # Deal with the centroid based variables
296            #=======================================
297            quantities = set(['elevation_c', 'friction_c', 'stage_c', 'xmomentum_c',
298                              'ymomentum_c', 'xvelocity_c', 'yvelocity_c', 'height_c'])
299            variables = set(fid.variables.keys())
300
301            quantities = list(quantities & variables)
302           
303            static_c_quantities = []
304            dynamic_c_quantities = []
305
306            for quantity in quantities:
307                # Test if quantity is static
308                if n_steps == fid.variables[quantity].shape[0]:
309                    dynamic_c_quantities.append(quantity)
310                else:
311                    static_c_quantities.append(quantity)
312               
313            for quantity in static_c_quantities:
314                out_s_c_quantities[quantity] = num.zeros((number_of_global_triangles,),num.float32)
315
316            # Quantities are stored as a 2D array of timesteps x data.
317            for quantity in dynamic_c_quantities:
318                out_d_c_quantities[quantity] = \
319                      num.zeros((n_steps,number_of_global_triangles),num.float32)
320                 
321            description = 'merged:' + getattr(fid, 'description')         
322            first_file = False
323
324
325        # Read in from files and add to global arrays
326
327        tri_l2g  = fid.variables['tri_l2g'][:]
328        node_l2g = fid.variables['node_l2g'][:]
329        tri_full_flag = fid.variables['tri_full_flag'][:]
330        volumes = num.array(fid.variables['volumes'][:],dtype=num.int)
331        l_volumes = num.zeros_like(volumes)
332        l_old_volumes = num.zeros_like(volumes)
333
334
335        # Change the local node ids to global id in the
336        # volume array
337
338        # FIXME SR: Surely we can knock up a numpy way of doing this
339        #for i in range(len(l_volumes)):
340        #    g_n0 = node_l2g[volumes[i,0]]
341        #    g_n1 = node_l2g[volumes[i,1]]
342        #    g_n2 = node_l2g[volumes[i,2]]
343        #
344        #    l_old_volumes[i,:] = [g_n0,g_n1,g_n2]
345
346        g_n0 = node_l2g[volumes[:,0]].reshape(-1,1)
347        g_n1 = node_l2g[volumes[:,1]].reshape(-1,1)
348        g_n2 = node_l2g[volumes[:,2]].reshape(-1,1)
349
350        #print g_n0.shape
351        l_volumes = num.hstack((g_n0,g_n1,g_n2))
352
353        #assert num.allclose(l_volumes, l_old_volumes)
354
355        # Just pick out the full triangles
356        ftri_ids = num.where(tri_full_flag>0)
357        ftri_l2g = num.compress(tri_full_flag, tri_l2g)
358
359        #print l_volumes
360        #print tri_full_flag
361        #print tri_l2g
362        #print ftri_l2g
363
364        fg_volumes = num.compress(tri_full_flag,l_volumes,axis=0)
365        g_volumes[ftri_l2g] = fg_volumes
366
367
368
369
370        #g_x[node_l2g] = fid.variables['x']
371        #g_y[node_l2g] = fid.variables['y']
372
373        g_points[node_l2g,0] = fid.variables['x']
374        g_points[node_l2g,1] = fid.variables['y']
375       
376
377        #print number_of_timesteps
378
379
380        # FIXME SR: It seems that some of the "ghost" node quantity values
381        # are being storded. We should only store those nodes which are associated with
382        # full triangles. So we need an index array of "full" nodes, ie those in
383        # full triangles
384
385        #use numpy.compress and numpy.unique to get "full nodes
386
387        f_volumes = num.compress(tri_full_flag,volumes,axis=0)
388        fl_nodes = num.unique(f_volumes)
389        f_node_l2g = node_l2g[fl_nodes]
390
391        #print len(node_l2g)
392        #print len(fl_nodes)
393
394        # Read in static quantities
395        for quantity in static_quantities:
396            #out_s_quantities[quantity][node_l2g] = \
397            #             num.array(fid.variables[quantity],dtype=num.float32)
398            q = fid.variables[quantity]
399            #print quantity, q.shape
400            out_s_quantities[quantity][f_node_l2g] = \
401                         num.array(q,dtype=num.float32)[fl_nodes]
402
403       
404        #Collate all dynamic quantities according to their timestep
405        for quantity in dynamic_quantities:
406            q = fid.variables[quantity]
407            #print q.shape
408            for i in range(n_steps):
409                #out_d_quantities[quantity][i][node_l2g] = \
410                #           num.array(q[i],dtype=num.float32)
411                out_d_quantities[quantity][i][f_node_l2g] = \
412                           num.array(q[i],dtype=num.float32)[fl_nodes]
413
414
415        # Read in static c quantities
416        for quantity in static_c_quantities:
417            #out_s_quantities[quantity][node_l2g] = \
418            #             num.array(fid.variables[quantity],dtype=num.float32)
419            q = fid.variables[quantity]
420            out_s_c_quantities[quantity][ftri_l2g] = \
421                         num.array(q,dtype=num.float32)[ftri_ids]
422
423       
424        #Collate all dynamic c quantities according to their timestep
425        for quantity in dynamic_c_quantities:
426            q = fid.variables[quantity]
427            #print q.shape
428            for i in range(n_steps):
429                out_d_c_quantities[quantity][i][ftri_l2g] = \
430                           num.array(q[i],dtype=num.float32)[ftri_ids]
431
432
433        fid.close()
434
435
436    #---------------------------
437    # Write out the SWW file
438    #---------------------------
439    #print g_points.shape
440
441    #print number_of_global_triangles
442    #print number_of_global_nodes
443
444
445    if verbose:
446            print 'Writing file ', output, ':'
447    fido = NetCDFFile(output, netcdf_mode_w)
448
449    sww = Write_sww(static_quantities, dynamic_quantities, static_c_quantities, dynamic_c_quantities)
450    sww.store_header(fido, starttime,
451                             number_of_global_triangles,
452                             number_of_global_nodes,
453                             description=description,
454                             sww_precision=netcdf_float32)
455
456
457
458    from anuga.coordinate_transforms.geo_reference import Geo_reference
459    geo_reference = Geo_reference()
460   
461    sww.store_triangulation(fido, g_points, g_volumes, points_georeference=geo_reference)
462
463    fido.order      = order
464    fido.xllcorner  = xllcorner;
465    fido.yllcorner  = yllcorner ;
466    fido.zone       = zone;
467    fido.false_easting  = false_easting;
468    fido.false_northing = false_northing;
469    fido.datum      = datum;
470    fido.projection = projection;
471       
472    sww.store_static_quantities(fido, verbose=verbose, **out_s_quantities)
473    sww.store_static_quantities_centroid(fido, verbose=verbose, **out_s_c_quantities)
474
475    # Write out all the dynamic quantities for each timestep
476
477    for i in range(n_steps):
478        fido.variables['time'][i] = times[i]
479       
480    for q in dynamic_quantities:
481        q_values = out_d_quantities[q]
482        for i in range(n_steps):
483            fido.variables[q][i] = q_values[i]
484       
485        # This updates the _range values
486        q_range = fido.variables[q + Write_sww.RANGE][:]
487        q_values_min = num.min(q_values)
488        if q_values_min < q_range[0]:
489            fido.variables[q + Write_sww.RANGE][0] = q_values_min
490        q_values_max = num.max(q_values)
491        if q_values_max > q_range[1]:
492            fido.variables[q + Write_sww.RANGE][1] = q_values_max       
493
494    for q in dynamic_c_quantities:
495        q_values = out_d_c_quantities[q]
496        for i in range(n_steps):
497            fido.variables[q][i] = q_values[i]
498
499                                       
500    #print out_s_quantities
501    #print out_d_quantities
502   
503    #print g_x
504    #print g_y
505
506    #print g_volumes
507   
508    fido.close()
509   
510    if delete_old:
511        import os
512        for filename in swwfiles:
513
514            if verbose:
515                print 'Deleting file ', filename, ':'
516            os.remove(filename)
517
518def _sww_merge_parallel_non_smooth(swwfiles, output,  verbose=False, delete_old=False):
519    """
520        Merge a list of sww files into a single file.
521
522        Used to merge files created by parallel runs.
523
524        The sww files to be merged must have exactly the same timesteps.
525
526        It is assumed that the separate sww files have been stored in non_smooth
527        format.
528
529        Note that some advanced information and custom quantities may not be
530        exported.
531
532        swwfiles is a list of .sww files to merge.
533        output is the output filename, including .sww extension.
534        verbose True to log output information
535    """
536
537    if verbose:
538        print "MERGING SWW Files"
539
540
541    first_file = True
542    tri_offset = 0
543    for filename in swwfiles:
544        if verbose:
545            print 'Reading file ', filename, ':'
546
547        fid = NetCDFFile(filename, netcdf_mode_r)
548
549        if first_file:
550
551            times    = fid.variables['time'][:]
552            n_steps = len(times)
553            number_of_timesteps = fid.dimensions['number_of_timesteps']
554            #print n_steps, number_of_timesteps
555            starttime = int(fid.starttime)
556
557            out_s_quantities = {}
558            out_d_quantities = {}
559
560
561            xllcorner = fid.xllcorner
562            yllcorner = fid.yllcorner
563
564            number_of_global_triangles = int(fid.number_of_global_triangles)
565            number_of_global_nodes     = int(fid.number_of_global_nodes)
566            number_of_global_triangle_vertices = 3*number_of_global_triangles
567
568
569            order      = fid.order
570            xllcorner  = fid.xllcorner;
571            yllcorner  = fid.yllcorner ;
572            zone       = fid.zone;
573            false_easting  = fid.false_easting;
574            false_northing = fid.false_northing;
575            datum      = fid.datum;
576            projection = fid.projection;
577
578            g_volumes = num.arange(number_of_global_triangles*3).reshape(-1,3)
579
580
581
582            g_x = num.zeros((number_of_global_triangle_vertices,),num.float32)
583            g_y = num.zeros((number_of_global_triangle_vertices,),num.float32)
584
585            g_points = num.zeros((number_of_global_triangle_vertices,2),num.float32)
586
587
588            quantities = set(['elevation', 'friction', 'stage', 'xmomentum',
589                              'ymomentum', 'xvelocity', 'yvelocity', 'height'])
590            variables = set(fid.variables.keys())
591
592            quantities = list(quantities & variables)
593
594            static_quantities = []
595            dynamic_quantities = []
596
597            for quantity in quantities:
598                # Test if elevation is static
599                if n_steps == fid.variables[quantity].shape[0]:
600                    dynamic_quantities.append(quantity)
601                else:
602                    static_quantities.append(quantity)
603
604            for quantity in static_quantities:
605                out_s_quantities[quantity] = num.zeros((3*number_of_global_triangles,),num.float32)
606
607            # Quantities are stored as a 2D array of timesteps x data.
608            for quantity in dynamic_quantities:
609                out_d_quantities[quantity] = \
610                      num.zeros((n_steps,3*number_of_global_triangles),num.float32)
611
612            description = 'merged:' + getattr(fid, 'description')
613            first_file = False
614
615
616        # Read in from files and add to global arrays
617
618        tri_l2g  = fid.variables['tri_l2g'][:]
619        node_l2g = fid.variables['node_l2g'][:]
620        tri_full_flag = fid.variables['tri_full_flag'][:]
621
622
623
624
625        f_ids = num.argwhere(tri_full_flag==1).reshape(-1,)
626
627        f_gids = tri_l2g[f_ids]
628
629
630        g_vids = (3*f_gids.reshape(-1,1) + num.array([0,1,2])).reshape(-1,)
631        l_vids = (3*f_ids.reshape(-1,1) + num.array([0,1,2])).reshape(-1,)
632
633
634        l_x = num.array(fid.variables['x'][:],dtype=num.float32)
635        l_y = num.array(fid.variables['y'][:],dtype=num.float32)
636
637       
638        g_x[g_vids] = l_x[l_vids]
639        g_y[g_vids] = l_y[l_vids]
640
641        g_points[g_vids,0] = g_x[g_vids]
642        g_points[g_vids,1] = g_y[g_vids]
643
644
645        # Read in static quantities
646        for quantity in static_quantities:
647            #out_s_quantities[quantity][node_l2g] = \
648            #             num.array(fid.variables[quantity],dtype=num.float32)
649            q = fid.variables[quantity]
650            #print quantity, q.shape
651            out_s_quantities[quantity][g_vids] = \
652                         num.array(q,dtype=num.float32)[l_vids]
653
654
655        #Collate all dynamic quantities according to their timestep
656        for quantity in dynamic_quantities:
657            q = fid.variables[quantity]
658            #print q.shape
659            for i in range(n_steps):
660                #out_d_quantities[quantity][i][node_l2g] = \
661                #           num.array(q[i],dtype=num.float32)
662                out_d_quantities[quantity][i][g_vids] = \
663                           num.array(q[i],dtype=num.float32)[l_vids]
664
665
666        fid.close()
667
668    #---------------------------
669    # Write out the SWW file
670    #---------------------------
671
672    if verbose:
673            print 'Writing file ', output, ':'
674
675    fido = NetCDFFile(output, netcdf_mode_w)
676    sww = Write_sww(static_quantities, dynamic_quantities)
677    sww.store_header(fido, starttime,
678                             number_of_global_triangles,
679                             number_of_global_triangles*3,
680                             description=description,
681                             sww_precision=netcdf_float32)
682
683
684    from anuga.coordinate_transforms.geo_reference import Geo_reference
685    geo_reference = Geo_reference()
686
687    sww.store_triangulation(fido, g_points, g_volumes, points_georeference=geo_reference)
688
689    fido.order      = order
690    fido.xllcorner  = xllcorner;
691    fido.yllcorner  = yllcorner ;
692    fido.zone       = zone;
693    fido.false_easting  = false_easting;
694    fido.false_northing = false_northing;
695    fido.datum      = datum;
696    fido.projection = projection;
697
698    sww.store_static_quantities(fido, verbose=verbose, **out_s_quantities)
699
700    # Write out all the dynamic quantities for each timestep
701
702    for i in range(n_steps):
703        fido.variables['time'][i] = times[i]
704
705    for q in dynamic_quantities:
706        q_values = out_d_quantities[q]
707        for i in range(n_steps):
708            fido.variables[q][i] = q_values[i]
709
710        # This updates the _range values
711        q_range = fido.variables[q + Write_sww.RANGE][:]
712        q_values_min = num.min(q_values)
713        if q_values_min < q_range[0]:
714            fido.variables[q + Write_sww.RANGE][0] = q_values_min
715        q_values_max = num.max(q_values)
716        if q_values_max > q_range[1]:
717            fido.variables[q + Write_sww.RANGE][1] = q_values_max
718
719
720    fido.close()
721
722    if delete_old:
723        import os
724        for filename in swwfiles:
725
726            if verbose:
727                print 'Deleting file ', filename, ':'
728            os.remove(filename)
729
730
731
732
733
734
735if __name__ == "__main__":
736
737    import argparse
738    from anuga.anuga_exceptions import ANUGAError
739
740
741    parser = argparse.ArgumentParser(description='Merge sww files created from parallel run')
742    parser.add_argument('-np', type=int, default = 4,
743                   help='number of processors used to produce sww files')
744    parser.add_argument('-f', type=str, default="domain",
745                   help='domain global name')
746    parser.add_argument('-v', nargs='?', type=bool, const=True, default=False,
747                   help='verbosity')
748    parser.add_argument('-delete_old', nargs='?', type=bool, const=True, default=False,
749                   help='Flag to delete the input files')
750    args = parser.parse_args()
751
752    np = args.np
753    domain_global_name = args.f
754    verbose = args.v
755    delete_old = args.delete_old
756
757
758    try:
759        sww_merge_parallel(domain_global_name, np, verbose, delete_old)
760    except:
761        msg = 'ERROR: When merging sww files %s '% domain_global_name
762        print msg
763        raise
Note: See TracBrowser for help on using the repository browser.