Added the ability to set linear multi-point equation constraints.

This commit is contained in:
Ryan Latture 2016-10-03 22:03:02 -07:00
parent f9974df76f
commit 692e8878e6
15 changed files with 293 additions and 36 deletions

View File

@ -21,11 +21,11 @@ Alternatively, I can package pre-built binaries if there is interest.
## Introduction ##
This contains a C++ implementation of 3D Euler-Bernoulli beam element formulation.
An analysis can be formulated in C++, through a command line interface via a config file, or using the graphical user interface.
An analysis can be formulated in C++, through a command line interface via a configuration file (in JSON format), or using the graphical user interface.
### Method 1: Using C++ ###
An analysis consists of the `fea::Job` as well as any boundary conditions (`fea::BC`), prescribed nodal forces (`fea::Force`), and ties (`fea::Tie`).
The latter of which ties to nodes together via a linear springs between all translational and rotational degrees of freedom.
An analysis consists of the `fea::Job` as well as any boundary conditions (`fea::BC`), prescribed nodal forces (`fea::Force`), ties (`fea::Tie`) and equation constraints (`fea::Equation`).
Ties to nodes together via a linear springs between all translational and rotational degrees of freedom, and equation constraints allow linear multi-point constraints to be applied to the model.
The `fea::Options` struct can be used to request results of the analysis be written to disk as well as modify various aspect of the analysis.
#### Forming the job ####
@ -156,7 +156,10 @@ The `fea::Summary` object can return a report of the analysis in the form of a s
// form an empty vector of ties since none were prescribed
std::vector<fea::Tie> tie_list;
fea::Summary summary = fea::solve(job, node_list, elem_list, bc_list, force_list, tie_list, opts);
// also create an empty list of equations as none were prescribed
std::vector<fea::Equation> eqn_list;
fea::Summary summary = fea::solve(job, node_list, elem_list, bc_list, force_list, tie_list, eqn_list, opts);
// print a report of the analysis
std::cout << summary.fullReport() << std::endl;
@ -173,6 +176,7 @@ Model parameters
BCs : 6
Ties : 0
Forces : 0
Equations : 0
Total time 0ms
Assembly time : 0ms
@ -191,6 +195,7 @@ Nodal Forces
Maximum : Node 1 DOF 1 Value 1.000
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
#### Ties ####
Ties are enforced by placing linear springs between all degrees of freedom for 2 nodes.
To form a tie specify the 2 nodes that will be linked as well as the spring constants for translational and rotational degrees of freedom.
@ -242,6 +247,24 @@ fea::Tie tie1(nn1, nn2, lmult, rmult);
std::vector<fea::Tie> tie_list = {tie1};
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
#### Equations ####
Equations are linear multi-point constraints that are applied to nodal degrees of freedom.
Each equation is composed of a list of terms that sum to zero, e.g. `t1 + t2 + t3 ... = 0`, where `tn` is the `n`th term.
Each term specifies the node number, degree of freedom and coefficient.
The node number and degree of freedom specify which nodal variable (either nodal displacement or rotation) is involved with the equation constraint,
and coefficient is multiplied by the specified nodal variable when forming the equation.
Note, the equation sums to zero, so in order to specify that 2 nodal degrees of freedom are equal their coefficients should be equal and opposite.
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~{.cpp}
// Create an empty equation
fea::Equation eqn;
// Stipulate that the x and y displacement for the first node must be equal
unsigned int node_number = 0;
eqn.terms.push_back(fea::Equation::Term(node_number, fea::DOF::DISPLACEMENT_X, 1.0));
eqn.terms.push_back(fea::Equation::Term(node_number, fea::DOF::DISPLACEMENT_Y, -1.0);
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
### Method 2: Using the command line interface ###
After using CMake to build the targets, an executable will be created that provide a command line interface (CLI) to the beam element code.
Once in the build directory, navigate to the `bin` folder containing fea_cmd. running `./fea_cmd -h` from the terminal will show the help
@ -251,12 +274,13 @@ An example is shown below.
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~{.txt}
{
"nodes" : "path/to/nodes.csv",
"elems" : "path/to/elems.csv",
"props" : "path/to/props.csv",
"bcs" : "path/to/bcs.csv",
"forces" : "path/to/forces.csv",
"ties" : "path/to/ties.csv",
"nodes" : "path/to/nodes.csv",
"elems" : "path/to/elems.csv",
"props" : "path/to/props.csv",
"bcs" : "path/to/bcs.csv",
"forces" : "path/to/forces.csv",
"ties" : "path/to/ties.csv",
"equations" : "path/to/equations.csv",
"options" : {
"epsilon" : 1.0E-14,
"csv_delimiter" : ",",
@ -275,7 +299,7 @@ An example is shown below.
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
The use of a JSON document avoids the need to set each of these options using command line options, which can become tedious when running multiple jobs.
The "nodes", "elems", and "props" keys are required. Keys "bcs", "forces", and "ties" are optional--if not provided the analysis will assume none were prescribed.
The "nodes", "elems", and "props" keys are required. Keys "bcs", "forces", "ties" and "equations" are optional--if not provided the analysis will assume none were prescribed.
If the "options" key is not provided the analysis will run with the default options.
Any of all of the "options" keys presented above can be used to customize the analysis.
If a key is not provided the default value is used in its place.
@ -330,7 +354,7 @@ elN_EA,elN_EIz,elN_EIy,elN_GJ,elN_nvec_x_comp,elN_nvec_y_comp,elN_nvec_z_comp
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
where each entry is a double and each line has 7 entries.
The "bcs" and "forces" CSV files have the same format.
The "bcs" and "forces" CSV files have the same format as each other.
Each line specifies the node number, degree of freedom, and value:
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~{.txt}
@ -358,6 +382,24 @@ tieN_node_num1,tieN_node_num2,tieN_lmult,tieN_rmult
where `lmult` is the spring constant for the translational degrees of freedom
and `rmult` is the spring constant for the rotational degrees of freedom.
Equation constraints are specified by a series of 3 items (representing a single term) repeated until the desired number of terms are created.
Each term is defined by the node index, degree of freedom for the specified node, and the coefficient that will multiply the nodal degree of freedom.
For example a single equation constraint that specifies that the x and y displacements for the first node must remain equal is given by:
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~{.txt}
0,0,1,0,1,-1
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
in general the equations CSV file will be resemble:
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~{.txt}
eq1_term1_node,eq1_term1_dof,eq1_term1_coeff,...
eq2_term1_node,eq2_term1_dof,eq2_term1_coeff,...
...
...
...
eqN_term1_node,eqN_term1_dof,eqN_term1_coeff,...
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
## Contact ##

View File

@ -68,6 +68,9 @@ int main(int argc, char *argv[])
// initialize empty vector of ties
std::vector<Tie> ties;
// initialize empty vector of equation constraints
std::vector<Equation> equations;
// initialize vector of prescribed forces
std::vector<Force> forces;
@ -75,9 +78,9 @@ int main(int argc, char *argv[])
Options opts;
// solve for nodal displacements
Summary summary = solve(job, bcs, forces, ties, opts);
Summary summary = solve(job, bcs, forces, ties, equations, opts);
// write report to terminal
std::cout << summary.FullReport() << std::endl;
return 0;
}
}

View File

@ -76,13 +76,16 @@ int main(int argc, char *argv[])
// apply force on node at (2,0,0)
std::vector<Force> forces= {Force(3, DOF::DISPLACEMENT_Y, 0.01)};
// initialize empty vector of equation constraints
std::vector<Equation> equations;
// use default options
Options opts;
// solve for nodal displacements
Summary summary = solve(job, bcs, forces, ties, opts);
Summary summary = solve(job, bcs, forces, ties, equations, opts);
// write report to terminal
std::cout << summary.FullReport() << std::endl;
return 0;
}
}

