DEMO_febio_0016_cube_viscoelastic_ramp_hold
Below is a demonstration for:
- Building geometry for a cube with hexahedral elements
- Defining the boundary conditions
- Coding the febio structure
- Running the model
- Importing and visualizing the displacement and stress results
Contents
Keywords
- febio_spec version 4.0
- febio, FEBio
- uniaxial loading
- compression, tension, compressive, tensile
- displacement control, displacement boundary condition
- hexahedral elements, hex8
- cube, box, rectangular
- static, solid
- hyperelastic, Ogden
- viscoelastic
- uncoupled, coupled
- ramp hold
- displacement logfile
- stress logfile
clear; close all; clc;
Plot settings
fontSize=20; faceAlpha1=0.8; markerSize=40; markerSize2=20; lineWidth=3;
Control parameters
% Path names defaultFolder = fileparts(fileparts(mfilename('fullpath'))); savePath=fullfile(defaultFolder,'data','temp'); % Defining file names febioFebFileNamePart='tempModel'; febioFebFileName=fullfile(savePath,[febioFebFileNamePart,'.feb']); %FEB file name febioLogFileName=[febioFebFileNamePart,'.txt']; %FEBio log file name febioLogFileName_disp=[febioFebFileNamePart,'_disp_out.txt']; %Log file name for exporting displacement febioLogFileName_force=[febioFebFileNamePart,'_force_out.txt']; %Log file name for exporting force febioLogFileName_stress=[febioFebFileNamePart,'_stress_out.txt']; %Log file name for exporting stress %Specifying dimensions and number of elements cubeSize=10; sampleWidth=cubeSize; %Width sampleThickness=cubeSize; %Thickness sampleHeight=cubeSize; %Height pointSpacings=5*ones(1,3); %Desired point spacing between nodes numElementsWidth=round(sampleWidth/pointSpacings(1)); %Number of elemens in dir 1 numElementsThickness=round(sampleThickness/pointSpacings(2)); %Number of elemens in dir 2 numElementsHeight=round(sampleHeight/pointSpacings(3)); %Number of elemens in dir 3 %Define applied displacement appliedStrain=0.3; %Linear strain (Only used to compute applied stretch) loadingOption='compression'; % or 'tension' switch loadingOption case 'compression' stretchLoad=1-appliedStrain; %The applied stretch for uniaxial loading case 'tension' stretchLoad=1+appliedStrain; %The applied stretch for uniaxial loading end displacementMagnitude=(stretchLoad*sampleHeight)-sampleHeight; %The displacement magnitude %Material parameter set c1=1e-3; %Shear-modulus-like parameter m1=2; %Material parameter setting degree of non-linearity k_factor=100; %Bulk modulus factor k=c1*k_factor; %Bulk modulus g1=1/4; %Viscoelastic QLV proportional coefficient t1=12; %Viscoelastic QLV time coefficient d=1e-9; %Density (not required for static analysis) uncoupledLaw=1; %1=uncoupled, 2=coupled analysisType='STATIC'; % FEA control settings t_load=5; %Time from start to max load t_step_ini1=t_load/20; %Initial desired step size numTimeSteps1=round(t_load/t_step_ini1); %Number of time steps desired t_step1=t_load/numTimeSteps1; %Step size dtmin1=t_step1/100; %Smallest allowed step size dtmax1=t_step1; %Largest allowed step size t_hold=40; t_step_ini2=t_step_ini1; %Initial desired step size numTimeSteps2=round(t_hold/t_step_ini2); %Number of time steps desired t_step2=t_hold/numTimeSteps2; %Step size dtmin2=t_step2/100; %Smallest allowed step size dtmax2=1; %Largest allowed step size max_refs=25; %Max reforms max_ups=0; %Set to zero to use full-Newton iterations opt_iter=6; %Optimum number of iterations max_retries=5; %Maximum number of retires runMode='external';%'internal';
Creating model geometry and mesh
A box is created with tri-linear hexahedral (hex8) elements using the hexMeshBox function. The function offers the boundary faces with seperate labels for the top, bottom, left, right, front, and back sides. As such these can be used to define boundary conditions on the exterior.
% Create a box with hexahedral elements cubeDimensions=[sampleWidth sampleThickness sampleHeight]; %Dimensions cubeElementNumbers=[numElementsWidth numElementsThickness numElementsHeight]; %Number of elements outputStructType=2; %A structure compatible with mesh view [meshStruct]=hexMeshBox(cubeDimensions,cubeElementNumbers,outputStructType); %Access elements, nodes, and faces from the structure E=meshStruct.elements; %The elements V=meshStruct.nodes; %The nodes (vertices) Fb=meshStruct.facesBoundary; %The boundary faces Cb=meshStruct.boundaryMarker; %The "colors" or labels for the boundary faces elementMaterialIndices=ones(size(E,1),1); %Element material indices
Plotting model boundary surfaces and a cut view
hFig=cFigure; subplot(1,2,1); hold on; title('Model boundary surfaces and labels','FontSize',fontSize); gpatch(Fb,V,Cb,'k',faceAlpha1); colormap(gjet(6)); icolorbar; axisGeom(gca,fontSize); hs=subplot(1,2,2); hold on; title('Cut view of solid mesh','FontSize',fontSize); optionStruct.hFig=[hFig hs]; meshView(meshStruct,optionStruct); axisGeom(gca,fontSize); drawnow;
Defining the boundary conditions
The visualization of the model boundary shows colors for each side of the cube. These labels can be used to define boundary conditions.
%Define supported node sets logicFace=Cb==1; %Logic for current face set Fr=Fb(logicFace,:); %The current face set bcSupportList_X=unique(Fr(:)); %Node set part of selected face logicFace=Cb==3; %Logic for current face set Fr=Fb(logicFace,:); %The current face set bcSupportList_Y=unique(Fr(:)); %Node set part of selected face logicFace=Cb==5; %Logic for current face set Fr=Fb(logicFace,:); %The current face set bcSupportList_Z=unique(Fr(:)); %Node set part of selected face %Prescribed displacement nodes logicPrescribe=Cb==6; %Logic for current face set Fr=Fb(logicPrescribe,:); %The current face set bcPrescribeList=unique(Fr(:)); %Node set part of selected face
Visualizing boundary conditions. Markers plotted on the semi-transparent model denote the nodes in the various boundary condition lists.
hf=cFigure; title('Boundary conditions','FontSize',fontSize); xlabel('X','FontSize',fontSize); ylabel('Y','FontSize',fontSize); zlabel('Z','FontSize',fontSize); hold on; gpatch(Fb,V,'kw','k',0.5); hl(1)=plotV(V(bcSupportList_X,:),'r.','MarkerSize',markerSize); hl(2)=plotV(V(bcSupportList_Y,:),'g.','MarkerSize',markerSize); hl(3)=plotV(V(bcSupportList_Z,:),'b.','MarkerSize',markerSize); hl(4)=plotV(V(bcPrescribeList,:),'k.','MarkerSize',markerSize); legend(hl,{'BC x support','BC y support','BC z support','BC z prescribe'}); axisGeom(gca,fontSize); camlight headlight; drawnow;
Defining the FEBio input structure
See also febioStructTemplate and febioStruct2xml and the FEBio user manual.
%Get a template with default settings [febio_spec]=febioStructTemplate; %febio_spec version febio_spec.ATTR.version='4.0'; %Module section febio_spec.Module.ATTR.type='solid'; %Control sections for each step febio_spec.Step.step{1}.Control=febio_spec.Control; %Copy from template febio_spec.Step.step{1}.ATTR.id=1; febio_spec.Step.step{1}.Control.analysis=analysisType; febio_spec.Step.step{1}.Control.time_steps=numTimeSteps1; febio_spec.Step.step{1}.Control.step_size=t_step1; febio_spec.Step.step{1}.Control.solver.max_refs=max_refs; febio_spec.Step.step{1}.Control.time_stepper.dtmin=dtmin1; febio_spec.Step.step{1}.Control.time_stepper.dtmax=dtmax1; febio_spec.Step.step{1}.Control.time_stepper.max_retries=max_retries; febio_spec.Step.step{1}.Control.time_stepper.opt_iter=opt_iter; febio_spec.Step.step{2}.Control=febio_spec.Control; %Copy from template febio_spec.Step.step{2}.ATTR.id=2; febio_spec.Step.step{2}.Control.analysis=analysisType; febio_spec.Step.step{2}.Control.time_steps=numTimeSteps2; febio_spec.Step.step{2}.Control.step_size=t_step2; febio_spec.Step.step{2}.Control.solver.max_refs=max_refs; febio_spec.Step.step{2}.Control.time_stepper.dtmin=dtmin2; febio_spec.Step.step{2}.Control.time_stepper.dtmax=dtmax2; febio_spec.Step.step{2}.Control.time_stepper.max_retries=max_retries; febio_spec.Step.step{2}.Control.time_stepper.opt_iter=opt_iter; %Remove control field (part of template) since step specific control sections are used febio_spec=rmfield(febio_spec,'Control'); %Material section materialName1='Material1'; febio_spec.Material.material{1}.ATTR.name=materialName1; switch uncoupledLaw case 1 %Viscoelastic part febio_spec.Material.material{1}.ATTR.type='uncoupled viscoelastic'; febio_spec.Material.material{1}.ATTR.id=1; febio_spec.Material.material{1}.g1=g1; febio_spec.Material.material{1}.t1=t1; febio_spec.Material.material{1}.k=k; febio_spec.Material.material{1}.density=d; %Elastic part febio_spec.Material.material{1}.elastic{1}.ATTR.type='Ogden'; febio_spec.Material.material{1}.elastic{1}.c1=c1; febio_spec.Material.material{1}.elastic{1}.m1=m1; febio_spec.Material.material{1}.elastic{1}.c2=c1; febio_spec.Material.material{1}.elastic{1}.m2=-m1; febio_spec.Material.material{1}.elastic{1}.density=d; case 2 %Viscoelastic part febio_spec.Material.material{1}.ATTR.type='viscoelastic'; febio_spec.Material.material{1}.ATTR.id=1; febio_spec.Material.material{1}.g1=g1; febio_spec.Material.material{1}.t1=t1; febio_spec.Material.material{1}.density=d; %Elastic part febio_spec.Material.material{1}.elastic{1}.ATTR.type='Ogden unconstrained'; febio_spec.Material.material{1}.elastic{1}.c1=c1; febio_spec.Material.material{1}.elastic{1}.m1=m1; febio_spec.Material.material{1}.elastic{1}.c2=c1; febio_spec.Material.material{1}.elastic{1}.m2=-m1; febio_spec.Material.material{1}.elastic{1}.cp=k; febio_spec.Material.material{1}.elastic{1}.density=d; end % Mesh section % -> Nodes febio_spec.Mesh.Nodes{1}.ATTR.name='Object1'; %The node set name febio_spec.Mesh.Nodes{1}.node.ATTR.id=(1:size(V,1))'; %The node id's febio_spec.Mesh.Nodes{1}.node.VAL=V; %The nodel coordinates % -> Elements partName1='Part1'; febio_spec.Mesh.Elements{1}.ATTR.name=partName1; %Name of this part febio_spec.Mesh.Elements{1}.ATTR.type='hex8'; %Element type febio_spec.Mesh.Elements{1}.elem.ATTR.id=(1:1:size(E,1))'; %Element id's febio_spec.Mesh.Elements{1}.elem.VAL=E; %The element matrix % -> NodeSets nodeSetName1='bcSupportList_X'; nodeSetName2='bcSupportList_Y'; nodeSetName3='bcSupportList_Z'; nodeSetName4='bcPrescribeList'; febio_spec.Mesh.NodeSet{1}.ATTR.name=nodeSetName1; febio_spec.Mesh.NodeSet{1}.VAL=mrow(bcSupportList_X); febio_spec.Mesh.NodeSet{2}.ATTR.name=nodeSetName2; febio_spec.Mesh.NodeSet{2}.VAL=mrow(bcSupportList_Y); febio_spec.Mesh.NodeSet{3}.ATTR.name=nodeSetName3; febio_spec.Mesh.NodeSet{3}.VAL=mrow(bcSupportList_Z); febio_spec.Mesh.NodeSet{4}.ATTR.name=nodeSetName4; febio_spec.Mesh.NodeSet{4}.VAL=mrow(bcPrescribeList); %MeshDomains section febio_spec.MeshDomains.SolidDomain.ATTR.name=partName1; febio_spec.MeshDomains.SolidDomain.ATTR.mat=materialName1; %Boundary condition section % -> Fix boundary conditions febio_spec.Boundary.bc{1}.ATTR.name='FixedDisplacement01'; febio_spec.Boundary.bc{1}.ATTR.type='zero displacement'; febio_spec.Boundary.bc{1}.ATTR.node_set=nodeSetName1; febio_spec.Boundary.bc{1}.x_dof=1; febio_spec.Boundary.bc{1}.y_dof=0; febio_spec.Boundary.bc{1}.z_dof=0; febio_spec.Boundary.bc{2}.ATTR.name='FixedDisplacement02'; febio_spec.Boundary.bc{2}.ATTR.type='zero displacement'; febio_spec.Boundary.bc{2}.ATTR.node_set=nodeSetName2; febio_spec.Boundary.bc{2}.x_dof=0; febio_spec.Boundary.bc{2}.y_dof=1; febio_spec.Boundary.bc{2}.z_dof=0; febio_spec.Boundary.bc{3}.ATTR.name='FixedDisplacement03'; febio_spec.Boundary.bc{3}.ATTR.type='zero displacement'; febio_spec.Boundary.bc{3}.ATTR.node_set=nodeSetName3; febio_spec.Boundary.bc{3}.x_dof=0; febio_spec.Boundary.bc{3}.y_dof=0; febio_spec.Boundary.bc{3}.z_dof=1; febio_spec.Boundary.bc{4}.ATTR.name='bcPrescribeList'; febio_spec.Boundary.bc{4}.ATTR.type='prescribed displacement'; febio_spec.Boundary.bc{4}.ATTR.node_set=nodeSetName4; febio_spec.Boundary.bc{4}.dof='z'; febio_spec.Boundary.bc{4}.value.ATTR.lc=1; febio_spec.Boundary.bc{4}.value.VAL=displacementMagnitude; febio_spec.Boundary.bc{4}.relative=0; %LoadData section % -> load_controller febio_spec.LoadData.load_controller{1}.ATTR.name='LC1'; febio_spec.LoadData.load_controller{1}.ATTR.id=1; febio_spec.LoadData.load_controller{1}.ATTR.type='loadcurve'; febio_spec.LoadData.load_controller{1}.interpolate='LINEAR'; febio_spec.LoadData.load_controller{1}.extend='CONSTANT'; febio_spec.LoadData.load_controller{1}.points.pt.VAL=[0 0;t_load 1;(t_load+t_hold) 1]; %Output section % -> log file febio_spec.Output.logfile.ATTR.file=febioLogFileName; febio_spec.Output.logfile.node_data{1}.ATTR.file=febioLogFileName_disp; febio_spec.Output.logfile.node_data{1}.ATTR.data='ux;uy;uz'; febio_spec.Output.logfile.node_data{1}.ATTR.delim=','; febio_spec.Output.logfile.node_data{2}.ATTR.file=febioLogFileName_force; febio_spec.Output.logfile.node_data{2}.ATTR.data='Rx;Ry;Rz'; febio_spec.Output.logfile.node_data{2}.ATTR.delim=','; febio_spec.Output.logfile.element_data{1}.ATTR.file=febioLogFileName_stress; febio_spec.Output.logfile.element_data{1}.ATTR.data='sz'; febio_spec.Output.logfile.element_data{1}.ATTR.delim=','; febio_spec.Output.plotfile.compression=0;
Quick viewing of the FEBio input file structure
The febView function can be used to view the xml structure in a MATLAB figure window.
febView(febio_spec); %Viewing the febio file
Exporting the FEBio input file
Exporting the febio_spec structure to an FEBio input file is done using the febioStruct2xml function.
febioStruct2xml(febio_spec,febioFebFileName); %Exporting to file and domNode
Running the FEBio analysis
To run the analysis defined by the created FEBio input file the runMonitorFEBio function is used. The input for this function is a structure defining job settings e.g. the FEBio input file name. The optional output runFlag informs the user if the analysis was run succesfully.
febioAnalysis.run_filename=febioFebFileName; %The input file name febioAnalysis.run_logname=febioLogFileName; %The name for the log file febioAnalysis.disp_on=1; %Display information on the command window febioAnalysis.runMode=runMode; [runFlag]=runMonitorFEBio(febioAnalysis);%START FEBio NOW!!!!!!!!
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% --------> RUNNING/MONITORING FEBIO JOB <-------- 01-Sep-2023 11:27:19 FEBio path: /home/kevin/FEBioStudio/bin/febio4 # Attempt removal of existing log files 01-Sep-2023 11:27:19 * Removal succesful 01-Sep-2023 11:27:19 # Attempt removal of existing .xplt files 01-Sep-2023 11:27:19 * Removal succesful 01-Sep-2023 11:27:19 # Starting FEBio... 01-Sep-2023 11:27:19 Max. total analysis time is: Inf s * Waiting for log file creation 01-Sep-2023 11:27:19 Max. wait time: 30 s * Log file found. 01-Sep-2023 11:27:19 # Parsing log file... 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 0.25 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 0.5 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 0.75 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 1 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 1.25 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 1.5 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 1.75 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 2 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 2.25 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 2.5 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 2.75 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 3 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 3.25 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 3.5 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 3.75 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 4 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 4.25 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 4.5 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 4.75 01-Sep-2023 11:27:19 number of iterations : 3 01-Sep-2023 11:27:19 number of reformations : 3 01-Sep-2023 11:27:19 ------- converged at time : 5 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 5.25 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 5.65 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 6.17 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 6.786 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 7.4788 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 8.23304 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 9.03643 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 9.87915 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 10.7533 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 11.6527 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 12.5721 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 13.5077 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 14.4562 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 15.4149 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 16.3819 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 17.3556 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 18.3344 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 19.3176 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 20.304 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 21.2932 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 22.2846 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 23.2777 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 24.2721 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 25.2677 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 26.2642 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 27.2613 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 28.2591 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 29.2573 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 30.2558 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 31.2546 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 32.2537 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 33.253 01-Sep-2023 11:27:19 number of iterations : 2 01-Sep-2023 11:27:19 number of reformations : 2 01-Sep-2023 11:27:19 ------- converged at time : 34.2524 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 35.2519 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 36.2515 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 37.2512 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 38.251 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 39.2508 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 40.2506 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 41.2505 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 42.2504 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 43.2503 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 44.2503 01-Sep-2023 11:27:19 number of iterations : 1 01-Sep-2023 11:27:19 number of reformations : 1 01-Sep-2023 11:27:19 ------- converged at time : 45 01-Sep-2023 11:27:19 Elapsed time : 0:00:00 01-Sep-2023 11:27:19 N O R M A L T E R M I N A T I O N # Done 01-Sep-2023 11:27:19 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
Import FEBio results
if runFlag==1 %i.e. a succesful run
Importing nodal displacements from a log file
dataStruct=importFEBio_logfile(fullfile(savePath,febioLogFileName_disp),0,1); %Access data N_disp_mat=dataStruct.data; %Displacement timeVec=dataStruct.time; %Time %Create deformed coordinate set V_DEF=N_disp_mat+repmat(V,[1 1 size(N_disp_mat,3)]);
Plotting the simulated results using anim8 to visualize and animate deformations
DN_magnitude=sqrt(sum(N_disp_mat(:,:,end).^2,2)); %Current displacement magnitude % Create basic view and store graphics handle to initiate animation hf=cFigure; %Open figure gtitle([febioFebFileNamePart,': Press play to animate']); title('Displacement magnitude [mm]','Interpreter','Latex') hp=gpatch(Fb,V_DEF(:,:,end),DN_magnitude,'k',1); %Add graphics object to animate hp.Marker='.'; hp.MarkerSize=markerSize2; hp.FaceColor='interp'; gpatch(Fb,V,0.5*ones(1,3),'k',0.25); %A static graphics object axisGeom(gca,fontSize); colormap(gjet(250)); colorbar; caxis([0 max(DN_magnitude)]); axis(axisLim(V_DEF)); %Set axis limits statically camlight headlight; % Set up animation features animStruct.Time=timeVec; %The time vector for qt=1:1:size(N_disp_mat,3) %Loop over time increments DN_magnitude=sqrt(sum(N_disp_mat(:,:,qt).^2,2)); %Current displacement magnitude %Set entries in animation structure animStruct.Handles{qt}=[hp hp]; %Handles of objects to animate animStruct.Props{qt}={'Vertices','CData'}; %Properties of objects to animate animStruct.Set{qt}={V_DEF(:,:,qt),DN_magnitude}; %Property values for to set in order to animate end anim8(hf,animStruct); %Initiate animation feature drawnow;
Importing element stress from a log file
dataStruct=importFEBio_logfile(fullfile(savePath,febioLogFileName_stress),0,1);
%Access data
E_stress_mat=dataStruct.data;
Plotting the simulated results using anim8 to visualize and animate deformations
[CV]=faceToVertexMeasure(E,V,E_stress_mat(:,:,end)); % Create basic view and store graphics handle to initiate animation hf=cFigure; %Open figure gtitle([febioFebFileNamePart,': Press play to animate']); title('$\sigma_{zz}$ [MPa]','Interpreter','Latex') hp=gpatch(Fb,V_DEF(:,:,end),CV,'k',1); %Add graphics object to animate hp.Marker='.'; hp.MarkerSize=markerSize2; hp.FaceColor='interp'; gpatch(Fb,V,0.5*ones(1,3),'k',0.25); %A static graphics object axisGeom(gca,fontSize); colormap(gjet(250)); colorbar; caxis([min(E_stress_mat(:)) max(E_stress_mat(:))]); axis(axisLim(V_DEF)); %Set axis limits statically camlight headlight; % Set up animation features animStruct.Time=timeVec; %The time vector for qt=1:1:size(N_disp_mat,3) %Loop over time increments [CV]=faceToVertexMeasure(E,V,E_stress_mat(:,:,qt)); %Set entries in animation structure animStruct.Handles{qt}=[hp hp]; %Handles of objects to animate animStruct.Props{qt}={'Vertices','CData'}; %Properties of objects to animate animStruct.Set{qt}={V_DEF(:,:,qt),CV}; %Property values for to set in order to animate end anim8(hf,animStruct); %Initiate animation feature drawnow;
Calculate metrics to visualize time-stress curve
stress_cauchy_sim=squeeze(mean(E_stress_mat,1)); %Cauchy stress sigma_z
Visualize stress-stretch curve
cFigure; hold on; title('Uniaxial stress-time curve','FontSize',fontSize); xlabel('Time [s]','FontSize',fontSize,'Interpreter','Latex'); ylabel('$\sigma_{zz}$ [MPa]','FontSize',fontSize,'Interpreter','Latex'); plot(timeVec(:),stress_cauchy_sim(:),'r-','lineWidth',lineWidth); view(2); axis tight; grid on; axis square; box on; set(gca,'FontSize',fontSize); drawnow;
Importing element stress from a log file
dataStruct=importFEBio_logfile(fullfile(savePath,febioLogFileName_force),0,1);
%Access data
N_force_mat=dataStruct.data;
force_sim=-squeeze(sum(N_force_mat(:,3,:),1));
Visualize stress-stretch curve
cFigure; hold on; title('Uniaxial stress-time curve','FontSize',fontSize); xlabel('Time [s]','FontSize',fontSize,'Interpreter','Latex'); ylabel('$F_{z}$ [N]','FontSize',fontSize,'Interpreter','Latex'); plot(timeVec(:),force_sim(:),'r-','lineWidth',lineWidth); view(2); axis tight; grid on; axis square; box on; set(gca,'FontSize',fontSize); drawnow;
end
GIBBON www.gibboncode.org
Kevin Mattheus Moerman, [email protected]
GIBBON footer text
License: https://github.com/gibbonCode/GIBBON/blob/master/LICENSE
GIBBON: The Geometry and Image-based Bioengineering add-On. A toolbox for image segmentation, image-based modeling, meshing, and finite element analysis.
Copyright (C) 2006-2023 Kevin Mattheus Moerman and the GIBBON contributors
This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/.