suanPan
Factory.hpp
Go to the documentation of this file.
1 /*******************************************************************************
2  * Copyright (C) 2017-2022 Theodore Chang
3  *
4  * This program is free software: you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation, either version 3 of the License, or
7  * (at your option) any later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with this program. If not, see <http://www.gnu.org/licenses/>.
16  ******************************************************************************/
31 #ifndef FACTORY_HPP
32 #define FACTORY_HPP
33 
34 #include <future>
35 
37 #include <Toolbox/container.h>
38 
39 enum class AnalysisType {
40  NONE,
41  DISP,
42  EIGEN,
43  BUCKLE,
44  STATICS,
45  DYNAMICS
46 };
47 
48 enum class StorageScheme {
49  FULL,
50  BAND,
51  BANDSYMM,
52  SYMMPACK,
53  SPARSE,
55 };
56 
57 enum class SolverType {
58  LAPACK,
59  SPIKE,
60  SUPERLU,
61  MUMPS,
62  CUDA,
63  PARDISO,
64  FGMRES
65 };
66 
67 template<sp_d T> class Factory final {
68  unsigned n_size = 0; // number of degrees of freedom
69  unsigned n_lobw = 0; // low bandwidth
70  unsigned n_upbw = 0; // up bandwidth
71  unsigned n_sfbw = n_lobw + n_upbw; // matrix storage offset
72  unsigned n_rfld = 0; // reference load size
73  unsigned n_mpc = 0; // multipoint constraint size
74  uword n_elem = 0;
75 
76  AnalysisType analysis_type = AnalysisType::NONE; // type of analysis
77  StorageScheme storage_type = StorageScheme::FULL; // type of analysis
78 
79  bool nlgeom = false;
80 
82  SolverSetting<T> setting{};
83 
84  T error = 0.; // error produced by certain solvers
85 
86  Col<T> ninja; // the result from A*X=B
87  Col<T> sushi; // modified right-hand side B
88 
89  suanpan::set<uword> reference_dof;
90  SpMat<T> reference_load;
91 
92  uvec auxiliary_encoding; // for constraints using multiplier method
93  Col<T> auxiliary_lambda; // for constraints using multiplier method
94  Col<T> auxiliary_resistance; // for constraints using multiplier method
95  Col<T> auxiliary_load; // for constraints using multiplier method
96  SpMat<T> auxiliary_stiffness; // for constraints using multiplier method
97 
98  SpCol<T> trial_constraint_resistance;
99  SpCol<T> current_constraint_resistance;
100 
101  T trial_time = 0.; // global trial (pseudo) time
102  T incre_time = 0.; // global incremental (pseudo) time
103  T current_time = 0.; // global current (pseudo) time
104  T pre_time = 0.; // global previous (pseudo) time
105 
106  T strain_energy = 0.;
107  T kinetic_energy = 0.;
108  T viscous_energy = 0.;
109  T complementary_energy = 0.;
110  Col<T> momentum;
111 
112  Col<T> trial_load_factor; // global trial load factor
113  Col<T> trial_load; // global trial load vector
114  Col<T> trial_settlement; // global trial displacement load vector
115  Col<T> trial_resistance; // global trial resistance vector
116  Col<T> trial_damping_force; // global trial damping force vector
117  Col<T> trial_inertial_force; // global trial inertial force vector
118  Col<T> trial_displacement; // global trial displacement vector
119  Col<T> trial_velocity; // global trial velocity vector
120  Col<T> trial_acceleration; // global trial acceleration vector
121  Col<T> trial_temperature; // global trial temperature vector
122 
123  Col<T> incre_load_factor; // global incremental load vector
124  Col<T> incre_load; // global incremental load vector
125  Col<T> incre_settlement; // global incremental displacement load vector
126  Col<T> incre_resistance; // global incremental resistance vector
127  Col<T> incre_damping_force; // global incremental damping force vector
128  Col<T> incre_inertial_force; // global incremental inertial force vector
129  Col<T> incre_displacement; // global incremental displacement vector
130  Col<T> incre_velocity; // global incremental velocity vector
131  Col<T> incre_acceleration; // global incremental acceleration vector
132  Col<T> incre_temperature; // global incremental temperature vector
133 
134  Col<T> current_load_factor; // global current load vector
135  Col<T> current_load; // global current load vector
136  Col<T> current_settlement; // global current displacement load vector
137  Col<T> current_resistance; // global current resistance vector
138  Col<T> current_damping_force; // global current damping force vector
139  Col<T> current_inertial_force; // global current inertial force vector
140  Col<T> current_displacement; // global current displacement vector
141  Col<T> current_velocity; // global current velocity vector
142  Col<T> current_acceleration; // global current acceleration vector
143  Col<T> current_temperature; // global current temperature vector
144 
145  Col<T> pre_load_factor; // global previous load vector
146  Col<T> pre_load; // global previous load vector
147  Col<T> pre_settlement; // global previous displacement load vector
148  Col<T> pre_resistance; // global previous resistance vector
149  Col<T> pre_damping_force; // global previous damping force vector
150  Col<T> pre_inertial_force; // global previous inertial force vector
151  Col<T> pre_displacement; // global previous displacement vector
152  Col<T> pre_velocity; // global previous velocity vector
153  Col<T> pre_acceleration; // global previous acceleration vector
154  Col<T> pre_temperature; // global previous temperature vector
155 
156  shared_ptr<MetaMat<T>> global_mass = nullptr; // global mass matrix
157  shared_ptr<MetaMat<T>> global_damping = nullptr; // global damping matrix
158  shared_ptr<MetaMat<T>> global_stiffness = nullptr; // global stiffness matrix
159  shared_ptr<MetaMat<T>> global_geometry = nullptr; // global geometry matrix
160 
161  std::vector<std::mutex> global_mutex = std::vector<std::mutex>(20);
162 
163  Col<T> eigenvalue; // eigenvalues
164 
165  Mat<T> eigenvector; // eigenvectors
166 
167  template<sp_d T1> friend unique_ptr<MetaMat<T1>> get_basic_container(const Factory<T1>*);
168  template<sp_d T1> friend unique_ptr<MetaMat<T1>> get_matrix_container(const Factory<T1>*);
169 
170 public:
171  const bool initialized = false;
172 
174 
175  void set_size(unsigned);
176  [[nodiscard]] unsigned get_size() const;
177 
178  void set_entry(uword);
179  [[nodiscard]] uword get_entry() const;
180 
181  void set_nlgeom(bool);
182  [[nodiscard]] bool is_nlgeom() const;
183 
185  [[nodiscard]] SolverType get_solver_type() const;
186 
188  [[nodiscard]] const SolverSetting<double>& get_solver_setting() const;
189 
191  [[nodiscard]] AnalysisType get_analysis_type() const;
192 
194  [[nodiscard]] StorageScheme get_storage_scheme() const;
195 
196  [[nodiscard]] bool is_sparse() const;
197 
198  void set_bandwidth(unsigned, unsigned);
199  void get_bandwidth(unsigned&, unsigned&) const;
200 
201  void update_reference_size();
202  void set_reference_size(unsigned);
203  [[nodiscard]] unsigned get_reference_size() const;
204 
205  void update_reference_dof(const uvec&);
207  [[nodiscard]] const suanpan::set<uword>& get_reference_dof() const;
208 
209  void set_error(const T&);
210  const T& get_error() const;
211 
212  /*************************INITIALIZER*************************/
213 
214  int initialize();
215 
216  void initialize_load_factor();
217  void initialize_load();
218  void initialize_settlement();
219  void initialize_resistance();
223  void initialize_velocity();
225  void initialize_temperature();
227 
228  void initialize_mass();
229  void initialize_damping();
230  void initialize_stiffness();
231  void initialize_geometry();
232  void initialize_eigen();
233 
234  /*************************SETTER*************************/
235 
236  void set_ninja(const Col<T>&);
237  void set_sushi(const Col<T>&);
238 
239  void set_mpc(unsigned);
240 
241  void set_reference_load(const SpMat<T>&);
242 
243  void set_trial_time(const T&);
244  void set_trial_load_factor(const Col<T>&);
245  void set_trial_load(const Col<T>&);
246  void set_trial_settlement(const Col<T>&);
247  void set_trial_resistance(const Col<T>&);
248  void set_trial_damping_force(const Col<T>&);
249  void set_trial_inertial_force(const Col<T>&);
250  void set_trial_displacement(const Col<T>&);
251  void set_trial_velocity(const Col<T>&);
252  void set_trial_acceleration(const Col<T>&);
253  void set_trial_temperature(const Col<T>&);
254 
255  void set_incre_time(const T&);
256  void set_incre_load_factor(const Col<T>&);
257  void set_incre_load(const Col<T>&);
258  void set_incre_settlement(const Col<T>&);
259  void set_incre_resistance(const Col<T>&);
260  void set_incre_damping_force(const Col<T>&);
261  void set_incre_inertial_force(const Col<T>&);
262  void set_incre_displacement(const Col<T>&);
263  void set_incre_velocity(const Col<T>&);
264  void set_incre_acceleration(const Col<T>&);
265  void set_incre_temperature(const Col<T>&);
266 
267  void set_current_time(const T&);
268  void set_current_load_factor(const Col<T>&);
269  void set_current_load(const Col<T>&);
270  void set_current_settlement(const Col<T>&);
271  void set_current_resistance(const Col<T>&);
272  void set_current_damping_force(const Col<T>&);
273  void set_current_inertial_force(const Col<T>&);
274  void set_current_displacement(const Col<T>&);
275  void set_current_velocity(const Col<T>&);
276  void set_current_acceleration(const Col<T>&);
277  void set_current_temperature(const Col<T>&);
278 
279  void set_pre_time(const T&);
280  void set_pre_load_factor(const Col<T>&);
281  void set_pre_load(const Col<T>&);
282  void set_pre_settlement(const Col<T>&);
283  void set_pre_resistance(const Col<T>&);
284  void set_pre_damping_force(const Col<T>&);
285  void set_pre_inertial_force(const Col<T>&);
286  void set_pre_displacement(const Col<T>&);
287  void set_pre_velocity(const Col<T>&);
288  void set_pre_acceleration(const Col<T>&);
289  void set_pre_temperature(const Col<T>&);
290 
291  void set_mass(const shared_ptr<MetaMat<T>>&);
292  void set_damping(const shared_ptr<MetaMat<T>>&);
293  void set_stiffness(const shared_ptr<MetaMat<T>>&);
294  void set_geometry(const shared_ptr<MetaMat<T>>&);
295 
296  void set_eigenvalue(const Col<T>&);
297  void set_eigenvector(const Mat<T>&);
298 
299  /*************************GETTER*************************/
300 
301  const Col<T>& get_ninja() const;
302  const Col<T>& get_sushi() const;
303 
304  [[nodiscard]] unsigned get_mpc() const;
305 
306  const SpMat<T>& get_reference_load() const;
307 
308  [[nodiscard]] const uvec& get_auxiliary_encoding() const;
309  const Col<T>& get_auxiliary_lambda() const;
310  const Col<T>& get_auxiliary_resistance() const;
311  const Col<T>& get_auxiliary_load() const;
312  const SpMat<T>& get_auxiliary_stiffness() const;
313 
314  const SpCol<T>& get_trial_constraint_resistance() const;
315  const SpCol<T>& get_current_constraint_resistance() const;
316 
321  const Col<T>& get_momentum();
322 
323  const T& get_trial_time() const;
324  const Col<T>& get_trial_load_factor() const;
325  const Col<T>& get_trial_load() const;
326  const Col<T>& get_trial_settlement() const;
327  const Col<T>& get_trial_resistance() const;
328  const Col<T>& get_trial_damping_force() const;
329  const Col<T>& get_trial_inertial_force() const;
330  const Col<T>& get_trial_displacement() const;
331  const Col<T>& get_trial_velocity() const;
332  const Col<T>& get_trial_acceleration() const;
333  const Col<T>& get_trial_temperature() const;
334 
335  const T& get_incre_time() const;
336  const Col<T>& get_incre_load_factor() const;
337  const Col<T>& get_incre_load() const;
338  const Col<T>& get_incre_settlement() const;
339  const Col<T>& get_incre_resistance() const;
340  const Col<T>& get_incre_damping_force() const;
341  const Col<T>& get_incre_inertial_force() const;
342  const Col<T>& get_incre_displacement() const;
343  const Col<T>& get_incre_velocity() const;
344  const Col<T>& get_incre_acceleration() const;
345  const Col<T>& get_incre_temperature() const;
346 
347  const T& get_current_time() const;
348  const Col<T>& get_current_load_factor() const;
349  const Col<T>& get_current_load() const;
350  const Col<T>& get_current_settlement() const;
351  const Col<T>& get_current_resistance() const;
352  const Col<T>& get_current_damping_force() const;
353  const Col<T>& get_current_inertial_force() const;
354  const Col<T>& get_current_displacement() const;
355  const Col<T>& get_current_velocity() const;
356  const Col<T>& get_current_acceleration() const;
357  const Col<T>& get_current_temperature() const;
358 
359  const T& get_pre_time() const;
360  const Col<T>& get_pre_load_factor() const;
361  const Col<T>& get_pre_load() const;
362  const Col<T>& get_pre_settlement() const;
363  const Col<T>& get_pre_resistance() const;
364  const Col<T>& get_pre_damping_force() const;
365  const Col<T>& get_pre_inertial_force() const;
366  const Col<T>& get_pre_displacement() const;
367  const Col<T>& get_pre_velocity() const;
368  const Col<T>& get_pre_acceleration() const;
369  const Col<T>& get_pre_temperature() const;
370 
371  const shared_ptr<MetaMat<T>>& get_mass() const;
372  const shared_ptr<MetaMat<T>>& get_damping() const;
373  const shared_ptr<MetaMat<T>>& get_stiffness() const;
374  const shared_ptr<MetaMat<T>>& get_geometry() const;
375 
376  std::mutex& get_auxiliary_encoding_mutex();
377  std::mutex& get_auxiliary_resistance_mutex();
378  std::mutex& get_auxiliary_load_mutex();
379  std::mutex& get_auxiliary_stiffness_mutex();
380 
382 
383  std::mutex& get_trial_load_mutex();
384  std::mutex& get_trial_settlement_mutex();
385 
386  std::mutex& get_mass_mutex();
387  std::mutex& get_damping_mutex();
388  std::mutex& get_stiffness_mutex();
389  std::mutex& get_geometry_mutex();
390 
391  const Col<T>& get_eigenvalue() const;
392  const Mat<T>& get_eigenvector() const;
393 
394  /*************************UPDATER*************************/
395 
396  void update_trial_time(const T&);
397  void update_trial_load_factor(const Col<T>&);
398  void update_trial_load(const Col<T>&);
399  void update_trial_settlement(const Col<T>&);
400  void update_trial_resistance(const Col<T>&);
401  void update_trial_damping_force(const Col<T>&);
402  void update_trial_inertial_force(const Col<T>&);
403  void update_trial_displacement(const Col<T>&);
404  void update_trial_velocity(const Col<T>&);
405  void update_trial_acceleration(const Col<T>&);
406  void update_trial_temperature(const Col<T>&);
407 
408  void update_incre_time(const T&);
409  void update_incre_load_factor(const Col<T>&);
410  void update_incre_load(const Col<T>&);
411  void update_incre_settlement(const Col<T>&);
412  void update_incre_resistance(const Col<T>&);
413  void update_incre_damping_force(const Col<T>&);
414  void update_incre_inertial_force(const Col<T>&);
415  void update_incre_displacement(const Col<T>&);
416  void update_incre_velocity(const Col<T>&);
417  void update_incre_acceleration(const Col<T>&);
418  void update_incre_temperature(const Col<T>&);
419 
420  void update_current_time(const T&);
421  void update_current_load_factor(const Col<T>&);
422  void update_current_load(const Col<T>&);
423  void update_current_settlement(const Col<T>&);
424  void update_current_resistance(const Col<T>&);
425  void update_current_damping_force(const Col<T>&);
426  void update_current_inertial_force(const Col<T>&);
427  void update_current_displacement(const Col<T>&);
428  void update_current_velocity(const Col<T>&);
429  void update_current_acceleration(const Col<T>&);
430  void update_current_temperature(const Col<T>&);
431 
432  void update_trial_time_by(const T&);
433  void update_trial_load_factor_by(const Col<T>&);
434  void update_trial_load_by(const Col<T>&);
435  void update_trial_settlement_by(const Col<T>&);
436  void update_trial_resistance_by(const Col<T>&);
437  void update_trial_damping_force_by(const Col<T>&);
438  void update_trial_inertial_force_by(const Col<T>&);
439  void update_trial_displacement_by(const Col<T>&);
440  void update_trial_velocity_by(const Col<T>&);
441  void update_trial_acceleration_by(const Col<T>&);
442  void update_trial_temperature_by(const Col<T>&);
443 
444  void update_incre_time_by(const T&);
445  void update_incre_load_factor_by(const Col<T>&);
446  void update_incre_load_by(const Col<T>&);
447  void update_incre_settlement_by(const Col<T>&);
448  void update_incre_resistance_by(const Col<T>&);
449  void update_incre_damping_force_by(const Col<T>&);
450  void update_incre_inertial_force_by(const Col<T>&);
451  void update_incre_displacement_by(const Col<T>&);
452  void update_incre_velocity_by(const Col<T>&);
453  void update_incre_acceleration_by(const Col<T>&);
454  void update_incre_temperature_by(const Col<T>&);
455 
456  void update_current_time_by(const T&);
457  void update_current_load_factor_by(const Col<T>&);
458  void update_current_load_by(const Col<T>&);
459  void update_current_settlement_by(const Col<T>&);
460  void update_current_resistance_by(const Col<T>&);
461  void update_current_damping_force_by(const Col<T>&);
462  void update_current_inertial_force_by(const Col<T>&);
463  void update_current_displacement_by(const Col<T>&);
464  void update_current_velocity_by(const Col<T>&);
465  void update_current_acceleration_by(const Col<T>&);
466  void update_current_temperature_by(const Col<T>&);
467 
468  /*************************FRIEND*************************/
469 
470  template<sp_d T1> friend Col<T1>& get_ninja(const shared_ptr<Factory<T1>>&);
471  template<sp_d T1> friend Col<T1>& get_sushi(const shared_ptr<Factory<T1>>&);
472 
473  template<sp_d T1> friend suanpan::set<uword>& get_reference_dof(const shared_ptr<Factory<T1>>&);
474  template<sp_d T1> friend SpMat<T1>& get_reference_load(const shared_ptr<Factory<T1>>&);
475 
476  template<sp_d T1> friend uvec& get_auxiliary_encoding(const shared_ptr<Factory<T1>>&);
477  template<sp_d T1> friend Col<T1>& get_auxiliary_lambda(const shared_ptr<Factory<T1>>&);
478  template<sp_d T1> friend Col<T1>& get_auxiliary_resistance(const shared_ptr<Factory<T1>>&);
479  template<sp_d T1> friend Col<T1>& get_auxiliary_load(const shared_ptr<Factory<T1>>&);
480  template<sp_d T1> friend SpMat<T1>& get_auxiliary_stiffness(const shared_ptr<Factory<T1>>&);
481 
482  template<sp_d T1> friend SpCol<T1>& get_trial_constraint_resistance(const shared_ptr<Factory<T1>>&);
483  template<sp_d T1> friend SpCol<T1>& get_current_constraint_resistance(const shared_ptr<Factory<T1>>&);
484 
485  template<sp_d T1> friend T& get_trial_time(const shared_ptr<Factory<T1>>&);
486  template<sp_d T1> friend Col<T1>& get_trial_load_factor(const shared_ptr<Factory<T1>>&);
487  template<sp_d T1> friend Col<T1>& get_trial_load(const shared_ptr<Factory<T1>>&);
488  template<sp_d T1> friend Col<T1>& get_trial_settlement(const shared_ptr<Factory<T1>>&);
489  template<sp_d T1> friend Col<T1>& get_trial_resistance(const shared_ptr<Factory<T1>>&);
490  template<sp_d T1> friend Col<T1>& get_trial_damping_force(const shared_ptr<Factory<T1>>&);
491  template<sp_d T1> friend Col<T1>& get_trial_inertial_force(const shared_ptr<Factory<T1>>&);
492  template<sp_d T1> friend Col<T1>& get_trial_displacement(const shared_ptr<Factory<T1>>&);
493  template<sp_d T1> friend Col<T1>& get_trial_velocity(const shared_ptr<Factory<T1>>&);
494  template<sp_d T1> friend Col<T1>& get_trial_acceleration(const shared_ptr<Factory<T1>>&);
495  template<sp_d T1> friend Col<T1>& get_trial_temperature(const shared_ptr<Factory<T1>>&);
496 
497  template<sp_d T1> friend T& get_incre_time(const shared_ptr<Factory<T1>>&);
498  template<sp_d T1> friend Col<T1>& get_incre_load_factor(const shared_ptr<Factory<T1>>&);
499  template<sp_d T1> friend Col<T1>& get_incre_load(const shared_ptr<Factory<T1>>&);
500  template<sp_d T1> friend Col<T1>& get_incre_settlement(const shared_ptr<Factory<T1>>&);
501  template<sp_d T1> friend Col<T1>& get_incre_resistance(const shared_ptr<Factory<T1>>&);
502  template<sp_d T1> friend Col<T1>& get_incre_damping_force(const shared_ptr<Factory<T1>>&);
503  template<sp_d T1> friend Col<T1>& get_incre_inertial_force(const shared_ptr<Factory<T1>>&);
504  template<sp_d T1> friend Col<T1>& get_incre_displacement(const shared_ptr<Factory<T1>>&);
505  template<sp_d T1> friend Col<T1>& get_incre_velocity(const shared_ptr<Factory<T1>>&);
506  template<sp_d T1> friend Col<T1>& get_incre_acceleration(const shared_ptr<Factory<T1>>&);
507  template<sp_d T1> friend Col<T1>& get_incre_temperature(const shared_ptr<Factory<T1>>&);
508 
509  template<sp_d T1> friend T& get_current_time(const shared_ptr<Factory<T1>>&);
510  template<sp_d T1> friend Col<T1>& get_current_load_factor(const shared_ptr<Factory<T1>>&);
511  template<sp_d T1> friend Col<T1>& get_current_load(const shared_ptr<Factory<T1>>&);
512  template<sp_d T1> friend Col<T1>& get_current_settlement(const shared_ptr<Factory<T1>>&);
513  template<sp_d T1> friend Col<T1>& get_current_resistance(const shared_ptr<Factory<T1>>&);
514  template<sp_d T1> friend Col<T1>& get_current_damping_force(const shared_ptr<Factory<T1>>&);
515  template<sp_d T1> friend Col<T1>& get_current_inertial_force(const shared_ptr<Factory<T1>>&);
516  template<sp_d T1> friend Col<T1>& get_current_displacement(const shared_ptr<Factory<T1>>&);
517  template<sp_d T1> friend Col<T1>& get_current_velocity(const shared_ptr<Factory<T1>>&);
518  template<sp_d T1> friend Col<T1>& get_current_acceleration(const shared_ptr<Factory<T1>>&);
519  template<sp_d T1> friend Col<T1>& get_current_temperature(const shared_ptr<Factory<T1>>&);
520 
521  template<sp_d T1> friend T& get_pre_time(const shared_ptr<Factory<T1>>&);
522  template<sp_d T1> friend Col<T1>& get_pre_load_factor(const shared_ptr<Factory<T1>>&);
523  template<sp_d T1> friend Col<T1>& get_pre_load(const shared_ptr<Factory<T1>>&);
524  template<sp_d T1> friend Col<T1>& get_pre_settlement(const shared_ptr<Factory<T1>>&);
525  template<sp_d T1> friend Col<T1>& get_pre_resistance(const shared_ptr<Factory<T1>>&);
526  template<sp_d T1> friend Col<T1>& get_pre_damping_force(const shared_ptr<Factory<T1>>&);
527  template<sp_d T1> friend Col<T1>& get_pre_inertial_force(const shared_ptr<Factory<T1>>&);
528  template<sp_d T1> friend Col<T1>& get_pre_displacement(const shared_ptr<Factory<T1>>&);
529  template<sp_d T1> friend Col<T1>& get_pre_velocity(const shared_ptr<Factory<T1>>&);
530  template<sp_d T1> friend Col<T1>& get_pre_acceleration(const shared_ptr<Factory<T1>>&);
531  template<sp_d T1> friend Col<T1>& get_pre_temperature(const shared_ptr<Factory<T1>>&);
532 
533  template<sp_d T1> friend shared_ptr<MetaMat<T1>>& get_mass(const shared_ptr<Factory<T1>>&);
534  template<sp_d T1> friend shared_ptr<MetaMat<T1>>& get_damping(const shared_ptr<Factory<T1>>&);
535  template<sp_d T1> friend shared_ptr<MetaMat<T1>>& get_stiffness(const shared_ptr<Factory<T1>>&);
536  template<sp_d T1> friend shared_ptr<MetaMat<T1>>& get_geometry(const shared_ptr<Factory<T1>>&);
537 
538  template<sp_d T1> friend Col<T1>& get_eigenvalue(const shared_ptr<Factory<T1>>&);
539  template<sp_d T1> friend Mat<T1>& get_eigenvector(const shared_ptr<Factory<T1>>&);
540 
541  /*************************STATUS*************************/
542 
543  void commit_energy();
544  void clear_energy();
545 
546  void commit_status();
547  void commit_time();
548  void commit_load_factor();
549  void commit_load();
550  void commit_settlement();
551  void commit_resistance();
552  void commit_damping_force();
553  void commit_inertial_force();
554  void commit_displacement();
555  void commit_velocity();
556  void commit_acceleration();
557  void commit_temperature();
559 
560  void commit_pre_status();
561  void commit_pre_time();
562  void commit_pre_load_factor();
563  void commit_pre_load();
564  void commit_pre_settlement();
565  void commit_pre_resistance();
569  void commit_pre_velocity();
571  void commit_pre_temperature();
572 
573  void clear_status();
574  void clear_time();
575  void clear_load_factor();
576  void clear_load();
577  void clear_settlement();
578  void clear_resistance();
579  void clear_damping_force();
580  void clear_inertial_force();
581  void clear_displacement();
582  void clear_velocity();
583  void clear_acceleration();
584  void clear_temperature();
586 
587  void reset_status();
588  void reset_time();
589  void reset_load_factor();
590  void reset_load();
591  void reset_settlement();
592  void reset_resistance();
593  void reset_damping_force();
594  void reset_inertial_force();
595  void reset_displacement();
596  void reset_velocity();
597  void reset_acceleration();
598  void reset_temperature();
600 
601  void clear_eigen();
602  void clear_mass();
603  void clear_damping();
604  void clear_stiffness();
605  void clear_geometry();
606  void clear_auxiliary();
607 
608  void reset();
609 
610  /*************************ASSEMBLER*************************/
611 
612  void assemble_resistance(const Mat<T>&, const uvec&);
613  void assemble_damping_force(const Mat<T>&, const uvec&);
614  void assemble_inertial_force(const Mat<T>&, const uvec&);
615  void assemble_mass(const Mat<T>&, const uvec&);
616  void assemble_damping(const Mat<T>&, const uvec&);
617  void assemble_stiffness(const Mat<T>&, const uvec&);
618  void assemble_geometry(const Mat<T>&, const uvec&);
619  void assemble_stiffness(const SpMat<T>&, const uvec&);
620 
621  /*************************UTILITY*************************/
622 
623  void print() const;
624 };
625 
626 template<sp_d T> Factory<T>::Factory(const unsigned D, const AnalysisType AT, const StorageScheme SS)
627  : n_size(D)
628  , analysis_type(AT)
629  , storage_type(SS) {}
630 
631 template<sp_d T> void Factory<T>::set_size(const unsigned D) {
632  if(D == n_size) return;
633  n_size = D;
634  access::rw(initialized) = false;
635 }
636 
637 template<sp_d T> unsigned Factory<T>::get_size() const { return n_size; }
638 
639 template<sp_d T> void Factory<T>::set_entry(const uword N) { n_elem = N; }
640 
641 template<sp_d T> uword Factory<T>::get_entry() const { return n_elem; }
642 
643 template<sp_d T> void Factory<T>::set_nlgeom(const bool B) {
644  if(B == nlgeom) return;
645  nlgeom = B;
646  access::rw(initialized) = false;
647 }
648 
649 template<sp_d T> bool Factory<T>::is_nlgeom() const { return nlgeom; }
650 
651 template<sp_d T> void Factory<T>::set_solver_type(const SolverType E) { solver = E; }
652 
653 template<sp_d T> SolverType Factory<T>::get_solver_type() const { return solver; }
654 
655 template<sp_d T> void Factory<T>::set_solver_setting(const SolverSetting<double>& SS) { setting = SS; }
656 
657 template<sp_d T> const SolverSetting<double>& Factory<T>::get_solver_setting() const { return setting; }
658 
659 template<sp_d T> void Factory<T>::set_analysis_type(const AnalysisType AT) {
660  if(AT == analysis_type) return;
661  analysis_type = AT;
662  access::rw(initialized) = false;
663 }
664 
665 template<sp_d T> AnalysisType Factory<T>::get_analysis_type() const { return analysis_type; }
666 
667 template<sp_d T> void Factory<T>::set_storage_scheme(const StorageScheme SS) {
668  if(SS == storage_type) return;
669  storage_type = SS;
670  access::rw(initialized) = false;
671 }
672 
673 template<sp_d T> StorageScheme Factory<T>::get_storage_scheme() const { return storage_type; }
674 
675 template<sp_d T> bool Factory<T>::is_sparse() const { return StorageScheme::SPARSE == storage_type || StorageScheme::SPARSESYMM == storage_type; }
676 
677 template<sp_d T> void Factory<T>::set_bandwidth(const unsigned L, const unsigned U) {
678  if(L == n_lobw && U == n_upbw) return;
679  n_lobw = L;
680  n_upbw = U;
681  n_sfbw = L + U;
682  access::rw(initialized) = false;
683 }
684 
685 template<sp_d T> void Factory<T>::get_bandwidth(unsigned& L, unsigned& U) const {
686  L = n_lobw;
687  U = n_upbw;
688 }
689 
690 template<sp_d T> void Factory<T>::update_reference_size() { n_rfld = static_cast<unsigned>(reference_dof.size()); }
691 
692 template<sp_d T> void Factory<T>::set_reference_size(const unsigned S) {
693  if(S == n_rfld) return;
694  n_rfld = S;
695 }
696 
697 template<sp_d T> unsigned Factory<T>::get_reference_size() const { return n_rfld; }
698 
699 template<sp_d T> void Factory<T>::update_reference_dof(const uvec& S) { reference_dof.insert(S.cbegin(), S.cend()); }
700 
701 template<sp_d T> void Factory<T>::set_reference_dof(const suanpan::set<uword>& D) { reference_dof = D; }
702 
703 template<sp_d T> const suanpan::set<uword>& Factory<T>::get_reference_dof() const { return reference_dof; }
704 
705 template<sp_d T> void Factory<T>::set_error(const T& E) { error = E; }
706 
707 template<sp_d T> const T& Factory<T>::get_error() const { return error; }
708 
709 template<sp_d T> int Factory<T>::initialize() {
710  reference_dof.clear(); // clear reference dof vector in every step
711 
712  if(initialized || n_size == 0) return 0;
713 
714  ninja.zeros(n_size);
715  sushi.zeros(n_size);
716 
717  reset();
718 
719  switch(analysis_type) {
720  case AnalysisType::DISP:
721  initialize_displacement();
722  break;
723  case AnalysisType::EIGEN:
724  initialize_mass();
725  initialize_stiffness();
726  initialize_eigen();
727  break;
730  initialize_load();
731  initialize_resistance();
732  initialize_displacement();
733  initialize_stiffness();
734  initialize_geometry();
735  break;
737  initialize_load();
738  initialize_resistance();
739  initialize_damping_force();
740  initialize_inertial_force();
741  initialize_displacement();
742  initialize_velocity();
743  initialize_acceleration();
744  initialize_mass();
745  initialize_damping();
746  initialize_stiffness();
747  initialize_geometry();
748  break;
749  case AnalysisType::NONE:
750  break;
751  }
752 
753  initialize_auxiliary_resistance();
754 
755  access::rw(initialized) = true;
756 
757  return 0;
758 }
759 
760 template<sp_d T> void Factory<T>::initialize_load_factor() {
761  if(n_rfld == 0) return;
762 
763  trial_load_factor.zeros(n_rfld);
764  incre_load_factor.zeros(n_rfld);
765  current_load_factor.zeros(n_rfld);
766 
767  reference_load.zeros(n_size, n_rfld);
768 }
769 
770 template<sp_d T> void Factory<T>::initialize_load() {
771  trial_load.zeros(n_size);
772  incre_load.zeros(n_size);
773  current_load.zeros(n_size);
774 }
775 
776 template<sp_d T> void Factory<T>::initialize_settlement() {
777  trial_settlement.zeros(n_size);
778  incre_settlement.zeros(n_size);
779  current_settlement.zeros(n_size);
780 }
781 
782 template<sp_d T> void Factory<T>::initialize_resistance() {
783  trial_resistance.zeros(n_size);
784  incre_resistance.zeros(n_size);
785  current_resistance.zeros(n_size);
786 }
787 
788 template<sp_d T> void Factory<T>::initialize_damping_force() {
789  trial_damping_force.zeros(n_size);
790  incre_damping_force.zeros(n_size);
791  current_damping_force.zeros(n_size);
792 }
793 
794 template<sp_d T> void Factory<T>::initialize_inertial_force() {
795  trial_inertial_force.zeros(n_size);
796  incre_inertial_force.zeros(n_size);
797  current_inertial_force.zeros(n_size);
798 }
799 
800 template<sp_d T> void Factory<T>::initialize_displacement() {
801  trial_displacement.zeros(n_size);
802  incre_displacement.zeros(n_size);
803  current_displacement.zeros(n_size);
804 }
805 
806 template<sp_d T> void Factory<T>::initialize_velocity() {
807  trial_velocity.zeros(n_size);
808  incre_velocity.zeros(n_size);
809  current_velocity.zeros(n_size);
810 }
811 
812 template<sp_d T> void Factory<T>::initialize_acceleration() {
813  trial_acceleration.zeros(n_size);
814  incre_acceleration.zeros(n_size);
815  current_acceleration.zeros(n_size);
816 }
817 
818 template<sp_d T> void Factory<T>::initialize_temperature() {
819  trial_temperature.zeros(n_size);
820  incre_temperature.zeros(n_size);
821  current_temperature.zeros(n_size);
822 }
823 
825  trial_constraint_resistance.zeros(n_size);
826  current_constraint_resistance.zeros(n_size);
827 }
828 
829 template<sp_d T> void Factory<T>::initialize_mass() { global_mass = get_matrix_container(this); }
830 
831 template<sp_d T> void Factory<T>::initialize_damping() { global_damping = get_matrix_container(this); }
832 
833 template<sp_d T> void Factory<T>::initialize_stiffness() { global_stiffness = get_matrix_container(this); }
834 
835 template<sp_d T> void Factory<T>::initialize_geometry() {
836  if(!nlgeom) return;
837 
838  global_geometry = get_matrix_container(this);
839 }
840 
841 template<sp_d T> void Factory<T>::initialize_eigen() {
842  eigenvalue.zeros(n_size);
843  eigenvector.zeros(n_size, n_size);
844 }
845 
846 template<sp_d T> void Factory<T>::set_ninja(const Col<T>& N) { ninja = N; }
847 
848 template<sp_d T> void Factory<T>::set_sushi(const Col<T>& S) { sushi = S; }
849 
850 template<sp_d T> void Factory<T>::set_mpc(const unsigned S) {
851  n_mpc = S;
852  auxiliary_encoding.zeros(n_mpc);
853  auxiliary_resistance.zeros(n_mpc);
854  auxiliary_load.zeros(n_mpc);
855  auxiliary_stiffness.zeros(n_size, n_mpc);
856 }
857 
858 template<sp_d T> void Factory<T>::set_reference_load(const SpMat<T>& L) { reference_load = L; }
859 
860 template<sp_d T> void Factory<T>::set_trial_time(const T& M) { trial_time = M; }
861 
862 template<sp_d T> void Factory<T>::set_trial_load_factor(const Col<T>& L) { trial_load_factor = L; }
863 
864 template<sp_d T> void Factory<T>::set_trial_load(const Col<T>& L) { trial_load = L; }
865 
866 template<sp_d T> void Factory<T>::set_trial_settlement(const Col<T>& S) { trial_settlement = S; }
867 
868 template<sp_d T> void Factory<T>::set_trial_resistance(const Col<T>& R) { trial_resistance = R; }
869 
870 template<sp_d T> void Factory<T>::set_trial_damping_force(const Col<T>& R) { trial_damping_force = R; }
871 
872 template<sp_d T> void Factory<T>::set_trial_inertial_force(const Col<T>& R) { trial_inertial_force = R; }
873 
874 template<sp_d T> void Factory<T>::set_trial_displacement(const Col<T>& D) { trial_displacement = D; }
875 
876 template<sp_d T> void Factory<T>::set_trial_velocity(const Col<T>& V) { trial_velocity = V; }
877 
878 template<sp_d T> void Factory<T>::set_trial_acceleration(const Col<T>& A) { trial_acceleration = A; }
879 
880 template<sp_d T> void Factory<T>::set_trial_temperature(const Col<T>& M) { trial_temperature = M; }
881 
882 template<sp_d T> void Factory<T>::set_incre_time(const T& M) { incre_time = M; }
883 
884 template<sp_d T> void Factory<T>::set_incre_load_factor(const Col<T>& L) { incre_load_factor = L; }
885 
886 template<sp_d T> void Factory<T>::set_incre_load(const Col<T>& L) { incre_load = L; }
887 
888 template<sp_d T> void Factory<T>::set_incre_settlement(const Col<T>& S) { incre_settlement = S; }
889 
890 template<sp_d T> void Factory<T>::set_incre_resistance(const Col<T>& R) { incre_resistance = R; }
891 
892 template<sp_d T> void Factory<T>::set_incre_damping_force(const Col<T>& R) { incre_damping_force = R; }
893 
894 template<sp_d T> void Factory<T>::set_incre_inertial_force(const Col<T>& R) { incre_inertial_force = R; }
895 
896 template<sp_d T> void Factory<T>::set_incre_displacement(const Col<T>& D) { incre_displacement = D; }
897 
898 template<sp_d T> void Factory<T>::set_incre_velocity(const Col<T>& V) { incre_velocity = V; }
899 
900 template<sp_d T> void Factory<T>::set_incre_acceleration(const Col<T>& A) { incre_acceleration = A; }
901 
902 template<sp_d T> void Factory<T>::set_incre_temperature(const Col<T>& M) { incre_temperature = M; }
903 
904 template<sp_d T> void Factory<T>::set_current_time(const T& M) { current_time = M; }
905 
906 template<sp_d T> void Factory<T>::set_current_load_factor(const Col<T>& L) { current_load_factor = L; }
907 
908 template<sp_d T> void Factory<T>::set_current_load(const Col<T>& L) { current_load = L; }
909 
910 template<sp_d T> void Factory<T>::set_current_settlement(const Col<T>& S) { current_settlement = S; }
911 
912 template<sp_d T> void Factory<T>::set_current_resistance(const Col<T>& R) { current_resistance = R; }
913 
914 template<sp_d T> void Factory<T>::set_current_damping_force(const Col<T>& R) { current_damping_force = R; }
915 
916 template<sp_d T> void Factory<T>::set_current_inertial_force(const Col<T>& R) { current_inertial_force = R; }
917 
918 template<sp_d T> void Factory<T>::set_current_displacement(const Col<T>& D) { current_displacement = D; }
919 
920 template<sp_d T> void Factory<T>::set_current_velocity(const Col<T>& V) { current_velocity = V; }
921 
922 template<sp_d T> void Factory<T>::set_current_acceleration(const Col<T>& A) { current_acceleration = A; }
923 
924 template<sp_d T> void Factory<T>::set_current_temperature(const Col<T>& M) { current_temperature = M; }
925 
926 template<sp_d T> void Factory<T>::set_pre_time(const T& M) { pre_time = M; }
927 
928 template<sp_d T> void Factory<T>::set_pre_load_factor(const Col<T>& L) { pre_load_factor = L; }
929 
930 template<sp_d T> void Factory<T>::set_pre_load(const Col<T>& L) { pre_load = L; }
931 
932 template<sp_d T> void Factory<T>::set_pre_settlement(const Col<T>& S) { pre_settlement = S; }
933 
934 template<sp_d T> void Factory<T>::set_pre_resistance(const Col<T>& R) { pre_resistance = R; }
935 
936 template<sp_d T> void Factory<T>::set_pre_damping_force(const Col<T>& R) { pre_damping_force = R; }
937 
938 template<sp_d T> void Factory<T>::set_pre_inertial_force(const Col<T>& R) { pre_inertial_force = R; }
939 
940 template<sp_d T> void Factory<T>::set_pre_displacement(const Col<T>& D) { pre_displacement = D; }
941 
942 template<sp_d T> void Factory<T>::set_pre_velocity(const Col<T>& V) { pre_velocity = V; }
943 
944 template<sp_d T> void Factory<T>::set_pre_acceleration(const Col<T>& A) { pre_acceleration = A; }
945 
946 template<sp_d T> void Factory<T>::set_pre_temperature(const Col<T>& M) { pre_temperature = M; }
947 
948 template<sp_d T> void Factory<T>::set_mass(const shared_ptr<MetaMat<T>>& M) { global_mass = M; }
949 
950 template<sp_d T> void Factory<T>::set_damping(const shared_ptr<MetaMat<T>>& C) { global_damping = C; }
951 
952 template<sp_d T> void Factory<T>::set_stiffness(const shared_ptr<MetaMat<T>>& K) { global_stiffness = K; }
953 
954 template<sp_d T> void Factory<T>::set_geometry(const shared_ptr<MetaMat<T>>& G) { global_geometry = G; }
955 
956 template<sp_d T> void Factory<T>::set_eigenvalue(const Col<T>& L) { eigenvalue = L; }
957 
958 template<sp_d T> void Factory<T>::set_eigenvector(const Mat<T>& V) { eigenvector = V; }
959 
960 template<sp_d T> const Col<T>& Factory<T>::get_ninja() const { return ninja; }
961 
962 template<sp_d T> const Col<T>& Factory<T>::get_sushi() const { return sushi; }
963 
964 template<sp_d T> unsigned Factory<T>::get_mpc() const { return n_mpc; }
965 
966 template<sp_d T> const SpMat<T>& Factory<T>::get_reference_load() const { return reference_load; }
967 
968 template<sp_d T> const uvec& Factory<T>::get_auxiliary_encoding() const { return auxiliary_encoding; }
969 
970 template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_lambda() const { return auxiliary_lambda; }
971 
972 template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_resistance() const { return auxiliary_resistance; }
973 
974 template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_load() const { return auxiliary_load; }
975 
976 template<sp_d T> const SpMat<T>& Factory<T>::get_auxiliary_stiffness() const { return auxiliary_stiffness; }
977 
978 template<sp_d T> const SpCol<T>& Factory<T>::get_trial_constraint_resistance() const { return trial_constraint_resistance; }
979 
980 template<sp_d T> const SpCol<T>& Factory<T>::get_current_constraint_resistance() const { return current_constraint_resistance; }
981 
982 template<sp_d T> T Factory<T>::get_strain_energy() { return strain_energy; }
983 
984 template<sp_d T> T Factory<T>::get_kinetic_energy() { return kinetic_energy; }
985 
986 template<sp_d T> T Factory<T>::get_viscous_energy() { return viscous_energy; }
987 
988 template<sp_d T> T Factory<T>::get_complementary_energy() { return complementary_energy; }
989 
990 template<sp_d T> const Col<T>& Factory<T>::get_momentum() { return momentum; }
991 
992 template<sp_d T> const T& Factory<T>::get_trial_time() const { return trial_time; }
993 
994 template<sp_d T> const Col<T>& Factory<T>::get_trial_load_factor() const { return trial_load_factor; }
995 
996 template<sp_d T> const Col<T>& Factory<T>::get_trial_load() const { return trial_load; }
997 
998 template<sp_d T> const Col<T>& Factory<T>::get_trial_settlement() const { return trial_settlement; }
999 
1000 template<sp_d T> const Col<T>& Factory<T>::get_trial_resistance() const { return trial_resistance; }
1001 
1002 template<sp_d T> const Col<T>& Factory<T>::get_trial_damping_force() const { return trial_damping_force; }
1003 
1004 template<sp_d T> const Col<T>& Factory<T>::get_trial_inertial_force() const { return trial_inertial_force; }
1005 
1006 template<sp_d T> const Col<T>& Factory<T>::get_trial_displacement() const { return trial_displacement; }
1007 
1008 template<sp_d T> const Col<T>& Factory<T>::get_trial_velocity() const { return trial_velocity; }
1009 
1010 template<sp_d T> const Col<T>& Factory<T>::get_trial_acceleration() const { return trial_acceleration; }
1011 
1012 template<sp_d T> const Col<T>& Factory<T>::get_trial_temperature() const { return trial_temperature; }
1013 
1014 template<sp_d T> const T& Factory<T>::get_incre_time() const { return incre_time; }
1015 
1016 template<sp_d T> const Col<T>& Factory<T>::get_incre_load_factor() const { return incre_load_factor; }
1017 
1018 template<sp_d T> const Col<T>& Factory<T>::get_incre_load() const { return incre_load; }
1019 
1020 template<sp_d T> const Col<T>& Factory<T>::get_incre_settlement() const { return incre_settlement; }
1021 
1022 template<sp_d T> const Col<T>& Factory<T>::get_incre_resistance() const { return incre_resistance; }
1023 
1024 template<sp_d T> const Col<T>& Factory<T>::get_incre_damping_force() const { return incre_damping_force; }
1025 
1026 template<sp_d T> const Col<T>& Factory<T>::get_incre_inertial_force() const { return incre_inertial_force; }
1027 
1028 template<sp_d T> const Col<T>& Factory<T>::get_incre_displacement() const { return incre_displacement; }
1029 
1030 template<sp_d T> const Col<T>& Factory<T>::get_incre_velocity() const { return incre_velocity; }
1031 
1032 template<sp_d T> const Col<T>& Factory<T>::get_incre_acceleration() const { return incre_acceleration; }
1033 
1034 template<sp_d T> const Col<T>& Factory<T>::get_incre_temperature() const { return incre_temperature; }
1035 
1036 template<sp_d T> const T& Factory<T>::get_current_time() const { return current_time; }
1037 
1038 template<sp_d T> const Col<T>& Factory<T>::get_current_load_factor() const { return current_load_factor; }
1039 
1040 template<sp_d T> const Col<T>& Factory<T>::get_current_load() const { return current_load; }
1041 
1042 template<sp_d T> const Col<T>& Factory<T>::get_current_settlement() const { return current_settlement; }
1043 
1044 template<sp_d T> const Col<T>& Factory<T>::get_current_resistance() const { return current_resistance; }
1045 
1046 template<sp_d T> const Col<T>& Factory<T>::get_current_damping_force() const { return current_damping_force; }
1047 
1048 template<sp_d T> const Col<T>& Factory<T>::get_current_inertial_force() const { return current_inertial_force; }
1049 
1050 template<sp_d T> const Col<T>& Factory<T>::get_current_displacement() const { return current_displacement; }
1051 
1052 template<sp_d T> const Col<T>& Factory<T>::get_current_velocity() const { return current_velocity; }
1053 
1054 template<sp_d T> const Col<T>& Factory<T>::get_current_acceleration() const { return current_acceleration; }
1055 
1056 template<sp_d T> const Col<T>& Factory<T>::get_current_temperature() const { return current_temperature; }
1057 
1058 template<sp_d T> const T& Factory<T>::get_pre_time() const { return pre_time; }
1059 
1060 template<sp_d T> const Col<T>& Factory<T>::get_pre_load_factor() const { return pre_load_factor; }
1061 
1062 template<sp_d T> const Col<T>& Factory<T>::get_pre_load() const { return pre_load; }
1063 
1064 template<sp_d T> const Col<T>& Factory<T>::get_pre_settlement() const { return pre_settlement; }
1065 
1066 template<sp_d T> const Col<T>& Factory<T>::get_pre_resistance() const { return pre_resistance; }
1067 
1068 template<sp_d T> const Col<T>& Factory<T>::get_pre_damping_force() const { return pre_damping_force; }
1069 
1070 template<sp_d T> const Col<T>& Factory<T>::get_pre_inertial_force() const { return pre_inertial_force; }
1071 
1072 template<sp_d T> const Col<T>& Factory<T>::get_pre_displacement() const { return pre_displacement; }
1073 
1074 template<sp_d T> const Col<T>& Factory<T>::get_pre_velocity() const { return pre_velocity; }
1075 
1076 template<sp_d T> const Col<T>& Factory<T>::get_pre_acceleration() const { return pre_acceleration; }
1077 
1078 template<sp_d T> const Col<T>& Factory<T>::get_pre_temperature() const { return pre_temperature; }
1079 
1080 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_mass() const { return global_mass; }
1081 
1082 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_damping() const { return global_damping; }
1083 
1084 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_stiffness() const { return global_stiffness; }
1085 
1086 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_geometry() const { return global_geometry; }
1087 
1088 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_encoding_mutex() { return global_mutex.at(0); }
1089 
1090 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_resistance_mutex() { return global_mutex.at(1); }
1091 
1092 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_load_mutex() { return global_mutex.at(2); }
1093 
1094 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_stiffness_mutex() { return global_mutex.at(3); }
1095 
1096 template<sp_d T> std::mutex& Factory<T>::get_trial_constraint_resistance_mutex() { return global_mutex.at(4); }
1097 
1098 template<sp_d T> std::mutex& Factory<T>::get_trial_load_mutex() { return global_mutex.at(5); }
1099 
1100 template<sp_d T> std::mutex& Factory<T>::get_trial_settlement_mutex() { return global_mutex.at(6); }
1101 
1102 template<sp_d T> std::mutex& Factory<T>::get_mass_mutex() { return global_mutex.at(7); }
1103 
1104 template<sp_d T> std::mutex& Factory<T>::get_damping_mutex() { return global_mutex.at(8); }
1105 
1106 template<sp_d T> std::mutex& Factory<T>::get_stiffness_mutex() { return global_mutex.at(9); }
1107 
1108 template<sp_d T> std::mutex& Factory<T>::get_geometry_mutex() { return global_mutex.at(10); }
1109 
1110 template<sp_d T> const Col<T>& Factory<T>::get_eigenvalue() const { return eigenvalue; }
1111 
1112 template<sp_d T> const Mat<T>& Factory<T>::get_eigenvector() const { return eigenvector; }
1113 
1114 template<sp_d T> void Factory<T>::update_trial_time(const T& M) {
1115  trial_time = M;
1116  incre_time = trial_time - current_time;
1117 }
1118 
1119 template<sp_d T> void Factory<T>::update_trial_load_factor(const Col<T>& L) {
1120  trial_load_factor = L;
1121  incre_load_factor = trial_load_factor - current_load_factor;
1122 }
1123 
1124 template<sp_d T> void Factory<T>::update_trial_load(const Col<T>& L) {
1125  trial_load = L;
1126  incre_load = trial_load - current_load;
1127 }
1128 
1129 template<sp_d T> void Factory<T>::update_trial_settlement(const Col<T>& S) {
1130  trial_settlement = S;
1131  incre_settlement = trial_settlement - current_settlement;
1132 }
1133 
1134 template<sp_d T> void Factory<T>::update_trial_resistance(const Col<T>& R) {
1135  trial_resistance = R;
1136  incre_resistance = trial_resistance - current_resistance;
1137 }
1138 
1139 template<sp_d T> void Factory<T>::update_trial_damping_force(const Col<T>& R) {
1140  trial_damping_force = R;
1141  incre_damping_force = trial_damping_force - current_damping_force;
1142 }
1143 
1144 template<sp_d T> void Factory<T>::update_trial_inertial_force(const Col<T>& R) {
1145  trial_inertial_force = R;
1146  incre_inertial_force = trial_inertial_force - current_inertial_force;
1147 }
1148 
1149 template<sp_d T> void Factory<T>::update_trial_displacement(const Col<T>& D) {
1150  trial_displacement = D;
1151  incre_displacement = trial_displacement - current_displacement;
1152 }
1153 
1154 template<sp_d T> void Factory<T>::update_trial_velocity(const Col<T>& V) {
1155  trial_velocity = V;
1156  incre_velocity = trial_velocity - current_velocity;
1157 }
1158 
1159 template<sp_d T> void Factory<T>::update_trial_acceleration(const Col<T>& A) {
1160  trial_acceleration = A;
1161  incre_acceleration = trial_acceleration - current_acceleration;
1162 }
1163 
1164 template<sp_d T> void Factory<T>::update_trial_temperature(const Col<T>& M) {
1165  trial_temperature = M;
1166  incre_temperature = trial_temperature - current_temperature;
1167 }
1168 
1169 template<sp_d T> void Factory<T>::update_incre_time(const T& M) {
1170  incre_time = M;
1171  trial_time = current_time + incre_time;
1172 }
1173 
1174 template<sp_d T> void Factory<T>::update_incre_load_factor(const Col<T>& L) {
1175  incre_load_factor = L;
1176  trial_load_factor = current_load_factor + incre_load_factor;
1177 }
1178 
1179 template<sp_d T> void Factory<T>::update_incre_load(const Col<T>& L) {
1180  incre_load = L;
1181  trial_load = current_load + incre_load;
1182 }
1183 
1184 template<sp_d T> void Factory<T>::update_incre_settlement(const Col<T>& S) {
1185  incre_settlement = S;
1186  trial_settlement = current_settlement + incre_settlement;
1187 }
1188 
1189 template<sp_d T> void Factory<T>::update_incre_resistance(const Col<T>& R) {
1190  incre_resistance = R;
1191  trial_resistance = current_resistance + incre_resistance;
1192 }
1193 
1194 template<sp_d T> void Factory<T>::update_incre_damping_force(const Col<T>& R) {
1195  incre_damping_force = R;
1196  trial_damping_force = current_damping_force + incre_damping_force;
1197 }
1198 
1199 template<sp_d T> void Factory<T>::update_incre_inertial_force(const Col<T>& R) {
1200  incre_inertial_force = R;
1201  trial_inertial_force = current_inertial_force + incre_inertial_force;
1202 }
1203 
1204 template<sp_d T> void Factory<T>::update_incre_displacement(const Col<T>& D) {
1205  incre_displacement = D;
1206  trial_displacement = current_displacement + incre_displacement;
1207 }
1208 
1209 template<sp_d T> void Factory<T>::update_incre_velocity(const Col<T>& V) {
1210  incre_velocity = V;
1211  trial_velocity = current_velocity + incre_velocity;
1212 }
1213 
1214 template<sp_d T> void Factory<T>::update_incre_acceleration(const Col<T>& A) {
1215  incre_acceleration = A;
1216  trial_acceleration = current_acceleration + incre_acceleration;
1217 }
1218 
1219 template<sp_d T> void Factory<T>::update_incre_temperature(const Col<T>& M) {
1220  incre_temperature = M;
1221  trial_temperature = current_temperature + incre_temperature;
1222 }
1223 
1224 template<sp_d T> void Factory<T>::update_current_time(const T& M) {
1225  trial_time = current_time = M;
1226  incre_time = 0.;
1227 }
1228 
1229 template<sp_d T> void Factory<T>::update_current_load_factor(const Col<T>& L) {
1230  trial_load_factor = current_load_factor = L;
1231  incre_load_factor.zeros();
1232 }
1233 
1234 template<sp_d T> void Factory<T>::update_current_load(const Col<T>& L) {
1235  trial_load = current_load = L;
1236  incre_load.zeros();
1237 }
1238 
1239 template<sp_d T> void Factory<T>::update_current_settlement(const Col<T>& S) {
1240  trial_settlement = current_settlement = S;
1241  incre_settlement.zeros();
1242 }
1243 
1244 template<sp_d T> void Factory<T>::update_current_resistance(const Col<T>& R) {
1245  trial_resistance = current_resistance = R;
1246  incre_resistance.zeros();
1247 }
1248 
1249 template<sp_d T> void Factory<T>::update_current_damping_force(const Col<T>& R) {
1250  trial_damping_force = current_damping_force = R;
1251  incre_damping_force.zeros();
1252 }
1253 
1254 template<sp_d T> void Factory<T>::update_current_inertial_force(const Col<T>& R) {
1255  trial_inertial_force = current_inertial_force = R;
1256  incre_inertial_force.zeros();
1257 }
1258 
1259 template<sp_d T> void Factory<T>::update_current_displacement(const Col<T>& D) {
1260  trial_displacement = current_displacement = D;
1261  incre_displacement.zeros();
1262 }
1263 
1264 template<sp_d T> void Factory<T>::update_current_velocity(const Col<T>& V) {
1265  trial_velocity = current_velocity = V;
1266  incre_velocity.zeros();
1267 }
1268 
1269 template<sp_d T> void Factory<T>::update_current_acceleration(const Col<T>& A) {
1270  trial_acceleration = current_acceleration = A;
1271  incre_acceleration.zeros();
1272 }
1273 
1274 template<sp_d T> void Factory<T>::update_current_temperature(const Col<T>& M) {
1275  trial_temperature = current_temperature = M;
1276  incre_temperature.zeros();
1277 }
1278 
1279 template<sp_d T> void Factory<T>::update_trial_time_by(const T& M) {
1280  trial_time += M;
1281  incre_time = trial_time - current_time;
1282 }
1283 
1284 template<sp_d T> void Factory<T>::update_trial_load_factor_by(const Col<T>& L) {
1285  trial_load_factor += L;
1286  incre_load_factor = trial_load_factor - current_load_factor;
1287 }
1288 
1289 template<sp_d T> void Factory<T>::update_trial_load_by(const Col<T>& L) {
1290  trial_load += L;
1291  incre_load = trial_load - current_load;
1292 }
1293 
1294 template<sp_d T> void Factory<T>::update_trial_settlement_by(const Col<T>& S) {
1295  trial_settlement += S;
1296  incre_settlement = trial_settlement - current_settlement;
1297 }
1298 
1299 template<sp_d T> void Factory<T>::update_trial_resistance_by(const Col<T>& R) {
1300  trial_resistance += R;
1301  incre_resistance = trial_resistance - current_resistance;
1302 }
1303 
1304 template<sp_d T> void Factory<T>::update_trial_damping_force_by(const Col<T>& R) {
1305  trial_damping_force += R;
1306  incre_damping_force = trial_damping_force - current_damping_force;
1307 }
1308 
1309 template<sp_d T> void Factory<T>::update_trial_inertial_force_by(const Col<T>& R) {
1310  trial_inertial_force += R;
1311  incre_inertial_force = trial_inertial_force - current_inertial_force;
1312 }
1313 
1314 template<sp_d T> void Factory<T>::update_trial_displacement_by(const Col<T>& D) {
1315  trial_displacement += D;
1316  incre_displacement = trial_displacement - current_displacement;
1317 }
1318 
1319 template<sp_d T> void Factory<T>::update_trial_velocity_by(const Col<T>& V) {
1320  trial_velocity += V;
1321  incre_velocity = trial_velocity - current_velocity;
1322 }
1323 
1324 template<sp_d T> void Factory<T>::update_trial_acceleration_by(const Col<T>& A) {
1325  trial_acceleration += A;
1326  incre_acceleration = trial_acceleration - current_acceleration;
1327 }
1328 
1329 template<sp_d T> void Factory<T>::update_trial_temperature_by(const Col<T>& M) {
1330  trial_temperature += M;
1331  incre_temperature = trial_temperature - current_temperature;
1332 }
1333 
1334 template<sp_d T> void Factory<T>::update_incre_time_by(const T& M) {
1335  incre_time += M;
1336  trial_time = current_time + incre_time;
1337 }
1338 
1339 template<sp_d T> void Factory<T>::update_incre_load_factor_by(const Col<T>& L) {
1340  incre_load_factor += L;
1341  trial_load_factor = current_load_factor + incre_load_factor;
1342 }
1343 
1344 template<sp_d T> void Factory<T>::update_incre_load_by(const Col<T>& L) {
1345  incre_load += L;
1346  trial_load = current_load + incre_load;
1347 }
1348 
1349 template<sp_d T> void Factory<T>::update_incre_settlement_by(const Col<T>& S) {
1350  incre_settlement += S;
1351  trial_settlement = current_settlement + incre_settlement;
1352 }
1353 
1354 template<sp_d T> void Factory<T>::update_incre_resistance_by(const Col<T>& R) {
1355  incre_resistance += R;
1356  trial_resistance = current_resistance + incre_resistance;
1357 }
1358 
1359 template<sp_d T> void Factory<T>::update_incre_damping_force_by(const Col<T>& R) {
1360  incre_damping_force += R;
1361  trial_damping_force = current_damping_force + incre_damping_force;
1362 }
1363 
1364 template<sp_d T> void Factory<T>::update_incre_inertial_force_by(const Col<T>& R) {
1365  incre_inertial_force += R;
1366  trial_inertial_force = current_inertial_force + incre_inertial_force;
1367 }
1368 
1369 template<sp_d T> void Factory<T>::update_incre_displacement_by(const Col<T>& D) {
1370  incre_displacement += D;
1371  trial_displacement = current_displacement + incre_displacement;
1372 }
1373 
1374 template<sp_d T> void Factory<T>::update_incre_velocity_by(const Col<T>& V) {
1375  incre_velocity += V;
1376  trial_velocity = current_velocity + incre_velocity;
1377 }
1378 
1379 template<sp_d T> void Factory<T>::update_incre_acceleration_by(const Col<T>& A) {
1380  incre_acceleration += A;
1381  trial_acceleration = current_acceleration + incre_acceleration;
1382 }
1383 
1384 template<sp_d T> void Factory<T>::update_incre_temperature_by(const Col<T>& M) {
1385  incre_temperature += M;
1386  trial_temperature = current_temperature + incre_temperature;
1387 }
1388 
1389 template<sp_d T> void Factory<T>::update_current_time_by(const T& M) {
1390  trial_time = current_time += M;
1391  incre_time = 0.;
1392 }
1393 
1394 template<sp_d T> void Factory<T>::update_current_load_factor_by(const Col<T>& L) {
1395  trial_load_factor = current_load_factor += L;
1396  incre_load_factor.zeros();
1397 }
1398 
1399 template<sp_d T> void Factory<T>::update_current_load_by(const Col<T>& L) {
1400  trial_load = current_load += L;
1401  incre_load.zeros();
1402 }
1403 
1404 template<sp_d T> void Factory<T>::update_current_settlement_by(const Col<T>& S) {
1405  trial_settlement = current_settlement += S;
1406  incre_settlement.zeros();
1407 }
1408 
1409 template<sp_d T> void Factory<T>::update_current_resistance_by(const Col<T>& R) {
1410  trial_resistance = current_resistance += R;
1411  incre_resistance.zeros();
1412 }
1413 
1414 template<sp_d T> void Factory<T>::update_current_damping_force_by(const Col<T>& R) {
1415  trial_damping_force = current_damping_force += R;
1416  incre_damping_force.zeros();
1417 }
1418 
1419 template<sp_d T> void Factory<T>::update_current_inertial_force_by(const Col<T>& R) {
1420  trial_inertial_force = current_inertial_force += R;
1421  incre_inertial_force.zeros();
1422 }
1423 
1424 template<sp_d T> void Factory<T>::update_current_displacement_by(const Col<T>& D) {
1425  trial_displacement = current_displacement += D;
1426  incre_displacement.zeros();
1427 }
1428 
1429 template<sp_d T> void Factory<T>::update_current_velocity_by(const Col<T>& V) {
1430  trial_velocity = current_velocity += V;
1431  incre_velocity.zeros();
1432 }
1433 
1434 template<sp_d T> void Factory<T>::update_current_acceleration_by(const Col<T>& A) {
1435  trial_acceleration = current_acceleration += A;
1436  incre_acceleration.zeros();
1437 }
1438 
1439 template<sp_d T> void Factory<T>::update_current_temperature_by(const Col<T>& M) {
1440  trial_temperature = current_temperature += M;
1441  incre_temperature.zeros();
1442 }
1443 
1444 template<sp_d T> void Factory<T>::commit_energy() {
1445  auto se = std::async([&] { if(!trial_resistance.empty() && !incre_displacement.empty()) strain_energy += .5 * dot(trial_resistance + current_resistance, incre_displacement); });
1446  auto ke = std::async([&] { if(!trial_inertial_force.empty() && !trial_velocity.empty()) kinetic_energy = .5 * dot(global_mass * trial_velocity, trial_velocity); });
1447  auto ve = std::async([&] { if(!trial_damping_force.empty() && !incre_displacement.empty()) viscous_energy += .5 * dot(trial_damping_force + current_damping_force, incre_displacement); });
1448  auto ce = std::async([&] { if(!trial_displacement.empty() && !incre_resistance.empty()) complementary_energy += .5 * dot(trial_displacement + current_displacement, incre_resistance); });
1449  auto mm = std::async([&] { if(!trial_inertial_force.empty() && !trial_velocity.empty()) momentum = global_mass * trial_velocity; });
1450 
1451  se.get();
1452  ke.get();
1453  ve.get();
1454  ce.get();
1455  mm.get();
1456 }
1457 
1458 template<sp_d T> void Factory<T>::clear_energy() {
1459  strain_energy = 0.;
1460  kinetic_energy = 0.;
1461  viscous_energy = 0.;
1462  complementary_energy = 0.;
1463  momentum.zeros();
1464 }
1465 
1466 template<sp_d T> void Factory<T>::commit_status() {
1467  commit_energy();
1468 
1469  commit_time();
1470  commit_load_factor();
1471  commit_load();
1472  commit_settlement();
1473  commit_resistance();
1474  commit_damping_force();
1475  commit_inertial_force();
1476  commit_displacement();
1477  commit_velocity();
1478  commit_acceleration();
1479  commit_temperature();
1480  commit_auxiliary_resistance();
1481 }
1482 
1483 template<sp_d T> void Factory<T>::commit_time() {
1484  current_time = trial_time;
1485  incre_time = 0.;
1486 }
1487 
1488 template<sp_d T> void Factory<T>::commit_load_factor() {
1489  if(trial_load_factor.is_empty()) return;
1490  current_load_factor = trial_load_factor;
1491  incre_load_factor.zeros();
1492 }
1493 
1494 template<sp_d T> void Factory<T>::commit_load() {
1495  if(trial_load.is_empty()) return;
1496  current_load = trial_load;
1497  incre_load.zeros();
1498 }
1499 
1500 template<sp_d T> void Factory<T>::commit_settlement() {
1501  if(trial_settlement.is_empty()) return;
1502  current_settlement = trial_settlement;
1503  incre_settlement.zeros();
1504 }
1505 
1506 template<sp_d T> void Factory<T>::commit_resistance() {
1507  if(trial_resistance.is_empty()) return;
1508  current_resistance = trial_resistance;
1509  incre_resistance.zeros();
1510 }
1511 
1512 template<sp_d T> void Factory<T>::commit_damping_force() {
1513  if(trial_damping_force.is_empty()) return;
1514  current_damping_force = trial_damping_force;
1515  incre_damping_force.zeros();
1516 }
1517 
1518 template<sp_d T> void Factory<T>::commit_inertial_force() {
1519  if(trial_inertial_force.is_empty()) return;
1520  current_inertial_force = trial_inertial_force;
1521  incre_inertial_force.zeros();
1522 }
1523 
1524 template<sp_d T> void Factory<T>::commit_displacement() {
1525  if(trial_displacement.is_empty()) return;
1526  current_displacement = trial_displacement;
1527  incre_displacement.zeros();
1528 }
1529 
1530 template<sp_d T> void Factory<T>::commit_velocity() {
1531  if(trial_velocity.is_empty()) return;
1532  current_velocity = trial_velocity;
1533  incre_velocity.zeros();
1534 }
1535 
1536 template<sp_d T> void Factory<T>::commit_acceleration() {
1537  if(trial_acceleration.is_empty()) return;
1538  current_acceleration = trial_acceleration;
1539  incre_acceleration.zeros();
1540 }
1541 
1542 template<sp_d T> void Factory<T>::commit_temperature() {
1543  if(trial_temperature.is_empty()) return;
1544  current_temperature = trial_temperature;
1545  incre_temperature.zeros();
1546 }
1547 
1549  if(trial_constraint_resistance.is_empty()) return;
1550  current_constraint_resistance = trial_constraint_resistance;
1551 }
1552 
1553 template<sp_d T> void Factory<T>::commit_pre_status() {
1554  commit_pre_time();
1555  commit_pre_load_factor();
1556  commit_pre_load();
1557  commit_pre_settlement();
1558  commit_pre_resistance();
1559  commit_pre_damping_force();
1560  commit_pre_inertial_force();
1561  commit_pre_displacement();
1562  commit_pre_velocity();
1563  commit_pre_acceleration();
1564  commit_pre_temperature();
1565 }
1566 
1567 template<sp_d T> void Factory<T>::commit_pre_time() { pre_time = current_time; }
1568 
1569 template<sp_d T> void Factory<T>::commit_pre_load_factor() { if(!current_load_factor.is_empty()) pre_load_factor = current_load_factor; }
1570 
1571 template<sp_d T> void Factory<T>::commit_pre_load() { if(!current_load.is_empty()) pre_load = current_load; }
1572 
1573 template<sp_d T> void Factory<T>::commit_pre_settlement() { if(!current_settlement.is_empty()) pre_settlement = current_settlement; }
1574 
1575 template<sp_d T> void Factory<T>::commit_pre_resistance() { if(!current_resistance.is_empty()) pre_resistance = current_resistance; }
1576 
1577 template<sp_d T> void Factory<T>::commit_pre_damping_force() { if(!current_damping_force.is_empty()) pre_damping_force = current_damping_force; }
1578 
1579 template<sp_d T> void Factory<T>::commit_pre_inertial_force() { if(!current_inertial_force.is_empty()) pre_inertial_force = current_inertial_force; }
1580 
1581 template<sp_d T> void Factory<T>::commit_pre_displacement() { if(!current_displacement.is_empty()) pre_displacement = current_displacement; }
1582 
1583 template<sp_d T> void Factory<T>::commit_pre_velocity() { if(!current_velocity.is_empty()) pre_velocity = current_velocity; }
1584 
1585 template<sp_d T> void Factory<T>::commit_pre_acceleration() { if(!current_acceleration.is_empty()) pre_acceleration = current_acceleration; }
1586 
1587 template<sp_d T> void Factory<T>::commit_pre_temperature() { if(!current_temperature.is_empty()) pre_temperature = current_temperature; }
1588 
1589 template<sp_d T> void Factory<T>::clear_status() {
1590  access::rw(initialized) = false;
1591 
1592  ninja.zeros();
1593 
1594  clear_energy();
1595 
1596  clear_time();
1597  clear_load_factor();
1598  clear_load();
1599  clear_settlement();
1600  clear_resistance();
1601  clear_damping_force();
1602  clear_inertial_force();
1603  clear_displacement();
1604  clear_velocity();
1605  clear_acceleration();
1606  clear_temperature();
1607  clear_auxiliary_resistance();
1608 }
1609 
1610 template<sp_d T> void Factory<T>::clear_time() { trial_time = incre_time = current_time = 0.; }
1611 
1612 template<sp_d T> void Factory<T>::clear_load_factor() {
1613  if(!pre_load_factor.is_empty()) pre_load_factor.zeros();
1614  if(!trial_load_factor.is_empty()) trial_load_factor.zeros();
1615  if(!incre_load_factor.is_empty()) incre_load_factor.zeros();
1616  if(!current_load_factor.is_empty()) current_load_factor.zeros();
1617 }
1618 
1619 template<sp_d T> void Factory<T>::clear_load() {
1620  if(!pre_load.is_empty()) pre_load.zeros();
1621  if(!trial_load.is_empty()) trial_load.zeros();
1622  if(!incre_load.is_empty()) incre_load.zeros();
1623  if(!current_load.is_empty()) current_load.zeros();
1624 }
1625 
1626 template<sp_d T> void Factory<T>::clear_settlement() {
1627  if(!pre_settlement.is_empty()) pre_settlement.zeros();
1628  if(!trial_settlement.is_empty()) trial_settlement.zeros();
1629  if(!incre_settlement.is_empty()) incre_settlement.zeros();
1630  if(!current_settlement.is_empty()) current_settlement.zeros();
1631 }
1632 
1633 template<sp_d T> void Factory<T>::clear_resistance() {
1634  if(!pre_resistance.is_empty()) pre_resistance.zeros();
1635  if(!trial_resistance.is_empty()) trial_resistance.zeros();
1636  if(!incre_resistance.is_empty()) incre_resistance.zeros();
1637  if(!current_resistance.is_empty()) current_resistance.zeros();
1638 }
1639 
1640 template<sp_d T> void Factory<T>::clear_damping_force() {
1641  if(!pre_damping_force.is_empty()) pre_damping_force.zeros();
1642  if(!trial_damping_force.is_empty()) trial_damping_force.zeros();
1643  if(!incre_damping_force.is_empty()) incre_damping_force.zeros();
1644  if(!current_damping_force.is_empty()) current_damping_force.zeros();
1645 }
1646 
1647 template<sp_d T> void Factory<T>::clear_inertial_force() {
1648  if(!pre_inertial_force.is_empty()) pre_inertial_force.zeros();
1649  if(!trial_inertial_force.is_empty()) trial_inertial_force.zeros();
1650  if(!incre_inertial_force.is_empty()) incre_inertial_force.zeros();
1651  if(!current_inertial_force.is_empty()) current_inertial_force.zeros();
1652 }
1653 
1654 template<sp_d T> void Factory<T>::clear_displacement() {
1655  if(!pre_displacement.is_empty()) pre_displacement.zeros();
1656  if(!trial_displacement.is_empty()) trial_displacement.zeros();
1657  if(!incre_displacement.is_empty()) incre_displacement.zeros();
1658  if(!current_displacement.is_empty()) current_displacement.zeros();
1659 }
1660 
1661 template<sp_d T> void Factory<T>::clear_velocity() {
1662  if(!pre_velocity.is_empty()) pre_velocity.zeros();
1663  if(!trial_velocity.is_empty()) trial_velocity.zeros();
1664  if(!incre_velocity.is_empty()) incre_velocity.zeros();
1665  if(!current_velocity.is_empty()) current_velocity.zeros();
1666 }
1667 
1668 template<sp_d T> void Factory<T>::clear_acceleration() {
1669  if(!pre_acceleration.is_empty()) pre_acceleration.zeros();
1670  if(!trial_acceleration.is_empty()) trial_acceleration.zeros();
1671  if(!incre_acceleration.is_empty()) incre_acceleration.zeros();
1672  if(!current_acceleration.is_empty()) current_acceleration.zeros();
1673 }
1674 
1675 template<sp_d T> void Factory<T>::clear_temperature() {
1676  if(!pre_temperature.is_empty()) pre_temperature.zeros();
1677  if(!trial_temperature.is_empty()) trial_temperature.zeros();
1678  if(!incre_temperature.is_empty()) incre_temperature.zeros();
1679  if(!current_temperature.is_empty()) current_temperature.zeros();
1680 }
1681 
1683  if(!trial_constraint_resistance.is_empty()) trial_constraint_resistance.zeros();
1684  if(!current_constraint_resistance.is_empty()) current_constraint_resistance.zeros();
1685 }
1686 
1687 template<sp_d T> void Factory<T>::reset_status() {
1688  ninja.zeros();
1689 
1690  reset_time();
1691  reset_load_factor();
1692  reset_load();
1693  reset_settlement();
1694  reset_resistance();
1695  reset_damping_force();
1696  reset_inertial_force();
1697  reset_displacement();
1698  reset_velocity();
1699  reset_acceleration();
1700  reset_temperature();
1701  reset_auxiliary_resistance();
1702 }
1703 
1704 template<sp_d T> void Factory<T>::reset_time() {
1705  trial_time = current_time;
1706  incre_time = 0.;
1707 }
1708 
1709 template<sp_d T> void Factory<T>::reset_load_factor() {
1710  if(trial_load_factor.is_empty()) return;
1711  trial_load_factor = current_load_factor;
1712  incre_load_factor.zeros();
1713 }
1714 
1715 template<sp_d T> void Factory<T>::reset_load() {
1716  if(trial_load.is_empty()) return;
1717  trial_load = current_load;
1718  incre_load.zeros();
1719 }
1720 
1721 template<sp_d T> void Factory<T>::reset_settlement() {
1722  if(trial_settlement.is_empty()) return;
1723  trial_settlement = current_settlement;
1724  incre_settlement.zeros();
1725 }
1726 
1727 template<sp_d T> void Factory<T>::reset_resistance() {
1728  if(trial_resistance.is_empty()) return;
1729  trial_resistance = current_resistance;
1730  incre_resistance.zeros();
1731 }
1732 
1733 template<sp_d T> void Factory<T>::reset_damping_force() {
1734  if(trial_damping_force.is_empty()) return;
1735  trial_damping_force = current_damping_force;
1736  incre_damping_force.zeros();
1737 }
1738 
1739 template<sp_d T> void Factory<T>::reset_inertial_force() {
1740  if(trial_inertial_force.is_empty()) return;
1741  trial_inertial_force = current_inertial_force;
1742  incre_inertial_force.zeros();
1743 }
1744 
1745 template<sp_d T> void Factory<T>::reset_displacement() {
1746  if(trial_displacement.is_empty()) return;
1747  trial_displacement = current_displacement;
1748  incre_displacement.zeros();
1749 }
1750 
1751 template<sp_d T> void Factory<T>::reset_velocity() {
1752  if(trial_velocity.is_empty()) return;
1753  trial_velocity = current_velocity;
1754  incre_velocity.zeros();
1755 }
1756 
1757 template<sp_d T> void Factory<T>::reset_acceleration() {
1758  if(trial_acceleration.is_empty()) return;
1759  trial_acceleration = current_acceleration;
1760  incre_acceleration.zeros();
1761 }
1762 
1763 template<sp_d T> void Factory<T>::reset_temperature() {
1764  if(trial_temperature.is_empty()) return;
1765  trial_temperature = current_temperature;
1766  incre_temperature.zeros();
1767 }
1768 
1770  if(trial_constraint_resistance.is_empty()) return;
1771  trial_constraint_resistance = current_constraint_resistance;
1772 }
1773 
1774 template<sp_d T> void Factory<T>::clear_eigen() {
1775  if(!eigenvalue.is_empty()) eigenvalue.zeros();
1776  if(!eigenvector.is_empty()) eigenvector.zeros();
1777 }
1778 
1779 template<sp_d T> void Factory<T>::clear_mass() { if(global_mass != nullptr) global_mass->zeros(); }
1780 
1781 template<sp_d T> void Factory<T>::clear_damping() { if(global_damping != nullptr) global_damping->zeros(); }
1782 
1783 template<sp_d T> void Factory<T>::clear_stiffness() { if(global_stiffness != nullptr) global_stiffness->zeros(); }
1784 
1785 template<sp_d T> void Factory<T>::clear_geometry() { if(global_geometry != nullptr) global_geometry->zeros(); }
1786 
1787 template<sp_d T> void Factory<T>::clear_auxiliary() {
1788  n_mpc = 0;
1789  auxiliary_load.reset();
1790  auxiliary_stiffness.set_size(n_size, 0);
1791  auxiliary_resistance.reset();
1792  auxiliary_encoding.reset();
1793 }
1794 
1795 template<sp_d T> void Factory<T>::reset() {
1796  global_mass = nullptr;
1797  global_damping = nullptr;
1798  global_stiffness = nullptr;
1799  global_geometry = nullptr;
1800 }
1801 
1802 template<sp_d T> void Factory<T>::assemble_resistance(const Mat<T>& ER, const uvec& EI) {
1803  if(ER.is_empty()) return;
1804  for(unsigned I = 0; I < EI.n_elem; ++I) trial_resistance(EI(I)) += ER(I);
1805 }
1806 
1807 template<sp_d T> void Factory<T>::assemble_damping_force(const Mat<T>& ER, const uvec& EI) {
1808  if(ER.is_empty()) return;
1809  for(unsigned I = 0; I < EI.n_elem; ++I) trial_damping_force(EI(I)) += ER(I);
1810 }
1811 
1812 template<sp_d T> void Factory<T>::assemble_inertial_force(const Mat<T>& ER, const uvec& EI) {
1813  if(ER.is_empty()) return;
1814  for(unsigned I = 0; I < EI.n_elem; ++I) trial_inertial_force(EI(I)) += ER(I);
1815 }
1816 
1817 template<sp_d T> void Factory<T>::assemble_mass(const Mat<T>& EM, const uvec& EI) {
1818  if(EM.is_empty()) return;
1819  for(unsigned I = 0; I < EI.n_elem; ++I) for(unsigned J = 0; J < EI.n_elem; ++J) global_mass->at(EI(J), EI(I)) += EM(J, I);
1820 }
1821 
1822 template<sp_d T> void Factory<T>::assemble_damping(const Mat<T>& EC, const uvec& EI) {
1823  if(EC.is_empty()) return;
1824  for(unsigned I = 0; I < EI.n_elem; ++I) for(unsigned J = 0; J < EI.n_elem; ++J) global_damping->at(EI(J), EI(I)) += EC(J, I);
1825 }
1826 
1827 template<sp_d T> void Factory<T>::assemble_stiffness(const Mat<T>& EK, const uvec& EI) {
1828  if(EK.is_empty()) return;
1829  for(unsigned I = 0; I < EI.n_elem; ++I) for(unsigned J = 0; J < EI.n_elem; ++J) global_stiffness->at(EI(J), EI(I)) += EK(J, I);
1830 }
1831 
1832 template<sp_d T> void Factory<T>::assemble_geometry(const Mat<T>& EG, const uvec& EI) {
1833  if(EG.is_empty() || !nlgeom) return;
1834  for(unsigned I = 0; I < EI.n_elem; ++I) for(unsigned J = 0; J < EI.n_elem; ++J) global_geometry->at(EI(J), EI(I)) += EG(J, I);
1835 }
1836 
1837 template<sp_d T> void Factory<T>::assemble_stiffness(const SpMat<T>& EK, const uvec& EI) {
1838  if(EK.is_empty()) return;
1839  for(auto I = EK.begin(); I != EK.end(); ++I) global_stiffness->at(EI(I.row()), EI(I.col())) += *I;
1840 }
1841 
1842 template<sp_d T> void Factory<T>::print() const { suanpan_info("This is a Factory object with size of %u.\n", n_size); }
1843 
1844 #endif // FACTORY_HPP
1845 
void reset(ExternalMaterialData *data, int *info)
Definition: ElasticExternal.cpp:74
unique_ptr< MetaMat< T > > get_matrix_container(const Factory< T > *const W)
Definition: FactoryHelper.hpp:61
A Factory class.
Definition: Factory.hpp:67
friend Col< T1 > & get_pre_inertial_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_resistance(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_acceleration(const shared_ptr< Factory< T1 >> &)
friend unique_ptr< MetaMat< T1 > > get_matrix_container(const Factory< T1 > *)
friend Col< T1 > & get_trial_temperature(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_inertial_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_damping_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_pre_load(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_acceleration(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_temperature(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_sushi(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_acceleration(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_inertial_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_auxiliary_load(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_velocity(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_pre_velocity(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_velocity(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_inertial_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_load(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_settlement(const shared_ptr< Factory< T1 >> &)
friend SpMat< T1 > & get_reference_load(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_pre_settlement(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_ninja(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_pre_acceleration(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_resistance(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_damping_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_damping_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_settlement(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_auxiliary_resistance(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_resistance(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_current_settlement(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_displacement(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_load(const shared_ptr< Factory< T1 >> &)
friend Mat< T1 > & get_eigenvector(const shared_ptr< Factory< T1 >> &)
friend unique_ptr< MetaMat< T1 > > get_basic_container(const Factory< T1 > *)
friend Col< T1 > & get_pre_displacement(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_pre_temperature(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_pre_damping_force(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_pre_resistance(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_velocity(const shared_ptr< Factory< T1 >> &)
friend SpMat< T1 > & get_auxiliary_stiffness(const shared_ptr< Factory< T1 >> &)
friend T & get_pre_time(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_displacement(const shared_ptr< Factory< T1 >> &)
friend T & get_current_time(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_incre_temperature(const shared_ptr< Factory< T1 >> &)
friend T & get_incre_time(const shared_ptr< Factory< T1 >> &)
const bool initialized
Definition: Factory.hpp:171
friend Col< T1 > & get_current_displacement(const shared_ptr< Factory< T1 >> &)
friend suanpan::set< uword > & get_reference_dof(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_trial_load(const shared_ptr< Factory< T1 >> &)
friend Col< T1 > & get_eigenvalue(const shared_ptr< Factory< T1 >> &)
friend T & get_trial_time(const shared_ptr< Factory< T1 >> &)
A MetaMat class that holds matrices.
Definition: MetaMat.hpp:39
void set_pre_inertial_force(const Col< T > &)
Definition: Factory.hpp:938
void set_eigenvalue(const Col< T > &)
Definition: Factory.hpp:956
void clear_acceleration()
Definition: Factory.hpp:1668
void initialize_load()
Definition: Factory.hpp:770
void set_mpc(unsigned)
Definition: Factory.hpp:850
void commit_temperature()
Definition: Factory.hpp:1542
void update_trial_time_by(const T &)
Definition: Factory.hpp:1279
bool is_sparse() const
Definition: Factory.hpp:675
void set_incre_velocity(const Col< T > &)
Definition: Factory.hpp:898
uword get_entry() const
Definition: Factory.hpp:641
void set_incre_temperature(const Col< T > &)
Definition: Factory.hpp:902
void print() const
Definition: Factory.hpp:1842
void clear_resistance()
Definition: Factory.hpp:1633
void initialize_auxiliary_resistance()
Definition: Factory.hpp:824
void clear_load_factor()
Definition: Factory.hpp:1612
void set_incre_acceleration(const Col< T > &)
Definition: Factory.hpp:900
int initialize()
Definition: Factory.hpp:709
void set_incre_time(const T &)
Definition: Factory.hpp:882
const Col< T > & get_current_acceleration() const
Definition: Factory.hpp:1054
void update_incre_damping_force_by(const Col< T > &)
Definition: Factory.hpp:1359
const Col< T > & get_pre_damping_force() const
Definition: Factory.hpp:1068
void reset_displacement()
Definition: Factory.hpp:1745
void update_incre_temperature_by(const Col< T > &)
Definition: Factory.hpp:1384
void set_trial_time(const T &)
Definition: Factory.hpp:860
void set_ninja(const Col< T > &)
Definition: Factory.hpp:846
void update_trial_temperature_by(const Col< T > &)
Definition: Factory.hpp:1329
bool is_nlgeom() const
Definition: Factory.hpp:649
void update_incre_load_factor_by(const Col< T > &)
Definition: Factory.hpp:1339
const Mat< T > & get_eigenvector() const
Definition: Factory.hpp:1112
const Col< T > & get_trial_velocity() const
Definition: Factory.hpp:1008
const Col< T > & get_trial_resistance() const
Definition: Factory.hpp:1000
void update_trial_load_factor_by(const Col< T > &)
Definition: Factory.hpp:1284
const Col< T > & get_ninja() const
Definition: Factory.hpp:960
void set_pre_velocity(const Col< T > &)
Definition: Factory.hpp:942
const Col< T > & get_current_load() const
Definition: Factory.hpp:1040
void set_trial_settlement(const Col< T > &)
Definition: Factory.hpp:866
void set_incre_load_factor(const Col< T > &)
Definition: Factory.hpp:884
void update_incre_load_by(const Col< T > &)
Definition: Factory.hpp:1344
void set_solver_type(SolverType)
Definition: Factory.hpp:651
const Col< T > & get_pre_load() const
Definition: Factory.hpp:1062
void commit_acceleration()
Definition: Factory.hpp:1536
AnalysisType get_analysis_type() const
Definition: Factory.hpp:665
void set_current_resistance(const Col< T > &)
Definition: Factory.hpp:912
void update_incre_displacement(const Col< T > &)
Definition: Factory.hpp:1204
const T & get_current_time() const
Definition: Factory.hpp:1036
void set_incre_displacement(const Col< T > &)
Definition: Factory.hpp:896
void update_current_acceleration_by(const Col< T > &)
Definition: Factory.hpp:1434
const uvec & get_auxiliary_encoding() const
Definition: Factory.hpp:968
void reset_auxiliary_resistance()
Definition: Factory.hpp:1769
void commit_resistance()
Definition: Factory.hpp:1506
const Col< T > & get_incre_inertial_force() const
Definition: Factory.hpp:1026
void set_storage_scheme(StorageScheme)
Definition: Factory.hpp:667
void initialize_damping()
Definition: Factory.hpp:831
void update_incre_velocity(const Col< T > &)
Definition: Factory.hpp:1209
void initialize_geometry()
Definition: Factory.hpp:835
void update_incre_displacement_by(const Col< T > &)
Definition: Factory.hpp:1369
void update_current_damping_force_by(const Col< T > &)
Definition: Factory.hpp:1414
void set_trial_velocity(const Col< T > &)
Definition: Factory.hpp:876
void update_current_resistance(const Col< T > &)
Definition: Factory.hpp:1244
void clear_inertial_force()
Definition: Factory.hpp:1647
void update_trial_damping_force_by(const Col< T > &)
Definition: Factory.hpp:1304
unsigned get_mpc() const
Definition: Factory.hpp:964
const Col< T > & get_current_settlement() const
Definition: Factory.hpp:1042
void assemble_mass(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1817
const SolverSetting< double > & get_solver_setting() const
Definition: Factory.hpp:657
void commit_load_factor()
Definition: Factory.hpp:1488
std::mutex & get_damping_mutex()
Definition: Factory.hpp:1104
void update_trial_velocity_by(const Col< T > &)
Definition: Factory.hpp:1319
void update_current_load_by(const Col< T > &)
Definition: Factory.hpp:1399
const shared_ptr< MetaMat< T > > & get_damping() const
Definition: Factory.hpp:1082
void set_pre_settlement(const Col< T > &)
Definition: Factory.hpp:932
void set_nlgeom(bool)
Definition: Factory.hpp:643
StorageScheme get_storage_scheme() const
Definition: Factory.hpp:673
void commit_inertial_force()
Definition: Factory.hpp:1518
void initialize_resistance()
Definition: Factory.hpp:782
const Col< T > & get_current_damping_force() const
Definition: Factory.hpp:1046
void set_solver_setting(const SolverSetting< double > &)
Definition: Factory.hpp:655
void update_trial_displacement_by(const Col< T > &)
Definition: Factory.hpp:1314
void reset_velocity()
Definition: Factory.hpp:1751
const Col< T > & get_trial_load_factor() const
Definition: Factory.hpp:994
void set_current_load(const Col< T > &)
Definition: Factory.hpp:908
std::mutex & get_trial_load_mutex()
Definition: Factory.hpp:1098
void initialize_stiffness()
Definition: Factory.hpp:833
void initialize_damping_force()
Definition: Factory.hpp:788
const Col< T > & get_pre_settlement() const
Definition: Factory.hpp:1064
void update_trial_temperature(const Col< T > &)
Definition: Factory.hpp:1164
void commit_velocity()
Definition: Factory.hpp:1530
T get_complementary_energy()
Definition: Factory.hpp:988
void reset()
Definition: Factory.hpp:1795
const Col< T > & get_incre_displacement() const
Definition: Factory.hpp:1028
void set_current_damping_force(const Col< T > &)
Definition: Factory.hpp:914
void update_trial_acceleration(const Col< T > &)
Definition: Factory.hpp:1159
const Col< T > & get_auxiliary_load() const
Definition: Factory.hpp:974
void commit_pre_settlement()
Definition: Factory.hpp:1573
const shared_ptr< MetaMat< T > > & get_mass() const
Definition: Factory.hpp:1080
void commit_status()
Definition: Factory.hpp:1466
void commit_pre_load_factor()
Definition: Factory.hpp:1569
const suanpan::set< uword > & get_reference_dof() const
Definition: Factory.hpp:703
void update_trial_inertial_force(const Col< T > &)
Definition: Factory.hpp:1144
void set_current_velocity(const Col< T > &)
Definition: Factory.hpp:920
void commit_displacement()
Definition: Factory.hpp:1524
void update_current_velocity_by(const Col< T > &)
Definition: Factory.hpp:1429
const Col< T > & get_incre_damping_force() const
Definition: Factory.hpp:1024
void update_trial_load_factor(const Col< T > &)
Definition: Factory.hpp:1119
void reset_load_factor()
Definition: Factory.hpp:1709
void clear_velocity()
Definition: Factory.hpp:1661
void update_current_load(const Col< T > &)
Definition: Factory.hpp:1234
void update_incre_settlement_by(const Col< T > &)
Definition: Factory.hpp:1349
unsigned get_reference_size() const
Definition: Factory.hpp:697
void update_current_acceleration(const Col< T > &)
Definition: Factory.hpp:1269
StorageScheme
Definition: Factory.hpp:48
void commit_pre_temperature()
Definition: Factory.hpp:1587
void commit_pre_status()
Definition: Factory.hpp:1553
const T & get_error() const
Definition: Factory.hpp:707
void set_trial_acceleration(const Col< T > &)
Definition: Factory.hpp:878
void set_error(const T &)
Definition: Factory.hpp:705
void update_trial_damping_force(const Col< T > &)
Definition: Factory.hpp:1139
const Col< T > & get_incre_settlement() const
Definition: Factory.hpp:1020
void reset_inertial_force()
Definition: Factory.hpp:1739
void set_reference_dof(const suanpan::set< uword > &)
Definition: Factory.hpp:701
const Col< T > & get_momentum()
Definition: Factory.hpp:990
void update_current_displacement(const Col< T > &)
Definition: Factory.hpp:1259
void clear_settlement()
Definition: Factory.hpp:1626
void set_pre_load(const Col< T > &)
Definition: Factory.hpp:930
std::mutex & get_auxiliary_stiffness_mutex()
Definition: Factory.hpp:1094
void commit_pre_damping_force()
Definition: Factory.hpp:1577
T get_strain_energy()
Definition: Factory.hpp:982
void commit_pre_inertial_force()
Definition: Factory.hpp:1579
void update_trial_settlement(const Col< T > &)
Definition: Factory.hpp:1129
void set_trial_inertial_force(const Col< T > &)
Definition: Factory.hpp:872
const Col< T > & get_eigenvalue() const
Definition: Factory.hpp:1110
void commit_pre_displacement()
Definition: Factory.hpp:1581
std::mutex & get_geometry_mutex()
Definition: Factory.hpp:1108
void commit_pre_acceleration()
Definition: Factory.hpp:1585
const Col< T > & get_incre_velocity() const
Definition: Factory.hpp:1030
void update_current_temperature_by(const Col< T > &)
Definition: Factory.hpp:1439
void set_incre_load(const Col< T > &)
Definition: Factory.hpp:886
void update_reference_size()
Definition: Factory.hpp:690
void set_pre_acceleration(const Col< T > &)
Definition: Factory.hpp:944
void set_incre_inertial_force(const Col< T > &)
Definition: Factory.hpp:894
T get_kinetic_energy()
Definition: Factory.hpp:984
void set_current_settlement(const Col< T > &)
Definition: Factory.hpp:910
void update_trial_settlement_by(const Col< T > &)
Definition: Factory.hpp:1294
void initialize_inertial_force()
Definition: Factory.hpp:794
void set_pre_resistance(const Col< T > &)
Definition: Factory.hpp:934
void commit_load()
Definition: Factory.hpp:1494
void update_trial_acceleration_by(const Col< T > &)
Definition: Factory.hpp:1324
const Col< T > & get_pre_acceleration() const
Definition: Factory.hpp:1076
void update_incre_acceleration_by(const Col< T > &)
Definition: Factory.hpp:1379
void commit_pre_load()
Definition: Factory.hpp:1571
void commit_pre_velocity()
Definition: Factory.hpp:1583
void update_trial_resistance_by(const Col< T > &)
Definition: Factory.hpp:1299
AnalysisType
Definition: Factory.hpp:39
void commit_pre_time()
Definition: Factory.hpp:1567
void commit_auxiliary_resistance()
Definition: Factory.hpp:1548
void update_current_settlement(const Col< T > &)
Definition: Factory.hpp:1239
void set_stiffness(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:952
void update_trial_load(const Col< T > &)
Definition: Factory.hpp:1124
void update_incre_load(const Col< T > &)
Definition: Factory.hpp:1179
const shared_ptr< MetaMat< T > > & get_stiffness() const
Definition: Factory.hpp:1084
void commit_settlement()
Definition: Factory.hpp:1500
const Col< T > & get_auxiliary_lambda() const
Definition: Factory.hpp:970
void reset_load()
Definition: Factory.hpp:1715
void reset_damping_force()
Definition: Factory.hpp:1733
void update_current_inertial_force_by(const Col< T > &)
Definition: Factory.hpp:1419
void set_trial_temperature(const Col< T > &)
Definition: Factory.hpp:880
std::mutex & get_auxiliary_encoding_mutex()
Definition: Factory.hpp:1088
SolverType
Definition: Factory.hpp:57
void assemble_damping_force(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1807
const Col< T > & get_pre_displacement() const
Definition: Factory.hpp:1072
void assemble_damping(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1822
const Col< T > & get_current_inertial_force() const
Definition: Factory.hpp:1048
std::mutex & get_mass_mutex()
Definition: Factory.hpp:1102
const Col< T > & get_pre_load_factor() const
Definition: Factory.hpp:1060
void update_incre_time(const T &)
Definition: Factory.hpp:1169
const Col< T > & get_trial_displacement() const
Definition: Factory.hpp:1006
void update_trial_velocity(const Col< T > &)
Definition: Factory.hpp:1154
void update_current_temperature(const Col< T > &)
Definition: Factory.hpp:1274
void reset_temperature()
Definition: Factory.hpp:1763
void assemble_stiffness(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1827
std::mutex & get_auxiliary_resistance_mutex()
Definition: Factory.hpp:1090
void set_current_acceleration(const Col< T > &)
Definition: Factory.hpp:922
void update_current_settlement_by(const Col< T > &)
Definition: Factory.hpp:1404
void initialize_settlement()
Definition: Factory.hpp:776
void clear_damping()
Definition: Factory.hpp:1781
const Col< T > & get_incre_load_factor() const
Definition: Factory.hpp:1016
void clear_energy()
Definition: Factory.hpp:1458
const SpMat< T > & get_reference_load() const
Definition: Factory.hpp:966
void update_incre_time_by(const T &)
Definition: Factory.hpp:1334
void initialize_displacement()
Definition: Factory.hpp:800
void reset_time()
Definition: Factory.hpp:1704
const Col< T > & get_current_displacement() const
Definition: Factory.hpp:1050
unsigned get_size() const
Definition: Factory.hpp:637
void reset_settlement()
Definition: Factory.hpp:1721
void clear_damping_force()
Definition: Factory.hpp:1640
void set_analysis_type(AnalysisType)
Definition: Factory.hpp:659
const Col< T > & get_incre_load() const
Definition: Factory.hpp:1018
void set_entry(uword)
Definition: Factory.hpp:639
const Col< T > & get_incre_acceleration() const
Definition: Factory.hpp:1032
void clear_geometry()
Definition: Factory.hpp:1785
T get_viscous_energy()
Definition: Factory.hpp:986
std::mutex & get_trial_constraint_resistance_mutex()
Definition: Factory.hpp:1096
const Col< T > & get_current_temperature() const
Definition: Factory.hpp:1056
void set_reference_size(unsigned)
Definition: Factory.hpp:692
void initialize_temperature()
Definition: Factory.hpp:818
void set_current_displacement(const Col< T > &)
Definition: Factory.hpp:918
void update_reference_dof(const uvec &)
Definition: Factory.hpp:699
const Col< T > & get_trial_inertial_force() const
Definition: Factory.hpp:1004
void set_current_inertial_force(const Col< T > &)
Definition: Factory.hpp:916
void set_trial_load(const Col< T > &)
Definition: Factory.hpp:864
void set_pre_time(const T &)
Definition: Factory.hpp:926
void update_incre_temperature(const Col< T > &)
Definition: Factory.hpp:1219
void update_current_velocity(const Col< T > &)
Definition: Factory.hpp:1264
const Col< T > & get_trial_settlement() const
Definition: Factory.hpp:998
void set_eigenvector(const Mat< T > &)
Definition: Factory.hpp:958
void initialize_load_factor()
Definition: Factory.hpp:760
const Col< T > & get_trial_load() const
Definition: Factory.hpp:996
void update_current_displacement_by(const Col< T > &)
Definition: Factory.hpp:1424
void set_trial_resistance(const Col< T > &)
Definition: Factory.hpp:868
const T & get_incre_time() const
Definition: Factory.hpp:1014
void set_pre_load_factor(const Col< T > &)
Definition: Factory.hpp:928
void update_current_damping_force(const Col< T > &)
Definition: Factory.hpp:1249
void update_current_load_factor(const Col< T > &)
Definition: Factory.hpp:1229
void commit_energy()
Definition: Factory.hpp:1444
void commit_damping_force()
Definition: Factory.hpp:1512
void update_incre_acceleration(const Col< T > &)
Definition: Factory.hpp:1214
const SpMat< T > & get_auxiliary_stiffness() const
Definition: Factory.hpp:976
void update_current_time_by(const T &)
Definition: Factory.hpp:1389
void assemble_inertial_force(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1812
void clear_status()
Definition: Factory.hpp:1589
void set_bandwidth(unsigned, unsigned)
Definition: Factory.hpp:677
void update_incre_velocity_by(const Col< T > &)
Definition: Factory.hpp:1374
const Col< T > & get_current_load_factor() const
Definition: Factory.hpp:1038
void commit_time()
Definition: Factory.hpp:1483
void clear_mass()
Definition: Factory.hpp:1779
void get_bandwidth(unsigned &, unsigned &) const
Definition: Factory.hpp:685
const Col< T > & get_pre_inertial_force() const
Definition: Factory.hpp:1070
void commit_pre_resistance()
Definition: Factory.hpp:1575
void reset_resistance()
Definition: Factory.hpp:1727
SolverType get_solver_type() const
Definition: Factory.hpp:653
std::mutex & get_auxiliary_load_mutex()
Definition: Factory.hpp:1092
const Col< T > & get_trial_acceleration() const
Definition: Factory.hpp:1010
void assemble_resistance(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1802
void set_incre_damping_force(const Col< T > &)
Definition: Factory.hpp:892
const SpCol< T > & get_trial_constraint_resistance() const
Definition: Factory.hpp:978
void update_incre_damping_force(const Col< T > &)
Definition: Factory.hpp:1194
const Col< T > & get_pre_velocity() const
Definition: Factory.hpp:1074
const Col< T > & get_current_resistance() const
Definition: Factory.hpp:1044
void initialize_velocity()
Definition: Factory.hpp:806
std::mutex & get_stiffness_mutex()
Definition: Factory.hpp:1106
void update_incre_inertial_force_by(const Col< T > &)
Definition: Factory.hpp:1364
const Col< T > & get_current_velocity() const
Definition: Factory.hpp:1052
void update_incre_load_factor(const Col< T > &)
Definition: Factory.hpp:1174
void set_incre_resistance(const Col< T > &)
Definition: Factory.hpp:890
void update_current_time(const T &)
Definition: Factory.hpp:1224
const Col< T > & get_pre_temperature() const
Definition: Factory.hpp:1078
void clear_load()
Definition: Factory.hpp:1619
const Col< T > & get_incre_temperature() const
Definition: Factory.hpp:1034
Factory(unsigned=0, AnalysisType=AnalysisType::NONE, StorageScheme=StorageScheme::FULL)
Definition: Factory.hpp:626
void set_size(unsigned)
Definition: Factory.hpp:631
void update_trial_displacement(const Col< T > &)
Definition: Factory.hpp:1149
const Col< T > & get_incre_resistance() const
Definition: Factory.hpp:1022
void assemble_geometry(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1832
const Col< T > & get_sushi() const
Definition: Factory.hpp:962
const Col< T > & get_auxiliary_resistance() const
Definition: Factory.hpp:972
const Col< T > & get_trial_damping_force() const
Definition: Factory.hpp:1002
void update_current_resistance_by(const Col< T > &)
Definition: Factory.hpp:1409
void set_damping(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:950
void set_pre_temperature(const Col< T > &)
Definition: Factory.hpp:946
void update_trial_load_by(const Col< T > &)
Definition: Factory.hpp:1289
const SpCol< T > & get_current_constraint_resistance() const
Definition: Factory.hpp:980
void clear_auxiliary_resistance()
Definition: Factory.hpp:1682
void update_incre_resistance_by(const Col< T > &)
Definition: Factory.hpp:1354
void set_incre_settlement(const Col< T > &)
Definition: Factory.hpp:888
void set_sushi(const Col< T > &)
Definition: Factory.hpp:848
void update_incre_settlement(const Col< T > &)
Definition: Factory.hpp:1184
void set_pre_displacement(const Col< T > &)
Definition: Factory.hpp:940
void set_current_load_factor(const Col< T > &)
Definition: Factory.hpp:906
void set_trial_displacement(const Col< T > &)
Definition: Factory.hpp:874
void update_trial_resistance(const Col< T > &)
Definition: Factory.hpp:1134
void clear_eigen()
Definition: Factory.hpp:1774
void clear_auxiliary()
Definition: Factory.hpp:1787
void set_geometry(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:954
void set_reference_load(const SpMat< T > &)
Definition: Factory.hpp:858
const T & get_pre_time() const
Definition: Factory.hpp:1058
void update_incre_resistance(const Col< T > &)
Definition: Factory.hpp:1189
void update_trial_time(const T &)
Definition: Factory.hpp:1114
void update_current_load_factor_by(const Col< T > &)
Definition: Factory.hpp:1394
void set_trial_load_factor(const Col< T > &)
Definition: Factory.hpp:862
void initialize_eigen()
Definition: Factory.hpp:841
void update_incre_inertial_force(const Col< T > &)
Definition: Factory.hpp:1199
void clear_time()
Definition: Factory.hpp:1610
void clear_displacement()
Definition: Factory.hpp:1654
const T & get_trial_time() const
Definition: Factory.hpp:992
void set_mass(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:948
const shared_ptr< MetaMat< T > > & get_geometry() const
Definition: Factory.hpp:1086
void reset_status()
Definition: Factory.hpp:1687
void update_current_inertial_force(const Col< T > &)
Definition: Factory.hpp:1254
const Col< T > & get_trial_temperature() const
Definition: Factory.hpp:1012
void set_trial_damping_force(const Col< T > &)
Definition: Factory.hpp:870
std::mutex & get_trial_settlement_mutex()
Definition: Factory.hpp:1100
void initialize_mass()
Definition: Factory.hpp:829
void set_pre_damping_force(const Col< T > &)
Definition: Factory.hpp:936
void clear_temperature()
Definition: Factory.hpp:1675
void initialize_acceleration()
Definition: Factory.hpp:812
void update_trial_inertial_force_by(const Col< T > &)
Definition: Factory.hpp:1309
void set_current_temperature(const Col< T > &)
Definition: Factory.hpp:924
void reset_acceleration()
Definition: Factory.hpp:1757
void set_current_time(const T &)
Definition: Factory.hpp:904
const Col< T > & get_pre_resistance() const
Definition: Factory.hpp:1066
void clear_stiffness()
Definition: Factory.hpp:1783
std::set< T > set
Definition: container.h:54
void suanpan_info(const char *M,...)
Definition: print.cpp:47