View File

@ -173,6 +173,13 @@ void MainWindow::setTiesText() {
}
}
void MainWindow::setEquationsText() {
QString filename = QFileDialog::getOpenFileName(this, tr("Select equations"), QDir::currentPath(), tr("Data files (*.txt *.csv);; All files (*)"));
if (!filename.isEmpty()) {
equationsLineEdit->setText(filename);
}
}
void MainWindow::updateProgressText() {
if (progress) {
progress->setLabelText(QString(feaProcess->readAllStandardOutput()));
@ -271,6 +278,15 @@ void MainWindow::createChooseFilesGroupBox()
row_counter++);
connect(loadTiesButton, &QPushButton::clicked, this, &MainWindow::setTiesText);
equationsLineEdit = new QLineEdit();
loadEquationsButton = new QPushButton("Equations");
loadEquationsButton->setStatusTip(tr("Choose file containing equation constraints"));
initializeChooseFilesRow(glayout,
equationsLineEdit,
loadEquationsButton,
row_counter++);
connect(loadEquationsButton, &QPushButton::clicked, this, &MainWindow::setEquationsText);
chooseFilesGroupBox->setLayout(glayout);
}

View File

@ -42,6 +42,7 @@ private slots:
void setBCsText();
void setForcesText();
void setTiesText();
void setEquationsText();
void updateProgressText();
private:
@ -89,6 +90,7 @@ private:
QPushButton *loadBCsButton;
QPushButton *loadForcesButton;
QPushButton *loadTiesButton;
QPushButton *loadEquationsButton;
QPushButton *submitButton;
QLineEdit *nodesLineEdit;
@ -97,6 +99,7 @@ private:
QLineEdit *bcsLineEdit;
QLineEdit *forcesLineEdit;
QLineEdit *tiesLineEdit;
QLineEdit *equationsLineEdit;
QCheckBox *nodalDispCheckBox;
QCheckBox *nodalForcesCheckBox;

