A First Vfx Program

This tutorial demonstrates how to create a complete structural finite element analysis using VfxTools from scratch. The example solves a simple 3D solid mechanics problem with two hexahedral elements, applying loads and boundary conditions, running the analysis, and extracting results.

Problem Description

The example models a simple cantilever beam structure consisting of two hexahedral (brick) elements:

                fixed at nodes 1,2,3,4,5,6 in all directions

             4------5------6
            /.     /.     /|
  y        / .    / .    / |
  |       /  .   /  .   /  |
  --x    /   1../...2../...3
 /     10-----11-----12   /
z       |  .   |  .   |  /
        | .    | .    | /
        |.     |.     |/
        7------8------9

        force applied on node 7 in z direction

The structure is fixed at the top face (nodes 1-6) and loaded with a force in the z-direction at the bottom node (node 7).

Step-by-Step Walkthrough

1. Create the Mesh

Create a Connect object to define the mesh topology and geometry:

vis_Connect* connect = vis_ConnectBegin();
vis_ConnectPre(connect, SYS_DOUBLE);
vis_ConnectDef(connect, MAX_NODE, MAX_ELEM);

Define element connectivity and coordinates:

for (Vint i = 0; i < MAX_ELEM; i++) {
    vis_ConnectSetTopology(connect, i + 1, SYS_SHAPEHEX, 2, 0, 0);
    vis_ConnectSetElemNode(connect, i + 1, connectivity[i]);
    vis_ConnectSetElemAssoc(connect, VIS_PROPID, i + 1, PROPERTY_ID);
}

Set the nodal coordinates in the Connect object:

for (Vint i = 0; i < MAX_NODE; i++) {
    vis_ConnectSetCoordsdv(connect, i + 1, coordinates[i]);
}

2. Create the Model and Associate the Mesh

Create a Model object to encapsulate all model data:

vis_Model* model = vis_ModelBegin();
vis_ModelSetObject(model, VIS_CONNECT, connect);

3. Define Element Properties

Create element properties (EProp) and associate them with material IDs:

vis_EProp* elementProperty = vis_EPropBegin();
vis_EPropDef(elementProperty, SYS_ELEM_SOLID);
vis_EPropSetValuei(elementProperty, EPROP_MID, MATERIAL_ID);

vsy_HashTable* elementPropertyHastTable = vsy_HashTableBegin();
vsy_HashTableInsert(elementPropertyHastTable, PROPERTY_ID, (Vobject*)elementProperty);
vis_ModelSetHashTable(model, VIS_EPROP, elementPropertyHastTable);

4. Define Material Properties

Define isotropic material properties (MProp):

vis_MProp* materialProperty = vis_MPropBegin();
vis_MPropDef(materialProperty, SYS_MAT_ISOTROPIC);
vis_MPropSetValuef(materialProperty, MPROP_E, 2.e+7f);      // Young's modulus
vis_MPropSetValuef(materialProperty, MPROP_NU, 0.3f);       // Poisson's ratio

vsy_HashTable* materialPropertyHashTable = vsy_HashTableBegin();
vsy_HashTableInsert(materialPropertyHashTable, MATERIAL_ID, (Vobject*)materialProperty);
vis_ModelSetHashTable(model, VIS_MPROP, materialPropertyHashTable);

5. Define Solution Properties

Configure the solution type (SProp) for static structural analysis:

vis_SProp* solutionProperty = vis_SPropBegin();
vis_SPropDef(solutionProperty, SYS_SOL_STATIC);
vis_SPropSetValuec(solutionProperty, SPROP_TITLE, (Vchar*)"Structural Analysis");
vis_SPropSetValuei(solutionProperty, SPROP_ANALYSIS, SYS_ANALYSIS_STRUCTURAL);
vis_SPropSetValuei(solutionProperty, SPROP_LCASE_NUM, 1);
vis_SPropSetValuei(solutionProperty, SPROP_LCASE, 1);
vis_SPropSetValuei(solutionProperty, SPROP_RCASE, 1);

// Specify which result types to save
Vint resultsToDump[3] = {SYS_RES_D, SYS_RES_S, SYS_RES_E};
vis_SPropSetValueiv(solutionProperty, SPROP_RESFILE, resultsToDump);

6. Define Boundary Conditions

Create restraint case (RCase) to fix nodes 1-6 in all directions:

vis_RCase* restraintCase = vis_RCaseBegin();
for (Vint node = 1; node <= 6; node++) {
    vis_RCaseSetSPC(restraintCase, node, SYS_DOF_TX, RCASE_FIXED, NULL, 0);
    vis_RCaseSetSPC(restraintCase, node, SYS_DOF_TY, RCASE_FIXED, NULL, 0);
    vis_RCaseSetSPC(restraintCase, node, SYS_DOF_TZ, RCASE_FIXED, NULL, 0);
}

