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