COOLFluiD  Release kernel
COOLFluiD is a Collaborative Simulation Environment (CSE) focused on complex MultiPhysics simulations.
atest-ufem-navier-stokes-mlaux.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 configuration
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 def run_case(modeltype):
16  implicit = modeltype == 'implicit'
17  # setup a model
18  model = root.create_component('NavierStokes'+modeltype, 'cf3.solver.ModelUnsteady')
19  domain = model.create_domain()
20  physics = model.create_physics('cf3.UFEM.NavierStokesPhysics')
21  solver = model.create_solver('cf3.UFEM.Solver')
22 
23  # Add the Navier-Stokes solver as an unsteady solver
24  if implicit:
25  ns_solver = solver.add_unsteady_solver('cf3.UFEM.NavierStokes')
26  else:
27  ns_solver = solver.add_unsteady_solver('cf3.UFEM.NavierStokesSemiImplicit')
28  ns_solver.options.pressure_rcg_solve = True
29  ns_solver.enable_body_force = True
30 
31  # Generate mesh
32  blocks = domain.create_component('blocks', 'cf3.mesh.BlockMesh.BlockArrays')
33  points = blocks.create_points(dimensions = 2, nb_points = 6)
34  points[0] = [0, 0.]
35  points[1] = [10., 0.]
36  points[2] = [0., 0.5]
37  points[3] = [10., 0.5]
38  points[4] = [0.,1.]
39  points[5] = [10., 1.]
40 
41  block_nodes = blocks.create_blocks(2)
42  block_nodes[0] = [0, 1, 3, 2]
43  block_nodes[1] = [2, 3, 5, 4]
44 
45  block_subdivs = blocks.create_block_subdivisions()
46  block_subdivs[0] = [40, 20]
47  block_subdivs[1] = [40, 20]
48 
49  gradings = blocks.create_block_gradings()
50  gradings[0] = [1., 1., 10., 10.]
51  gradings[1] = [1., 1., 0.1, 0.1]
52 
53  left_patch = blocks.create_patch_nb_faces(name = 'left', nb_faces = 2)
54  left_patch[0] = [2, 0]
55  left_patch[1] = [4, 2]
56 
57  bottom_patch = blocks.create_patch_nb_faces(name = 'bottom', nb_faces = 1)
58  bottom_patch[0] = [0, 1]
59 
60  top_patch = blocks.create_patch_nb_faces(name = 'top', nb_faces = 1)
61  top_patch[0] = [5, 4]
62 
63  right_patch = blocks.create_patch_nb_faces(name = 'right', nb_faces = 2)
64  right_patch[0] = [1, 3]
65  right_patch[1] = [3, 5]
66 
67  blocks.partition_blocks(nb_partitions = cf.Core.nb_procs(), direction = 1)
68 
69  mesh = domain.create_component('Mesh', 'cf3.mesh.Mesh')
70  blocks.create_mesh(mesh.uri())
71 
72  if not implicit:
73  link_horizontal = domain.create_component('LinkHorizontal', 'cf3.mesh.actions.LinkPeriodicNodes')
74  link_horizontal.mesh = mesh
75  link_horizontal.source_region = mesh.topology.right
76  link_horizontal.destination_region = mesh.topology.left
77  link_horizontal.translation_vector = [-10., 0.]
78  link_horizontal.execute()
79 
80  ns_solver.regions = [mesh.topology.uri()]
81 
82  if implicit:
83  lss = ns_solver.LSS
84  lss.SolutionStrategy.Parameters.preconditioner_type = 'ML'
85  lss.SolutionStrategy.Parameters.PreconditionerTypes.ML.MLSettings.add_parameter(name = 'aggregation: aux: enable', value = True)
86  lss.SolutionStrategy.Parameters.PreconditionerTypes.ML.MLSettings.add_parameter(name = 'aggregation: aux: threshold', value = 0.0001)
87  else:
88  for strat in [ns_solver.children.FirstPressureStrategy, ns_solver.children.SecondPressureStrategy]:
89  strat.MLParameters.aggregation_type = 'Uncoupled'
90  strat.MLParameters.max_levels = 4
91  strat.MLParameters.smoother_sweeps = 3
92  strat.MLParameters.coarse_type = 'Amesos-KLU'
93  strat.MLParameters.add_parameter(name = 'repartition: start level', value = 0)
94  strat.MLParameters.add_parameter(name = 'repartition: max min ratio', value = 1.1)
95  strat.MLParameters.add_parameter(name = 'aggregation: aux: enable', value = True)
96  strat.MLParameters.add_parameter(name = 'aggregation: aux: threshold', value = 0.05)
97  lss = ns_solver.VelocityLSS.LSS
98  lss.SolutionStrategy.Parameters.LinearSolverTypes.Belos.solver_type = 'Block GMRES'
99  lss.SolutionStrategy.Parameters.LinearSolverTypes.Belos.SolverTypes.BlockGMRES.convergence_tolerance = 1e-6
100  lss.SolutionStrategy.Parameters.LinearSolverTypes.Belos.SolverTypes.BlockGMRES.maximum_iterations = 300
101  lss.SolutionStrategy.Parameters.LinearSolverTypes.Belos.SolverTypes.BlockGMRES.num_blocks = 100
102 
103  u_in = [1., 0.]
104 
105  solver.create_fields()
106 
107  #initial condition for the velocity. Unset variables (i.e. the pressure) default to zero
108  if implicit:
109  solver.InitialConditions.navier_stokes_solution.Velocity = u_in
110  else:
111  ic_u = solver.InitialConditions.NavierStokes.create_initial_condition(builder_name = 'cf3.UFEM.InitialConditionFunction', field_tag = 'navier_stokes_u_solution')
112  ic_u.variable_name = 'Velocity'
113  ic_u.regions = [mesh.topology.uri()]
114  ic_u.value = ['0', '0']
115  ic_g = solver.InitialConditions.NavierStokes.create_initial_condition(builder_name = 'cf3.UFEM.InitialConditionFunction', field_tag = 'body_force')
116  ic_g.variable_name = 'Force'
117  ic_g.regions = [mesh.topology.uri()]
118  ic_g.value = ['2', '0']
119 
120  # Physical constants
121  physics.options().set('density', 1000.)
122  physics.options().set('dynamic_viscosity', 10.)
123 
124  # Boundary conditions
125  if implicit:
126  bc_u = ns_solver.BoundaryConditions
127  bc_p = bc_u
128  bc_u.add_constant_bc(region_name = 'left', variable_name = 'Velocity').value = u_in
129  bc_p.add_constant_bc(region_name = 'right', variable_name = 'Pressure').value = 0.
130  else:
131  bc_u = ns_solver.VelocityLSS.BC
132  bc_p = ns_solver.PressureLSS.BC
133 
134  bc_u.add_constant_bc(region_name = 'bottom', variable_name = 'Velocity').value = [0., 0.]
135  bc_u.add_constant_bc(region_name = 'top', variable_name = 'Velocity').value = [0., 0.]
136 
137  # Time setup
138  time = model.create_time()
139  time.time_step = 0.1
140  time.end_time = 2.*time.time_step
141 
142  model.simulate()
143  domain.write_mesh(cf.URI('mlaux-{type}.pvtu'.format(type=modeltype)))
144 
145  model.print_timing_tree()
146 
147 run_case('implicit')
148 run_case('semi')
Send comments to:
COOLFluiD Web Admin