COOLFluiD  Release kernel
COOLFluiD is a Collaborative Simulation Environment (CSE) focused on complex MultiPhysics simulations.
atest-ufem-navier-stokes-laminar-channel-2d.py
Go to the documentation of this file.
1 import sys
2 import coolfluid as cf
3 
4 # Some shortcuts
5 root = cf.Core.root()
6 env = cf.Core.environment()
7 
8 # Global confifuration
9 env.assertion_throws = False
10 env.assertion_backtrace = False
11 env.exception_backtrace = False
12 env.regist_signal_handlers = False
13 env.log_level = 4
14 
15 # setup a model
16 model = root.create_component('NavierStokes', 'cf3.solver.ModelUnsteady')
17 domain = model.create_domain()
18 physics = model.create_physics('cf3.UFEM.NavierStokesPhysics')
19 solver = model.create_solver('cf3.UFEM.Solver')
20 
21 # Add the Navier-Stokes solver as an unsteady solver
22 ns_solver = solver.add_unsteady_solver('cf3.UFEM.NavierStokes')
23 
24 ic_visc = solver.InitialConditions.create_initial_condition(builder_name = 'cf3.UFEM.InitialConditionFunction', field_tag = 'navier_stokes_viscosity')
25 ic_visc.variable_name = 'EffectiveViscosity'
26 
27 # Generate mesh
28 blocks = domain.create_component('blocks', 'cf3.mesh.BlockMesh.BlockArrays')
29 points = blocks.create_points(dimensions = 2, nb_points = 6)
30 points[0] = [0, 0.]
31 points[1] = [10., 0.]
32 points[2] = [0., 0.5]
33 points[3] = [10., 0.5]
34 points[4] = [0.,1.]
35 points[5] = [10., 1.]
36 
37 block_nodes = blocks.create_blocks(2)
38 block_nodes[0] = [0, 1, 3, 2]
39 block_nodes[1] = [2, 3, 5, 4]
40 
41 block_subdivs = blocks.create_block_subdivisions()
42 block_subdivs[0] = [40, 20]
43 block_subdivs[1] = [40, 20]
44 
45 gradings = blocks.create_block_gradings()
46 gradings[0] = [1., 1., 1., 1.]
47 gradings[1] = [1., 1., 1., 1.]
48 
49 left_patch = blocks.create_patch_nb_faces(name = 'left', nb_faces = 2)
50 left_patch[0] = [2, 0]
51 left_patch[1] = [4, 2]
52 
53 bottom_patch = blocks.create_patch_nb_faces(name = 'bottom', nb_faces = 1)
54 bottom_patch[0] = [0, 1]
55 
56 top_patch = blocks.create_patch_nb_faces(name = 'top', nb_faces = 1)
57 top_patch[0] = [5, 4]
58 
59 right_patch = blocks.create_patch_nb_faces(name = 'right', nb_faces = 2)
60 right_patch[0] = [1, 3]
61 right_patch[1] = [3, 5]
62 
63 blocks.partition_blocks(nb_partitions = cf.Core.nb_procs(), direction = 0)
64 
65 mesh = domain.create_component('Mesh', 'cf3.mesh.Mesh')
66 blocks.create_mesh(mesh.uri())
67 
68 ns_solver.regions = [mesh.topology.uri()]
69 
70 u_in = [2., 0.]
71 
72 solver.create_fields()
73 
74 #initial condition for the velocity. Unset variables (i.e. the pressure) default to zero
75 solver.InitialConditions.navier_stokes_solution.Velocity = u_in
76 ic_visc.value = ['10. + 2*sin(2/pi*x)']
77 ic_visc.regions = [mesh.topology.uri()]
78 ic_visc.execute()
79 domain.write_mesh(cf.URI('laminar-channel-2d_output-init.pvtu'))
80 
81 # Physical constants
82 physics.options().set('density', 1000.)
83 physics.options().set('dynamic_viscosity', 10.)
84 
85 # Boundary conditions
86 bc = ns_solver.BoundaryConditions
87 bc.add_constant_bc(region_name = 'left', variable_name = 'Velocity').value = u_in
88 bc.add_constant_bc(region_name = 'bottom', variable_name = 'Velocity').value = [0., 0.]
89 bc.add_constant_bc(region_name = 'top', variable_name = 'Velocity').value = [0., 0.]
90 bc.add_constant_bc(region_name = 'right', variable_name = 'Pressure').value = 0.
91 
92 pressure_integral = solver.add_unsteady_solver('cf3.UFEM.SurfaceIntegral')
93 pressure_integral.set_field(variable_name = 'Pressure', field_tag = 'navier_stokes_solution')
94 pressure_integral.regions = [mesh.topology.access_component('bottom').uri()]
95 pressure_integral.history = solver.create_component('ForceHistory', 'cf3.solver.History')
96 pressure_integral.history.file = cf.URI('force-implicit.tsv')
97 pressure_integral.history.dimension = 2
98 
99 # Time setup
100 time = model.create_time()
101 time.options().set('time_step', 0.1)
102 
103 #ns_solver.options().set('disabled_actions', ['SolveLSS'])
104 
105 # Setup a time series write
106 final_end_time = 10.
107 save_interval = 1.
108 current_end_time = 0.
109 iteration = 0
110 
111 writer = root.create_component('Writer', 'cf3.mesh.gmsh.Writer')
112 writer.fields = [mesh.geometry.navier_stokes_solution.uri()]
113 writer.enable_overlap = True
114 writer.mesh = mesh
115 
116 while current_end_time < final_end_time:
117  current_end_time += save_interval
118  time.options().set('end_time', current_end_time)
119  model.simulate()
120  writer.file = cf.URI('laminar-channel-2d_output-' +str(iteration) + '.msh')
121  writer.execute()
122  domain.write_mesh(cf.URI('laminar-channel-2d_output-' +str(iteration) + '.pvtu'))
123  iteration += 1
124  if iteration == 1:
125  solver.options().set('disabled_actions', ['InitialConditions'])
126 
127 # print timings
128 model.print_timing_tree()
common::URI uri(ComponentWrapper &self)
Send comments to:
COOLFluiD Web Admin