COOLFluiD  Release kernel
COOLFluiD is a Collaborative Simulation Environment (CSE) focused on complex MultiPhysics simulations.
atest-ufem-demo-chorin.py
Go to the documentation of this file.
1 import coolfluid as cf
2 import xml.etree.ElementTree as ET
3 
4 # Global configuration
5 cf.env.log_level = 4
6 cf.env.assertion_throws = False
7 cf.env.assertion_backtrace = False
8 cf.env.exception_backtrace = False
9 cf.env.exception_outputs = False
10 cf.env.regist_signal_handlers = False
11 
12 n = 10
13 
14 # Print out the problem size to the testing dashboard
15 measurement = ET.Element('DartMeasurement', name = 'Problem size', type = 'numeric/integer')
16 measurement.text = str(n)
17 print ET.tostring(measurement)
18 
19 # Setup a model
20 model = cf.Core.root().create_component('Model', 'cf3.solver.ModelUnsteady')
21 # The domain holds the mesh
22 domain = model.create_domain()
23 # Physical model keeps track of all the variables in the problem and any material properties (absent here)
24 physics = model.create_physics('cf3.UFEM.NavierStokesPhysics')
25 physics.density = 1.
26 physics.dynamic_viscosity = 1.
27 # Manager for finite element solvers
28 solver = model.create_solver('cf3.UFEM.Solver')
29 ns_solver = solver.add_unsteady_solver('cf3.UFEM.demo.NavierStokesChorin')
30 
31 # Generate a 2D channel mesh
32 mesh = domain.create_component('Mesh', 'cf3.mesh.Mesh')
33 mesh_generator = domain.create_component("MeshGenerator","cf3.mesh.SimpleMeshGenerator")
34 mesh_generator.mesh = mesh.uri()
35 mesh_generator.nb_cells = [n,n]
36 mesh_generator.lengths = [10.,2.]
37 mesh_generator.offsets = [0.,0.]
38 mesh_generator.execute()
39 
40 # Triangulate it
41 triangulator = domain.create_component('triangulator', 'cf3.mesh.MeshTriangulator')
42 triangulator.mesh = mesh
43 triangulator.execute()
44 
45 # Set the region over which the solver operates
46 ns_solver.regions = [mesh.topology.uri()]
47 
48 # Initial conditions
49 ic_u = solver.InitialConditions.create_initial_condition(builder_name = 'cf3.UFEM.InitialConditionFunction', field_tag = 'navier_stokes_u_solution')
50 ic_u.variable_name = 'Velocity'
51 ic_u.options.field_space_name = 'cf3.mesh.LagrangeP2'
52 ic_u.value = ['0', '0']
53 ic_u.regions = ns_solver.regions
54 
55 ic_p = solver.InitialConditions.create_initial_condition(builder_name = 'cf3.UFEM.InitialConditionFunction', field_tag = 'navier_stokes_p_solution')
56 ic_p.variable_name = 'Pressure'
57 ic_p.value = ['20-2*x']
58 ic_p.regions = ns_solver.regions
59 
60 # Boundary conditions
61 ns_solver.VelocityBC.regions = ns_solver.regions
62 ns_solver.PressureBC.regions = ns_solver.regions
63 # Set no-slip at the walls
64 bc_wall = ns_solver.VelocityBC.create_bc_action(region_name = 'top', builder_name = 'cf3.UFEM.BCDirichletFunction')
65 bc_wall.solving_for_difference = False
66 bc_wall.variable_name = 'Velocity'
67 bc_wall.field_tag = 'navier_stokes_u_solution'
68 bc_wall.space = 'cf3.mesh.LagrangeP2'
69 bc_wall.value = ['0', '0']
70 bc_wall.regions = [mesh.topology.top.uri(), mesh.topology.bottom.uri()]
71 # Imposed pressure difference
72 bc_p_out = ns_solver.PressureBC.create_bc_action(region_name = 'right', builder_name = 'cf3.UFEM.BCDirichletFunction')
73 bc_p_out.solving_for_difference = False
74 bc_p_out.variable_name = 'Pressure'
75 bc_p_out.field_tag = 'navier_stokes_p_solution'
76 bc_p_out.value = ['0']
77 bc_p_in = ns_solver.PressureBC.create_bc_action(region_name = 'left', builder_name = 'cf3.UFEM.BCDirichletFunction')
78 bc_p_in.solving_for_difference = False
79 bc_p_in.variable_name = 'Pressure'
80 bc_p_in.field_tag = 'navier_stokes_p_solution'
81 bc_p_in.value = ['20']
82 
83 # Create the fields, to ensure LSS creation
84 solver.create_fields()
85 # LSS options
86 lss = ns_solver.children.AuxiliaryLSS.create_lss()
87 lss.SolutionStrategy.Parameters.preconditioner_type = 'ML'
88 lss.SolutionStrategy.Parameters.PreconditionerTypes.ML.MLSettings.default_values = 'SA'
89 lss.SolutionStrategy.Parameters.PreconditionerTypes.ML.MLSettings.aggregation_type = 'MIS'
90 lss.SolutionStrategy.Parameters.PreconditionerTypes.ML.MLSettings.smoother_type = 'symmetric block Gauss-Seidel'
91 lss.SolutionStrategy.Parameters.PreconditionerTypes.ML.MLSettings.smoother_sweeps = 2
92 lss.SolutionStrategy.Parameters.PreconditionerTypes.ML.MLSettings.smoother_pre_or_post = 'both'
93 
94 # Timestepping
95 time = model.create_time()
96 time.time_step = 0.1
97 time.end_time = 10.*time.time_step
98 
99 # Run the simulation
100 model.simulate()
101 
102 # Print timings
103 model.print_timing_tree()
104 
105 # Save result
106 domain.write_mesh(cf.URI('chorin.msh'))
boost::python::object create_component(ComponentWrapper &self, const std::string &name, const std::string &builder_name)
Send comments to:
COOLFluiD Web Admin