7. Define Loads

Create load case (LCase) to apply forces at node 7:

Vdouble forceApplied[3] = {0., 0., 1.e+6};  // 1 MN in z-direction
vis_LCase* loadCase = vis_LCaseBegin();

vis_LCaseAddConcdv(loadCase, 7, LCASE_FORCE, forceApplied);

8. Run the Solver

Create a ProSolve object and execute the analysis:

vfx_ProSolve* prosolve = vfx_ProSolveBegin();
vfx_ProSolveSetParami(prosolve, PROSOLVE_SAVERES, SYS_ON);
vfx_ProSolveSetParami(prosolve, PROSOLVE_PRINTMODEL, SYS_ON); // Print model summary
vfx_ProSolveSetObject(prosolve, VIS_MODEL, (Vobject*)model);
vfx_ProSolveInit(prosolve);
vfx_ProSolveExec(prosolve);

9. Extract Results

Extract displacement results using RProp and State:

vis_RProp* resultProperty = vis_RPropBegin();
vis_RPropDef(resultProperty, SYS_NODE, SYS_NONE);
vis_RPropSetType(resultProperty, SYS_RES_D);  // Displacement results

vis_State* resultState = vis_StateBegin();
vis_StateDef(resultState, MAX_NODE, SYS_NODE, SYS_NONE, SYS_SIXDOF);
vfx_ProSolveState(prosolve, resultProperty, resultState);

// Retrieve and print displacements
Vint maxResultData, maxLoc, maxSection;
vis_StateDataMax(resultState, &maxResultData, &maxLoc, &maxSection);
vsy::vector<Vfloat> result(maxResultData);

for (Vint n = 1; n <= MAX_NODE; ++n) {
    vis_StateData(resultState, 1, &n, result.data());
    printf("%5d %12.4e %12.4e %12.4e\n", n, result[0], result[1], result[2]);
}

Similarly, extract stress and strain results by modifying the result property:

// For stresses
vis_RPropDef(resultProperty, SYS_ELEM, SYS_NODE);
vis_RPropSetType(resultProperty, SYS_RES_S);

// For strains
vis_RPropSetType(resultProperty, SYS_RES_E);

10. Clean Up

Release all allocated objects:

vis_RPropEnd(resultProperty);
vis_StateEnd(resultState);
vis_ModelDelete(model);
vis_ModelEnd(model);
vfx_ProSolveEnd(prosolve);

Complete Source Code

#include <stdio.h>
#include <vector>
#include "sam/vfx/vfx.h"
#include "sam/base/license.h"
#include "sam/hoops_license.h"
/*
                    fixed at nodes 1,2,3,4,5,6 in all directions

                 4------5------6
                /.     /.     /|
  y            / .    / .    / |
  |           /  .   /  .   /  |
  --x        /   1../...2../...3
 /         10-----11-----12   /
z           |  .   |  .   |  /
            | .    | .    | /
            |.     |.     |/
            7------8------9

            force applied on node 7 in z direction

*/
#define MAX_ELEM    2
#define MAX_NODE    12
#define PROPERTY_ID 7
#define MATERIAL_ID 3

static Vint connectivity[MAX_ELEM][8] = {{1, 2, 5, 4, 7, 8, 11, 10}, {2, 3, 6, 5, 8, 9, 12, 11}};

static Vdouble coordinates[MAX_NODE][3] = {{0., 0., 0.}, {1., 0., 0.}, {2., 0., 0.}, {0., 1., 0.}, {1., 1., 0.}, {2., 1., 0.},
                                           {0., 0., 1.}, {1., 0., 1.}, {2., 0., 1.}, {0., 1., 1.}, {1., 1., 1.}, {2., 1., 1.}};