View File

@ -214,6 +214,69 @@ namespace fea {
node_number_1(_node_number_1), node_number_2(_node_number_2), lmult(_lmult), rmult(_rmult) { };
};
/**
* @brief A linear multipoint constraint.
* @details Equation constraints are defined by a series of terms that
* sum to zero, e.g. `t1 + t2 + t3 ... = 0`, where `tn` is the `n`th term.
* Each term specifies the node number, degree of freedom and coefficient.
* The node number and degree of freedom specify which nodal variable
* (either nodal displacement or rotation) is involved with the equation
* constraint, and coefficient is multiplied by the specified nodal variable
* when forming the equation. Note, the equation sums to zero, so in order
* to specify that 2 nodal degrees of freedom are equal their coefficients
* should be equal and opposite.
*
* @code
* // Create an empty equation
* fea::Equation eqn;
*
* // Stipulate that the x and y displacement for the first node must be equal
* unsigned int node_number = 0;
* eqn.terms.push_back(fea::Equation::Term(node_number, fea::DOF::DISPLACEMENT_X, 1.0));
* eqn.terms.push_back(fea::Equation::Term(node_number, fea::DOF::DISPLACEMENT_Y, -1.0);
* @endcode
*/
struct Equation {
/**
* @brief A single term in the equation constraint.
* @details Each term defines the node number, degree of freedom and
* coefficient.
*/
struct Term {
unsigned int node_number;/**<Index of the node in the node list.*/
unsigned int dof;/**<Degree of freedom. @sa fea::DOF*/
double coefficient;/**<Coefficient to multiply the nodal variable by.*/
/**
* Default constructor
*/
Term() : node_number(0), dof(0), coefficient(0) {};
/**
* @brief Constructor
* @param node_number. `unsigned int`. Index of the node within the node list.
* @param dof. `unsigned int`. Degree of freedom for the specified node.
* @param coefficient. `double`. coefficient to multiply the corresponding nodal variable by.
*/
Term(unsigned int _node_number, unsigned int _dof, double _coefficient) :
node_number(_node_number), dof(_dof), coefficient(_coefficient) {};
};
std::vector<Term> terms;/**<A list of terms that sum to zero.*/
/**
* Default constructor
*/
Equation() : terms(0) {};
/**
* @brief Constructor
* @param terms. `std::vector<Term>`. A list of terms that sum to zero.k
*/
Equation(const std::vector<Term> &_terms) : terms(_terms) {};
};
/**
* @brief An element of the mesh. Contains the indices of the two `fea::Node`'s that form the element as well
* as the properties of the element given by the `fea::Props` struct.

View File

@ -85,6 +85,14 @@ namespace fea {
*/
std::vector<Tie> createTieVecFromJSON(const rapidjson::Document &config_doc);
/**
* Parses the file indicated by the "equations" key in `config_doc` into a vector of `fea::Equation`'s.
*
* @param config_doc `rapidjson::Document`. Document storing the file name containing the prescribed forces.
* @return Equation constraints. `std::vector<Equation>`.
*/
std::vector<Equation> createEquationVecFromJSON(const rapidjson::Document &config_doc);
/**
* Creates vectors of `fea::Node`'s and `fea::Elem`'s from the files specified in `config_doc`. A
* `fea::Job` is created from the node and element vectors and returned.

View File

@ -114,6 +114,11 @@ namespace fea {
*/
unsigned long num_ties;
/**
* The number of equation constraints in the analysis.
*/
unsigned long num_eqns;
/**
* The resultant nodal displacement from the FE analysis.
* `nodal_displacements` is a 2D vector where each row

View File

@ -179,6 +179,8 @@ namespace fea {
*/
void loadBCs(SparseMat &Kg, ForceVector &force_vec, const std::vector<BC> &BCs, unsigned int num_nodes);
void loadEquations(SparseMat &Kg, const std::vector<Equation> &equations, unsigned int num_nodes, unsigned int num_bcs);
/**
* @brief Loads any tie constraints into the set of triplets that will become the global stiffness matrix.
* @details Tie constraints are enforced via linear springs between the 2 specified
@ -228,6 +230,7 @@ namespace fea {
const std::vector<BC> &BCs,
const std::vector<Force> &forces,
const std::vector<Tie> &ties,
const std::vector<Equation> &equations,
const Options &options);
} // namespace fea

View File

@ -46,9 +46,14 @@ fea::Summary runAnalysis(const rapidjson::Document &config_doc) {
forces = fea::createForceVecFromJSON(config_doc);
}
std::vector<fea::Equation> equations;
if (config_doc.HasMember("equations")) {
equations = fea::createEquationVecFromJSON(config_doc);
}
fea::Options options = fea::createOptionsFromJSON(config_doc);
return fea::solve(job, bcs, forces, ties, options);
return fea::solve(job, bcs, forces, ties, equations, options);
}
int main(int argc, char *argv[]) {

View File

@ -181,15 +181,36 @@ namespace fea {
return ties_out;
}
std::vector<Equation> createEquationVecFromJSON(const rapidjson::Document &config_doc) {
std::vector< std::vector<double> > eqns_vec;
fea::createVectorFromJSON(config_doc, "equations", eqns_vec);
std::vector<Equation> eqns_out(eqns_vec.size());
for (size_t i = 0; i < eqns_vec.size(); ++i) {
if (eqns_vec[i].size() % 3 != 0) {
throw std::runtime_error(
(boost::format("Row %d in equations does not specify [node number,dof,coefficient,...] for each term.") %
i).str()
);
}
Equation eqn;
for (size_t j = 0; j < eqns_vec[i].size() / 3; ++j) {
eqn.terms.push_back(Equation::Term(
(unsigned int) eqns_vec[i][3 * j],
(unsigned int) eqns_vec[i][3 * j + 1],
eqns_vec[i][3 * j + 2]));
}
eqns_out[i] = eqn;
}
return eqns_out;
}
Job createJobFromJSON(const rapidjson::Document &config_doc) {
try {
std::vector<Node> nodes = createNodeVecFromJSON(config_doc);
std::vector<Elem> elems = createElemVecFromJSON(config_doc);
return Job(nodes, elems);
}
catch (std::runtime_error &e) {
throw;
}
std::vector<Node> nodes = createNodeVecFromJSON(config_doc);
std::vector<Elem> elems = createElemVecFromJSON(config_doc);
return Job(nodes, elems);
}
Options createOptionsFromJSON(const rapidjson::Document &config_doc) {

View File

@ -96,6 +96,7 @@ namespace fea {
num_bcs(0),
num_forces(0),
num_ties(0),
num_eqns(0),
nodal_displacements(0),
nodal_forces(0),
tie_forces(0) {
@ -112,6 +113,7 @@ namespace fea {
fe_params[2] = fe_param_pair("BCs", num_bcs);
fe_params[3] = fe_param_pair("Ties", num_ties);
fe_params[4] = fe_param_pair("Forces ", num_forces);
fe_params[5] = fe_param_pair("Equations ", num_eqns);
// get the maximum number of digits in the model parameters for formatting
int max_digits = 1;

View File

@ -311,6 +311,21 @@ namespace fea {
}
};
void loadEquations(SparseMat &Kg, const std::vector<Equation> &equations, unsigned int num_nodes, unsigned int num_bcs) {
size_t row_idx, col_idx;
const unsigned int dofs_per_elem = DOF::NUM_DOFS;
const unsigned int global_add_idx = dofs_per_elem * num_nodes + num_bcs;
for (size_t i = 0; i < equations.size(); ++i) {
row_idx = global_add_idx + i;
for (size_t j = 0; j < equations[i].terms.size(); ++j) {
col_idx = dofs_per_elem * equations[i].terms[j].node_number + equations[i].terms[j].dof;
Kg.insert(row_idx, col_idx) = equations[i].terms[j].coefficient;
Kg.insert(col_idx, row_idx) = equations[i].terms[j].coefficient;
}
}
};
void loadTies(std::vector<Eigen::Triplet<double> > &triplets, const std::vector<Tie> &ties) {
const unsigned int dofs_per_elem = DOF::NUM_DOFS;
unsigned int nn1, nn2;
@ -384,6 +399,7 @@ namespace fea {
const std::vector<BC> &BCs,
const std::vector<Force> &forces,
const std::vector<Tie> &ties,
const std::vector<Equation> &equations,
const Options &options) {
auto initial_start_time = std::chrono::high_resolution_clock::now();
@ -396,7 +412,7 @@ namespace fea {
const unsigned int dofs_per_elem = DOF::NUM_DOFS;
// calculate size of global stiffness matrix and force vector
const int size = dofs_per_elem * job.nodes.size() + BCs.size();
const unsigned long size = dofs_per_elem * job.nodes.size() + BCs.size() + equations.size();
// construct global stiffness matrix and force vector
SparseMat Kg(size, size);
@ -419,11 +435,15 @@ namespace fea {
// load prescribed boundary conditions into stiffness matrix and force vector
loadBCs(Kg, force_vec, BCs, job.nodes.size());
if (equations.size() > 0) {
loadEquations(Kg, equations, job.nodes.size(), BCs.size());
}
// load prescribed forces into force vector
if (forces.size() > 0) {
loadForces(force_vec, forces);
}
// compress global stiffness matrix since all non-zero values have been added.
Kg.prune(1.e-14);
Kg.makeCompressed();

View File

@ -185,7 +185,8 @@ TEST_F(beamFEATest, AssemblesGlobalStiffness) {
TEST_F(beamFEATest, CorrectNodalDisplacementsNoTies) {
std::vector<Tie> ties;
Summary summary = solve(JOB_L_BRACKET, BCS_L_BRACKET, FORCES_L_BRACKET, ties, Options());
std::vector<Equation> equations;
Summary summary = solve(JOB_L_BRACKET, BCS_L_BRACKET, FORCES_L_BRACKET, ties, equations, Options());
// the first 4 rows check nodal displacements
// the last row is associated with the reaction
@ -236,7 +237,9 @@ TEST_F(beamFEATest, CorrectNodalDisplacementsWithStiffTies) {
std::vector<Force> forces;
Summary summary = solve(job_tie, bcs, forces, ties, Options());
std::vector<Equation> equations;
Summary summary = solve(job_tie, bcs, forces, ties, equations, Options());
// The verification program used to generate expected values outputs data as floats
std::vector<std::vector<double> > expected = {{0., 0., 0., 0., 0., 0.},
@ -255,13 +258,41 @@ TEST_F(beamFEATest, CorrectNodalDisplacementsWithStiffTies) {
}
}
TEST_F(beamFEATest, CorrectDisplacementWithEquationsCantileverBeam) {
std::vector<BC> bcs = {BC(0, 0, 0.1),
BC(0, 1, 0.0),
BC(0, 2, 0.0),
BC(0, 3, 0.0),
BC(0, 4, 0.0),
BC(0, 5, 0.0)};
std::vector<Tie> ties;
std::vector<Force> forces;
Equation eqn;
eqn.terms.push_back(Equation::Term(0, 0, 1));
eqn.terms.push_back(Equation::Term(1, 0, 1));
std::vector<Equation> equations = {eqn};
Summary summary = solve(JOB_CANTILEVER, bcs, forces, ties, equations, Options());
std::vector<std::vector<double> > expected = {{0.1, 0.,0., 0., 0., 0.},
{-0.1, 0., 0., 0., 0., 0.}};
for (size_t i = 0; i < summary.nodal_displacements.size(); ++i) {
for (size_t j = 0; j < summary.nodal_displacements[i].size(); ++j)
EXPECT_DOUBLE_EQ(expected[i][j], summary.nodal_displacements[i][j]);
}
}
// This is simple cantilever beam with a point load at the
// free end. This checks the end displacements match the analytical results.
TEST_F(beamFEATest, CorrectTipDisplacementCantileverBeam) {
std::vector<Tie> ties;
std::vector<Equation> equations;
Summary summary = solve(JOB_CANTILEVER, BCS_CANTILEVER, FORCES_CANTILEVER, ties, Options());
Summary summary = solve(JOB_CANTILEVER, BCS_CANTILEVER, FORCES_CANTILEVER, ties, equations, Options());
// the first row checks nodal displacements of fixed node are zero
// the second row checks the analytical result for tip displacement
@ -282,17 +313,15 @@ TEST_F(beamFEATest, CorrectTipDisplacementCantileverBeam) {
TEST_F(beamFEATest, CorrectTipForcesCantileverBeam) {
std::vector<Tie> ties;
std::vector<Equation> equations;
std::vector<Force> forces;
std::vector<BC> bcs = BCS_CANTILEVER;
bcs.push_back(BC(1, 0, 0.1));
bcs.push_back(BC(1, 1, 0.1));
Options opts;
opts.save_report = true;
opts.save_nodal_forces = true;
opts.save_nodal_displacements = true;
Summary summary = solve(JOB_CANTILEVER, bcs, forces, ties, opts);
Summary summary = solve(JOB_CANTILEVER, bcs, forces, ties, equations, opts);
std::vector<std::vector<double> > expected = {{-0.1, -0.3, 0., 0., 0., -0.3},
{0.1, 0.3, 0., 0.0, 0.0, 0.0}};
@ -337,11 +366,12 @@ TEST_F(beamFEATest, CorrectDisplacementWeakTies) {
std::vector<Tie> ties = {Tie(1, 2, 0.01, 0.01)};
std::vector<Force> forces;
std::vector<Equation> equations;
Options opts;
opts.epsilon = 1e-10;
Summary summary = solve(job_tie, bcs, forces, ties, opts);
Summary summary = solve(job_tie, bcs, forces, ties, equations, opts);
std::vector<std::vector<double> > expected = {{0.0, 0.0, 0.0, 0.0, 0.0, 0.0},
{0.0, 0.0, 0.0, 0.0, 0.0, 0.0},
@ -385,11 +415,12 @@ TEST_F(beamFEATest, CorrectForcesWeakTies) {
std::vector<Tie> ties = {Tie(1, 2, 0.01, 0.01)};
std::vector<Force> forces;
std::vector<Equation> equations;
Options opts;
opts.epsilon = 1e-10;
Summary summary = solve(job_tie, bcs, forces, ties, opts);
Summary summary = solve(job_tie, bcs, forces, ties, equations, opts);
std::vector<std::vector<double> > expected = {{0.005, 0.0, 0.0, 0.005, 0.0, 0.0}};

View File

@ -230,6 +230,38 @@ TEST(SetupTest, CreatesCorrectTiesFromJSON) {
}
}
TEST(SetupTest, CreatesCorrectEquationsFromJSON) {
std::string equations_file = "CreatesCorrectEquations.csv";
std::string json = "{\"equations\":\"" + equations_file + "\"}\n";
std::string filename = "CreatesCorrectEquations.json";
writeStringToTxt(filename, json);
rapidjson::Document doc = parseJSONConfig(filename);
std::vector<std::vector<double> > expected = {{1, 2, 3, 4, 5, 6},
{7, 8, 9, 10, 11, 12}};
CSVParser csv;
csv.write(equations_file, expected, 1, ",");
std::vector<Equation> equations = createEquationVecFromJSON(doc);
for (size_t i = 0; i < equations.size(); ++i) {
for (size_t j = 0; j < equations[i].terms.size() / 3; ++j) {
EXPECT_EQ((unsigned int) expected[i][3 * j], equations[i].terms[j].node_number);
EXPECT_EQ((unsigned int) expected[i][3 * j + 1], equations[i].terms[j].dof);
EXPECT_DOUBLE_EQ(expected[i][3 * j + 2], equations[i].terms[j].coefficient);
}
}
if (std::remove(filename.c_str()) != 0) {
std::cerr << "Error removing test csv file " << filename << ".\n";
}
if (std::remove(equations_file.c_str()) != 0) {
std::cerr << "Error removing test csv file " << equations_file << ".\n";
}
}
TEST(SetupTest, CreatesCorrectJobFromJSON) {
std::string elems_file = "CreatesCorrectJob_elems.csv";
std::string props_file = "CreatesCorrectJob_props.csv";