6.4  Superelement import#

6.4.1  Generic representations in SDT#

For interfacing with external finite element software the well documented superelement formalism is used. This formalism is largely used by Multibody Dynamic Software (Simpack, Adams, Excite, ...) and thus widely documented.

A superelement representation of the model is of the form

 
    (6.97)

In general, the reduction is performed so that the DOFs retained {qR} are related to the original DOFs of a larger model by a Rayleigh Ritz reduction basis T using

 
    (6.98)

This representation is fairly standard. The data structure representation within SDT is described in section 7.6. SDT/FEMLink supports import from various FEM codes and more details are given in section 6.4.2 for NASTRAN and section 6.4.4 for Abaqus.

For Craig-Bampton type reduction (enforced displacement on a set of interface DOF I and fixed interface modes on other DOF C), the reduction basis has the form

 
    (6.99)

which verifies the constraints on basis columns that TII=I and TIC=0.

For the free mode variant (McNeal) (see sdtweb('fe_reduc','free')), the form is

 
    (6.100)

with no TI columns.

The observation formalism of SDT which is applicable to both test/analysis sensors (see sdtweb('sensor')) and strain observations used for non-linearities .

 
    (6.101)

Standard notions that have an equivalent in other code are

  • qI interface DOF are directly comparable in SDT and other software when using a DofSet entry that contains an identity enforced motion matrix on a set of DOF idof, typically known as
    • MASTER degree of freedom. That can be initialized with a call of the form model=fe_case(model, 'DofSet', 'In',struct('DOF',idof, 'def', speye(length(idof)))).
    • NASTRAN calls this the BSET, Abaqus uses *RETAINED NODAL DOF cards, ANSYS uses M commands.
  • TC columns may correspond to generalized or internal DOFs. External software will often require that a DOF number be associated to these interfaces. NASTRAN uses a QSET. Abaqus uses xxx. ANSYS uses additional node numbers, with no clearly defined rule, which may cause issues when coupling multiple superelements. xxx
  • bres the independent vectors used to generate residual loads are found as columns of a DofLoad entry.
    • For point loads on a set of DOF adof the entry would be struct('DOF',adof, 'def',speye(length(adof))).
    • NASTRAN allows uses of static loads, point loads defined using USET,U6 entries, viscous damping forces, ...
  • [c] observation matrices do not exist in other environments. The strategy usually retained for export is to add additional nodes of a set ObsNode to correspond to the observations of interest and define the observation equation using a multiple point constraint. This is achieved by modifying the model using fe_sens('MeshSensAsMPC'.
  • rigid assuming that the 3D motion of all nodes in the set depend rigidly on the center node motion. Can be used to define motion of sensor nodes. This tends to over-stiffen the area of connected nodes.
  • rbe3 assuming that the center node moves as the mean motion of the set of nodes is often considered to observe motion of nodes. The case dependent observations of SDT are more general, but may correspond to rbe3.

Implementations of these are discussed for Abaqus section 6.4.5, NASTRAN section 6.4.2, ANSYS section 6.4.6.

6.4.2  NASTRAN Craig-Bampton example#

For a demo, see d_cms('Tuto'). The superelement generation by NASTRAN is saved to an .op2 file that is automatically transformed to the SDT superelement format by FEMLink. A sample file is given in ubeamse.dat.

ASSIGN OUTPUT2='./ubeam_se.op2',UNIT=30
$
ID DFR
SOL 101
GEOMCHECK NONE
TIME 100
$
CEND
TITLE=Generic computation of mode shapes
METHOD=1
DISP(PLOT) = ALL
SPCFORCES(PLOT)=ALL
MPCFORCES(PLOT)=ALL
$ Now extract stresses on base
SET 101=1 THRU 16
STRESS(SORT)=101
$
MPC=1
SPC=1
$ RESVEC(NOINRL)= YES
EXTSEOUT(ASMBULK,EXTBULK,EXTID=100,DMIGOP2=30)
PARAM,POST,-2
PARAM,BAILOUT,-1
$
BEGIN BULK   
$EIGRL,SID,V1,V2,ND,MSGLV,MAXSET,SHFSCL,NORM
EIGRL,1,,,20
$ DOF and nodes to support modal DOF
QSET1,0,1000001,THRU,1000050 
SPOINT,1000001,THRU,1000050
$ Master DOF 4 base corners
BSET1,123,1,5,8,12
$ 
$ Residual on 3 DOF of input node 104
USET,U6,104,123
$
$ Residual associated with CAMP1 vector
CDAMP1  161     2       114     1       244     1
PDAMP*  2               1.              
$
include 'ubeam_include.bdf'
ENDDATA

The resulting basis has the following form

 
    (6.102)

The op2 file contains nodes and superelement definition. It is advised to read the bulk file to obtain a model containing elements and material properties.

  • The interface DOF are defined in NASTRAN usingBset cards. These are stored in SDT as a DofSet entry to the model.
  • QSET correspond to modal/generalized DOFs. These require the definition of a QSET card (to declare existing DOFs), SPOINT grids (to have node numbers to support these QSET DOFs). Note also that the SPOINT numbers should be distinct from other NodeId. The number of modes defined in the EIRGL card should be lower than the the number of SPOINT and the QSET card.
  • Fixed interface modes φc are computed by specifying the EIRGL card.
  • Residual loads bres are defined as follows and lead to additional shapes using the residual vector procedures of NASTRAN
    • point loads simply declared using the USET,U6 card
    • relative loads simply obtained by declaring a CDAMP element that generates a relative viscous load between its two nodes.

6.4.3  Free mode using NASTRAN#

When using a free mode computation, NASTRAN provides mechanisms to compute residual vectors, you should just insert the RESVEC=YES card. An example is given in the ubeamfr.dat file. The resulting basis has the following form

 
    (6.103)

The main mechanisms to generate residual vectors are

  • Free interface modes φc are computed by specifying the EIRGL card.
  • Residual loads bres are defined as follows and lead to additional shapes using the residual vector procedures of NASTRAN
    • point loads simply declared using the USET,U6,NodeId,DofList card. Alternatively RVDOF (MSC but possibly not NX-NASTRAN) can given a list of up to four NodeId,Dof per card.
    • relative loads simply obtained by declaring a CDAMP element that generates a relative viscous load between its two nodes.

6.4.4  Abaqus Craig-Bampton example#

Superelement generation in Abaqus requires at least two consecutive steps, a frequency step *FREQUENCY to compute the reduction basis and a substructure step *SUBSTRUCTURE GENERATE to assemble the model and perform projections. Preloading is possible with preliminary steps. ABAQUS gives the possibility to compute residual vectors with option ,RESIDUAL MODES in the *FREQUENCY card. Depending on the mode resolution strategy residual vectors computation capabilities and input differ.

  • EIGENSOLVER=LANCZOS: this is the standard strategy, that has limitations on very large models. This implementation allows computation of residual vectors from load fields. The associated load vectors must be resolved beforehand in a *STATIC, PERTURBATION step in which *LOAD CASE entries will be defined. A residual vector associated to each load case will then be computed during the frequency step.
  • EIGENSOLVER=AMS: this is the multi-level strategy that is usually required for very large models. This implementation has more limitations as only single DOF residual vectors are computed. No load case needs to be defined, the list of DOF for which residual vectors will be computed is given in the *FREQUENCY card input lines. To overcome the single DOF limitation, one can use MPC with control points to recover a residual vector associated to a distributed force field.

For a Craig-Bampton implementation, one needs to clamp the interface DOF in the *FREQUENCY step with *BOUNDARY and specify them in the *SUBSTRUCTURE GENERATE step with *RETAINED NODAL DOF. There is no specific DOF sets requirements. If the retained nodal DOF are not clamped in the frequency step, one will then get a MacNeal basis, and hybrid formulations with partially clamped sets is also possible.

Substructure output must be done to the RESULTS FILE (.fil extension). A direct import in SDT is then possible. sdtm.nodeLoadSE('FileName.fil') is the code independent entry point. d_cms('TutoAbqExp') illustrates a modeshape expansion case. This is the optimal solution as it is a documented binary file. The use of *INSTANCE for mesh definition is not recommended due to renumbering complexity. In particular, mapping between internal and global numbering is not written in the .fil file. SDT then exploits the .dat file in which the mapping is printed.

See SeGenResidual.inp

6.4.5  Free mode using ABAQUS#

ABAQUS does not natively support generation of substructures with free mode reduction, an interface is always required. The workaround is to use a disconnected node as the interface to recover the free modes reduction. The generation is then exactly similar to the Craig-Bampton version, but with the handling of a spurious node bearing the interface DOF. Implementation is integrated in abaqusJobWrite functionality to ease usage. A free mode reduction job generation is provided in the following example.

% ABAQUS free mode reduction job generation
% demonstration model
mo1=demosdt('demoubeam-noplot');
% Export mesh in input format
finp=fullfile(pwd,'ubeam_mesh.inp'); % mesh file name
abaqus(['write' finp],mo1); % export command

% Prepare job generation % setup context: data storage RT=struct('nmap',vhandle.nmap); % declare mesh file for input handling RT.nmap('MeshFile')=finp; % declare job name and working directory, store RJ=struct('Code','abaqus','Job','ubeam_job','RelDir',pwd); RT.nmap('CurJob')=RJ;

% Declare job: working model, spurious node and steps li={['MeshCfg{MeshFile};'... 'SimuCfg{Inc:SpurNode{opt,-bcstore}:'... 'EIG{opt,-bc}:SEGen{opt,AMS}};']};

% Write Job abaqus('JobWrite',RT,li)

The superelement import with spurious node cleanup can then be obtained with FEMLinkRead

% SE import
RT=struct('In',{{'se','ubeam_job.fil',pwd}});
SE=femlink('Read',RT);

6.4.6  ANSYS Craig-Bampton example#

The cards typically used for superelement generation in ANSYS are

  • antype,substr specify FE substructure generation in after /SOLU
  • cmsopt,fix,Nshapes,,,,,tcms card generates .sub. Use ans2sdt('subSE','file.sub') to import matrices from .sub, restitution from .tcms and mesh from .cdb files using the same file root.
  • resvec,on to use residual vectors in the basis
  • seopt,name,MatType,1 with MatType=2 for mass and stiffness, and 3 for stiffness, mass, viscous damping, name must be defined with card /FILENAME before the analysis.
  • m,NodeId,all master DOF definition repeat card for the various interface nodes. You will have to replace all with UX,,UY,UZ,ROTX,ROTY,ROTZ if the DOF is used by an element that supports multi-physics.
/FILENAME,ubeam_se  ! name must be the one used in SeOpt command
/PREP7
!...
! Use command F to apply loads that will define the residual vector
/SOLU

antype,substr ! substructure analysis

CmsOpt,Fix,20,,,,,TCMS RESVEC, ON SeOpt,ubeamse_ans,3,1,0,, ! 3(all matrices,2 for m and k), 1 to print

! Define list of master DOF, you cannot use ALL if the elements support multi-physics M,1,UX,,,UY,UZ,ROTX,ROTY,ROTZ M,5,UX,,,UY,UZ,ROTX,ROTY,ROTZ M,8,UX,,,UY,UZ,ROTX,ROTY,ROTZ M,12,UX,,,UY,UZ,ROTX,ROTY,ROTZ

SAVE ! save .db file SOLVE ! generate the matrices FINISH