fe_case#
UI function to handle FEM computation cases
Syntax
Case = fe_case(Case,'EntryType','Entry Name',Data) fe_case(model,'command' ...)
Description
FEM computation cases contains information other than nodes and elements used to describe a FEM computation. Currently supported entries in the case stack are
cyclic | (SDT) used to support cyclic symmetry conditions |
DofLoad | loads defined on DOFs (handled by fe_load) |
DofSet | (SDT) imposed displacements on DOFs |
FixDof | used to eliminated DOFs specified by the stack data |
FSurf | surface load defined on element faces (handled by fe_load). This will be phased out since surface load elements associated with volume loads entries are more general. |
FVol | volume loads defined on elements (handled by fe_load) |
info | used to stored non standard entries |
KeepDof | (obsolete) used to eliminated DOFs not specified by the stack data. These entries are less general than FixDof and should be avoided. |
map | field of normals at nodes |
mpc | multiple point constraints |
rbe3 | a flavor of MPC that enforce motion of a node a weighted average |
par | are used to define physical parameters (see upcom Par commands) |
rigid | linear constraints associated with rigid links |
SensDof | (SDT) Sensor definitions |
fe_case is called by the user to initialize (when Case is not provided as first argument) or modify cases (Case is provided).
Accepted commands are
Get, T, Set, Remove, Reset ...#
- [Case,CaseName]=fe_case(model,'GetCase') returns the current case.
GetCasei returns case number i (order in the model stack). GetCaseName returns a case with name Name and creates it if it does not exist. Note that the Case name cannot start with Case. - data=fe_case(model,'GetData EntryName') returns data associated with the case entry EntryName.
- model=fe_case(model,'SetData EntryName',data) sets data associated with the case entry EntryName.
- [Case,NNode,ModelDOF]=fe_case(model,'GetT'); returns a congruent transformation matrix which verifies constraints. Details are given in section 7.14. CaseDof=fe_case(model,'GetTDOF') returns the case DOF (for model DOF use feutil('getdof',model)). If fields Case.T and Case.DOF are already defined, they will be reused. Use command option new to force a reset of these fields.
- model=fe_case(model,'Remove','EntryName') removes the entry with name EntryName.
- Reset empties all information in the case stored in a model structure model = fe_case(model,'reset')
- fe_case SetCurve has a load reference a curve in model Stack. For example model=fe_case(model,'SetCurve','Point load 1','input'); associates Point load 1 to curve input. See section 7.9 for more details on curves format and fe_case SetCurve for details on the input syntax.
- stack_get applies the command to the case rather than the model. For example des = fe_case(model,'stack_get','par')
- stack_set applies the command to the case rather than the model. For example model = fe_case(model,'stack_set','info','Value',1)
- stack_rm applies the command to the case rather than the model. For example model = fe_case(model,'stack_rm','par')
Commands for advanced constraint generation
AutoSPC#
Analyses the rank of the stiffness matrix at each node and generates a fixdof case entry for DOFs found to be singular:
model = fe_case(model,'autospc')
Assemble#
Calls used to assemble the matrices of a model. See fe_mknl Assemble and section 4.10.7 for optimized assembly strategies.
Build Sec epsl d#
model = fe_cyclic('build (N) epsl (d)',model,LeftNodeSelect) is used to append a cyclic constraint entry in the current case.
ConnectionEqualDOF#
fe_caseg('Connection EqualDOF',model,'name',DOF1,DOF2) generates a set of MPC connecting each DOF of the vector DOF1 (slaves) to corresponding DOF in DOF2 (masters). DOF1 and DOF2 can be a list of NodeId, in that case all corresponding DOF are connected, or only DOF given as a -dof DOFs command option.
Following example defines 2 disjointed cubes and connects them with a set of MPC between DOFs along x and y of the given nodes,
% Build a Multiple Point Constraint (MPC) with DOF equalization % Generate a cube model cf=feplot; cf.model=femesh('testhexa8'); % duplicate the cube and translate cf.mdl=feutil('repeatsel 2 0.0 0.0 1.5',cf.mdl); % build the connection cf.mdl=fe_caseg('Connection EqualDOF -id7 -dof 1 2',cf.mdl, ... 'link1',[5:8]',[9:12]'); % display the result in feplot cf.sel='reset'; % reset feplot display % open feplot pro and view the built connection fecom(cf,'promodelviewon');fecom(cf,'curtab Cases','link1');
The option -id i can be added to the command to specify a MPC ID i for export to other software. Silent mode is obtained by adding ; at the end of the command.
By default a DOF input mismatch will generate an error. Command option -safe allows DOF mismatch in the input by applying the constraint only to DOF existing in both lists. If no such DOF exists the constraint is not created.
ConnectionPivot#
This command generates a set of MPC defining a pivot connection between two sets of nodes. It is meant for use with volume or shell models with no common nodes. For beams the pin flags (columns 9:10 of the element row) are typically more appropriate, see beam1for more details.
The command specifies the DOFs constraint at the pivot (in the example DOF 6 is free), the local z direction and the location of the pivot node. One then gives the model, the connection name, and node selections for the two sets of nodes.
% Build a pivot connection between plates model=demosdt('demoTwoPlate'); model=fe_caseg('Connection Pivot 12345 0 0 1 .5 .5 -3 -id 1111', ... model,'pivot','group1','group2'); def=fe_eig(model);feplot(model,def)
The option -id i can be added to the command to specify a MPC ID i for export to other software. Silent mode is obtained by adding ; at the end of the command.
ConnectionSurface#
ConnectionSurface implements node to surface connections trough constraints or elasticity. fe_caseg('ConnectionSurface DOFs',model,'name',NodeSel1,Eltsel2) generates a set of MPC connecting of DOFs of a set of nodes selected by NodeSel1 (this is a node selection string) to a surface selected by EltSel2 (this is an element selection string). ConnectionSurface performs a match between two selections using feutilb Match and exploits the result with feutilb MpcFromMatch.
The following example links x and z translations of two plates
% Build a surface connection between two plates model=demosdt('demoTwoPlate'); model=fe_caseg('Connection surface 13 -MaxDist0.1',model,'surface', ... 'z==0', ... % Selection of nodes to connect 'withnode {z==.1 & y<0.5 & x<0.5}'); % Selection of elements for matching def=fe_eig(model);feplot(model,def)
Accepted command options are
- Auto will run an automated refinement of then provided element selections element selection to locate areas of possible interactions.
- -aTol provides a custom tolerance in Auto mode to detect intersecting volume extensions where the match will be performed. By default one will consider 10 times the mesh characteristic length.
- -id i can be added to the command to specify a MPC ID i for export to other software.
- -Radius val can be used to increase the search radius for the feutilb Match operation.
- -radEstval can be used to exploit a radius based on the average mesh edge length of the elements selected for matching multiplied by val (0.1 to get 10% of the average mesh edge length). This command is exclusive with -Radius, the priority is on -Radius
- -MaxDist val eliminates matched node with distance to the matched point within the element higher than val. This is typically useful for matches on surfaces where the node can often be external. Using a -MaxDist is required for -Dof.
- -kp val is used to give the stiffness (force/length) for a penalty based implementation of the constraint. The stiffness matrix of the penalized bilateral connection is stored in a superelement with the constraint name.
- -KpAutoval is used is -kp is not present to ask for an automated estimation of the penalization stiffness based on mesh size and flange materials. The objective is to get a saturated stiffness not altering numerical conditionning. val is optionnal. It wiil be used as a correction factor to the default computed stiffness. To get 10% of the automated stiffness use 0.1.
- -dens uses a slave surface. In conjunction with -kp the coefficient provided is used as a surface stiffness density. With this option, the first selection must rethrow a face selection.
- -Dof val can be used to build surface connections of non structural DOFs (thermal fields, ...).
- -MatchS uses a surface based matching strategy that may be significantly faster.
- -disjCut will attempt at splitting the generated connection by disjointed connected areas of the surface (second selection), the result is either a series of mpc or a model with multiple SE depending on the mode.
- Silent mode is obtained by adding ; at the end of the command.
It is also possible to define the ConnectionSurface implicitly, to let the constraint resolution be performed after full model assembly. The ConnectionSurface is then defined as an MPC, which data structure features fields .type equal to ConnectionSurface with possible command options, and field .sel giving in a cell array a sequence {NodeSel1, EltSel2}, as defined in the explicit definition. The following example presents the implicit ConnectionSurface definition equivalent to the above explicit one.
% Build a surface connection between two plates % using implicit selections model=demosdt('demoTwoPlate'); model=fe_case(model,'mpc','surface',... struct('type','Connection surface 13 -MaxDist0.1',... 'sel',{{'z==0','withnode {z==.1 & y<0.5 & x<0.5}'}})); def=fe_eig(model);feplot(model,def)
% Build a penalized surface connection % with a given sitffness density between two plates model=demosdt('demoTwoPlate'); model=fe_caseg('Connection surface 123 -MaxDist 0.1 -kp1e8 -dens',model,... 'surface',... 'withnode{z==0}&selface',... 'withnode {z==.1 & y<0.5 & x<0.5}') def=fe_eig(model);cf=feplot(model,def); fecom(cf,'promodelinit'); fecom(cf,'curtabStack','surface'); fecom(cf,'proviewon');
Warning volume matching requires that nodes are within the element. To allow exterior nodes, you should add a & selface at the end of the element selection string for matching.
ConnectionScrew#
fe_caseg('Connection Screw',model,'name',data)
This command generates a set of RBE3 defining a screw connection. Nodes to be connected are defined in planes from their distance to the axis of the screw. The connected nodes define a master set enforcing the motion of a node taken on the axis of the screw with a set of RBE3 (plane type 1) or rigid links (plane type 0) ring for each plane.
In the case where rigid links are defined, the command appends a group of rigid elements to the model case.
Real screws can be represented by beams connecting all the axis slave nodes, this option is activated by adding the field MatProId in the data structure.
data defining the screw is a data structure with following fields:
Origin | a vector [x0 y0 z0] defining the origin of the screw. |
axis | a vector [nx ny nz] defining the direction of the screw axis. |
radius | defines the radius of the screw. |
planes | a matrix with as many lines as link rings. Each row is of the form [z0 type ProId zTol rad stype zTol2] where |
z0 is the plane distance to the origin along the axis of the screw | |
type is the type of link: 0 for rigid and 1 for rbe3 | |
ProId is the ProId of the elements containing nodes to connect. This limits the plane search to the elements of given ProId. By default, a zero value can be used, in which case all elements will be considered for the search | |
zTol is the plane position tolerance, nodes within z0-zTol to z0+zTol will be detected | |
rad is the radius considered for this plane detection, if a zero value is given the base radius is used | |
stype defines the node search type. A value of 0 (default) will use a spherical search of radius rad aorund the origin (only practical for perfectly planar definitions). A value of 1 will use a cylindrical node search along the screw axis from the origin, with symmetric distance from the origin defined by zTol. A value of 2 implements a cylindrical node search with non-symmetric height tolerances from origin, using from zTol to zTol2 | |
zTol2 second side height tolerance for stype=2 (non-symmetric height cylinder based node search) | |
MatProId | Optional. If present beams are added to connect slave nodes at the center of each link ring. It is a vector [MatId ProId] defining the MatId and the ProId of the beams. For new MatId, default material is steel and for new ProId, default beam section is a circle with provided radius. |
MasterCelas | Optional. It defines the celas element which is added if this field is present. It is of the form [0 0 -DofID1 DofID2 ProID EltID Kv Mv Cv Bv]. The first node of the celas is the slave node of the rbe3 ring and the second is added at the same location. This can be useful to reduce a superelement keeping the center of the rings in the interface. |
NewNode | Optional. If it is omitted or equal to 1 then a new slave node is added to the model at the centers of the link rings. If it equals to 0, existent model node can be kept. |
Nnode | Optional. Gives the number of points to retain in each plane. |
For each plane, nodes are searched following the stype strategy. The found nodes are then connected to the center node which is strictly defined at height z0 on the axis provided. The heights provided as z0, zTol and zTol2 must be understood along the axis provided and not as function of the main frame coordinates.
In the case of a rigid connection, nodes detection should be non intersecting to avoid multiple slaves. Overlapping slave node selection is avoided by sequentially eliminating used nodes in the following detections. Selection priority is thus performed following the plane order sequence.
One can also define more generally planes as a cell array whose each row defines a plane and is of the form {z0 type st} where z0 and type are defined above and st is a FindNode string. st can contain $FieldName tokens that will be replaced by corresponding data.FieldName value (for example 'cyl<= $radius o $Origin $axis & inElt{ProId $ProId}' will select nodes in cylinder of radius data.radius, origin data.Origin and axis data.axis, and in elements of ProId data.ProId).
Silent mode is obtained by adding ; at the end of the command.
Following example creates a test model, and adds 2 rbe3 rings in 2 planes.
% Sample connection builds commands for screws using rigid or RBE3 model=demosdt('demoscrew layer 0 40 20 3 3 layer 0 40 20 4'); % create model r1=struct('Origin',[20 10 0],'axis',[0 0 1],'radius',3, ... 'planes',[1.5 1 111 1 3.1; 5.0 1 112 1 4;], ... 'MasterCelas',[0 0 -123456 123456 10 0 1e14], ... 'NewNode',0); model=fe_caseg('ConnectionScrew',model,'screw1',r1); cf=feplot(model); % show model fecom('promodelviewon');fecom('curtab Cases','screw1');% alternative definintion using a beam model=demosdt('demoscrew layer 0 40 20 3 3 layer 0 40 20 4'); % create model r1=struct('Origin',[20 10 0],'axis',[0 0 1],'radius',3, ... 'planes',[1.5 1 111 1 3.1; 5.0 1 112 1 4;], ... 'MasterCelas',[0 0 -123456 123456 10 0 1e14], ... 'MatProId',[110 1001],... 'NewNode',0); model=fe_caseg('ConnectionScrew',model,'screw1',r1); cf=feplot(model); % show model fecom('promodelviewon');fecom('curtab Cases','screw1');
% alternative definition with a load, two beam elements are created model=demosdt('demoscrew layer 0 40 20 3 3 layer 0 40 20 4'); % create model model=fe_caseg('ConnectionScrew -load1e5;',model,'screw1',r1); def=fe_eig(model,[5 15 1e3]);
% alternative definition with a load, two beam elements are created % and a pin flag is added to release the beam compression model=demosdt('demoscrew layer 0 40 20 3 3 layer 0 40 20 4'); % create model model=fe_caseg('ConnectionScrew -load1e5 -pin1;',model,'screw1',r1); def1=fe_eig(model,[5 15 1e3]);
% a new rigid body mode has been added due to the pin flag addition [def.data(7) def1.data(7)]
Command option -loadval allows defining a loading force of amplitude val to the screw in the case where a beam is added to model the screw (through the MatId optional field). To this mean the last beam element (in the order defined by the planes entry) is split in two at a tenth of its length and a compression force is added to the larger element that is exclusively inside the beam. In complement, command option -pinpdof allows defining pin flags with identifiers pdof to the compressed beam1element.
Entries
The following paragraphs list available entries not handled by fe_load or upcom.
cyclic (SDT)#
cyclic entries are used to define sector edges for cyclic symmetry computations. They are generated using the fe_cyclic Build command.
FixDof#
FixDof entries correspond to rows of the Case.Stack cell array giving {'FixDof', Name, Data}. Name is a string identifying the entry. data is a column DOF definition vector (see section 7.10) or a string defining a node selection command. You can also use
data=struct('data',DataStringOrDof,'ID',ID) to specify a identifier.
You can now add DOF and ID specifications to the findnode command. For example 'x==0 -dof 1 2 -ID 101' fixes DOFs x and y on the x==0 plane and generates an data.ID field equal to 101 (for use in other software).
The following command gives syntax examples. An example is given at the end of the fe_case documentation.
% Declare a clamping constraint with fixdof model = fe_case(model,'FixDof','clamped dofs','z==0', ... 'FixDof','SimpleSupport','x==1 & y==1 -DOF 3', ... 'FixDof','DofList',[1.01;2.01;2.02], ... 'FixDof','AllDofAtNode',[5;6], ... 'FixDof','DofAtAllNode',[.05]);
map#
map entries are used to define maps for normals at nodes. These entries are typically used by shell elements or by meshing tools. Data is a structure with fields
- .normal a N by 3 matrix giving the normal at each node or element
- .ID a N by 1 vector giving identifiers. For normals at integration points, element coordinates can be given as two or three additional columns.
- .opt an option vector. opt(1) gives the type of map (1 for normals at element centers, 2 for normals at nodes, 3 normals at integration points specified as additional columns of Data.ID).
- .vertex an optional N by 3 matrix giving the location of each vector specified in .normal. This can be used for plotting.
MPC#
MPC (multiple point constraint) entries are rows of the Case.Stack cell array giving {'MPC', Name, Data}. Name is a string identifying the entry. Data is a structure with fields Data.ID positive integer for identification. Data.c is a sparse matrix whose columns correspond to DOFs in Data.DOF. c is the constraint matrix such that [c] {q} = {0} for q defined on DOF.
Data.slave is an optional vector of slave DOFs in Data.DOF. If the vector does not exist, it is filled by feutil FixMpcMaster.
Note that the current implementation has no provision for using local coordinates in the definition of MPC (they are assumed to be defined using global coordinates).
par (SDT)#
par entries are used to define variable coefficients in element selections. It is nominally used through upcom Par commands but other routines may also use it [38].
High level calls to define model parameters are packaged in fe_caseg Par.
At a lower level, command ParAdd allows quickly defining a parameter:
model=fe_case(model,'ParAddtype nom min max scale',pname,par)
Inputs defines the parameter as
-
type the types are associted to the MatType with tokens that refer to an internal numeric value
- k for stiffness (1).
- m for mass (2).
- c for viscous damping (3.1).
- t for shell thickness (3).
- ki for hyteretic damping (4).
- kg for non-linear geometry stiffness (5).
- 0 for no matrix association (0)
- -2 internal only (-2), to force a superelement matrix not to undergo low-level assembly checks for parametric assemblies (MatType option -1
- nom the parameter nominal value.
- min the parameter minimal value.
- max the parameter maximal value.
- scale the parameter varying scale,
- 1 linear variation
- 2 or lo for logarithmic variation
- more types are available in generic parameters, see fe_range.
- pname defines the parameter name.
- par provides the parameter defintion as structure, string inputs will override the original values. At low level, the following fields are admissible
- .coef a 5 value row vector providing the parameter values as [type nom min max scale], defined above. type is used for assembly. The following values are overriden by fe_rangeif a more advanced definition is used.
- .sel provides the element selection on which the parameter is applied. During assembly, a submodel based on the selection is assembled with the required type to provide the associated matrix.
- .zCoef (optionnal) defines a non standard weighting coefficient rule for the matrix.
- More entries can be added conforming to fe_rangeparameter definition.
RBE3 (SDT)#
rbe3 constraints enforce the motion of a slave node as a weighted average of master nodes. Two definition strategies are supported in SDT, either direct or implicit. There are known robustness problems with the current implementation of this constraint.
The direct definition explicitly declares each node with coupled DOFs and weighting in a data field. Several rbe3 constrains can be declared in data.data. Each row of data.data codes a set of constraints following the format
Rbe3ID NodeIdSlave DofSlave Weight1 DofMaster1 NodeId1 Weight2 ...
DofMaster and DofSlave code which DOFs are used (123 for translations, 123456 for both translations and rotations). You can obtain the expression of the RBE3 as a MPC constraint using data=fe_mpc('rbe3c',model,'CaseEntryName').
When reading NASTRAN models an alternate definition
Rbe3ID NodeIdSlave DofSlave Weight DofMaster NodeId1 NodeId2 ... may exist. If the automated attempt to detect this format fails you can fix the entry using model=fe_mpc('FixRbe3 Alt',model).
The implicit definition handles Node Selectors described in section 7.11 to define the rbe3. The input is then a structure:
% Define a RBE3 constraint data=struct('SlaveSel','NodeSel',... 'MasterSel','NodeSel',... 'DOF', DofSlave,... 'MasterDOF', DofMaster);
SlaveSel is the slave node selection (typically a single node), MasterSel is the master node selection, DOF is the declaration of the slave node coupling, MasterDOF is the declaration of the master nodes coupling (same for all master nodes).
Grounding or coupling the slave node movement is possible through the use of a celas, as shown in the example below featuring an implicit rbe3 definition. In a practical approach, the slave node is duplicated and a celas element is generated between the two, which allows the definition of global movement stiffness. Constraining the rotation of a drilled block around its bore axis is considered using a global rotation stiffness.
% Integrated generation of an RBE3 constraint in a model % Definition of a drilled block around y model=feutil('ObjectHoleInBlock 0 0 0 1 0 0 0 1 0 2 2 2 .5 4 4 4'); model=fe_mat('DefaultIl',model); % default material properties model=fe_mat('defaultPl',model); % default element integration properties % Generation of the bore surface node set [i1,r1]=feutil('Findnode cyl ==0.5 o 0 0 0 0 1 0',model); model=feutil('AddsetNodeId',model,'bolt',r1(:,1)); % Generation of the slave node driving the global bore movement model.Node(end+[1:2],1:7)=[242 0 0 0 0 0 0;244 0 0 0 0 0 0]; % Addition of the celas element between the slave node and its duplicate model.Elt(end+[1:2],1:7)=[inf abs('celas') 0;242 244 123456 0 0 0 1e11]; model=feutil('AddSetNodeId',model,'ref_rot',244); % Definition of the RBE3 constraint data=struct('SlaveSel','setname ref_rot',... 'MasterSel','setname bolt',... 'DOF',123456,... % Slave node constrained on 6 DOF 'MasterDOF',123); % Master only use translation model=fe_case(model,'rbe3','block_mov',data); % Grounding the global y rotation (leaving the celas stiffness work) model=fe_case(model,'fixdof','ClampBlockRot',242.05); % 5 rigid body modes model obtained def=fe_eig(model,[5 20 1e3]); cf=feplot(model,def);fecom('curtabCases','rbe3');fecom('ProViewOn');
rigid#
See details under rigid which also illustrates the RigidAppend command.
Sens ... (SDT)#
SensDof entries are detailed in section 4.7. Command options vel and acc can be used to specify that certain sensors should measure velocity or acceleration. They are stored as rows of the Case.Stack cell array giving {'SensDof', Name, data}..
To properly retrieve a unique SensDof from the model, command [wire,name]=fe_case('GetSensDof',model) looks in the model Case with this strategy :
- If only one SensDof is defined, return this SensDof and its name
- If several SensDofs are defined return SensDof Test if there, else return first SensDof in the list
- Empty return if no SensDof found
To get back the observation matrix, use the command Sens=fe_case(model,'sens','SensName') as detailed in Sens for both full and reduced models.
R1=fe_case('sensobserve',model,'SensEntryName',def); iiplot(R1) can be used to extract observations at sensors associated with a given response. The SensEntryName can be omitted if a single sensor set exist. Sens=fe_case(model,'sens','SensName');R1=fe_case('sensobserve',Sens,def); is also acceptable
un=0#
model=fe_case(model,'un=0','Normal motion',map); where map gives normals at nodes generates an mpc case entry that enforces the condition {u}T{n}=0 at each node of the map.
SetCurve#
To associate a time variation to a compatible case entry, one adds a field curve to the case entry structure. This field is a cell array that is of the same length as the number of solicitation contained in the case entry.
Each curve definition in the cell array can be defined as either
- a string referring to the name of a curve stacked in the model (recommended)
- a curve structure
- a string that will be interpreted on the fly by fe_curvewhen the load is assembled, see fe_curve('TestList') to get the corresponding strings
The assignation is performed using
model = fe_case(model,'SetCurve',EntryName,CurveName,Curve,ind);
with
- EntryName the case entry to which the curve will be assigned. Use ? to find name automatically if only one exists.
- CurveName a string or a cell array of string with the name of the curves to assign
- Curve (optional) a curve or a cell array of curves that will be assigned (if not in model stack), they will be set in the model stack and only their names will be mentioned in the case entry
- ind (optional) the index of the curves to assign in the curve field, if several solicitation are present in the case entry considered. If ind is omitted the whole field curve of the case entry will be replaced by CurveName.
In practice, a variant call is supported for retro-compatibility but is not recommended for use,
model = fe_case(model,'SetCurve',EntryName,Curve,ind);
allows a direct assignation of non stacked curves to the case entry with the same behavior than for the classical way.
Multiple curve assignation at once to a specific EntryName is supported with the following rules
- CurveName, Curve (optional) and ind (mandatory) have the same sizes. In this case, all given curves will be assigned to the case entry with their provided index
- A singleCurveName and Curve is provided with a vector of indices. In this case, all indexed curves will be assigned to the new provided one
To remove a curve assignation to a case entry. Command
model = fe_case(model,'SetCurve',EntryName,'remove');
will remove the field curve from case entry EntryName.
The flexibility of the command imposes some restriction to the curve names. Name remove and TestVal with Val begin a keyword used by fe_curve Test cannot be used.
The following example illustrate the use of SetCurve to assign curves to case entries
% Sample calls to assign curves to load cases % generate a sample cube model model=femesh('testhexa8'); % clamp the cube bottom model=fe_case(model,'FixDof','clamped dofs','z==0'); % load a DOF of the cube base model=fe_case(model,'DofLoad','in',struct('def',1,'DOF',5.02)); % generate a curve loading transient pattern R1=fe_curve('testramp t1.005 yf1'); % assign the curve to the load case model=fe_case(model,'SetCurve','in','tramp',R1);% add a new load case with two sollicitations model=fe_case(model,'DofLoad','in2',... struct('def',[1 0;0 1],'DOF',[6.02;6.03])); % assign a new transient variation to both directions model=fe_case(model,'SetCurve','in2','tramp1',... fe_curve('testramp t0.5 yf1'),1:2); % modify the first direction only to tramp instead of tramp1 model=fe_case(model,'SetCurve','in2','tramp',1);
% remove the curve assigned to input in model=fe_case(model,'SetCurve','in','remove')
Examples#
Here is an example combining various fe_case commands
% Sample fe_case commands for boundary conditions, connections, and loads femesh('reset'); model = femesh('test ubeam plot'); % specifying clamped dofs (FixDof) model = fe_case(model,'FixDof','clamped dofs','z==0'); % creating a volume load data = struct('sel','GroupAll','dir',[1 0 0]); model = fe_case(model,'FVol','Volumic load',data); % assemble active DOFs and matrices model=fe_mknl(model); % assemble RHS (volumic load) Load = fe_load(model,'Case1'); % compute static response kd=ofact(model.K{2});def.def= kd\Load.def; ofact('clear',kd) Case=fe_case(model,'gett'); def.DOF=Case.DOF; % plot displacements feplot('initdef',def); fecom(';undef;triax;showpatch;promodelinit');
See also