From a59e96105ecf2882b6e7c05d9b39970c11ee3b9b Mon Sep 17 00:00:00 2001 From: krande Date: Mon, 5 Jul 2021 16:48:25 +0200 Subject: [PATCH] Minor bugfixes in code aster eigenvalue analysis step and related tests --- images/tests/test_code_aster_fem.py | 8 +++++++- src/ada/fem/io/code_aster/writer.py | 24 ++++++++++++++++++------ 2 files changed, 25 insertions(+), 7 deletions(-) diff --git a/images/tests/test_code_aster_fem.py b/images/tests/test_code_aster_fem.py index 5052c36a5..b21851354 100644 --- a/images/tests/test_code_aster_fem.py +++ b/images/tests/test_code_aster_fem.py @@ -36,10 +36,13 @@ def test_bm_eigenfrequency_beam(self): p.gmsh.mesh(0.1) fix_set = p.fem.add_set(FemSet("bc_nodes", get_beam_end_nodes(bm), "nset")) - a.fem.add_bc(Bc("Fixed", fix_set, [1, 2, 3, 4, 5, 6])) + p.fem.add_bc(Bc("Fixed", fix_set, [1, 2, 3, 4, 5, 6])) a.fem.add_step(Step("Eigen", "eigenfrequency")) + a.to_fem("Cantilever_aba_EIG_bm", "abaqus", overwrite=True) + a.to_fem("Cantilever_ses_EIG_bm", "sesam", overwrite=True) + res = a.to_fem("Cantilever_CA_EIG_bm", "code_aster", overwrite=True, execute=True) f = h5py.File(res.results_file_path) @@ -63,6 +66,9 @@ def test_bm_eigenfrequency_shell(self): a.fem.add_step(Step("Eigen", "eigenfrequency")) + a.to_fem("Cantilever_CA_EIG_sh", "abaqus", overwrite=True) + a.to_fem("Cantilever_CA_EIG_sh", "sesam", overwrite=True) + res = a.to_fem("Cantilever_CA_EIG_sh", "code_aster", overwrite=True, execute=True) f = h5py.File(res.results_file_path) diff --git a/src/ada/fem/io/code_aster/writer.py b/src/ada/fem/io/code_aster/writer.py index 70b35abf7..d6bd6dd7e 100644 --- a/src/ada/fem/io/code_aster/writer.py +++ b/src/ada/fem/io/code_aster/writer.py @@ -5,7 +5,7 @@ import numpy as np from ada.config import Settings as _Settings -from ada.fem import ElemShapes, FemSection +from ada.fem import ElemShapes, FemSection, Step from ada.fem.containers import FemSections from ..utils import _folder_prep, get_fem_model_from_assembly @@ -90,7 +90,7 @@ def write_to_comm(assembly, part): materials_str = "\n".join([write_material(mat) for mat in part.materials]) sections_str = write_sections(part.fem.sections) - bc_str = "\n".join([write_boundary_condition(bc) for bc in assembly.fem.bcs]) + bc_str = "\n".join([write_boundary_condition(bc) for bc in assembly.fem.bcs + part.fem.bcs]) step_str = "\n".join([write_step(s, part) for s in assembly.fem.steps]) model_type_str = "" if len(part.fem.sections.beams) > 0: @@ -390,10 +390,18 @@ def step_static_str(step, part): )""" -def step_eig_str(step, part): - if len(part.get_assembly().fem.bcs) != 1: +def step_eig_str(step: Step, part): + if len(part.fem.bcs) == 1: + bcs = part.fem.bcs + elif len(part.get_assembly().fem.bcs) == 1: + bcs = part.get_assembly().fem.bcs + else: raise NotImplementedError("Number of BC sets is for now limited to 1 for eigenfrequency analysis") - bc = part.get_assembly().fem.bcs[0] + + eig_map = dict(sorensen="SORENSEN", lanczos="TRI_DIAG") + eig_type = step.metadata.get("eig_method", "sorensen") + eig_method = eig_map[eig_type] + bc = bcs[0] return f""" #modal analysis ASSEMBLAGE( @@ -407,9 +415,13 @@ def step_eig_str(step, part): _F(MATRICE=CO('mass'), OPTION ='MASS_MECA', ), ), ) +# Using Subspace Iteration method ('SORENSEN' AND 'PLUS_PETITE') +# See https://www.code-aster.org/V2/UPLOAD/DOC/Formations/01-modal-analysis.pdf for more information +# + modes = CALC_MODES( CALC_FREQ=_F(NMAX_FREQ={step.eigenmodes}, ) , - SOLVEUR_MODAL=_F(METHODE='SORENSEN'), + SOLVEUR_MODAL=_F(METHODE='{eig_method}'), MATR_MASS=mass, MATR_RIGI=stiff, OPTION='PLUS_PETITE',