/*----------------------------------------------------------------------
                      Structural Analysis
----------------------------------------------------------------------*/
int
main()
{
    vsy_LicenseValidate(HOOPS_LICENSE);

    /*  create connect/mesh object */
    vis_Connect* connect = vis_ConnectBegin();
    vis_ConnectPre(connect, SYS_DOUBLE);
    vis_ConnectDef(connect, MAX_NODE, MAX_ELEM);

    for (Vint i = 0; i < MAX_ELEM; i++) {
        /* set element connectivity */
        vis_ConnectSetTopology(connect, i + 1, SYS_SHAPEHEX, 2, 0, 0);
        vis_ConnectSetElemNode(connect, i + 1, connectivity[i]);

        /* assign property id */
        vis_ConnectSetElemAssoc(connect, VIS_PROPID, i + 1, PROPERTY_ID);
    }
    /* set nodal coordinates */
    for (Vint i = 0; i < MAX_NODE; i++) {
        vis_ConnectSetCoordsdv(connect, i + 1, coordinates[i]);
    }

    /* create and fill model object */
    vis_Model* model = vis_ModelBegin();

    /* saving mesh into the model */
    vis_ModelSetObject(model, VIS_CONNECT, connect);

    /* define element properties and save into the model*/
    vis_EProp* elementProperty = vis_EPropBegin();
    vis_EPropDef(elementProperty, SYS_ELEM_SOLID);
    vis_EPropSetValuei(elementProperty, EPROP_MID, MATERIAL_ID);
    vsy_HashTable* elementPropertyHastTable = vsy_HashTableBegin();
    vsy_HashTableInsert(elementPropertyHastTable, PROPERTY_ID, (Vobject*)elementProperty);
    vis_ModelSetHashTable(model, VIS_EPROP, elementPropertyHastTable);

    /* define material properties and save into the model*/
    vis_MProp* materialProperty = vis_MPropBegin();
    vis_MPropDef(materialProperty, SYS_MAT_ISOTROPIC);
    vis_MPropSetValuef(materialProperty, MPROP_E, 2.e+7f);
    vis_MPropSetValuef(materialProperty, MPROP_NU, 0.3f);
    vsy_HashTable* materialPropertyHashTable = vsy_HashTableBegin();
    vsy_HashTableInsert(materialPropertyHashTable, MATERIAL_ID, (Vobject*)materialProperty);
    vis_ModelSetHashTable(model, VIS_MPROP, materialPropertyHashTable);

    /* Setting up solution property and save into the model*/
    vis_SProp* solutionProperty = vis_SPropBegin();
    vis_SPropDef(solutionProperty, SYS_SOL_STATIC);
    vis_SPropSetValuec(solutionProperty, SPROP_TITLE, (Vchar*)"Structural Analysis");
    vis_SPropSetValuei(solutionProperty, SPROP_ANALYSIS, SYS_ANALYSIS_STRUCTURAL);
    vis_SPropSetValuei(solutionProperty, SPROP_LCASE_NUM, 1);
    vis_SPropSetValuei(solutionProperty, SPROP_LCASE, 1);
    vis_SPropSetValuei(solutionProperty, SPROP_RCASE, 1);
    vis_SPropSetValuei(solutionProperty, SPROP_RESFILE_NUM, 3);
    Vint resultsToDump[3] = {SYS_RES_D, SYS_RES_S, SYS_RES_E}; /* displacements, stresses, strains */
    vis_SPropSetValueiv(solutionProperty, SPROP_RESFILE, resultsToDump);
    vsy_List* solutionPropertyList = vsy_ListBegin();
    vsy_ListInsert(solutionPropertyList, 1, (Vobject*)solutionProperty);
    vis_ModelSetList(model, VIS_SPROP, solutionPropertyList);

    /* Setting up restraints and save into the model*/
    vis_RCase* restraintCase = vis_RCaseBegin();
    /*Fixing the top face nodes in all directions */
    for (Vint node = 1; node <= 6; node++) {
        vis_RCaseSetSPC(restraintCase, node, SYS_DOF_TX, RCASE_FIXED, NULL, 0);
        vis_RCaseSetSPC(restraintCase, node, SYS_DOF_TY, RCASE_FIXED, NULL, 0);
        vis_RCaseSetSPC(restraintCase, node, SYS_DOF_TZ, RCASE_FIXED, NULL, 0);
    }
    vsy_HashTable* restraintCaseHashTable = vsy_HashTableBegin();
    vsy_HashTableInsert(restraintCaseHashTable, 1, (Vobject*)restraintCase);
    vis_ModelSetHashTable(model, VIS_RCASE, restraintCaseHashTable);

    /* Setting up load case and save into the model*/
    Vdouble forceApplied[3] = {0., 0., 1.e+6};
    auto loadCase = vis_LCaseBegin();

    /* Applying force of 1 MN in z-direction at node 7 */
    vis_LCaseAddConcdv(loadCase, 7, LCASE_FORCE, forceApplied);

    vsy_HashTable* loadCaseHashTable = vsy_HashTableBegin();
    vsy_HashTableInsert(loadCaseHashTable, 1, (Vobject*)loadCase);
    vis_ModelSetHashTable(model, VIS_LCASE, loadCaseHashTable);

    /* Setting up ProSolve and feed model into it*/
    vfx_ProSolve* prosolve = vfx_ProSolveBegin();
    vfx_ProSolveSetParami(prosolve, PROSOLVE_SAVERES, SYS_ON);

    /*Print model summary*/
    vfx_ProSolveSetParami(prosolve, PROSOLVE_PRINTMODEL, SYS_ON);
    vfx_ProSolveSetObject(prosolve, VIS_MODEL, (Vobject*)model);
    vfx_ProSolveInit(prosolve);

    /* run ProSolve for structural problem */
    vfx_ProSolveExec(prosolve);

    /* creating result property and state to get displacements */
    vis_RProp* resultProperty = vis_RPropBegin();
    vis_RPropDef(resultProperty, SYS_NODE, SYS_NONE);
    vis_RPropSetType(resultProperty, SYS_RES_D);

    vis_State* resultState = vis_StateBegin();
    vis_StateDef(resultState, MAX_NODE, SYS_NODE, SYS_NONE, SYS_SIXDOF);
    vfx_ProSolveState(prosolve, resultProperty, resultState);
    printf("Displacements\n");

    Vint maxResultData, maxLoc, maxSection;
    vis_StateDataMax(resultState, &maxResultData, &maxLoc, &maxSection);

    std::vector<Vfloat> result(maxResultData);
    for (Vint n = 1; n <= MAX_NODE; ++n) {
        vis_StateData(resultState, 1, &n, result.data());
        printf("%5d %12.4e %12.4e %12.4e\n", n, result[0], result[1], result[2]);
    }
    printf("\n");

    /*modifying the result property and state to get stress results*/
    vis_RPropDef(resultProperty, SYS_ELEM, SYS_NODE);
    vis_RPropSetType(resultProperty, SYS_RES_S);
    vis_StateDef(resultState, MAX_ELEM, SYS_ELEM, SYS_NODE, SYS_TENSOR);
    vfx_ProSolveState(prosolve, resultProperty, resultState);

    vis_StateDataMax(resultState, &maxResultData, &maxLoc, &maxSection);
    result.resize(maxResultData);
    printf("Stresses\n");
    for (Vint n = 1; n <= MAX_ELEM; ++n) {
        vis_StateData(resultState, 1, &n, result.data());
        printf("%5d\n", n);
        for (Vint i = 0; i < 8; ++i) {
            printf("   %5d %12.4e %12.4e %12.4e %12.4e %12.4e %12.4e\n", i + 1, result[6 * i], result[6 * i + 1],
                   result[6 * i + 2], result[6 * i + 3], result[6 * i + 4], result[6 * i + 5]);
        }
    }
    printf("\n");

    /*modifying the result property and state to get strain results*/
    vis_RPropSetType(resultProperty, SYS_RES_E);
    vis_StateDef(resultState, MAX_ELEM, SYS_ELEM, SYS_NODE, SYS_TENSOR);
    vfx_ProSolveState(prosolve, resultProperty, resultState);
    vis_StateDataMax(resultState, &maxResultData, &maxLoc, &maxSection);
    result.resize(maxResultData);
    printf("Strains\n");
    for (Vint n = 1; n <= MAX_ELEM; ++n) {
        vis_StateData(resultState, 1, &n, result.data());
        printf("%5d\n", n);
        for (Vint i = 0; i < 8; ++i) {
            printf("   %5d %12.4e %12.4e %12.4e %12.4e %12.4e %12.4e\n", i + 1, result[6 * i], result[6 * i + 1],
                   result[6 * i + 2], result[6 * i + 3], result[6 * i + 4], result[6 * i + 5]);
        }
    }
    printf("\n");

    /* end objects */
    vis_RPropEnd(resultProperty);
    vis_StateEnd(resultState);
    vis_ModelDelete(model);
    vis_ModelEnd(model);
    vfx_ProSolveEnd(prosolve);
    return 0;
}

Key Concepts

  • Connect: Defines mesh topology and geometry
  • Model: Container for all finite element model data
  • EProp: Element properties linking elements to materials
  • MProp: Material properties (elastic constants)
  • SProp: Solution properties (analysis type, result requests)
  • RCase: Restraint case (boundary conditions)
  • LCase: Load case (applied forces)
  • ProSolve: Solver engine that performs the analysis
  • RProp/State: Result extraction and post-processing

Expected Output

The program will output: 1. The model summary 2. Displacements at all 12 nodes where nodes 1-6 should have zero displacement 3. Stresses at element level 4. Strains at element level 5. Log files with solution summary

This example demonstrates the complete workflow from model creation to result extraction, providing a foundation for more complex analyses.