Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Driver.hpp
Go to the documentation of this file.
1 //==============================================================================
2 //
3 // Seismo Virtual Laboratory
4 // Module for Serial and Parallel Analysis of seismic
5 // wave propagation and soil-structure interaction simulation
6 // Copyright (C) 2018-2021, The California Institute of Technology
7 // All Rights Reserved.
8 //
9 // Commercial use of this program without express permission of the California
10 // Institute of Technology, is strictly prohibited. See file "COPYRIGHT" in
11 // main directory for information on usage and redistribution, and for a
12 // DISCLAIMER OF ALL WARRANTIES.
13 //
14 //==============================================================================
15 //
16 // Written by:
17 // Danilo S. Kusanovic (dkusanov@caltech.edu)
18 // Elnaz E. Seylabi (elnaze@unr.edu)
19 //
20 // Supervised by:
21 // Domniki M. Asimaki (domniki@caltech.edu)
22 //
23 // References :
24 // [1]
25 //
26 // Description:
27 //This file contains the "parsing" functions to create, modify and
28 //analyze a finite element model provided in JSON format.
29 //------------------------------------------------------------------------------
30 
31 #ifndef _DRIVER_HPP_
32 #define _DRIVER_HPP_
33 
34 #include <stdio.h>
35 #include <stdlib.h>
36 #include <cstring>
37 #include <fstream>
38 #include <iostream>
39 #include <algorithm>
40 
41 #include "Viscous1DLinear.hpp"
42 #include "Elastic1DLinear.hpp"
43 #include "Hertzian1DLinear.hpp"
44 #include "Elastic2DPlaneStress.hpp"
45 #include "Elastic2DPlaneStrain.hpp"
46 #include "PlasticPlaneStrainJ2.hpp"
47 #include "PlasticPlaneStrainBA.hpp"
48 #include "Elastic3DLinear.hpp"
49 #include "Plastic1DJ2.hpp"
50 #include "Plastic3DJ2.hpp"
51 #include "Plastic3DBA.hpp"
52 
53 #include "Elastic1DGap.hpp"
54 #include "Plastic1DGap.hpp"
55 #include "Elastic1DFiber.hpp"
56 #include "Steel1DFiber.hpp"
57 #include "Concrete1DFiber.hpp"
58 
59 #include "Lin2DRectangular.hpp"
60 #include "Lin3DRectangular.hpp"
61 #include "Lin2DAngle.hpp"
62 #include "Lin3DAngle.hpp"
63 #include "Lin2DChannel.hpp"
64 #include "Lin3DChannel.hpp"
65 #include "Lin2DTee.hpp"
66 #include "Lin3DTee.hpp"
67 #include "Lin2DWideFlange.hpp"
68 #include "Lin3DWideFlange.hpp"
69 #include "Lin2DCircular.hpp"
70 #include "Lin3DCircular.hpp"
71 #include "Lin2DRectangularTube.hpp"
72 #include "Lin3DRectangularTube.hpp"
73 #include "Lin2DCircularTube.hpp"
74 #include "Lin3DCircularTube.hpp"
75 #include "Lin3DThinArea.hpp"
76 #include "Lin2DUserDefined.hpp"
77 #include "Lin3DUserDefined.hpp"
78 
79 #include "Fib3DLineSection.hpp"
80 
81 #include "ZeroLength1D.hpp"
82 #include "lin2DTruss2.hpp"
83 #include "kin2DTruss2.hpp"
84 #include "lin3DTruss2.hpp"
85 #include "kin3DTruss2.hpp"
86 #include "lin2DTruss3.hpp"
87 #include "lin3DTruss3.hpp"
88 #include "lin2DFrame2.hpp"
89 #include "kin2DFrame2.hpp"
90 #include "lin3DFrame2.hpp"
91 #include "lin2DTria3.hpp"
92 #include "lin2DTria6.hpp"
93 #include "lin2DQuad4.hpp"
94 #include "lin2DQuad8.hpp"
95 #include "PML2DQuad4.hpp"
96 #include "kin2DQuad4.hpp"
97 #include "PML2DQuad8.hpp"
98 #include "lin3DShell4.hpp"
99 #include "lin3DTetra4.hpp"
100 #include "lin3DTetra10.hpp"
101 #include "lin3DHexa8.hpp"
102 #include "kin3DHexa8.hpp"
103 #include "PML3DHexa8.hpp"
104 #include "lin3DHexa20.hpp"
105 #include "EQlin2DQuad4.hpp"
106 #include "TIEQlin2DQuad4.hpp"
107 #include "UnxBoucWen2DLink.hpp"
108 #include "UnxBoucWen3DLink.hpp"
109 #include "HDRBYamamoto2DLink.hpp"
110 #include "HDRBYamamoto3DLink.hpp"
111 #include "null2DFrame2.hpp"
112 #include "null3DFrame2.hpp"
113 
114 #include "StaticAnalysis.hpp"
115 #include "DynamicAnalysis.hpp"
116 
117 #include "Linear.hpp"
118 #include "NewtonRaphson.hpp"
119 //#include "GeneralDisplacementPath.hpp"
120 
121 #include "QuasiStatic.hpp"
122 #include "NewmarkBeta.hpp"
123 #include "CompositeBathe.hpp"
124 #include "CentralDifference.hpp"
125 #include "ExtendedNewmarkBeta.hpp"
126 
127 #include "EigenSolver.hpp"
128 #include "MumpsSolver.hpp"
129 #include "PetscSolver.hpp"
130 
131 #include "RSJparser.hpp"
132 #include "Definitions.hpp"
133 #include "Profiler.hpp"
134 
141 
144 template<typename T>
145 void
146 sortVector(std::vector<T>& v){
147  std::sort(v.begin(), v.end());
148 }
149 
155 template<typename T>
156 std::vector<T>
157 setOperation(std::vector<T>& v1, std::vector<T>& v2, std::string op){
158  //The resulting vector
159  std::vector<T> v;
160 
161  //The possible options for set operations
162  if( strcasecmp(op.c_str(),"DIFFERENCE") == 0) {
163  std::set_difference(v1.begin(), v1.end(), v2.begin(), v2.end(), std::inserter(v, v.begin()));
164  }
165  else if(strcasecmp(op.c_str(),"INTERSECTION") == 0){
166  std::set_intersection(v1.begin(), v1.end(), v2.begin(), v2.end(), back_inserter(v));
167  }
168  else if(strcasecmp(op.c_str(),"UNION") == 0){
169  std::set_union(v1.begin(), v1.end(), v2.begin(), v2.end(), std::back_inserter(v));
170  }
171 
172  //Removes possible repeated values in result vector
173  if(v.size() != 0){
174  auto last = std::unique(v.begin(), v.end());
175  v.erase(last, v.end());
176  }
177 
178  return v;
179 }
180 
185 template<typename T>
186 std::vector<T>
187 GetIDsFromJSON(RSJresource& jsonFile, std::string Name){
188  //Number of Entities in Property:
189  unsigned int num = jsonFile[Name].as_object().size();
190  std::vector<T> Tags(num);
191 
192  //Goes Over the identifiers of the JSON Property
193  unsigned int k = 0;
194  for(auto it = jsonFile[Name].as_object().begin(); it != jsonFile[Name].as_object().end(); ++it){
195  //Property identifier
196  unsigned int Tag = std::stoi(it->first);
197  Tags[k] = Tag;
198  k++;
199  }
200 
201  //Sorts the identifier list
202  sortVector<T>(Tags);
203 
204  return Tags;
205 }
206 
211 template<typename T>
212 std::vector<T>
213 GetIDsFromMESH(std::shared_ptr<Mesh>& mesh, std::string Name){
214  //The resulting ID vector
215  std::vector<T> Tags = mesh->GetVectorIDs<T>(Name);
216 
217  //Sorts the identifier list
218  sortVector<T>(Tags);
219 
220  return Tags;
221 }
222 
228 template<typename T>
229 std::map<std::string, std::vector<T> >
230 Entities2Update(std::shared_ptr<Mesh>& mesh, RSJresource& jsonFile, std::string Name){
231  //Entity Objects To Update
232  std::vector<T> v1 = GetIDsFromJSON<T>(jsonFile, Name);
233  std::vector<T> v2 = GetIDsFromMESH<T>(mesh, Name);
234 
235  //Map with indexes for which Node Objects needs to be added, deleted, or modified
236  std::map<std::string, std::vector<T> > Tags;
237  Tags["add"] = setOperation<T>(v1,v2,"Difference");
238  Tags["del"] = setOperation<T>(v2,v1,"Difference");
239  Tags["mod"] = setOperation<T>(v1,v2,"Intersection");
240 
241  return Tags;
242 }
243 
249 std::string
250 GetPartitionName(std::string theFile, int k, bool cond){
251  //Auxiliary variables.
252  std::string auxName = theFile;
253  std::string subDomainMesh;
254 
255  //Conver process number into string.
256  std::stringstream process;
257  process << k;
258 
259  //Replace strings.
260  size_t pos;
261  pos = auxName.find("$");
262  if(pos !=std::string::npos)
263  auxName.replace(pos, std::string("$").length(), process.str());
264 
265  //Modified File Name.
266  if(cond){
267  subDomainMesh = filePath + "/" + auxName;
268  }
269  else{
270  subDomainMesh = auxName;
271  }
272 
273  return subDomainMesh;
274 }
275 
280 std::string
281 GetSpacedName(std::string theFile, std::string toReplace){
282  //Auxiliary variables.
283  std::string auxName = theFile;
284  std::string subDomainMesh;
285 
286  //Replace strings.
287  size_t pos;
288  pos = auxName.find("~");
289  while(pos != std::string::npos){
290  auxName.replace(pos, std::string("~").length(), toReplace);
291  pos = auxName.find("~");
292  }
293 
294  //Modified File Name.
295  subDomainMesh = auxName;
296 
297  return subDomainMesh;
298 }
299 
303 void
304 UpdateNodes(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
305  //Node Identifiers To Be Updated
306  std::map<std::string, std::vector<unsigned int> > Tags = Entities2Update<unsigned int>(theMesh, jsonFile, "Nodes");
307 
308  //Nodes to be Removed from Mesh Object
309  for(unsigned int k = 0; k < Tags["del"].size(); k++){
310  unsigned int Tag = Tags["del"][k];
311  theMesh->DelNode(Tag);
312  }
313 
314  //Nodes to be Added To Mesh Object
315  for(unsigned int k = 0; k < Tags["add"].size(); k++){
316  unsigned int Tag = Tags["add"][k];
317 
318  //Node identifier and number of degree-of-freedom:
319  std::string nTag = std::to_string(Tag);
320  unsigned int nDofs = jsonFile["Nodes"][nTag]["ndof"].as<int>();
321 
322  //Allocate memory for Free/Total degree of freedom lists.
323  std::vector<int> freeDofs(nDofs,0);
324  std::vector<int> totalDofs(nDofs,0);
325 
326  //Allocate memory for the coordinate vector.
327  Eigen::VectorXd coordinates(nDimensions);
328 
329  for(int k = 0; k < jsonFile["Nodes"][nTag]["freedof"].size(); ++k)
330  freeDofs[k] = jsonFile["Nodes"][nTag]["freedof"][k].as<int>();
331 
332  for(int k = 0; k < jsonFile["Nodes"][nTag]["totaldof"].size(); ++k)
333  totalDofs[k] = jsonFile["Nodes"][nTag]["totaldof"][k].as<int>();
334 
335  for(int k = 0; k < jsonFile["Nodes"][nTag]["coords"].size(); ++k)
336  coordinates(k) = jsonFile["Nodes"][nTag]["coords"][k].as<double>();
337 
338  //Checks if the node is free/fixed.
339  bool IsFixed = false;
340  if(std::find(freeDofs.begin(), freeDofs.end(), -1) != freeDofs.end())
341  IsFixed = true;
342 
343  //Creates a node object.
344  std::shared_ptr<Node> theNode = std::make_shared<Node>(nDofs, coordinates, IsFixed);
345 
346  //Sets preliminary degree-of-freedom numbering.
347  theNode->SetFreeDegreeOfFreedom(freeDofs);
348  theNode->SetTotalDegreeOfFreedom(totalDofs);
349 
350  //Stores the information in node container.
351  theMesh->AddNode(Tag, theNode);
352  }
353 
354  //Nodes to be Modified in Mesh Object
355  if( Tags["mod"].size() > 0){
356  std::map<unsigned int, std::shared_ptr<Node> > theNodes = theMesh->GetNodes();
357 
358  for(unsigned int k = 0; k < Tags["mod"].size(); k++){
359  unsigned int Tag = Tags["mod"][k];
360 
361  //Node identifier and number of degree-of-freedom:
362  std::string nTag = std::to_string(Tag);
363  unsigned int nDofs = jsonFile["Nodes"][nTag]["ndof"].as<int>();
364 
365  //Allocate memory and get information from JSON file
366  std::vector<int> freeDofs(nDofs,0);
367  for(int k = 0; k < jsonFile["Nodes"][nTag]["freedof"].size(); ++k)
368  freeDofs[k] = jsonFile["Nodes"][nTag]["freedof"][k].as<int>();
369 
370  std::vector<int> totalDofs(nDofs,0);
371  for(int k = 0; k < jsonFile["Nodes"][nTag]["totaldof"].size(); ++k)
372  totalDofs[k] = jsonFile["Nodes"][nTag]["totaldof"][k].as<int>();
373 
374  Eigen::VectorXd coordinates(nDimensions);
375  for(int k = 0; k < jsonFile["Nodes"][nTag]["coords"].size(); ++k)
376  coordinates(k) = jsonFile["Nodes"][nTag]["coords"][k].as<double>();
377 
378  //Checks if the node is free/fixed.
379  bool IsFixed = false;
380  if(std::find(freeDofs.begin(), freeDofs.end(), -1) != freeDofs.end())
381  IsFixed = true;
382 
383  //Modify the Node's properties according to JSON file.
384  theNodes[Tag]->SetAsFixed(IsFixed);
385  theNodes[Tag]->SetCoordinates(coordinates);
386  theNodes[Tag]->SetFreeDegreeOfFreedom(freeDofs);
387  theNodes[Tag]->SetTotalDegreeOfFreedom(totalDofs);
388  }
389  }
390 }
391 
395 void
396 UpdateMasses(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
397  //Mass Identifiers To Be Updated
398  std::map<std::string, std::vector<unsigned int> > Tags = Entities2Update<unsigned int>(theMesh, jsonFile, "Masses");
399 
400  //Masses to be Removed from Mesh Object
401  for(unsigned int k = 0; k < Tags["del"].size(); k++){
402  unsigned int Tag = Tags["del"][k];
403  theMesh->DelMass(Tag);
404  }
405 
406  //Masses to be Added To Mesh Object
407  for(unsigned int k = 0; k < Tags["add"].size(); k++){
408  //Mass identifier and number of degree-of-freedom:
409  unsigned int Tag = Tags["add"][k];
410  std::string nTag = std::to_string(Tag);
411  unsigned int nDofs = jsonFile["Masses"][nTag]["ndof"].as<int>();
412 
413  Eigen::VectorXd Mass(nDofs);
414  for(unsigned int k = 0; k < nDofs; ++k)
415  Mass(k) = jsonFile["Masses"][nTag]["mass"][k].as<double>();
416 
417  //Adds the nodal mass to the node.
418  theMesh->AddMass(Tag, Mass);
419  }
420 
421  //Masses to be Modified in Mesh Object
422  for(unsigned int k = 0; k < Tags["mod"].size(); k++){
423  unsigned int Tag = Tags["mod"][k];
424  std::string nTag = std::to_string(Tag);
425  unsigned int nDofs = jsonFile["Masses"][nTag]["ndof"].as<int>();
426 
427  Eigen::VectorXd Mass(nDofs);
428  for(unsigned int k = 0; k < nDofs; ++k)
429  Mass(k) = jsonFile["Masses"][nTag]["mass"][k].as<double>();
430 
431  //Sets (replace) the new mass vector.
432  theMesh->AddMass(Tag, Mass);
433  }
434 }
435 
439 void
440 UpdateConstraints(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
441  //Constraint Identifiers To Be Updated
442  std::map<std::string, std::vector<int> > Tags = Entities2Update<int>(theMesh, jsonFile, "Constraints");
443 
444  //Constraint to be Removed from Mesh Object
445  for(unsigned int k = 0; k < Tags["del"].size(); k++){
446  int Tag = Tags["del"][k];
447  theMesh->DelConstraint(Tag);
448  }
449 
450  //Constraint to be Added To Mesh Object
451  for(unsigned int k = 0; k < Tags["add"].size(); k++){
452  //Node identifier and number of degree-of-freedom:
453  int Tag = Tags["add"][k];
454  std::string cTag = std::to_string(Tag);
455 
456  //Constraint and slave identifiers
457  unsigned int stag = jsonFile["Constraints"][cTag]["stag"].as<int>();
458  unsigned int nCombs = jsonFile["Constraints"][cTag]["mtag"].size();
459 
460  //The master identifiers and combinational factors
461  std::vector<double> factors(nCombs);
462  std::vector<unsigned int> mtag(nCombs);
463 
464  for(unsigned int k = 0; k < nCombs; ++k){
465  mtag[k] = jsonFile["Constraints"][cTag]["mtag"][k].as<int>();
466  factors[k] = jsonFile["Constraints"][cTag]["factor"][k].as<double>();
467  }
468 
469  //Creates a constaint object.
470  std::shared_ptr<Constraint> theConstraint = std::make_shared<Constraint>(stag, mtag, factors);
471 
472  //Stores the information in constraint container.
473  theMesh->AddConstraint(Tag, theConstraint);
474  }
475 
476  //Constraint to be Modified in Mesh Object
477  if( Tags["mod"].size() > 0){
478  std::map<int, std::shared_ptr<Constraint> > theConstraints = theMesh->GetConstraints();
479 
480  for(unsigned int k = 0; k < Tags["mod"].size(); k++){
481  int Tag = Tags["mod"][k];
482  std::string cTag = std::to_string(Tag);
483 
484  //Constraint and slave identifiers
485  unsigned int stag = jsonFile["Constraints"][cTag]["stag"].as<int>();
486  unsigned int nCombs = jsonFile["Constraints"][cTag]["mtag"].size();
487 
488  //The master identifiers and combinational factors
489  std::vector<double> factors(nCombs);
490  std::vector<unsigned int> mtag(nCombs);
491 
492  for(unsigned int k = 0; k < nCombs; ++k){
493  mtag[k] = jsonFile["Constraints"][cTag]["mtag"][k].as<int>();
494  factors[k] = jsonFile["Constraints"][cTag]["factor"][k].as<double>();
495  }
496 
497  //Modify the Constraint's properties according to JSON file.
498  theConstraints[Tag]->SetSlaveInformation(stag);
499  theConstraints[Tag]->SetMasterInformation(mtag);
500  theConstraints[Tag]->SetCombinationFactors(factors);
501  }
502  }
503 }
504 
508 void
509 UpdateSupportMotion(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
510  //Support Motions Identifiers To Be Updated
511  std::map<std::string, std::vector<int> > Tags = Entities2Update<int>(theMesh, jsonFile, "Supports");
512 
513  //Support Motions to be Removed from Mesh Object
514  for(unsigned int k = 0; k < Tags["del"].size(); k++){
515  unsigned int Tag = Tags["del"][k];
516  theMesh->DelSupportMotion(Tag);
517  }
518 
519  //Support Motions to be Added To Mesh Object
520  for(unsigned int k = 0; k < Tags["add"].size(); k++){
521  //Support motion name and identifier.
522  unsigned int Tag = Tags["add"][k];
523  std::string nTag = std::to_string(Tag);
524  std::string Name = jsonFile["Supports"][nTag]["type"].as<std::string>();
525 
526  for(int n = 0; n < jsonFile["Supports"][nTag]["dof"].size(); ++n){
527  unsigned int dof;
528  std::vector<double> Xo;
529 
530  if(strcasecmp(Name.c_str(),"CONSTANT") == 0){
531  Xo.resize(1);
532  Xo[0] = jsonFile["Supports"][nTag]["value"][n].as<double>();
533  }
534  else if(strcasecmp(Name.c_str(),"TIMESERIES") == 0){
535  std::string File = jsonFile["Supports"][nTag]["file"][n].as<std::string>();
536  //Loads the time history values into memory.
537  std::ifstream motion(File.c_str());
538 
539  //The File is Opened and Ready to be Loaded.
540  if(motion.is_open()){
541  //Number of time steps.
542  unsigned int nt;
543  motion >> nt;
544 
545  //Time-history load vector.
546  Xo.resize(nt);
547  for(unsigned int j = 0; j < nt; j++)
548  motion >> Xo[j];
549  }
550 
551  motion.close();
552  }
553  dof = jsonFile["Supports"][nTag]["dof"][n].as<int>();
554 
555  //Assign support motion to the node.
556  theMesh->SetSupportMotion(Tag, dof, Xo);
557  }
558  }
559 
560  //No Support Motion Modification. Support should be erased from previous analysis
561 }
562 
566 void
567 UpdateMaterials(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
568  //List of Materials already defined
569  std::vector<unsigned int> Tags = theMesh->GetVectorIDs<unsigned int>("Materials");
570 
571  //Materials to be Added To Mesh Object
572  for(auto it = jsonFile["Materials"].as_object().begin(); it != jsonFile["Materials"].as_object().end(); ++it){
573  //Material name and identifier.
574  unsigned int Tag = std::stoi(it->first);
575 
576  //Whether the material should be created
577  if (std::find(Tags.begin(), Tags.end(), Tag) == Tags.end()) {
578  //Material name and identifier.
579  unsigned int Tag = std::stoi(it->first);
580  std::string Name = it->second["name"].as<std::string>();
581 
582  //Creates a material object.
583  std::unique_ptr<Material> theMaterial;
584 
585  if(strcasecmp(Name.c_str(),"Elastic1DLinear") == 0){
586  double E = it->second["attributes"]["E"].as<double>(0.0);
587  double nu = it->second["attributes"]["nu"].as<double>(0.0);
588  double rho = it->second["attributes"]["rho"].as<double>(0.0);
589 
590  //Instantiate the material object.
591  theMaterial = std::make_unique<Elastic1DLinear>(E, nu, rho);
592  }
593  else if(strcasecmp(Name.c_str(),"Hertzian1DLinear") == 0){
594  double k1 = it->second["attributes"]["k1"].as<double>(0.0);
595  double k2 = it->second["attributes"]["k2"].as<double>(0.0);
596  double k3 = it->second["attributes"]["k3"].as<double>(0.0);
597  double rho = it->second["attributes"]["rho"].as<double>(0.0);
598 
599  //Instantiate the material object.
600  theMaterial = std::make_unique<Hertzian1DLinear>(k1, k2, k3, rho);
601  }
602  else if(strcasecmp(Name.c_str(),"Viscous1DLinear") == 0){
603  double eta = it->second["attributes"]["eta"].as<double>(0.0);
604 
605  //Instantiate the material object.
606  theMaterial =std::make_unique<Viscous1DLinear>(eta);
607  }
608  else if(strcasecmp(Name.c_str(),"Elastic2DPlaneStrain") == 0){
609  double E = it->second["attributes"]["E"].as<double>(0.0);
610  double nu = it->second["attributes"]["nu"].as<double>(0.0);
611  double rho = it->second["attributes"]["rho"].as<double>(0.0);
612 
613  //Instantiate the material object.
614  theMaterial = std::make_unique<Elastic2DPlaneStrain>(E, nu, rho);
615  }
616  else if(strcasecmp(Name.c_str(),"Elastic2DPlaneStress") == 0){
617  double E = it->second["attributes"]["E"].as<double>(0.0);
618  double nu = it->second["attributes"]["nu"].as<double>(0.0);
619  double rho = it->second["attributes"]["rho"].as<double>(0.0);
620 
621  //Instantiate the material object.
622  theMaterial = std::make_unique<Elastic2DPlaneStress>(E, nu, rho);
623  }
624  else if(strcasecmp(Name.c_str(),"Elastic3DLinear") == 0){
625  double E = it->second["attributes"]["E"].as<double>(0.0);
626  double nu = it->second["attributes"]["nu"].as<double>(0.0);
627  double rho = it->second["attributes"]["rho"].as<double>(0.0);
628 
629  //Instantiate the material object.
630  theMaterial = std::make_unique<Elastic3DLinear>(E, nu, rho);
631  }
632  else if(strcasecmp(Name.c_str(),"Plastic1DJ2") == 0){
633  double E = it->second["attributes"]["E"].as<double>(0.0);
634  double nu = it->second["attributes"]["nu"].as<double>(0.0);
635  double rho = it->second["attributes"]["rho"].as<double>(0.0);
636  double K = it->second["attributes"]["k"].as<double>(0.0);
637  double H = it->second["attributes"]["h"].as<double>(0.0);
638  double Sy = it->second["attributes"]["Sy"].as<double>(0.0);
639 
640  //Instantiate the material object.
641  theMaterial = std::make_unique<Plastic1DJ2>(E, nu, rho, K, H, Sy);
642  }
643  else if(strcasecmp(Name.c_str(),"PlasticPlaneStrainJ2") == 0){
644  double K = it->second["attributes"]["K"].as<double>(0.0);
645  double G = it->second["attributes"]["G"].as<double>(0.0);
646  double rho = it->second["attributes"]["rho"].as<double>(0.0);
647  double H = it->second["attributes"]["h"].as<double>(0.0);
648  double Sy = it->second["attributes"]["Sy"].as<double>(0.0);
649  double beta = it->second["attributes"]["beta"].as<double>(0.0);
650 
651  //Instantiate the material object.
652  theMaterial = std::make_unique<PlasticPlaneStrainJ2>(K, G, rho, H, beta, Sy);
653  }
654  else if(strcasecmp(Name.c_str(),"PlasticPlaneStrainBA") == 0){
655  double K = it->second["attributes"]["K"].as<double>(0.0);
656  double G = it->second["attributes"]["G"].as<double>(0.0);
657  double rho = it->second["attributes"]["rho"].as<double>(0.0);
658  double h0 = it->second["attributes"]["h0"].as<double>(0.0);
659  double h = it->second["attributes"]["h"].as<double>(0.0);
660  double m = it->second["attributes"]["m"].as<double>(0.0);
661  double Su = it->second["attributes"]["Su"].as<double>(0.0);
662  double beta = it->second["attributes"]["beta"].as<double>(0.0);
663 
664  //Instantiate the material object.
665  theMaterial = std::make_unique<PlasticPlaneStrainBA>(K, G, rho, h0, h, m, Su, beta);
666  }
667  else if(strcasecmp(Name.c_str(),"Plastic3DJ2") == 0){
668  double K = it->second["attributes"]["K"].as<double>(0.0);
669  double G = it->second["attributes"]["G"].as<double>(0.0);
670  double rho = it->second["attributes"]["rho"].as<double>(0.0);
671  double h = it->second["attributes"]["h"].as<double>(0.0);
672  double Sy = it->second["attributes"]["Sy"].as<double>(0.0);
673  double beta = it->second["attributes"]["beta"].as<double>(0.0);
674 
675  //Instantiate the material object.
676  theMaterial = std::make_unique<Plastic3DJ2>(K, G, rho, h, beta, Sy);
677  }
678  else if(strcasecmp(Name.c_str(),"Plastic3DBA") == 0){
679  double K = it->second["attributes"]["K"].as<double>(0.0);
680  double G = it->second["attributes"]["G"].as<double>(0.0);
681  double rho = it->second["attributes"]["rho"].as<double>(0.0);
682  double h0 = it->second["attributes"]["h0"].as<double>(0.0);
683  double h = it->second["attributes"]["h"].as<double>(0.0);
684  double m = it->second["attributes"]["m"].as<double>(0.0);
685  double Su = it->second["attributes"]["Su"].as<double>(0.0);
686  double beta = it->second["attributes"]["beta"].as<double>(0.0);
687 
688  //Instantiate the material object.
689  theMaterial = std::make_unique<Plastic3DBA>(K, G, rho, h0, h, m, Su, beta);
690  }
691  else if(strcasecmp(Name.c_str(),"Elastic1DFiber") == 0){
692  double E = it->second["attributes"]["E"].as<double>(0.0);
693  double nu = it->second["attributes"]["nu"].as<double>(0.0);
694  double rho = it->second["attributes"]["rho"].as<double>(0.0);
695 
696  //Instantiate the material (fiber) object.
697  theMaterial = std::make_unique<Elastic1DFiber>(E, nu, rho);
698  }
699  else if(strcasecmp(Name.c_str(),"Elastic1DGap") == 0){
700  double E = it->second["attributes"]["E"].as<double>(0.0);
701  double gap = it->second["attributes"]["gap"].as<double>(0.0);
702  double behavior = it->second["attributes"]["behavior"].as<bool>(false);
703 
704  //Instantiate the material (fiber) object.
705  theMaterial = std::make_unique<Elastic1DGap>(E, gap, behavior);
706  }
707  else if(strcasecmp(Name.c_str(),"Plastic1DGap") == 0){
708  double E = it->second["attributes"]["E"].as<double>(0.0);
709  double fy = it->second["attributes"]["fy"].as<double>(0.0);
710  double eta = it->second["attributes"]["ratio"].as<double>(0.0);
711  double gap = it->second["attributes"]["gap"].as<double>(0.0);
712  double behavior = it->second["attributes"]["behavior"].as<bool>(false);
713 
714  //Instantiate the material (fiber) object.
715  theMaterial = std::make_unique<Plastic1DGap>(E, fy, gap, eta, behavior);
716  }
717  else if(strcasecmp(Name.c_str(),"Steel1DFiber") == 0){
718  double E = it->second["attributes"]["E"].as<double>(0.0);
719  double fy = it->second["attributes"]["fy"].as<double>(0.0);
720  double b = it->second["attributes"]["b"].as<double>(0.0);
721  double R0 = it->second["attributes"]["R0"].as<double>(15.00);
722  double cR1 = it->second["attributes"]["cR1"].as<double>(0.925);
723  double cR2 = it->second["attributes"]["cR2"].as<double>(0.150);
724  double a1 = it->second["attributes"]["a1"].as<double>(0.0);
725  double a2 = it->second["attributes"]["a2"].as<double>(1.0);
726  double a3 = it->second["attributes"]["a3"].as<double>(0.0);
727  double a4 = it->second["attributes"]["a4"].as<double>(1.0);
728  double nu = it->second["attributes"]["nu"].as<double>(0.33);
729  double rho = it->second["attributes"]["rho"].as<double>(0.0);
730 
731  //Instantiate the material (fiber) object.
732  theMaterial = std::make_unique<Steel1DFiber>(fy, E, b, R0, cR1, cR2, a1, a2, a3, a4, nu, rho);
733  }
734  else if(strcasecmp(Name.c_str(),"Concrete1DFiber") == 0){
735  double fc = it->second["attributes"]["fc"].as<double>();
736  double ecc = it->second["attributes"]["ecc"].as<double>(-0.002);
737  double fcu = it->second["attributes"]["fcu"].as<double>();
738  double ecu = it->second["attributes"]["ecu"].as<double>(-0.012);
739  double ratio = it->second["attributes"]["ratio"].as<double>(0.1);
740  double ft = it->second["attributes"]["ft"].as<double>(0.0);
741  double Et = it->second["attributes"]["Et"].as<double>(0.0);
742  double nu = it->second["attributes"]["nu"].as<double>(0.25);
743  double rho = it->second["attributes"]["rho"].as<double>(0.0);
744 
745  //Instantiate the material (fiber) object.
746  theMaterial = std::make_unique<Concrete1DFiber>(fc, ecc, fcu, ecu, ratio, ft, Et, nu, rho);
747  }
748 
749  //TODO: Add more material models here.
750 
751  //Stores the information in material container.
752  theMesh->AddMaterial(Tag, theMaterial);
753  }
754  }
755 }
756 
760 void
761 UpdateSections(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
762  //List of Sections already defined
763  std::vector<unsigned int> Tags = theMesh->GetVectorIDs<unsigned int>("Sections");
764 
765  //Materials to be Added To Mesh Object
766  for(auto it = jsonFile["Sections"].as_object().begin(); it != jsonFile["Sections"].as_object().end(); ++it){
767  //Material name and identifier.
768  unsigned int Tag = std::stoi(it->first);
769 
770  //Whether the material should be created
771  if (std::find(Tags.begin(), Tags.end(), Tag) == Tags.end()) {
772  //Section name and identifier.
773  unsigned int Tag = std::stoi(it->first);
774  std::string Name = it->second["name"].as<std::string>();
775  std::string Model = it->second["model"].as<std::string>();
776 
777  //Creates a section object.
778  std::unique_ptr<Section> theSection;
779 
780  if(strcasecmp(Model.c_str(),"Plain") == 0){
781  //Insertion point, section rotation and material identifier
782  double theta = it->second["attributes"]["theta"].as<double>(0.0);
783  unsigned int ip = it->second["attributes"]["ip"].as<int>(10);
784  unsigned int matTag = it->second["attributes"]["material"].as<int>();
785 
786  if(strcasecmp(Name.c_str(),"Lin2DRectangular") == 0){
787  double h = it->second["attributes"]["h"].as<double>();
788  double b = it->second["attributes"]["b"].as<double>();
789 
790  //Instantiate the section object.
791  theSection = std::make_unique<Lin2DRectangular>(h, b, theMesh->GetMaterial(matTag), theta, ip);
792  }
793  else if(strcasecmp(Name.c_str(),"Lin3DRectangular") == 0){
794  double h = it->second["attributes"]["h"].as<double>();
795  double b = it->second["attributes"]["b"].as<double>();
796 
797  //Instantiate the section object.
798  theSection = std::make_unique<Lin3DRectangular>(h, b, theMesh->GetMaterial(matTag), theta, ip);
799  }
800  else if(strcasecmp(Name.c_str(),"Lin2DRectangularTube") == 0){
801  double h = it->second["attributes"]["h"].as<double>();
802  double b = it->second["attributes"]["b"].as<double>();
803  double tw = it->second["attributes"]["tw"].as<double>();
804  double tf = it->second["attributes"]["tf"].as<double>();
805 
806  //Instantiate the section object.
807  theSection = std::make_unique<Lin2DRectangularTube>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
808  }
809  else if(strcasecmp(Name.c_str(),"Lin3DRectangularTube") == 0){
810  double h = it->second["attributes"]["h"].as<double>();
811  double b = it->second["attributes"]["b"].as<double>();
812  double tw = it->second["attributes"]["tw"].as<double>();
813  double tf = it->second["attributes"]["tf"].as<double>();
814 
815  //Instantiate the section object.
816  theSection = std::make_unique<Lin3DRectangularTube>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
817  }
818  else if(strcasecmp(Name.c_str(),"Lin2DCircular") == 0){
819  double r = it->second["attributes"]["r"].as<double>();
820 
821  //Instantiate the section object.
822  theSection = std::make_unique<Lin2DCircular>(r, theMesh->GetMaterial(matTag), theta, ip);
823  }
824  else if(strcasecmp(Name.c_str(),"Lin3DCircular") == 0){
825  double r = it->second["attributes"]["r"].as<double>();
826 
827  //Instantiate the section object.
828  theSection = std::make_unique<Lin3DCircular>(r, theMesh->GetMaterial(matTag), theta, ip);
829  }
830  else if(strcasecmp(Name.c_str(),"Lin2DCircularTube") == 0){
831  double re = it->second["attributes"]["re"].as<double>();
832  double ri = it->second["attributes"]["ri"].as<double>();
833 
834  //Instantiate the section object.
835  theSection = std::make_unique<Lin2DCircularTube>(re, ri, theMesh->GetMaterial(matTag), theta, ip);
836  }
837  else if(strcasecmp(Name.c_str(),"Lin3DCircularTube") == 0){
838  double re = it->second["attributes"]["re"].as<double>();
839  double ri = it->second["attributes"]["ri"].as<double>();
840 
841  //Instantiate the section object.
842  theSection = std::make_unique<Lin3DCircularTube>(re, ri, theMesh->GetMaterial(matTag), theta, ip);
843  }
844  else if(strcasecmp(Name.c_str(),"Lin2DAngle") == 0){
845  double h = it->second["attributes"]["h"].as<double>();
846  double b = it->second["attributes"]["b"].as<double>();
847  double tw = it->second["attributes"]["tw"].as<double>();
848  double tf = it->second["attributes"]["tf"].as<double>();
849 
850  //Instantiate the section object.
851  theSection = std::make_unique<Lin2DAngle>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
852  }
853  else if(strcasecmp(Name.c_str(),"Lin3DAngle") == 0){
854  double h = it->second["attributes"]["h"].as<double>();
855  double b = it->second["attributes"]["b"].as<double>();
856  double tw = it->second["attributes"]["tw"].as<double>();
857  double tf = it->second["attributes"]["tf"].as<double>();
858 
859  //Instantiate the section object.
860  theSection = std::make_unique<Lin3DAngle>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
861  }
862  else if(strcasecmp(Name.c_str(),"Lin2DChannel") == 0){
863  double h = it->second["attributes"]["h"].as<double>();
864  double b = it->second["attributes"]["b"].as<double>();
865  double tw = it->second["attributes"]["tw"].as<double>();
866  double tf = it->second["attributes"]["tf"].as<double>();
867 
868  //Instantiate the section object.
869  theSection = std::make_unique<Lin2DChannel>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
870  }
871  else if(strcasecmp(Name.c_str(),"Lin3DChannel") == 0){
872  double h = it->second["attributes"]["h"].as<double>();
873  double b = it->second["attributes"]["b"].as<double>();
874  double tw = it->second["attributes"]["tw"].as<double>();
875  double tf = it->second["attributes"]["tf"].as<double>();
876 
877  //Instantiate the section object.
878  theSection = std::make_unique<Lin3DChannel>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
879  }
880  if(strcasecmp(Name.c_str(),"Lin2DTee") == 0){
881  double h = it->second["attributes"]["h"].as<double>();
882  double b = it->second["attributes"]["b"].as<double>();
883  double tw = it->second["attributes"]["tw"].as<double>();
884  double tf = it->second["attributes"]["tf"].as<double>();
885 
886  //Instantiate the section object.
887  theSection = std::make_unique<Lin2DTee>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
888  }
889  if(strcasecmp(Name.c_str(),"Lin3DTee") == 0){
890  double h = it->second["attributes"]["h"].as<double>();
891  double b = it->second["attributes"]["b"].as<double>();
892  double tw = it->second["attributes"]["tw"].as<double>();
893  double tf = it->second["attributes"]["tf"].as<double>();
894 
895  //Instantiate the section object.
896  theSection = std::make_unique<Lin3DTee>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
897  }
898  else if(strcasecmp(Name.c_str(),"Lin2DWideFlange") == 0){
899  double h = it->second["attributes"]["h"].as<double>();
900  double b = it->second["attributes"]["b"].as<double>();
901  double tw = it->second["attributes"]["tw"].as<double>();
902  double tf = it->second["attributes"]["tf"].as<double>();
903 
904  //Instantiate the section object.
905  theSection = std::make_unique<Lin2DWideFlange>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
906  }
907  else if(strcasecmp(Name.c_str(),"Lin3DWideFlange") == 0){
908  double h = it->second["attributes"]["h"].as<double>();
909  double b = it->second["attributes"]["b"].as<double>();
910  double tw = it->second["attributes"]["tw"].as<double>();
911  double tf = it->second["attributes"]["tf"].as<double>();
912 
913  //Instantiate the section object.
914  theSection = std::make_unique<Lin3DWideFlange>(h, b, tw, tf, theMesh->GetMaterial(matTag), theta, ip);
915  }
916  else if(strcasecmp(Name.c_str(),"Lin2DUserDefined") == 0){
917  std::vector<double> prop(3, 0.0);
918  prop[0] = it->second["attributes"]["A"].as<double>();
919  prop[1] = it->second["attributes"]["As2"].as<double>(0.0);
920  prop[2] = it->second["attributes"]["I33"].as<double>();
921 
922  //Instantiate the section object.
923  theSection = std::make_unique<Lin2DUserDefined>(prop, theMesh->GetMaterial(matTag), theta);
924  }
925  else if(strcasecmp(Name.c_str(),"Lin3DUserDefined") == 0){
926  std::vector<double> prop(7, 0.0);
927  prop[0] = it->second["attributes"]["A"].as<double>();
928  prop[3] = it->second["attributes"]["J"].as<double>();
929  prop[1] = it->second["attributes"]["As2"].as<double>(0.0);
930  prop[2] = it->second["attributes"]["As3"].as<double>(0.0);
931  prop[4] = it->second["attributes"]["I22"].as<double>();
932  prop[5] = it->second["attributes"]["I33"].as<double>();
933  prop[6] = it->second["attributes"]["I23"].as<double>(0.0);
934 
935  //Instantiate the section object.
936  theSection = std::make_unique<Lin3DUserDefined>(prop, theMesh->GetMaterial(matTag), theta);
937  }
938  else if(strcasecmp(Name.c_str(),"Lin3DThinArea") == 0){
939  double th = it->second["attributes"]["th"].as<double>();
940 
941  //Instantiate the section object.
942  theSection = std::make_unique<Lin3DThinArea>(th, theMesh->GetMaterial(matTag));
943  }
944  }
945  else if(strcasecmp(Model.c_str(),"Fiber") == 0){
946  if(strcasecmp(Name.c_str(),"Fib3DLineSection") == 0){
947  //Initialize the vector with fiber coordinates, and area
948  unsigned int numOfFibers = it->second["attributes"]["fiber"].size();
949  std::vector<double> zi(numOfFibers);
950  std::vector<double> yi(numOfFibers);
951  std::vector<double> Ai(numOfFibers);
952  std::vector<std::unique_ptr<Material> > fibers(numOfFibers);
953 
954  //Generates the Fiber section object
955  for(unsigned int k =0; k < numOfFibers; k++){
956  unsigned int fibTag = it->second["attributes"]["fiber"][k].as<int>();
957  zi[k] = it->second["attributes"]["zi"][k].as<double>();
958  yi[k] = it->second["attributes"]["yi"][k].as<double>();
959  Ai[k] = it->second["attributes"]["Ai"][k].as<double>();
960 
961  theMesh->AddFiber(fibTag, fibers[k]);
962  }
963  //Cross section Shear factor (Ash/A): example 5/6 rectangle, 0.9 circular
964  double kappa2 = it->second["attributes"]["kappa2"].as<double>(1.00);
965  double kappa3 = it->second["attributes"]["kappa3"].as<double>(1.00);
966 
967  //Section Insertion Point bounding Box
968  double h = it->second["attributes"]["h"].as<double>();
969  double b = it->second["attributes"]["b"].as<double>();
970  unsigned int ip = it->second["attributes"]["ip"].as<int>(10);
971 
972  //Section angle of rotation
973  double theta = it->second["attributes"]["theta"].as<double>(0.0);
974 
975  //Instantiate the section object.
976  theSection = std::make_unique<Fib3DLineSection>(h, b, fibers, zi, yi, Ai, kappa2, kappa3, theta, ip);
977  }
978  else if(strcasecmp(Name.c_str(),"Fib3DAreaSection") == 0){
979  //Instantiate the section object.
980  //theSection = std::make_unique<Fib3DAreaSection>(...);
981  }
982  }
983  else if(strcasecmp(Model.c_str(),"General") == 0){
984  //TODO: Parse general sections components.
985  }
986 
987  //TODO: Add more section models here.
988 
989  //Stores the section in the mesh.
990  theMesh->AddSection(Tag, theSection);
991  }
992  }
993 }
994 
998 void
999 UpdateElements(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
1000  //Element Identifiers To Be Updated
1001  std::map<std::string, std::vector<unsigned int> > Tags = Entities2Update<unsigned int>(theMesh, jsonFile, "Elements");
1002 
1003  //Element to be Removed from Mesh Object
1004  for(unsigned int k = 0; k < Tags["del"].size(); k++){
1005  unsigned int Tag = Tags["del"][k];
1006  theMesh->DelElement(Tag);
1007  }
1008 
1009  //Element to be Added To Mesh Object
1010  for(unsigned int k = 0; k < Tags["add"].size(); k++){
1011  unsigned int Tag = Tags["add"][k];
1012 
1013  //Node identifier and number of degree-of-freedom:
1014  std::string eTag = std::to_string(Tag);
1015  std::string Name = jsonFile["Elements"][eTag]["name"].as<std::string>();
1016 
1017  //Node connectivity
1018  unsigned int n = jsonFile["Elements"][eTag]["conn"].size();
1019  std::vector<unsigned int> nodes(n);
1020  for(unsigned int k = 0; k < n; ++k)
1021  nodes[k] = jsonFile["Elements"][eTag]["conn"][k].as<int>();
1022 
1023  //Creates an element object.
1024  std::shared_ptr<Element> theElement;
1025 
1026  if(strcasecmp(Name.c_str(),"lin2DTruss2") == 0){
1027  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1028  double parameters = jsonFile["Elements"][eTag]["attributes"]["area"].as<double>();
1029 
1030  //Instantiate the lin2DTruss2 element.
1031  theElement = std::make_shared<lin2DTruss2>(nodes, theMesh->GetMaterial(matID), parameters);
1032  }
1033  else if(strcasecmp(Name.c_str(),"kin2DTruss2") == 0){
1034  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1035  double parameters = jsonFile["Elements"][eTag]["attributes"]["area"].as<double>();
1036 
1037  //Instantiate the kin2DTruss2 element.
1038  theElement = std::make_shared<kin2DTruss2>(nodes, theMesh->GetMaterial(matID), parameters);
1039  }
1040  else if(strcasecmp(Name.c_str(),"lin3DTruss2") == 0){
1041  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1042  double parameters = jsonFile["Elements"][eTag]["attributes"]["area"].as<double>();
1043 
1044  //Instantiate the lin3DTruss2 element.
1045  theElement = std::make_shared<lin3DTruss2>(nodes, theMesh->GetMaterial(matID), parameters);
1046  }
1047  else if(strcasecmp(Name.c_str(),"kin3DTruss2") == 0){
1048  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1049  double parameters = jsonFile["Elements"][eTag]["attributes"]["area"].as<double>();
1050 
1051  //Instantiate the kin3DTruss2 element.
1052  theElement = std::make_shared<kin3DTruss2>(nodes, theMesh->GetMaterial(matID), parameters);
1053  }
1054  else if(strcasecmp(Name.c_str(),"lin2DTruss3") == 0){
1055  double parameters = jsonFile["Elements"][eTag]["attributes"]["area"].as<double>();
1056  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(3);
1057  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1058  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1059 
1060  //Instantiate the lin2DTruss3 element.
1061  theElement = std::make_shared<lin2DTruss3>(nodes, theMesh->GetMaterial(matID), parameters, Quadrature, nGauss);
1062  }
1063  else if(strcasecmp(Name.c_str(),"lin3DTruss3") == 0){
1064  double parameters = jsonFile["Elements"][eTag]["attributes"]["area"].as<double>();
1065  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(3);
1066  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1067  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1068 
1069  //Instantiate the lin3DTruss3 element.
1070  theElement = std::make_shared<lin3DTruss3>(nodes, theMesh->GetMaterial(matID), parameters, Quadrature, nGauss);
1071  }
1072  else if(strcasecmp(Name.c_str(),"ZeroLength1D") == 0){
1073  unsigned int dir = jsonFile["Elements"][eTag]["attributes"]["dir"].as<int>();
1074  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1075 
1076  //Instantiate the zerolength1D element.
1077  theElement = std::make_shared<ZeroLength1D>(nodes, theMesh->GetMaterial(matID), dir);
1078  }
1079  else if(strcasecmp(Name.c_str(),"UnxBoucWen2DLink") == 0){
1080  std::vector<double> variables(2, 0.0);
1081  std::vector<double> parameters(4, 0.0);
1082  double tol = jsonFile["Elements"][eTag]["attributes"]["tol"].as<double>(1E-06);
1083  unsigned int nmax = jsonFile["Elements"][eTag]["attributes"]["nmax"].as<int>(50);
1084  unsigned int dir = jsonFile["Elements"][eTag]["attributes"]["dir"].as<int>();
1085  unsigned int dim = jsonFile["Elements"][eTag]["attributes"]["dim"].as<int>();
1086 
1087  variables[0] = jsonFile["Elements"][eTag]["attributes"]["Fy"].as<double>();
1088  variables[1] = jsonFile["Elements"][eTag]["attributes"]["k0"].as<double>();
1089  parameters[0] = jsonFile["Elements"][eTag]["attributes"]["alpha"].as<double>(1.0);
1090  parameters[1] = jsonFile["Elements"][eTag]["attributes"]["eta"].as<double>(1.0);
1091  parameters[2] = jsonFile["Elements"][eTag]["attributes"]["beta"].as<double>(0.5);
1092  parameters[3] = jsonFile["Elements"][eTag]["attributes"]["gamma"].as<double>(0.5);
1093 
1094  //Instantiate the UnxBoucWen2DLink element.
1095  theElement = std::make_shared<UnxBoucWen2DLink>(nodes, parameters, variables, dim, dir, tol, nmax);
1096  }
1097  else if(strcasecmp(Name.c_str(),"UnxBoucWen3DLink") == 0){
1098  std::vector<double> variables(2, 0.0);
1099  std::vector<double> parameters(4, 0.0);
1100  double tol = jsonFile["Elements"][eTag]["attributes"]["tol"].as<double>(1E-06);
1101  unsigned int nmax = jsonFile["Elements"][eTag]["attributes"]["nmax"].as<int>(50);
1102  unsigned int dir = jsonFile["Elements"][eTag]["attributes"]["dir"].as<int>();
1103  unsigned int dim = jsonFile["Elements"][eTag]["attributes"]["dim"].as<int>();
1104 
1105  variables[0] = jsonFile["Elements"][eTag]["attributes"]["Fy"].as<double>();
1106  variables[1] = jsonFile["Elements"][eTag]["attributes"]["k0"].as<double>();
1107  parameters[0] = jsonFile["Elements"][eTag]["attributes"]["alpha"].as<double>(1.0);
1108  parameters[1] = jsonFile["Elements"][eTag]["attributes"]["eta"].as<double>(1.0);
1109  parameters[2] = jsonFile["Elements"][eTag]["attributes"]["beta"].as<double>(0.5);
1110  parameters[3] = jsonFile["Elements"][eTag]["attributes"]["gamma"].as<double>(0.5);
1111 
1112  //Instantiate the UnxBoucWen3DLink element.
1113  theElement = std::make_shared<UnxBoucWen3DLink>(nodes, parameters, variables, dim, dir, tol, nmax);
1114  }
1115  else if(strcasecmp(Name.c_str(),"HDRBYamamoto2DLink") == 0){
1116  double De = jsonFile["Elements"][eTag]["attributes"]["De"].as<double>(1.3);
1117  double Di = jsonFile["Elements"][eTag]["attributes"]["Di"].as<double>(0.3);
1118  double Hr = jsonFile["Elements"][eTag]["attributes"]["Hr"].as<double>(0.261);
1119  unsigned int dim = jsonFile["Elements"][eTag]["attributes"]["dim"].as<int>();
1120 
1121  //Instantiate the HDRBYamamoto2DLink element.
1122  theElement = std::make_shared<HDRBYamamoto2DLink>(nodes, De, Di, Hr, dim);
1123  }
1124  else if(strcasecmp(Name.c_str(),"HDRBYamamoto3DLink") == 0){
1125  double De = jsonFile["Elements"][eTag]["attributes"]["De"].as<double>(1.3);
1126  double Di = jsonFile["Elements"][eTag]["attributes"]["Di"].as<double>(0.3);
1127  double Hr = jsonFile["Elements"][eTag]["attributes"]["Hr"].as<double>(0.261);
1128  unsigned int dim = jsonFile["Elements"][eTag]["attributes"]["dim"].as<int>();
1129 
1130  //Instantiate the HDRBYamamoto3DLink element.
1131  theElement = std::make_shared<HDRBYamamoto3DLink>(nodes, De, Di, Hr, dim);
1132  }
1133  else if(strcasecmp(Name.c_str(),"lin2DFrame2") == 0){
1134  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(3);
1135  unsigned int secID = jsonFile["Elements"][eTag]["attributes"]["section"].as<int>();
1136  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1137  std::string Formulation = jsonFile["Elements"][eTag]["attributes"]["formulation"].as<std::string>("Bernoulli");
1138 
1139  //Frame element formulation.
1140  bool Condition = false;
1141  if(strcasecmp(Formulation.c_str(),"Timoshenko") == 0)
1142  Condition = true;
1143 
1144  //Instantiate the lin2DFrame2 element.
1145  theElement = std::make_shared<lin2DFrame2>(nodes, theMesh->GetSection(secID), Condition, Quadrature, nGauss);
1146  }
1147  else if(strcasecmp(Name.c_str(),"kin2DFrame2") == 0){
1148  unsigned int secID = jsonFile["Elements"][eTag]["attributes"]["section"].as<int>();
1149 
1150  //Instantiate the kin2DFrame2 element.
1151  theElement = std::make_shared<kin2DFrame2>(nodes, theMesh->GetSection(secID));
1152  }
1153  else if(strcasecmp(Name.c_str(),"lin3DFrame2") == 0){
1154  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(3);
1155  unsigned int secID = jsonFile["Elements"][eTag]["attributes"]["section"].as<int>();
1156  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1157  std::string Formulation = jsonFile["Elements"][eTag]["attributes"]["formulation"].as<std::string>("Bernoulli");
1158 
1159  //Frame element formulation.
1160  bool Condition = false;
1161  if(strcasecmp(Formulation.c_str(),"Timoshenko") == 0)
1162  Condition = true;
1163 
1164  //Instantiate the lin3DFrame2 element.
1165  theElement = std::make_shared<lin3DFrame2>(nodes, theMesh->GetSection(secID), Condition, Quadrature, nGauss);
1166  }
1167  else if(strcasecmp(Name.c_str(),"lin2DTria3") == 0){
1168  double thickness = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1169  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(4);
1170  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1171  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1172 
1173  //Instantiate the lin2DTria3 element.
1174  theElement = std::make_shared<lin2DTria3>(nodes, theMesh->GetMaterial(matID), thickness, Quadrature, nGauss);
1175  }
1176  else if(strcasecmp(Name.c_str(),"lin2DTria6") == 0){
1177  double thickness = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1178  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(9);
1179  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1180  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1181 
1182  //Instantiate the lin2DTria6 element.
1183  theElement = std::make_shared<lin2DTria6>(nodes, theMesh->GetMaterial(matID), thickness, Quadrature, nGauss);
1184  }
1185  else if(strcasecmp(Name.c_str(),"lin2DQuad4") == 0){
1186  double thickness = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1187  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(4);
1188  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1189  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1190 
1191  //Instantiate the lin2DQuad4 element.
1192  theElement = std::make_shared<lin2DQuad4>(nodes, theMesh->GetMaterial(matID), thickness, Quadrature, nGauss);
1193  }
1194  else if(strcasecmp(Name.c_str(),"lin2DQuad8") == 0){
1195  double thickness = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1196  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(9);
1197  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1198  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1199 
1200  //Instantiate the lin2DQuad8 element.
1201  theElement = std::make_shared<lin2DQuad8>(nodes, theMesh->GetMaterial(matID), thickness, Quadrature, nGauss);
1202  }
1203  else if(strcasecmp(Name.c_str(),"PML2DQuad4") == 0){
1204  std::vector<double> parameters(8);
1205  parameters[0] = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1206  parameters[1] = jsonFile["Elements"][eTag]["attributes"]["n"].as<double>();
1207  parameters[2] = jsonFile["Elements"][eTag]["attributes"]["L"].as<double>();
1208  parameters[3] = jsonFile["Elements"][eTag]["attributes"]["R"].as<double>();
1209  parameters[4] = jsonFile["Elements"][eTag]["attributes"]["x0"][0].as<double>();
1210  parameters[5] = jsonFile["Elements"][eTag]["attributes"]["x0"][1].as<double>();
1211  parameters[6] = jsonFile["Elements"][eTag]["attributes"]["npml"][0].as<double>();
1212  parameters[7] = jsonFile["Elements"][eTag]["attributes"]["npml"][1].as<double>();
1213 
1214  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(4);
1215  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1216  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1217 
1218  //Instantiate the PML2DQuad4 element.
1219  theElement = std::make_shared<PML2DQuad4>(nodes, theMesh->GetMaterial(matID), parameters, Quadrature, nGauss);
1220  }
1221  else if(strcasecmp(Name.c_str(),"PML2DQuad8") == 0){
1222  std::vector<double> parameters(8);
1223  parameters[0] = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1224  parameters[1] = jsonFile["Elements"][eTag]["attributes"]["n"].as<double>();
1225  parameters[2] = jsonFile["Elements"][eTag]["attributes"]["L"].as<double>();
1226  parameters[3] = jsonFile["Elements"][eTag]["attributes"]["R"].as<double>();
1227  parameters[4] = jsonFile["Elements"][eTag]["attributes"]["x0"][0].as<double>();
1228  parameters[5] = jsonFile["Elements"][eTag]["attributes"]["x0"][1].as<double>();
1229  parameters[6] = jsonFile["Elements"][eTag]["attributes"]["npml"][0].as<double>();
1230  parameters[7] = jsonFile["Elements"][eTag]["attributes"]["npml"][1].as<double>();
1231 
1232  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(9);
1233  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1234  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1235 
1236  //Instantiate the PML2DQuad4 element.
1237  theElement = std::make_shared<PML2DQuad8>(nodes, theMesh->GetMaterial(matID), parameters, Quadrature, nGauss);
1238  }
1239  else if(strcasecmp(Name.c_str(),"kin2DQuad4") == 0){
1240  double thickness = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1241  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(4);
1242  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1243  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1244 
1245  //Instantiate the kin2DQuad4 element.
1246  theElement = std::make_shared<kin2DQuad4>(nodes, theMesh->GetMaterial(matID), thickness, Quadrature, nGauss);
1247  }
1248  else if(strcasecmp(Name.c_str(),"lin3DShell4") == 0){
1249  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(9);
1250  unsigned int secID = jsonFile["Elements"][eTag]["attributes"]["section"].as<int>();
1251  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1252 
1253  //Instantiate the lin3DShell4 element.
1254  theElement = std::make_shared<lin3DShell4>(nodes, theMesh->GetSection(secID), Quadrature, nGauss);
1255  }
1256  else if(strcasecmp(Name.c_str(),"lin3DTetra4") == 0){
1257  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(4);
1258  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1259  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1260 
1261  //Instantiate the lin3DTetra4 element.
1262  theElement = std::make_shared<lin3DTetra4>(nodes, theMesh->GetMaterial(matID), Quadrature, nGauss);
1263  }
1264  else if(strcasecmp(Name.c_str(),"lin3DTetra10") == 0){
1265  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(11);
1266  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1267  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1268 
1269  //Instantiate the lin3DTetra10 element.
1270  theElement = std::make_shared<lin3DTetra10>(nodes, theMesh->GetMaterial(matID), Quadrature, nGauss);
1271  }
1272  else if(strcasecmp(Name.c_str(),"lin3DHexa8") == 0){
1273  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(8);
1274  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1275  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1276 
1277  //Instantiate the lin3DHexa8 element.
1278  theElement = std::make_shared<lin3DHexa8>(nodes, theMesh->GetMaterial(matID), Quadrature, nGauss);
1279  }
1280  else if(strcasecmp(Name.c_str(),"kin3DHexa8") == 0){
1281  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(8);
1282  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1283  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1284 
1285  //Instantiate the kin3DHexa8 element.
1286  theElement = std::make_shared<kin3DHexa8>(nodes, theMesh->GetMaterial(matID), Quadrature, nGauss);
1287  }
1288  else if(strcasecmp(Name.c_str(),"PML3DHexa8") == 0){
1289  std::vector<double> parameters(9);
1290  parameters[0] = jsonFile["Elements"][eTag]["attributes"]["n"].as<double>();
1291  parameters[1] = jsonFile["Elements"][eTag]["attributes"]["L"].as<double>();
1292  parameters[2] = jsonFile["Elements"][eTag]["attributes"]["R"].as<double>();
1293  parameters[3] = jsonFile["Elements"][eTag]["attributes"]["x0"][0].as<double>();
1294  parameters[4] = jsonFile["Elements"][eTag]["attributes"]["x0"][1].as<double>();
1295  parameters[5] = jsonFile["Elements"][eTag]["attributes"]["x0"][2].as<double>();
1296  parameters[6] = jsonFile["Elements"][eTag]["attributes"]["npml"][0].as<double>();
1297  parameters[7] = jsonFile["Elements"][eTag]["attributes"]["npml"][1].as<double>();
1298  parameters[8] = jsonFile["Elements"][eTag]["attributes"]["npml"][2].as<double>();
1299 
1300  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(8);
1301  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1302  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1303 
1304  //Instantiate the PML3DHexa8 element.
1305  theElement = std::make_shared<PML3DHexa8>(nodes, theMesh->GetMaterial(matID), parameters, Quadrature, nGauss);
1306  }
1307  else if(strcasecmp(Name.c_str(),"lin3DHexa20") == 0){
1308  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(27);
1309  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1310  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1311 
1312  //Instantiate the lin3DHexa20 element.
1313  theElement = std::make_shared<lin3DHexa20>(nodes, theMesh->GetMaterial(matID), Quadrature, nGauss);
1314  }
1315  else if(strcasecmp(Name.c_str(),"TIEQlin2DQuad4") == 0){
1316  double cf1 = jsonFile["Elements"][eTag]["attributes"]["cf1"].as<double>();
1317  double cf2 = jsonFile["Elements"][eTag]["attributes"]["cf2"].as<double>();
1318  double zref = jsonFile["Elements"][eTag]["attributes"]["zref"].as<double>();
1319  double eref = jsonFile["Elements"][eTag]["attributes"]["eref"].as<double>();
1320  double thickness = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1321  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1322  std::string eType = jsonFile["Elements"][eTag]["attributes"]["type"].as<std::string>("Darandelli");
1323  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(4);
1324  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1325 
1326  //Instantiate the TIEQlin2DQuad4 element.
1327  theElement = std::make_shared<TIEQlin2DQuad4>(nodes, theMesh->GetMaterial(matID), thickness, Quadrature, nGauss, eType, zref, cf1, cf2, eref);
1328  }
1329  else if(strcasecmp(Name.c_str(),"EQlin2DQuad4") == 0){
1330  double cf1 = jsonFile["Elements"][eTag]["attributes"]["cf1"].as<double>();
1331  double cf2 = jsonFile["Elements"][eTag]["attributes"]["cf2"].as<double>();
1332  double zref = jsonFile["Elements"][eTag]["attributes"]["zref"].as<double>();
1333  double thickness = jsonFile["Elements"][eTag]["attributes"]["th"].as<double>();
1334  unsigned int matID = jsonFile["Elements"][eTag]["attributes"]["material"].as<int>();
1335  std::string eType = jsonFile["Elements"][eTag]["attributes"]["type"].as<std::string>("Darandelli");
1336  unsigned int nGauss = jsonFile["Elements"][eTag]["attributes"]["np"].as<int>(4);
1337  std::string Quadrature = jsonFile["Elements"][eTag]["attributes"]["rule"].as<std::string>("GAUSS");
1338 
1339  //Instantiate the EQlin2DQuad4 element.
1340  theElement = std::make_shared<EQlin2DQuad4>(nodes, theMesh->GetMaterial(matID), thickness, Quadrature, nGauss, eType, zref, cf1, cf2);
1341  }
1342  else if(strcasecmp(Name.c_str(),"null2DFrame2") == 0){
1343  //Instantiate the null2DFrame2 element.
1344  theElement = std::make_shared<null2DFrame2>(nodes);
1345  }
1346  else if(strcasecmp(Name.c_str(),"null3DFrame2") == 0){
1347  //Instantiate the null3DFrame2 element.
1348  theElement = std::make_shared<null3DFrame2>(nodes);
1349  }
1350 
1351  //TODO: Add more elements models here.
1352 
1353  //Stores the information in element container.
1354  theMesh->AddElement(Tag, theElement);
1355  }
1356 }
1357 
1361 void
1362 UpdateDampings(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
1363  //Damping Identifiers To Be Updated
1364  std::map<std::string, std::vector<unsigned int> > Tags = Entities2Update<unsigned int>(theMesh, jsonFile, "Dampings");
1365 
1366  //Damping to be Removed from Mesh Object
1367  for(unsigned int k = 0; k < Tags["add"].size(); k++){
1368  unsigned int Tag = Tags["add"][k];
1369 
1370  //Damping name and identifier.
1371  std::string dTag = std::to_string(Tag);
1372  std::string Name = jsonFile["Dampings"][dTag]["name"].as<std::string>();
1373 
1374  //Damping parameters
1375  std::vector<double> parameters;
1376  if(strcasecmp(Name.c_str(),"Free") == 0){
1377  parameters.resize(2);
1378  parameters[0] = 0.0;
1379  parameters[1] = 0.0;
1380  }
1381  else if(strcasecmp(Name.c_str(),"Rayleigh") == 0){
1382  parameters.resize(2);
1383  parameters[0] = jsonFile["Dampings"][dTag]["attributes"]["am"].as<double>(0.0);
1384  parameters[1] = jsonFile["Dampings"][dTag]["attributes"]["ak"].as<double>(0.0);
1385  }
1386  else if(strcasecmp(Name.c_str(),"Caughey") == 0){
1387  //TODO: Caughey-viscous damping is not implemented yet.
1388  //parameters.resize(4);
1389  }
1390  else if(strcasecmp(Name.c_str(),"Capped") == 0){
1391  //TODO: Capped-viscous damping is not implemented yet.
1392  }
1393 
1394  //Damping associated to elements
1395  unsigned int n = jsonFile["Dampings"][dTag]["attributes"]["list"].size();
1396  std::vector<unsigned int> eList(n);
1397 
1398  for(unsigned int k = 0; k < n; k++)
1399  eList[k] = jsonFile["Dampings"][dTag]["attributes"]["list"][k].as<int>();
1400 
1401  //Creates a damping object.
1402  std::shared_ptr<Damping> theDamping = std::make_shared<Damping>(Name, parameters);
1403 
1404  //Stores the damping in the mesh.
1405  theMesh->AddDamping(Tag, theDamping);
1406 
1407  //Assign damping to group of elements.
1408  theMesh->SetDamping(Tag, eList);
1409  }
1410 
1411  //Damping to be Modified in Mesh Object
1412  if( Tags["mod"].size() > 0){
1413  std::map<unsigned int, std::shared_ptr<Damping> > theDampings = theMesh->GetDampings();
1414 
1415  for(unsigned int k = 0; k < Tags["mod"].size(); k++){
1416  unsigned int Tag = Tags["mod"][k];
1417 
1418  //Damping name and identifier.
1419  std::string dTag = std::to_string(Tag);
1420  std::string Name = jsonFile["Dampings"][dTag]["name"].as<std::string>();
1421 
1422  //Damping parameters
1423  std::vector<double> parameters;
1424  if(strcasecmp(Name.c_str(),"Free") == 0){
1425  parameters.resize(2);
1426  parameters[0] = 0.0;
1427  parameters[1] = 0.0;
1428  }
1429  else if(strcasecmp(Name.c_str(),"Rayleigh") == 0){
1430  parameters.resize(2);
1431  parameters[0] = jsonFile["Dampings"][dTag]["attributes"]["am"].as<double>(0.0);
1432  parameters[1] = jsonFile["Dampings"][dTag]["attributes"]["ak"].as<double>(0.0);
1433  }
1434  else if(strcasecmp(Name.c_str(),"Caughey") == 0){
1435  //TODO: Caughey-viscous damping is not implemented yet.
1436  //parameters.resize(4);
1437  }
1438  else if(strcasecmp(Name.c_str(),"Capped") == 0){
1439  //TODO: Capped-viscous damping is not implemented yet.
1440  }
1441 
1442  //Update damping associated to elements
1443  unsigned int n = jsonFile["Dampings"][dTag]["attributes"]["list"].size();
1444  std::vector<unsigned int> eList(n);
1445 
1446  for(unsigned int k = 0; k < n; k++)
1447  eList[k] = jsonFile["Dampings"][dTag]["attributes"]["list"][k].as<int>();
1448 
1449  //Update at Damping and Element level.
1450  theDampings[Tag]->SetName(Name);
1451  theDampings[Tag]->SetParameters(parameters);
1452  theMesh->SetDamping(Tag, eList);
1453  }
1454  }
1455 }
1456 
1460 void
1461 UpdateLoads(std::shared_ptr<Mesh> &theMesh, RSJresource& jsonFile){
1462  //Load Identifiers To Be Updated
1463  std::map<std::string, std::vector<unsigned int> > Tags = Entities2Update<unsigned int>(theMesh, jsonFile, "Loads");
1464 
1465  //Load to be Removed from Mesh Object
1466  for(unsigned int k = 0; k < Tags["del"].size(); k++){
1467  int Tag = Tags["del"][k];
1468  theMesh->DelLoad(Tag);
1469  }
1470 
1471  //Load to be Added To Mesh Object
1472  for(unsigned int k = 0; k < Tags["add"].size(); k++){
1473  unsigned int Tag = Tags["add"][k];
1474 
1475  //Load name and identifier.
1476  std::string lTag = std::to_string(Tag);
1477  std::string Name = jsonFile["Loads"][lTag]["name"].as<std::string>();
1478 
1479  //Creates a load object.
1480  std::shared_ptr<Load> theLoad;
1481 
1482  if(strcasecmp(Name.c_str(),"PointLoad") == 0){
1483  //Node identifiers that share this load
1484  unsigned int n = jsonFile["Loads"][lTag]["attributes"]["list"].size();
1485  std::vector<unsigned int> nodes(n);
1486  for(unsigned int k = 0; k < n; k++)
1487  nodes[k] = jsonFile["Loads"][lTag]["attributes"]["list"][k].as<int>();
1488 
1489  //Direction of the applied load
1490  unsigned int ndirs = jsonFile["Loads"][lTag]["attributes"]["dir"].size();
1491  Eigen::VectorXd direction(ndirs);
1492  for(unsigned int k = 0; k < ndirs; k++)
1493  direction(k) = jsonFile["Loads"][lTag]["attributes"]["dir"][k].as<double>();
1494 
1495  //The assigned load pattern function
1496  std::vector<double> value;
1497  std::string pattern = jsonFile["Loads"][lTag]["attributes"]["type"].as<std::string>();
1498  std::string function = jsonFile["Loads"][lTag]["attributes"]["name"].as<std::string>();
1499 
1500  if(strcasecmp(pattern.c_str(),"CONCENTRATED") == 0){
1501  if (strcasecmp(function.c_str(),"CONSTANT") == 0){
1502  //The constant force magnitude
1503  value.resize(1);
1504  value[0] = jsonFile["Loads"][lTag]["attributes"]["mag"].as<double>();
1505 
1506  //Creates the constant point load.
1507  theLoad = std::make_shared<Load>(direction, value, POINTLOAD_CONCENTRATED_CONSTANT);
1508 
1509  //Add the nodes that shares this load.
1510  theLoad->AddNodes(nodes);
1511  }
1512  else if(strcasecmp(function.c_str(),"TIMESERIES") == 0){
1513  //Loads the time history values into memory
1514  std::string pathfile = jsonFile["Loads"][lTag]["attributes"]["file"].as<std::string>();
1515  std::ifstream load(pathfile.c_str());
1516 
1517  //The File is Opened and Ready to be Loaded.
1518  if(load.is_open()){
1519  //Number of time steps.
1520  unsigned int nt;
1521  load >> nt;
1522 
1523  //Time-history load vector.
1524  value.resize(nt);
1525  for(unsigned int j = 0; j < nt; j++)
1526  load >> value[j];
1527  }
1528  load.close();
1529 
1530  //Creates the timeseries point load.
1531  theLoad = std::make_shared<Load>(direction, value, POINTLOAD_CONCENTRATED_DYNAMIC);
1532 
1533  //Add the nodes that shares this load.
1534  theLoad->AddNodes(nodes);
1535  }
1536  }
1537  else if(strcasecmp(pattern.c_str(),"BODY") == 0){
1538  if (strcasecmp(function.c_str(),"CONSTANT") == 0){
1539  //The constant force magnitude
1540  value.resize(1);
1541  value[0] = jsonFile["Loads"][lTag]["attributes"]["mag"].as<double>();
1542 
1543  //Creates the constant point load.
1544  theLoad = std::make_shared<Load>(direction, value, POINTLOAD_BODY_CONSTANT);
1545 
1546  //Add the nodes that shares this load.
1547  theLoad->AddNodes(nodes);
1548  }
1549  else if(strcasecmp(function.c_str(),"TIMESERIES") == 0){
1550  //Loads the time history values into memory
1551  std::string pathfile = jsonFile["Loads"][lTag]["attributes"]["file"].as<std::string>();
1552  std::ifstream load(pathfile.c_str());
1553 
1554  //The File is Opened and Ready to be Loaded.
1555  if(load.is_open()){
1556  //Number of time steps.
1557  unsigned int nt;
1558  load >> nt;
1559 
1560  //Time-history load vector.
1561  value.resize(nt);
1562  for(unsigned int j = 0; j < nt; j++)
1563  load >> value[j];
1564  }
1565  load.close();
1566 
1567  //Creates the timeseries point load.
1568  theLoad = std::make_shared<Load>(direction, value, POINTLOAD_BODY_DYNAMIC);
1569 
1570  //Add the nodes that shares this load.
1571  theLoad->AddNodes(nodes);
1572  }
1573  }
1574  }
1575  else if(strcasecmp(Name.c_str(),"ElementLoad") == 0){
1576  //The function and pattern applied to these elements
1577  std::string pattern = jsonFile["Loads"][lTag]["attributes"]["type"].as<std::string>();
1578  std::string function = jsonFile["Loads"][lTag]["attributes"]["name"].as<std::string>();
1579 
1580  //Node identifiers that share this load
1581  unsigned int n = jsonFile["Loads"][lTag]["attributes"]["list"].size();
1582  std::vector<unsigned int> elements(n);
1583  for(unsigned int k = 0; k < n; k++)
1584  elements[k] = jsonFile["Loads"][lTag]["attributes"]["list"][k].as<int>();
1585 
1586  if(strcasecmp(pattern.c_str(),"SURFACE") == 0){
1587  //Direction of the applied surface load
1588  unsigned int ndirs = jsonFile["Loads"][lTag]["attributes"]["dir"].size();
1589  Eigen::VectorXd direction(ndirs);
1590  for(unsigned int k = 0; k < ndirs; k++)
1591  direction(k) = jsonFile["Loads"][lTag]["attributes"]["dir"][k].as<double>();
1592 
1593  //Face (surfaces) identifiers that share this load
1594  unsigned int m = jsonFile["Loads"][lTag]["attributes"]["list"].size();
1595  std::vector<unsigned int> faces(m);
1596 
1597  for(unsigned int k = 0; k < m; k++){
1598  std::string Tag = std::to_string(elements[k]);
1599  unsigned int sTag = jsonFile["Surfaces"][Tag]["face"].as<int>();
1600  unsigned int eTag = jsonFile["Surfaces"][Tag]["element"].as<int>();
1601 
1602  faces[k] = sTag;
1603  elements[k] = eTag;
1604  }
1605 
1606  //The assigned load pattern function
1607  std::vector<double> value;
1608  if (strcasecmp(function.c_str(),"CONSTANT") == 0){
1609  //The constant force magnitude
1610  value.resize(1);
1611  value[0] = jsonFile["Loads"][lTag]["attributes"]["mag"].as<double>();
1612 
1613  //Creates the constant surface load.
1614  theLoad = std::make_shared<Load>(direction, value, ELEMENTLOAD_SURFACE_CONSTANT);
1615 
1616  //Add the element and face list that shares this load.
1617  theLoad->AddFaces(faces);
1618  theLoad->AddElements(elements);
1619  }
1620  else if(strcasecmp(function.c_str(),"TIMESERIES") == 0){
1621  //TODO: Implement Dynamic Element Surface Load.
1622  }
1623  }
1624  else if(strcasecmp(pattern.c_str(),"BODY") == 0){
1625  //Direction of the applied body load
1626  unsigned int ndirs = jsonFile["Loads"][lTag]["attributes"]["dir"].size();
1627  Eigen::VectorXd direction(ndirs);
1628  for(unsigned int k = 0; k < ndirs; k++)
1629  direction(k) = jsonFile["Loads"][lTag]["attributes"]["dir"][k].as<double>();
1630 
1631  //The assigned load pattern function
1632  std::vector<double> value;
1633  if (strcasecmp(function.c_str(),"CONSTANT") == 0){
1634  //The constant force magnitude
1635  value.resize(1);
1636  value[0] = jsonFile["Loads"][lTag]["attributes"]["mag"].as<double>();
1637 
1638  //Creates the constant body load.
1639  theLoad = std::make_shared<Load>(direction, value, ELEMENTLOAD_BODY_CONSTANT);
1640 
1641  //Add the element list that shares this load.
1642  theLoad->AddElements(elements);
1643  }
1644  else if(strcasecmp(function.c_str(),"TIMESERIES") == 0){
1645  //Loads the time history values into memory
1646  std::string pathfile = jsonFile["Loads"][lTag]["attributes"]["file"].as<std::string>();
1647  std::ifstream load(pathfile.c_str());
1648 
1649  //The File is Opened and Ready to be Loaded.
1650  if(load.is_open()){
1651  //Number of time steps.
1652  unsigned int nt;
1653  load >> nt;
1654 
1655  //Time-history load vector.
1656  value.resize(nt);
1657  for(unsigned int j = 0; j < nt; j++)
1658  load >> value[j];
1659  }
1660  load.close();
1661 
1662  //Creates the variable body load.
1663  theLoad = std::make_shared<Load>(direction, value, ELEMENTLOAD_BODY_DYNAMIC);
1664 
1665  //Add the element list that shares this load.
1666  theLoad->AddElements(elements);
1667  }
1668  }
1669  else if(strcasecmp(pattern.c_str(),"GENERALWAVE") == 0){
1670  //Creates the domain reduction load.
1671  theLoad = std::make_shared<Load>(ELEMENTLOAD_DOMAIN_REDUCTION);
1672  theLoad->AddElements(elements);
1673 
1674  //Gets all nodes involved in elements.
1675  std::map<unsigned int, bool> nodes;
1676  std::vector<unsigned int> elemIDs = theLoad->GetElements();
1677 
1678  //Mesh node and element information.
1679  std::map<unsigned int, std::shared_ptr<Node> > theNodes = theMesh->GetNodes();
1680  std::map<unsigned int, std::shared_ptr<Element> > theElements = theMesh->GetElements();
1681 
1682  for(unsigned int k = 0; k < elemIDs.size(); k++){
1683  std::vector<unsigned int> nodeIDs = theElements[elemIDs[k]]->GetNodes();
1684 
1685  for(unsigned int j = 0; j < nodeIDs.size(); j++)
1686  nodes[nodeIDs[j]] = true;
1687  }
1688 
1689  //Reads/Copy the node information in specified folder.
1690  std::string pathfile = jsonFile["Loads"][lTag]["attributes"]["file"].as<std::string>();
1691 
1692  for(auto it : nodes){
1693  auto &ind = it.first;
1694  std::string LoadFile = GetPartitionName(pathfile, ind, false);
1695 
1696  //Loads the time history values into memory.
1697  std::ifstream load(LoadFile.c_str());
1698 
1699  //The File is Opened and Ready to be Loaded.
1700  if (load.is_open()){
1701  //Number of time-steps and field-components.
1702  bool cond;
1703  unsigned int nt, nFields;
1704  load >> nt >> nFields >> cond;
1705 
1706  //Time-history displacement field.
1707  Eigen::MatrixXd Signal(nt,nFields);
1708 
1709  for(unsigned int i = 0; i < nt; i++){
1710  for(unsigned int j = 0; j < nFields; j++)
1711  load >> Signal(i,j);
1712  }
1713 
1714  //The node is exterior (change the sign).
1715  if(cond)
1716  Signal = -1.00*Signal;
1717 
1718  theLoad->AddDRMCondition(ind,cond);
1719  theNodes[ind]->SetDomainReductionMotion(Signal);
1720  }
1721  load.close();
1722  }
1723  }
1724  }
1725  else if(strcasecmp(Name.c_str(),"SupportMotion") == 0){
1726  //Creates the support motion as point load.
1727  theLoad = std::make_shared<Load>(POINTLOAD_SUPPORT_MOTION);
1728 
1729  unsigned int n = jsonFile["Loads"][lTag]["attributes"]["list"].size();
1730  std::vector<unsigned int> nodes(n);
1731  for(unsigned int k = 0; k < n; k++)
1732  nodes[k] = jsonFile["Loads"][lTag]["attributes"]["list"][k].as<int>();
1733 
1734  theLoad->AddNodes(nodes);
1735  }
1736 
1737  //Stores the load in the mesh.
1738  theMesh->AddLoad(Tag, theLoad);
1739  }
1740 }
1741 
1747 bool
1748 UpdateAnalysis(std::shared_ptr<Mesh> &theMesh, std::unique_ptr<Analysis> &theAnalysis, std::vector<std::shared_ptr<Recorder> > &Recorders, std::map<unsigned int, std::shared_ptr<LoadCombo> > &LoadCombos, std::string InputFile){
1749  //Gets the partitioned mesh file name.
1750  std::string file2open = GetPartitionName(InputFile, rank, true);
1751  file2open = GetSpacedName(file2open, " ");
1752 
1753  //Opens the JSON file
1754  std::ifstream file2stream(file2open.c_str());
1755 
1756  if(file2stream.is_open()){
1757  RSJresource jsonFile(file2stream);
1758 
1759  //The pointes to define the simulation
1760  std::unique_ptr<LinearSystem> theSolver;
1761  std::shared_ptr<Algorithm> theAlgorithm;
1762  std::shared_ptr<Integrator> theIntegrator;
1763 
1764  //Creates the solver
1765  std::string solver = jsonFile["Simulations"]["attributes"]["solver"]["name"].as<std::string>();
1766  unsigned int Tag = jsonFile["Simulations"]["combo"].as<int>();
1767 
1768  if(strcasecmp(solver.c_str(),"EIGEN") == 0){
1769  bool update = jsonFile["Simulations"]["attributes"]["solver"]["update"].as<bool>(false);
1770  theSolver = std::make_unique<EigenSolver>(update);
1771  }
1772  else if(strcasecmp(solver.c_str(),"MUMPS") == 0){
1773  bool update = jsonFile["Simulations"]["attributes"]["solver"]["update"].as<bool>(false);
1774  unsigned int option = jsonFile["Simulations"]["attributes"]["solver"]["option"].as<int>(2);
1775  theSolver = std::make_unique<MumpsSolver>(option, update);
1776  }
1777  else if(strcasecmp(solver.c_str(),"PETSC") == 0){
1778  unsigned int dnz = jsonFile["Simulations"]["attributes"]["solver"]["d_nz"].as<int>();
1779  unsigned int onz = jsonFile["Simulations"]["attributes"]["solver"]["o_nz"].as<int>();
1780  unsigned int ksp = jsonFile["Simulations"]["attributes"]["solver"]["option"].as<int>(3);
1781  double tol = jsonFile["Simulations"]["attributes"]["solver"]["tol"].as<double>(1e-08);
1782  theSolver = std::make_unique<PetscSolver>(dnz, onz, tol, ksp);
1783  }
1784 
1785  //Creates the algorithm
1786  std::string algorithm = jsonFile["Simulations"]["attributes"]["algorithm"]["name"].as<std::string>();
1787 
1788  if(strcasecmp(algorithm.c_str(),"LINEAR") == 0){
1789  theAlgorithm = std::make_shared<Linear>(theSolver, theMesh);
1790  }
1791  else if(strcasecmp(algorithm.c_str(),"NEWTON") == 0){
1792  double ConvergenceTol = jsonFile["Simulations"]["attributes"]["algorithm"]["cnvgtol"].as<double>(1e-06);
1793  unsigned int nMaxIteration = jsonFile["Simulations"]["attributes"]["algorithm"]["nstep"].as<int>(50);
1794  unsigned int ConvergenceTest = jsonFile["Simulations"]["attributes"]["algorithm"]["cnvgtest"].as<int>(4);
1795  theAlgorithm = std::make_shared<NewtonRaphson>(theSolver, theMesh, ConvergenceTol, nMaxIteration, ConvergenceTest);
1796  }
1797  //else if(strcasecmp(algorithm.c_str(),"CONTROLPATH") == 0){
1798  // theAlgorithm = std::make_shared<GeneralDisplacementPath>(theSolver, theIntegrator, theMesh);
1799  //}
1800 
1801  //Creates the integrator
1802  std::string integrator = jsonFile["Simulations"]["attributes"]["integrator"]["name"].as<std::string>();
1803  double dt = jsonFile["Simulations"]["attributes"]["integrator"]["dt"].as<double>(0.0);
1804  double mtol = jsonFile["Simulations"]["attributes"]["integrator"]["mtol"].as<double>(1e-12);
1805  double ktol = jsonFile["Simulations"]["attributes"]["integrator"]["ktol"].as<double>(1e-12);
1806  double ftol = jsonFile["Simulations"]["attributes"]["integrator"]["ftol"].as<double>(1e-12);
1807 
1808  if(strcasecmp(integrator.c_str(),"STATIC") == 0){
1809  theIntegrator = std::make_shared<QuasiStatic>(theMesh, mtol, ktol, ftol);
1810  }
1811  else if(strcasecmp(integrator.c_str(),"NEWMARK") == 0){
1812  theIntegrator = std::make_shared<NewmarkBeta>(theMesh, dt, mtol, ktol, ftol);
1813  }
1814  else if(strcasecmp(integrator.c_str(),"BATHE") == 0){
1815  theIntegrator = std::make_shared<CompositeBathe>(theMesh, dt, mtol, ktol, ftol);
1816  }
1817  else if(strcasecmp(integrator.c_str(),"CENTRALDIFFERENCE") == 0){
1818  theIntegrator = std::make_shared<CentralDifference>(theMesh, dt, mtol, ktol, ftol);
1819  }
1820  else if(strcasecmp(integrator.c_str(),"EXTENDEDNEWMARK") == 0){
1821  theIntegrator = std::make_shared<ExtendedNewmarkBeta>(theMesh, dt, mtol, ktol, ftol);
1822  }
1823 
1824  //Creates Circular dependency between algorithm and integrator.
1825  theIntegrator->SetAlgorithm(theAlgorithm);
1826  theAlgorithm->SetIntegrator(theIntegrator);
1827 
1828  //Creates the analysis case.
1829  std::string analysis = jsonFile["Simulations"]["attributes"]["analysis"]["name"].as<std::string>();
1830  unsigned int nt = jsonFile["Simulations"]["attributes"]["analysis"]["nt"].as<int>(1);
1831 
1832  if(strcasecmp(analysis.c_str(),"STATIC") == 0){
1833  double Factor = 1.00 / (double)nt;
1834  theAlgorithm->SetLoadFactor(Factor);
1835  theAnalysis = std::make_unique<StaticAnalysis>(theMesh, theAlgorithm, theIntegrator, LoadCombos[Tag], nt);
1836  }
1837  else if(strcasecmp(analysis.c_str(),"DYNAMIC") == 0){
1838  theAlgorithm->SetLoadFactor(1.00);
1839  theAnalysis = std::make_unique<DynamicAnalysis>(theMesh, theAlgorithm, theIntegrator, LoadCombos[Tag], nt);
1840  }
1841 
1842  //Sets the Recorders for the analysis.
1843  for(unsigned int k = 0; k < Recorders.size(); k++)
1844  theAnalysis->SetRecorder(Recorders[k]);
1845 
1846  file2stream.close();
1847  }
1848  else{
1849  return true;
1850  }
1851 
1852  return false;
1853 }
1854 
1858 void
1859 UpdateRecorders(std::vector<std::shared_ptr<Recorder> > &Recorders, std::string InputFile){
1860  //Gets the partitioned mesh file name.
1861  std::string file2open = GetPartitionName(InputFile, rank, true);
1862  file2open = GetSpacedName(file2open, " ");
1863 
1864  //Opens the JSON file
1865  std::ifstream file2stream(file2open.c_str());
1866 
1867  if(file2stream.is_open()){
1868  RSJresource jsonFile(file2stream);
1869 
1870  //Recorder Objects
1871  for(auto it = jsonFile["Recorders"].as_object().begin(); it != jsonFile["Recorders"].as_object().end(); ++it){
1872  //Recorder name and identifier.
1873  std::string name = it->second["name"].as<std::string>();
1874  std::string file = it->second["file"].as<std::string>();
1875  unsigned int precision = it->second["ndps"].as<int>();
1876  unsigned int nsample = it->second["nsamp"].as<int>();
1877 
1878  //Creates a Recorder Object.
1879  std::shared_ptr<Recorder> theRecorder;
1880 
1881  if(strcasecmp(name.c_str(),"PARAVIEW") == 0){
1882  unsigned int features = it->second["features"].as<int>();
1883 
1884  //Creates the Paraview Recorder object.
1885  theRecorder = std::make_shared<Recorder>(file, name, features, nsample, precision);
1886  }
1887  else if(strcasecmp(name.c_str(),"SECTION") == 0){
1888  unsigned int nlist = it->second["list"].size();
1889  unsigned int ndims = it->second["coords"].size();
1890  std::string response = it->second["resp"].as<std::string>();
1891 
1892  std::vector<double> coord(ndims);
1893  for(unsigned int k = 0; k < ndims; k++)
1894  coord[k] = it->second["coords"][k].as<double>();
1895 
1896  std::vector<unsigned int> IDs(nlist);
1897  for(unsigned int k = 0; k < nlist; k++)
1898  IDs[k] = it->second["list"][k].as<int>();
1899 
1900  //Creates the Paraview Recorder object.
1901  theRecorder = std::make_shared<Recorder>(file, name, response, coord, IDs, nsample, precision);
1902  }
1903  else{
1904  std::string response = it->second["resp"].as<std::string>();
1905  unsigned int nlist = it->second["list"].size();
1906  std::vector<unsigned int> IDs(nlist);
1907 
1908  for(unsigned int k = 0; k < nlist; k++)
1909  IDs[k] = it->second["list"][k].as<int>();
1910 
1911  //Creates the Point/Element Recorder object.
1912  theRecorder = std::make_shared<Recorder>(file, name, response, IDs, nsample, precision);
1913  }
1914 
1915  //Adds the recorder to the analysis.
1916  Recorders.push_back(theRecorder);
1917  }
1918 
1919  file2stream.close();
1920  }
1921  else{
1922  std::cout << "\x1B[31m ERROR: \x1B[0mThe JSON file in \'Driver::UpdateRecorders()\' in Processor [" << rank << "] couldn't be opened. \n";
1923  }
1924 }
1925 
1929 void
1930 UpdateCombinations(std::map<unsigned int, std::shared_ptr<LoadCombo> > &LoadCombos, std::string InputFile){
1931  //Gets the partitioned mesh file name.
1932  std::string file2open = GetPartitionName(InputFile, rank, true);
1933  file2open = GetSpacedName(file2open, " ");
1934 
1935  //Opens the JSON file
1936  std::ifstream file2stream(file2open.c_str());
1937 
1938  if(file2stream.is_open()){
1939  RSJresource jsonFile(file2stream);
1940 
1941  //Combination Objects
1942  for(auto it = jsonFile["Combinations"].as_object().begin(); it != jsonFile["Combinations"].as_object().end(); ++it){
1943  //Combination name and identifier.
1944  unsigned int Tag = std::stoi(it->first);
1945  std::string Name = it->second["name"].as<std::string>();
1946 
1947  //Load identifiers and combinational factors
1948  std::vector<double> factors;
1949  std::vector<unsigned int> loads;
1950 
1951  if( it->second["attributes"]["load"].exists() ){
1952  unsigned int nCombs = it->second["attributes"]["load"].size();
1953  loads.resize(nCombs);
1954  factors.resize(nCombs);
1955 
1956  for(unsigned int k = 0; k < loads.size(); k++){
1957  loads[k] = it->second["attributes"]["load"][k].as<int>();
1958  factors[k] = it->second["attributes"]["factor"][k].as<double>();
1959  }
1960  }
1961 
1962  //Creates the load combination.
1963  std::shared_ptr<LoadCombo> theCombo = std::make_unique<LoadCombo>(Name, loads, factors);
1964 
1965  //Adds the combination to the analysis.
1966  LoadCombos[Tag] = theCombo;
1967  }
1968 
1969  file2stream.close();
1970  }
1971  else{
1972  std::cout << "\x1B[31m ERROR: \x1B[0mThe JSON file in \'Driver::UpdateCombinations()\' in Processor [" << rank << "] couldn't be opened. \n";
1973  }
1974 }
1975 
1980 bool
1981 UpdateMesh(std::shared_ptr<Mesh> &theMesh, std::string InputFile){
1982  //Gets the partitioned mesh file name.
1983  std::string file2open = GetPartitionName(InputFile, rank, true);
1984  file2open = GetSpacedName(file2open, " ");
1985 
1986  //Opens the JSON file
1987  std::ifstream file2stream(file2open.c_str());
1988 
1989  if(file2stream.is_open()){
1990  RSJresource jsonFile(file2stream);
1991 
1992  //Global Variables
1993  if( jsonFile["Global"].exists() ){
1994  nDimensions = jsonFile["Global"]["ndim"].as<int>();
1995  numberOfTotalDofs = jsonFile["Global"]["ntotal"].as<int>();
1996  numberOfFreeDofs = jsonFile["Global"]["nfree"].as<int>();
1997  UpdateOption = jsonFile["Global"]["update"].as<std::string>("RESTART");
1998  std::string MassForm = jsonFile["Global"]["massform"].as<std::string>("CONSISTENT");
1999 
2000  //Mass formulation for elements mass matrix
2001  if(strcasecmp(MassForm.c_str(),"LUMPED") == 0){
2002  MassFormulation = true;
2003  }
2004  else if(strcasecmp(MassForm.c_str(),"CONSISTENT") == 0){
2005  MassFormulation = false;
2006  }
2007  }
2008 
2009  //Node Objects
2010  UpdateNodes(theMesh, jsonFile);
2011 
2012  //Mass Objects
2013  UpdateMasses(theMesh, jsonFile);
2014 
2015  //Constraint Objects
2016  UpdateConstraints(theMesh, jsonFile);
2017 
2018  //Support Motion Objects
2019  UpdateSupportMotion(theMesh, jsonFile);
2020 
2021  //Material Objects
2022  UpdateMaterials(theMesh, jsonFile);
2023 
2024  //Section Objects
2025  UpdateSections(theMesh, jsonFile);
2026 
2027  //Element Objects
2028  UpdateElements(theMesh, jsonFile);
2029 
2030  //Damping Objects
2031  UpdateDampings(theMesh, jsonFile);
2032 
2033  //Load Objects
2034  UpdateLoads(theMesh, jsonFile);
2035 
2036  //TODO: Include File inside JSON
2037 
2038  file2stream.close();
2039  }
2040  else{
2041  return true;
2042  }
2043 
2044  return false;
2045 }
2046 
2048 void
2050  //Starts profiling this function.
2051  PROFILE_FUNCTION();
2052 
2053  if(driverFile){
2054  //Creates the Mesh subdomain Pointer.
2055  std::shared_ptr<Mesh> theMesh = std::make_shared<Mesh>();
2056 
2057  //Performs the one simulations after the other.
2058  for(unsigned int k = 0; k < fileName.size(); k++){
2059  //Vector of Recorders to store solution.
2060  std::vector<std::shared_ptr<Recorder> > Recorders;
2061 
2062  //Load combination cases.
2063  std::map<unsigned int, std::shared_ptr<LoadCombo> > LoadCombos;
2064 
2065  //Update the Mesh object from json file
2066  bool FailMesh = UpdateMesh(theMesh, fileName[k]);
2067 
2068  if(!FailMesh){
2069  //Sets up memory allocation
2070  theMesh->Initialize();
2071 
2072  //Update the Combination object from json file
2073  UpdateCombinations(LoadCombos, fileName[k]);
2074 
2075  //Update the Combination object from json file
2076  UpdateRecorders(Recorders, fileName[k]);
2077 
2078  //Update the Analysis object from json file
2079  std::unique_ptr<Analysis> theAnalysis;
2080  bool FailAnalysis = UpdateAnalysis(theMesh, theAnalysis, Recorders, LoadCombos, fileName[k]);
2081 
2082  if(!FailAnalysis){
2083  //Runs the current Analysis.
2084  bool StopSimulation = theAnalysis->Analyze();
2085 
2086  if(StopSimulation){
2087  if(rank == 0)
2088  std::cout << " **A PROBLEM WAS ENCOUNTERED. THE ANALYSIS WILL BE TERMINATED**\n\n";
2089  return;
2090  }
2091  }
2092  else{
2093  std::cout << "\x1B[31m ERROR: \x1B[0mThe JSON file \'" << fileName[k] << "\' in Driver::UpdateAnalysis() in Processor [" << rank << "] couldn't be opened. \n";
2094  return;
2095  }
2096  }
2097  else{
2098  std::cout << "\x1B[31m ERROR: \x1B[0mThe JSON file \'" << fileName[k] << "\' in Driver::UpdateMesh() in Processor [" << rank << "] couldn't be opened. \n";
2099  return;
2100  }
2101  }
2102  }
2103 }
2104 
2105 #endif
void UpdateCombinations(std::map< unsigned int, std::shared_ptr< LoadCombo > > &LoadCombos, std::string InputFile)
Updates the Combination Entities required for the Analysis.
Definition: Driver.hpp:1930
RSJobject & as_object(bool force=false)
Definition: RSJparser.hpp:653
This file contains the "Lin3DRectangular" section declarations, which defines a Rectangular geometry ...
#define POINTLOAD_CONCENTRATED_DYNAMIC
Define global load pattern time dependant concentrated load.
Definition: Definitions.hpp:122
This file contains the "lin2DQuad8" linearized eight-node element declarations, which defines an elem...
This file contains the "CentralDifference" integration method, which inncludes any damping and inerti...
This file contains the "Lin2DTee" section declarations, which defines a circular (T) geometry for 2D ...
This file contains the "Steel1DFiber" material declarations, which defines a steel material for 1D...
This file contains the "lin2DQuad4" linearized four-node element declarations, which defines an eleme...
This file contains the "kin2DQuad4" kinematic four-node element declarations, which defines an elemen...
void RunDriverFile()
Runs the User&#39;s JSON Input file.
Definition: Driver.hpp:2049
int rank
The processor number.
This file contains the "Lin3DTee" section declarations, which defines a Tee (T) geometry for 3D analy...
This file contains the "lin3DTetra10" linearized ten-node element declarations, which defines an elem...
#define POINTLOAD_CONCENTRATED_CONSTANT
Define global load pattern constant concentrated load.
Definition: Definitions.hpp:119
This file contains the "Lin2DAngle" section declarations, which defines a Angle geometry in 2D analys...
This file contains the "lin2DTria6" linearized six-node element declarations, which defines an elemen...
This file contains the "Plastic3DBA" material declarations, which defines an triaxial plastic boundin...
This file contains the "null2DFrame2" null frame element declarations, which defines a frame-like ele...
This file contains the "lin3DTruss3" linearized three-node element declarations, which defines an ele...
void UpdateMasses(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Masses in Mesh Object.
Definition: Driver.hpp:396
std::string GetPartitionName(std::string theFile, int k, bool cond)
Sets the partition subdomain tag number.
Definition: Driver.hpp:250
This file contains the "Lin3DCircularTube" section declarations, which defines a hollow circular geom...
This file contains the "Elastic1DLinear" material declarations, which defines an uniaxial isotropic l...
This file contains the "Lin3DChannel" section declarations, which defines a circular (C) geometry for...
unsigned int nDimensions
The problem dimension (1D, 2D, 3D).
This file contains the "Lin2DRectangularTube" section declarations, which defines a hollow rectangula...
This file contains the "Lin3DRectangularTube" section declarations, which defines a hollow rectangula...
#define POINTLOAD_BODY_DYNAMIC
Define global load pattern time dependant for a body load.
Definition: Definitions.hpp:128
void UpdateRecorders(std::vector< std::shared_ptr< Recorder > > &Recorders, std::string InputFile)
Updates the Recorder Entities required for the Analysis.
Definition: Driver.hpp:1859
bool UpdateAnalysis(std::shared_ptr< Mesh > &theMesh, std::unique_ptr< Analysis > &theAnalysis, std::vector< std::shared_ptr< Recorder > > &Recorders, std::map< unsigned int, std::shared_ptr< LoadCombo > > &LoadCombos, std::string InputFile)
Creates a new Analysis object from the json file.
Definition: Driver.hpp:1748
std::vector< T > GetIDsFromMESH(std::shared_ptr< Mesh > &mesh, std::string Name)
Query the MESH object and obtains the Identifier list.
Definition: Driver.hpp:213
This file contains the "Plastic3DJ2" material declarations, which defines an triaxial isotropic plast...
This file contains the "Plastic1DJ2" material declarations, which defines an uniaxial isotropic plast...
This file contains the "Plastic1DGap" material declarations, which defines a plastic gap material for...
This file contains the "Lin3DUserDefined" section declarations, which defines a general geometry for ...
This file contains the "Fib3DLineSection" section declarations, which defines a Rectangular geometry ...
void UpdateSections(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Section Entities in Mesh Object.
Definition: Driver.hpp:761
This file contains the "Static Analysis object" declarations, and updates information in the mesh...
This file contains the "PML3DHexa8" linearized four-node perfectly matched layer element declarations...
This file sets the global variables to be used during SeismoVLAB execution.
This file contains the "null2DFrame2" null frame element declarations, which defines a frame-like ele...
This file contains the "PlasticPlaneStrainBA" material declarations, which defines an biaxial plastic...
#define ELEMENTLOAD_SURFACE_CONSTANT
Define global load pattern for a constant surface load.
Definition: Definitions.hpp:131
This file contains the "lin3DHexa20" linearized twenty-node element declarations, which defines an el...
This file contains the "PML2DQuad8" linearized eight-node perfectly matched layer element declaration...
This file contains the "Lin3DCircular" section declarations, which defines a circular geometry for 3D...
This file contains the "kin2DFrame2" finite kinematic two-node element declarations, which defines a frame element in a finite element mesh.
This file contains the "Lin2DCircular" section declarations, which defines a Circular geometry to be ...
void UpdateLoads(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Loads Entities in Mesh Object.
Definition: Driver.hpp:1461
This file contains the "PML2DQuad4" linearized four-node perfectly matched layer element declarations...
std::string GetSpacedName(std::string theFile, std::string toReplace)
Fix blank spaces provided by user in path (if any).
Definition: Driver.hpp:281
This file contains the "lin2DFrame2" linearized two-node element declarations, which defines a frame ...
This file contains the "NewtonRaphson" algorithm declarations, which solves the non-linear system of ...
std::vector< T > GetIDsFromJSON(RSJresource &jsonFile, std::string Name)
Query the JSON file and obtains the Identifier list.
Definition: Driver.hpp:187
This file contains the implementation of "A Ridiculously Simple JSON Parser for C++." This great parsing class was taken from http://subhrajit.net/.
std::vector< T > setOperation(std::vector< T > &v1, std::vector< T > &v2, std::string op)
Performs union, difference or intersection operation between two vectors.
Definition: Driver.hpp:157
void UpdateElements(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Element Entities in Mesh Object.
Definition: Driver.hpp:999
This file contains the "PlasticPlaneStrainJ2" material declarations, which defines an biaxial plastic...
Definition: RSJparser.hpp:240
This file contains the "Lin2DCircularTube" section declarations, which defines a hollow circular geom...
This file contains the "Elastic2DPlaneStress" material declarations, which defines an biaxial isotrop...
This file contains the "Lin3DWideFlange" section declarations, which defines a WideFlange geometry (I...
This file contains the "EigenSolver" solver declaration; the object performs Cholesky-decomposition o...
This file contains the "lin2DTria3" linearized three-node element declarations, which defines an elem...
void SetRecorder(std::shared_ptr< Recorder > &recorder)
Sets the recorder for the analysis.
This file contains the "lin2DTruss2" linearized two-node element declarations, which defines an eleme...
std::vector< std::string > fileName
The file name to be loaded.
This file contains the "Lin2DWideFlange" section declarations, which defines a circular (I) geometry ...
This file contains the "Elastic3DLinear" material declarations, which defines a triaxial isotropic li...
This file contains the "Lin2DRectangular" section declarations, which defines a rectangular geometry ...
This file contains the "kin2DTruss2" kinematic two-node element declarations, which defines an elemen...
This file contains the "EQlin2DQuad4" equivalent linear four-node element declarations, which defines an element in a finite element mesh.
This file contains the "Hertzian1DLinear" material declarations, which defines an uniaxial isotropic ...
dataType as(const dataType &def=dataType())
Definition: RSJparser.hpp:309
void sortVector(std::vector< T > &v)
Sorts a vector of tag number form smaller to larger.
Definition: Driver.hpp:146
This file contains the "ExtendedNewmarkBeta" integration declaration, the routine includes any dampin...
This file contains the "CompositeBathe" integration declaration, the file integrator includes any dam...
#define ELEMENTLOAD_BODY_CONSTANT
Define global load pattern for a constant body load.
Definition: Definitions.hpp:137
int size
The number of partitions.
virtual bool Analyze()=0
Performs the required analysis on the domain.
std::string to_string(RSJresourceType rt)
Definition: RSJparser.hpp:65
This file contains the "lin3DTetra4" linearized eight-node element declarations, which defines an ele...
This file contains the "Lin2DChannel" section declarations, which defines a channel geometry in 2D an...
This file contains the "Linear" algorithm declarations, which solves the linear system of equations i...
bool driverFile
Whether the driver (JSON) file is provided.
This file contains the "Lin3DAngle" section declarations, which defines a circular (L) geometry for 3...
This file contains the "Elastic1DGap" material declarations, which defines an elastic gap material fo...
#define ELEMENTLOAD_DOMAIN_REDUCTION
Define global load pattern time dependant general wave load.
Definition: Definitions.hpp:143
This file contains the "ZeroLength1D" two-node element declarations, which defines an element of zero...
void UpdateSupportMotion(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Support Motions in Mesh Object.
Definition: Driver.hpp:509
#define POINTLOAD_SUPPORT_MOTION
Define global load pattern constant/time-dependant support motion load.
Definition: Definitions.hpp:146
This file contains the "Concrete1DFiber" material declarations, which defines a concrete material for...
void UpdateConstraints(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Constraint Entities in Mesh Object.
Definition: Driver.hpp:440
This file contains the "lin3DTruss2" linearized two-node element declarations, which defines an eleme...
bool MassFormulation
The element mass formulation.
This file contains the "Profiler object" declarations to compute how long it takes for a program to b...
#define POINTLOAD_BODY_CONSTANT
Define global load pattern for a node body load.
Definition: Definitions.hpp:125
This file contains the "kin3DTruss2" kinematic two-node element declarations, which defines an elemen...
This file contains the "kin3DHexa8" finite kinematic eight-node element declarations, which defines an element in a finite element mesh.
This file contains the "lin3DHexa8" linearized eight-node element declarations, which defines an elem...
This file contains the "NewmarkBeta" integration declaration, the routine include any damping and ine...
This file contains the "lin3DShell4" linearized four-node element declarations, which defines a frame...
void UpdateMaterials(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Material Entities in Mesh Object.
Definition: Driver.hpp:567
This file contains the "Dynamic Analysis object" declarations, and updates information in the mesh...
This file contains the "PetscSolver" solver declaration; this object is an iterative solver...
This file contains the "Lin2DUserDefined" section declarations, which defines a general geometry for ...
This file contains the "lin2DTruss3" linearized three-node element declarations, which defines an ele...
void UpdateNodes(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Node Entities in Mesh Object.
Definition: Driver.hpp:304
This file contains the "lin3DFrame2" linearized two-node element declarations, which defines a frame ...
void UpdateDampings(std::shared_ptr< Mesh > &theMesh, RSJresource &jsonFile)
Updates the Damping Entities in Mesh Object.
Definition: Driver.hpp:1362
#define PROFILE_FUNCTION()
Definition: Profiler.hpp:210
This file contains the "Steel1DFiber" material declarations, which defines a linear material for 1D...
std::string filePath
The folder path where the file is loaded.
This file contains the "EQlin2DQuad4" equivalent linear four-node element declarations, which defines an element in a finite element mesh.
bool UpdateMesh(std::shared_ptr< Mesh > &theMesh, std::string InputFile)
Populate the Mesh object with the json provided entities.
Definition: Driver.hpp:1981
unsigned int numberOfFreeDofs
Total number of free-degree-of-freedom.
std::map< std::string, std::vector< T > > Entities2Update(std::shared_ptr< Mesh > &mesh, RSJresource &jsonFile, std::string Name)
The Entities indexes to be updated.
Definition: Driver.hpp:230
std::string UpdateOption
The update option for member in Mesh.
This file contains the "Viscous1DLinear" material declarations, which defines an uniaxial viscous mat...
This file contains the "Lin3DThinArea" section declarations, which defines an area geometry for 3D an...
This file contains the "Elastic2DPlaneStrain" material declarations, which.
This file contains the "MumpsSolver" solver declaration; this object performs a (MU)ltifrontal (M)ass...
#define ELEMENTLOAD_BODY_DYNAMIC
Define global load pattern for a time dependant body load.
Definition: Definitions.hpp:140
This file contains the "QuasiStatic" integrator declaration, this integrator neglects any damping and...
unsigned int numberOfTotalDofs
Total number of total-degree-of-freedom.