suanPan
Factory.hpp
Go to the documentation of this file.
1 /*******************************************************************************
2  * Copyright (C) 2017-2024 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 #include <Toolbox/container.h>
36 #include <Element/MappingDOF.h>
37 #include <Domain/MetaMat/MetaMat>
38 
39 #ifdef SUANPAN_MAGMA
40 #include <magmasparse.h>
41 #endif
42 
43 enum class AnalysisType {
44  NONE,
45  DISP,
46  EIGEN,
47  BUCKLE,
48  STATICS,
49  DYNAMICS
50 };
51 
52 enum class StorageScheme {
53  FULL,
54  BAND,
55  BANDSYMM,
56  SYMMPACK,
57  SPARSE,
59 };
60 
61 enum class SolverType {
62  LAPACK,
63  SPIKE,
64  SUPERLU,
65  MUMPS,
66  CUDA,
67  PARDISO,
68  FGMRES,
69  MAGMA,
70  LIS
71 };
72 
73 template<sp_d T> class Factory final {
74  unsigned n_size = 0; // number of degrees of freedom
75  unsigned n_lobw = 0; // low bandwidth
76  unsigned n_upbw = 0; // up bandwidth
77  unsigned n_sfbw = n_lobw + n_upbw; // matrix storage offset
78  unsigned n_rfld = 0; // reference load size
79  unsigned n_mpc = 0; // multipoint constraint size
80  uword n_elem = 0;
81 
82  AnalysisType analysis_type = AnalysisType::NONE; // type of analysis
83  StorageScheme storage_type = StorageScheme::FULL; // type of analysis
84 
85 #ifdef SUANPAN_MAGMA
86  magma_dopts magma_setting{};
87 #endif
88 
89  bool nlgeom = false;
90  bool nonviscous = false;
91 
93  SolverSetting<T> setting{};
94 
95  T error = T(0); // error produced by certain solvers
96 
97  Col<T> ninja; // the result from A*X=B
98  Col<T> sushi; // modified right-hand side B
99 
100  suanpan::set<uword> reference_dof;
101  SpMat<T> reference_load;
102 
103  uvec auxiliary_encoding; // for constraints using multiplier method
104  Col<T> auxiliary_lambda; // for constraints using multiplier method
105  Col<T> auxiliary_resistance; // for constraints using multiplier method
106  Col<T> auxiliary_load; // for constraints using multiplier method
107  SpMat<T> auxiliary_stiffness; // for constraints using multiplier method
108 
109  SpCol<T> trial_constraint_resistance;
110  SpCol<T> current_constraint_resistance;
111 
112  T trial_time = T(0); // global trial (pseudo) time
113  T incre_time = T(0); // global incremental (pseudo) time
114  T current_time = T(0); // global current (pseudo) time
115  T pre_time = T(0); // global previous (pseudo) time
116 
117  T strain_energy = T(0);
118  T kinetic_energy = T(0);
119  T viscous_energy = T(0);
120  T nonviscous_energy = T(0);
121  T complementary_energy = T(0);
122  Col<T> momentum;
123 
124  Col<T> trial_load_factor; // global trial load factor
125  Col<T> trial_load; // global trial load vector
126  Col<T> trial_settlement; // global trial displacement load vector
127  Col<T> trial_resistance; // global trial resistance vector
128  Col<T> trial_damping_force; // global trial damping force vector
129  Col<T> trial_nonviscous_force; // global trial nonviscous damping force vector
130  Col<T> trial_inertial_force; // global trial inertial force vector
131  Col<T> trial_displacement; // global trial displacement vector
132  Col<T> trial_velocity; // global trial velocity vector
133  Col<T> trial_acceleration; // global trial acceleration vector
134  Col<T> trial_temperature; // global trial temperature vector
135 
136  Col<T> incre_load_factor; // global incremental load vector
137  Col<T> incre_load; // global incremental load vector
138  Col<T> incre_settlement; // global incremental displacement load vector
139  Col<T> incre_resistance; // global incremental resistance vector
140  Col<T> incre_damping_force; // global incremental damping force vector
141  Col<T> incre_nonviscous_force; // global incremental nonviscous damping force vector
142  Col<T> incre_inertial_force; // global incremental inertial force vector
143  Col<T> incre_displacement; // global incremental displacement vector
144  Col<T> incre_velocity; // global incremental velocity vector
145  Col<T> incre_acceleration; // global incremental acceleration vector
146  Col<T> incre_temperature; // global incremental temperature vector
147 
148  Col<T> current_load_factor; // global current load vector
149  Col<T> current_load; // global current load vector
150  Col<T> current_settlement; // global current displacement load vector
151  Col<T> current_resistance; // global current resistance vector
152  Col<T> current_damping_force; // global current damping force vector
153  Col<T> current_nonviscous_force; // global current nonviscous damping force vector
154  Col<T> current_inertial_force; // global current inertial force vector
155  Col<T> current_displacement; // global current displacement vector
156  Col<T> current_velocity; // global current velocity vector
157  Col<T> current_acceleration; // global current acceleration vector
158  Col<T> current_temperature; // global current temperature vector
159 
160  Col<T> pre_load_factor; // global previous load vector
161  Col<T> pre_load; // global previous load vector
162  Col<T> pre_settlement; // global previous displacement load vector
163  Col<T> pre_resistance; // global previous resistance vector
164  Col<T> pre_damping_force; // global previous damping force vector
165  Col<T> pre_nonviscous_force; // global previous nonviscous damping force vector
166  Col<T> pre_inertial_force; // global previous inertial force vector
167  Col<T> pre_displacement; // global previous displacement vector
168  Col<T> pre_velocity; // global previous velocity vector
169  Col<T> pre_acceleration; // global previous acceleration vector
170  Col<T> pre_temperature; // global previous temperature vector
171 
172  shared_ptr<MetaMat<T>> global_mass = nullptr; // global mass matrix
173  shared_ptr<MetaMat<T>> global_damping = nullptr; // global damping matrix
174  shared_ptr<MetaMat<T>> global_nonviscous = nullptr; // global nonviscous damping matrix
175  shared_ptr<MetaMat<T>> global_stiffness = nullptr; // global stiffness matrix
176  shared_ptr<MetaMat<T>> global_geometry = nullptr; // global geometry matrix
177 
178  std::vector<std::mutex> global_mutex{20};
179 
180  Col<T> eigenvalue; // eigenvalues
181 
182  Mat<T> eigenvector; // eigenvectors
183 
184  unique_ptr<MetaMat<T>> get_basic_container();
185  unique_ptr<MetaMat<T>> get_matrix_container();
186 
187  void assemble_matrix_helper(shared_ptr<MetaMat<T>>&, const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
188 
189 public:
190  const bool initialized = false;
191 
193 
194  void set_size(unsigned);
195  [[nodiscard]] unsigned get_size() const;
196 
197  void set_entry(uword);
198  [[nodiscard]] uword get_entry() const;
199 
200  void set_nlgeom(bool);
201  [[nodiscard]] bool is_nlgeom() const;
202 
203  void set_nonviscous(bool);
204  [[nodiscard]] bool is_nonviscous() const;
205 
207  [[nodiscard]] SolverType get_solver_type() const;
208 
210  [[nodiscard]] const SolverSetting<double>& get_solver_setting() const;
211 
212 #ifdef SUANPAN_MAGMA
213  void set_solver_setting(const magma_dopts& magma_opt) { magma_setting = magma_opt; }
214 
215  [[nodiscard]] const magma_dopts& get_magma_setting() const { return magma_setting; }
216 #endif
217 
219  [[nodiscard]] AnalysisType get_analysis_type() const;
220 
222  [[nodiscard]] StorageScheme get_storage_scheme() const;
223 
224  [[nodiscard]] bool is_sparse() const;
225 
226  void set_bandwidth(unsigned, unsigned);
227  [[nodiscard]] std::pair<unsigned, unsigned> get_bandwidth() const;
228 
229  void update_reference_size();
230  void set_reference_size(unsigned);
231  [[nodiscard]] unsigned get_reference_size() const;
232 
233  void update_reference_dof(const uvec&);
235  [[nodiscard]] const suanpan::set<uword>& get_reference_dof() const;
236 
237  void set_error(T);
238  T get_error() const;
239 
240  /*************************INITIALIZER*************************/
241 
242  int initialize();
243 
244  void initialize_load_factor();
245  void initialize_load();
246  void initialize_settlement();
247  void initialize_resistance();
252  void initialize_velocity();
254  void initialize_temperature();
256 
257  void initialize_mass();
258  void initialize_damping();
259  void initialize_nonviscous();
260  void initialize_stiffness();
261  void initialize_geometry();
262  void initialize_eigen();
263 
264  /*************************SETTER*************************/
265 
266  void set_ninja(const Col<T>&);
267  void set_sushi(const Col<T>&);
268 
269  void update_sushi_by(const Col<T>&);
270 
271  void set_mpc(unsigned);
272 
273  void set_reference_load(const SpMat<T>&);
274 
275  void set_trial_time(T);
276  void set_trial_load_factor(const Col<T>&);
277  void set_trial_load(const Col<T>&);
278  void set_trial_settlement(const Col<T>&);
279  void set_trial_resistance(const Col<T>&);
280  void set_trial_damping_force(const Col<T>&);
281  void set_trial_nonviscous_force(const Col<T>&);
282  void set_trial_inertial_force(const Col<T>&);
283  void set_trial_displacement(const Col<T>&);
284  void set_trial_velocity(const Col<T>&);
285  void set_trial_acceleration(const Col<T>&);
286  void set_trial_temperature(const Col<T>&);
287 
288  void set_incre_time(T);
289  void set_incre_load_factor(const Col<T>&);
290  void set_incre_load(const Col<T>&);
291  void set_incre_settlement(const Col<T>&);
292  void set_incre_resistance(const Col<T>&);
293  void set_incre_damping_force(const Col<T>&);
294  void set_incre_nonviscous_force(const Col<T>&);
295  void set_incre_inertial_force(const Col<T>&);
296  void set_incre_displacement(const Col<T>&);
297  void set_incre_velocity(const Col<T>&);
298  void set_incre_acceleration(const Col<T>&);
299  void set_incre_temperature(const Col<T>&);
300 
301  void set_current_time(T);
302  void set_current_load_factor(const Col<T>&);
303  void set_current_load(const Col<T>&);
304  void set_current_settlement(const Col<T>&);
305  void set_current_resistance(const Col<T>&);
306  void set_current_damping_force(const Col<T>&);
307  void set_current_nonviscous_force(const Col<T>&);
308  void set_current_inertial_force(const Col<T>&);
309  void set_current_displacement(const Col<T>&);
310  void set_current_velocity(const Col<T>&);
311  void set_current_acceleration(const Col<T>&);
312  void set_current_temperature(const Col<T>&);
313 
314  void set_pre_time(T);
315  void set_pre_load_factor(const Col<T>&);
316  void set_pre_load(const Col<T>&);
317  void set_pre_settlement(const Col<T>&);
318  void set_pre_resistance(const Col<T>&);
319  void set_pre_damping_force(const Col<T>&);
320  void set_pre_nonviscous_force(const Col<T>&);
321  void set_pre_inertial_force(const Col<T>&);
322  void set_pre_displacement(const Col<T>&);
323  void set_pre_velocity(const Col<T>&);
324  void set_pre_acceleration(const Col<T>&);
325  void set_pre_temperature(const Col<T>&);
326 
327  void set_mass(const shared_ptr<MetaMat<T>>&);
328  void set_damping(const shared_ptr<MetaMat<T>>&);
329  void set_nonviscous(const shared_ptr<MetaMat<T>>&);
330  void set_stiffness(const shared_ptr<MetaMat<T>>&);
331  void set_geometry(const shared_ptr<MetaMat<T>>&);
332 
333  void set_eigenvalue(const Col<T>&);
334  void set_eigenvector(const Mat<T>&);
335 
336  /*************************GETTER*************************/
337 
338  const Col<T>& get_ninja() const;
339  const Col<T>& get_sushi() const;
340 
341  [[nodiscard]] unsigned get_mpc() const;
342 
343  const SpMat<T>& get_reference_load() const;
344 
345  [[nodiscard]] const uvec& get_auxiliary_encoding() const;
346  const Col<T>& get_auxiliary_lambda() const;
347  const Col<T>& get_auxiliary_resistance() const;
348  const Col<T>& get_auxiliary_load() const;
349  const SpMat<T>& get_auxiliary_stiffness() const;
350 
351  const SpCol<T>& get_trial_constraint_resistance() const;
352  const SpCol<T>& get_current_constraint_resistance() const;
353 
359  const Col<T>& get_momentum();
360 
361  T get_trial_time() const;
362  const Col<T>& get_trial_load_factor() const;
363  const Col<T>& get_trial_load() const;
364  const Col<T>& get_trial_settlement() const;
365  const Col<T>& get_trial_resistance() const;
366  const Col<T>& get_trial_damping_force() const;
367  const Col<T>& get_trial_nonviscous_force() const;
368  const Col<T>& get_trial_inertial_force() const;
369  const Col<T>& get_trial_displacement() const;
370  const Col<T>& get_trial_velocity() const;
371  const Col<T>& get_trial_acceleration() const;
372  const Col<T>& get_trial_temperature() const;
373 
374  T get_incre_time() const;
375  const Col<T>& get_incre_load_factor() const;
376  const Col<T>& get_incre_load() const;
377  const Col<T>& get_incre_settlement() const;
378  const Col<T>& get_incre_resistance() const;
379  const Col<T>& get_incre_damping_force() const;
380  const Col<T>& get_incre_nonviscous_force() const;
381  const Col<T>& get_incre_inertial_force() const;
382  const Col<T>& get_incre_displacement() const;
383  const Col<T>& get_incre_velocity() const;
384  const Col<T>& get_incre_acceleration() const;
385  const Col<T>& get_incre_temperature() const;
386 
387  T get_current_time() const;
388  const Col<T>& get_current_load_factor() const;
389  const Col<T>& get_current_load() const;
390  const Col<T>& get_current_settlement() const;
391  const Col<T>& get_current_resistance() const;
392  const Col<T>& get_current_damping_force() const;
393  const Col<T>& get_current_nonviscous_force() const;
394  const Col<T>& get_current_inertial_force() const;
395  const Col<T>& get_current_displacement() const;
396  const Col<T>& get_current_velocity() const;
397  const Col<T>& get_current_acceleration() const;
398  const Col<T>& get_current_temperature() const;
399 
400  T get_pre_time() const;
401  const Col<T>& get_pre_load_factor() const;
402  const Col<T>& get_pre_load() const;
403  const Col<T>& get_pre_settlement() const;
404  const Col<T>& get_pre_resistance() const;
405  const Col<T>& get_pre_damping_force() const;
406  const Col<T>& get_pre_nonviscous_force() const;
407  const Col<T>& get_pre_inertial_force() const;
408  const Col<T>& get_pre_displacement() const;
409  const Col<T>& get_pre_velocity() const;
410  const Col<T>& get_pre_acceleration() const;
411  const Col<T>& get_pre_temperature() const;
412 
413  const shared_ptr<MetaMat<T>>& get_mass() const;
414  const shared_ptr<MetaMat<T>>& get_damping() const;
415  const shared_ptr<MetaMat<T>>& get_nonviscous() const;
416  const shared_ptr<MetaMat<T>>& get_stiffness() const;
417  const shared_ptr<MetaMat<T>>& get_geometry() const;
418 
419  std::mutex& get_auxiliary_encoding_mutex();
420  std::mutex& get_auxiliary_resistance_mutex();
421  std::mutex& get_auxiliary_load_mutex();
422  std::mutex& get_auxiliary_stiffness_mutex();
423 
425 
426  std::mutex& get_trial_load_mutex();
427  std::mutex& get_trial_settlement_mutex();
428  std::mutex& get_reference_load_mutex();
429 
430  std::mutex& get_mass_mutex();
431  std::mutex& get_damping_mutex();
432  std::mutex& get_nonviscous_mutex();
433  std::mutex& get_stiffness_mutex();
434  std::mutex& get_geometry_mutex();
435 
436  const Col<T>& get_eigenvalue() const;
437  const Mat<T>& get_eigenvector() const;
438 
439  /*************************UPDATER*************************/
440 
441  void update_trial_time(T);
442  void update_trial_load_factor(const Col<T>&);
443  void update_trial_load(const Col<T>&);
444  void update_trial_settlement(const Col<T>&);
445  void update_trial_resistance(const Col<T>&);
446  void update_trial_damping_force(const Col<T>&);
447  void update_trial_nonviscous_force(const Col<T>&);
448  void update_trial_inertial_force(const Col<T>&);
449  void update_trial_displacement(const Col<T>&);
450  void update_trial_velocity(const Col<T>&);
451  void update_trial_acceleration(const Col<T>&);
452  void update_trial_temperature(const Col<T>&);
453 
454  void update_incre_time(T);
455  void update_incre_load_factor(const Col<T>&);
456  void update_incre_load(const Col<T>&);
457  void update_incre_settlement(const Col<T>&);
458  void update_incre_resistance(const Col<T>&);
459  void update_incre_damping_force(const Col<T>&);
460  void update_incre_nonviscous_force(const Col<T>&);
461  void update_incre_inertial_force(const Col<T>&);
462  void update_incre_displacement(const Col<T>&);
463  void update_incre_velocity(const Col<T>&);
464  void update_incre_acceleration(const Col<T>&);
465  void update_incre_temperature(const Col<T>&);
466 
467  void update_current_time(T);
468  void update_current_load_factor(const Col<T>&);
469  void update_current_load(const Col<T>&);
470  void update_current_settlement(const Col<T>&);
471  void update_current_resistance(const Col<T>&);
472  void update_current_damping_force(const Col<T>&);
473  void update_current_nonviscous_force(const Col<T>&);
474  void update_current_inertial_force(const Col<T>&);
475  void update_current_displacement(const Col<T>&);
476  void update_current_velocity(const Col<T>&);
477  void update_current_acceleration(const Col<T>&);
478  void update_current_temperature(const Col<T>&);
479 
480  void update_trial_time_by(T);
481  void update_trial_load_factor_by(const Col<T>&);
482  void update_trial_load_by(const Col<T>&);
483  void update_trial_settlement_by(const Col<T>&);
484  void update_trial_resistance_by(const Col<T>&);
485  void update_trial_damping_force_by(const Col<T>&);
486  void update_trial_nonviscous_force_by(const Col<T>&);
487  void update_trial_inertial_force_by(const Col<T>&);
488  void update_trial_displacement_by(const Col<T>&);
489  void update_trial_velocity_by(const Col<T>&);
490  void update_trial_acceleration_by(const Col<T>&);
491  void update_trial_temperature_by(const Col<T>&);
492 
493  void update_incre_time_by(T);
494  void update_incre_load_factor_by(const Col<T>&);
495  void update_incre_load_by(const Col<T>&);
496  void update_incre_settlement_by(const Col<T>&);
497  void update_incre_resistance_by(const Col<T>&);
498  void update_incre_damping_force_by(const Col<T>&);
499  void update_incre_nonviscous_force_by(const Col<T>&);
500  void update_incre_inertial_force_by(const Col<T>&);
501  void update_incre_displacement_by(const Col<T>&);
502  void update_incre_velocity_by(const Col<T>&);
503  void update_incre_acceleration_by(const Col<T>&);
504  void update_incre_temperature_by(const Col<T>&);
505 
507  void update_current_load_factor_by(const Col<T>&);
508  void update_current_load_by(const Col<T>&);
509  void update_current_settlement_by(const Col<T>&);
510  void update_current_resistance_by(const Col<T>&);
511  void update_current_damping_force_by(const Col<T>&);
512  void update_current_nonviscous_force_by(const Col<T>&);
513  void update_current_inertial_force_by(const Col<T>&);
514  void update_current_displacement_by(const Col<T>&);
515  void update_current_velocity_by(const Col<T>&);
516  void update_current_acceleration_by(const Col<T>&);
517  void update_current_temperature_by(const Col<T>&);
518 
519  /*************************FRIEND*************************/
520 
521  Col<T>& modify_ninja();
522  Col<T>& modify_sushi();
523 
525  SpMat<T>& modify_reference_load();
526 
528  Col<T>& modify_auxiliary_lambda();
529  Col<T>& modify_auxiliary_resistance();
530  Col<T>& modify_auxiliary_load();
531  SpMat<T>& modify_auxiliary_stiffness();
532 
535 
536  T& modify_trial_time();
537  Col<T>& modify_trial_load_factor();
538  Col<T>& modify_trial_load();
539  Col<T>& modify_trial_settlement();
540  Col<T>& modify_trial_resistance();
541  Col<T>& modify_trial_damping_force();
543  Col<T>& modify_trial_inertial_force();
544  Col<T>& modify_trial_displacement();
545  Col<T>& modify_trial_velocity();
546  Col<T>& modify_trial_acceleration();
547  Col<T>& modify_trial_temperature();
548 
549  T& modify_incre_time();
550  Col<T>& modify_incre_load_factor();
551  Col<T>& modify_incre_load();
552  Col<T>& modify_incre_settlement();
553  Col<T>& modify_incre_resistance();
554  Col<T>& modify_incre_damping_force();
556  Col<T>& modify_incre_inertial_force();
557  Col<T>& modify_incre_displacement();
558  Col<T>& modify_incre_velocity();
559  Col<T>& modify_incre_acceleration();
560  Col<T>& modify_incre_temperature();
561 
563  Col<T>& modify_current_load_factor();
564  Col<T>& modify_current_load();
565  Col<T>& modify_current_settlement();
566  Col<T>& modify_current_resistance();
570  Col<T>& modify_current_displacement();
571  Col<T>& modify_current_velocity();
572  Col<T>& modify_current_acceleration();
573  Col<T>& modify_current_temperature();
574 
575  T& modify_pre_time();
576  Col<T>& modify_pre_load_factor();
577  Col<T>& modify_pre_load();
578  Col<T>& modify_pre_settlement();
579  Col<T>& modify_pre_resistance();
580  Col<T>& modify_pre_damping_force();
581  Col<T>& modify_pre_nonviscous_force();
582  Col<T>& modify_pre_inertial_force();
583  Col<T>& modify_pre_displacement();
584  Col<T>& modify_pre_velocity();
585  Col<T>& modify_pre_acceleration();
586  Col<T>& modify_pre_temperature();
587 
588  shared_ptr<MetaMat<T>>& modify_mass();
589  shared_ptr<MetaMat<T>>& modify_damping();
590  shared_ptr<MetaMat<T>>& modify_nonviscous();
591  shared_ptr<MetaMat<T>>& modify_stiffness();
592  shared_ptr<MetaMat<T>>& modify_geometry();
593 
594  Col<T>& modify_eigenvalue();
595  Mat<T>& modify_eigenvector();
596 
597  /*************************STATUS*************************/
598 
599  void commit_energy();
600  void clear_energy();
601 
602  void commit_status();
603  void commit_time();
604  void commit_load_factor();
605  void commit_load();
606  void commit_settlement();
607  void commit_resistance();
608  void commit_damping_force();
610  void commit_inertial_force();
611  void commit_displacement();
612  void commit_velocity();
613  void commit_acceleration();
614  void commit_temperature();
616 
617  void commit_pre_status();
618  void commit_pre_time();
619  void commit_pre_load_factor();
620  void commit_pre_load();
621  void commit_pre_settlement();
622  void commit_pre_resistance();
627  void commit_pre_velocity();
629  void commit_pre_temperature();
630 
631  void clear_status();
632  void clear_time();
633  void clear_load_factor();
634  void clear_load();
635  void clear_settlement();
636  void clear_resistance();
637  void clear_damping_force();
638  void clear_nonviscous_force();
639  void clear_inertial_force();
640  void clear_displacement();
641  void clear_velocity();
642  void clear_acceleration();
643  void clear_temperature();
645 
646  void reset_status();
647  void reset_time();
648  void reset_load_factor();
649  void reset_load();
650  void reset_settlement();
651  void reset_resistance();
652  void reset_damping_force();
653  void reset_nonviscous_force();
654  void reset_inertial_force();
655  void reset_displacement();
656  void reset_velocity();
657  void reset_acceleration();
658  void reset_temperature();
660 
661  void clear_eigen();
662  void clear_mass();
663  void clear_damping();
664  void clear_nonviscous();
665  void clear_stiffness();
666  void clear_geometry();
667  void clear_auxiliary();
668 
669  void reset();
670 
671  /*************************ASSEMBLER*************************/
672 
673  void assemble_resistance(const Mat<T>&, const uvec&);
674  void assemble_damping_force(const Mat<T>&, const uvec&);
675  void assemble_nonviscous_force(const Mat<T>&, const uvec&);
676  void assemble_inertial_force(const Mat<T>&, const uvec&);
677 
678  void assemble_mass(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
679  void assemble_damping(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
680  void assemble_nonviscous(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
681  void assemble_stiffness(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
682  void assemble_geometry(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
683 
684  void assemble_stiffness(const SpMat<T>&, const uvec&);
685 
686  /*************************UTILITY*************************/
687 
688  void print() const;
689 };
690 
691 template<sp_d T> Factory<T>::Factory(const unsigned D, const AnalysisType AT, const StorageScheme SS)
692  : n_size(D)
693  , analysis_type(AT)
694  , storage_type(SS) {}
695 
696 template<sp_d T> void Factory<T>::set_size(const unsigned D) {
697  if(D == n_size) return;
698  n_size = D;
699  access::rw(initialized) = false;
700 }
701 
702 template<sp_d T> unsigned Factory<T>::get_size() const { return n_size; }
703 
704 template<sp_d T> void Factory<T>::set_entry(const uword N) {
705  n_elem = N;
706  if(n_elem > std::numeric_limits<int>::max()) throw invalid_argument("too many elements");
707 }
708 
709 template<sp_d T> uword Factory<T>::get_entry() const { return n_elem; }
710 
711 template<sp_d T> void Factory<T>::set_nlgeom(const bool B) {
712  if(B == nlgeom) return;
713  nlgeom = B;
714  access::rw(initialized) = false;
715 }
716 
717 template<sp_d T> bool Factory<T>::is_nlgeom() const { return nlgeom; }
718 
719 template<sp_d T> void Factory<T>::set_nonviscous(const bool B) {
720  if(B == nonviscous) return;
721  nonviscous = B;
722  access::rw(initialized) = false;
723 }
724 
725 template<sp_d T> bool Factory<T>::is_nonviscous() const { return nonviscous; }
726 
727 template<sp_d T> void Factory<T>::set_solver_type(const SolverType E) { solver = E; }
728 
729 template<sp_d T> SolverType Factory<T>::get_solver_type() const { return solver; }
730 
731 template<sp_d T> void Factory<T>::set_solver_setting(const SolverSetting<double>& SS) { setting = SS; }
732 
733 template<sp_d T> const SolverSetting<double>& Factory<T>::get_solver_setting() const { return setting; }
734 
735 template<sp_d T> void Factory<T>::set_analysis_type(const AnalysisType AT) {
736  if(AT == analysis_type) return;
737  analysis_type = AT;
738  access::rw(initialized) = false;
739 }
740 
741 template<sp_d T> AnalysisType Factory<T>::get_analysis_type() const { return analysis_type; }
742 
743 template<sp_d T> void Factory<T>::set_storage_scheme(const StorageScheme SS) {
744  if(SS == storage_type) return;
745  storage_type = SS;
746  access::rw(initialized) = false;
747 }
748 
749 template<sp_d T> StorageScheme Factory<T>::get_storage_scheme() const { return storage_type; }
750 
751 template<sp_d T> bool Factory<T>::is_sparse() const { return StorageScheme::SPARSE == storage_type || StorageScheme::SPARSESYMM == storage_type; }
752 
753 template<sp_d T> void Factory<T>::set_bandwidth(const unsigned L, const unsigned U) {
754  if(L == n_lobw && U == n_upbw) return;
755  n_lobw = L;
756  n_upbw = U;
757  n_sfbw = L + U;
758  access::rw(initialized) = false;
759 }
760 
761 template<sp_d T> std::pair<unsigned, unsigned> Factory<T>::get_bandwidth() const { return {n_lobw, n_upbw}; }
762 
763 template<sp_d T> void Factory<T>::update_reference_size() { n_rfld = static_cast<unsigned>(reference_dof.size()); }
764 
765 template<sp_d T> void Factory<T>::set_reference_size(const unsigned S) {
766  if(S == n_rfld) return;
767  n_rfld = S;
768 }
769 
770 template<sp_d T> unsigned Factory<T>::get_reference_size() const { return n_rfld; }
771 
772 template<sp_d T> void Factory<T>::update_reference_dof(const uvec& S) { reference_dof.insert(S.cbegin(), S.cend()); }
773 
774 template<sp_d T> void Factory<T>::set_reference_dof(const suanpan::set<uword>& D) { reference_dof = D; }
775 
776 template<sp_d T> const suanpan::set<uword>& Factory<T>::get_reference_dof() const { return reference_dof; }
777 
778 template<sp_d T> void Factory<T>::set_error(const T E) { error = E; }
779 
780 template<sp_d T> T Factory<T>::get_error() const { return error; }
781 
782 template<sp_d T> int Factory<T>::initialize() {
783  reference_dof.clear(); // clear reference dof vector in every step
784 
785  if(initialized || n_size == 0) return 0;
786 
787  ninja.zeros(n_size);
788  sushi.zeros(n_size);
789 
790  reset();
791 
792  switch(analysis_type) {
793  case AnalysisType::DISP:
794  initialize_displacement();
795  break;
796  case AnalysisType::EIGEN:
797  initialize_mass();
798  initialize_stiffness();
799  initialize_eigen();
800  break;
803  initialize_load();
804  initialize_resistance();
805  initialize_displacement();
806  initialize_stiffness();
807  initialize_geometry();
808  break;
810  initialize_load();
811  initialize_resistance();
812  initialize_damping_force();
813  initialize_nonviscous_force();
814  initialize_inertial_force();
815  initialize_displacement();
816  initialize_velocity();
817  initialize_acceleration();
818  initialize_mass();
819  initialize_damping();
820  initialize_nonviscous();
821  initialize_stiffness();
822  initialize_geometry();
823  break;
824  case AnalysisType::NONE:
825  break;
826  }
827 
828  initialize_auxiliary_resistance();
829 
830  access::rw(initialized) = true;
831 
832  return 0;
833 }
834 
835 template<sp_d T> void Factory<T>::initialize_load_factor() {
836  if(n_rfld == 0) return;
837 
838  trial_load_factor.zeros(n_rfld);
839  incre_load_factor.zeros(n_rfld);
840  current_load_factor.zeros(n_rfld);
841 
842  reference_load.zeros(n_size, n_rfld);
843 }
844 
845 template<sp_d T> void Factory<T>::initialize_load() {
846  trial_load.zeros(n_size);
847  incre_load.zeros(n_size);
848  current_load.zeros(n_size);
849 }
850 
851 template<sp_d T> void Factory<T>::initialize_settlement() {
852  trial_settlement.zeros(n_size);
853  incre_settlement.zeros(n_size);
854  current_settlement.zeros(n_size);
855 }
856 
857 template<sp_d T> void Factory<T>::initialize_resistance() {
858  trial_resistance.zeros(n_size);
859  incre_resistance.zeros(n_size);
860  current_resistance.zeros(n_size);
861 }
862 
863 template<sp_d T> void Factory<T>::initialize_damping_force() {
864  trial_damping_force.zeros(n_size);
865  incre_damping_force.zeros(n_size);
866  current_damping_force.zeros(n_size);
867 }
868 
870  trial_nonviscous_force.zeros(n_size);
871  incre_nonviscous_force.zeros(n_size);
872  current_nonviscous_force.zeros(n_size);
873 }
874 
875 template<sp_d T> void Factory<T>::initialize_inertial_force() {
876  trial_inertial_force.zeros(n_size);
877  incre_inertial_force.zeros(n_size);
878  current_inertial_force.zeros(n_size);
879 }
880 
881 template<sp_d T> void Factory<T>::initialize_displacement() {
882  trial_displacement.zeros(n_size);
883  incre_displacement.zeros(n_size);
884  current_displacement.zeros(n_size);
885 }
886 
887 template<sp_d T> void Factory<T>::initialize_velocity() {
888  trial_velocity.zeros(n_size);
889  incre_velocity.zeros(n_size);
890  current_velocity.zeros(n_size);
891 }
892 
893 template<sp_d T> void Factory<T>::initialize_acceleration() {
894  trial_acceleration.zeros(n_size);
895  incre_acceleration.zeros(n_size);
896  current_acceleration.zeros(n_size);
897 }
898 
899 template<sp_d T> void Factory<T>::initialize_temperature() {
900  trial_temperature.zeros(n_size);
901  incre_temperature.zeros(n_size);
902  current_temperature.zeros(n_size);
903 }
904 
906  trial_constraint_resistance.zeros(n_size);
907  current_constraint_resistance.zeros(n_size);
908 }
909 
910 template<sp_d T> void Factory<T>::initialize_mass() { global_mass = get_matrix_container(); }
911 
912 template<sp_d T> void Factory<T>::initialize_damping() { global_damping = get_matrix_container(); }
913 
914 template<sp_d T> void Factory<T>::initialize_nonviscous() {
915  if(!nonviscous) return;
916 
917  global_nonviscous = get_matrix_container();
918 }
919 
920 template<sp_d T> void Factory<T>::initialize_stiffness() { global_stiffness = get_matrix_container(); }
921 
922 template<sp_d T> void Factory<T>::initialize_geometry() {
923  if(!nlgeom) return;
924 
925  global_geometry = get_matrix_container();
926 }
927 
928 template<sp_d T> void Factory<T>::initialize_eigen() {
929  eigenvalue.zeros(n_size);
930  eigenvector.zeros(n_size, n_size);
931 }
932 
933 template<sp_d T> void Factory<T>::set_ninja(const Col<T>& N) { ninja = N; }
934 
935 template<sp_d T> void Factory<T>::set_sushi(const Col<T>& S) { sushi = S; }
936 
937 template<sp_d T> void Factory<T>::update_sushi_by(const Col<T>& S) { sushi += S; }
938 
939 template<sp_d T> void Factory<T>::set_mpc(const unsigned S) {
940  n_mpc = S;
941  auxiliary_encoding.zeros(n_mpc);
942  auxiliary_resistance.zeros(n_mpc);
943  auxiliary_load.zeros(n_mpc);
944  auxiliary_stiffness.zeros(n_size, n_mpc);
945 }
946 
947 template<sp_d T> void Factory<T>::set_reference_load(const SpMat<T>& L) { reference_load = L; }
948 
949 template<sp_d T> void Factory<T>::set_mass(const shared_ptr<MetaMat<T>>& M) { global_mass = M; }
950 
951 template<sp_d T> void Factory<T>::set_damping(const shared_ptr<MetaMat<T>>& C) { global_damping = C; }
952 
953 template<sp_d T> void Factory<T>::set_nonviscous(const shared_ptr<MetaMat<T>>& C) { global_nonviscous = C; }
954 
955 template<sp_d T> void Factory<T>::set_stiffness(const shared_ptr<MetaMat<T>>& K) { global_stiffness = K; }
956 
957 template<sp_d T> void Factory<T>::set_geometry(const shared_ptr<MetaMat<T>>& G) { global_geometry = G; }
958 
959 template<sp_d T> void Factory<T>::set_eigenvalue(const Col<T>& L) { eigenvalue = L; }
960 
961 template<sp_d T> void Factory<T>::set_eigenvector(const Mat<T>& V) { eigenvector = V; }
962 
963 template<sp_d T> const Col<T>& Factory<T>::get_ninja() const { return ninja; }
964 
965 template<sp_d T> const Col<T>& Factory<T>::get_sushi() const { return sushi; }
966 
967 template<sp_d T> unsigned Factory<T>::get_mpc() const { return n_mpc; }
968 
969 template<sp_d T> const SpMat<T>& Factory<T>::get_reference_load() const { return reference_load; }
970 
971 template<sp_d T> const uvec& Factory<T>::get_auxiliary_encoding() const { return auxiliary_encoding; }
972 
973 template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_lambda() const { return auxiliary_lambda; }
974 
975 template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_resistance() const { return auxiliary_resistance; }
976 
977 template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_load() const { return auxiliary_load; }
978 
979 template<sp_d T> const SpMat<T>& Factory<T>::get_auxiliary_stiffness() const { return auxiliary_stiffness; }
980 
981 template<sp_d T> const SpCol<T>& Factory<T>::get_trial_constraint_resistance() const { return trial_constraint_resistance; }
982 
983 template<sp_d T> const SpCol<T>& Factory<T>::get_current_constraint_resistance() const { return current_constraint_resistance; }
984 
985 template<sp_d T> T Factory<T>::get_strain_energy() { return strain_energy; }
986 
987 template<sp_d T> T Factory<T>::get_kinetic_energy() { return kinetic_energy; }
988 
989 template<sp_d T> T Factory<T>::get_viscous_energy() { return viscous_energy; }
990 
991 template<sp_d T> T Factory<T>::get_nonviscous_energy() { return nonviscous_energy; }
992 
993 template<sp_d T> T Factory<T>::get_complementary_energy() { return complementary_energy; }
994 
995 template<sp_d T> const Col<T>& Factory<T>::get_momentum() { return momentum; }
996 
997 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_mass() const { return global_mass; }
998 
999 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_damping() const { return global_damping; }
1000 
1001 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_nonviscous() const { return global_nonviscous; }
1002 
1003 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_stiffness() const { return global_stiffness; }
1004 
1005 template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_geometry() const { return global_geometry; }
1006 
1007 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_encoding_mutex() { return global_mutex.at(0); }
1008 
1009 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_resistance_mutex() { return global_mutex.at(1); }
1010 
1011 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_load_mutex() { return global_mutex.at(2); }
1012 
1013 template<sp_d T> std::mutex& Factory<T>::get_auxiliary_stiffness_mutex() { return global_mutex.at(3); }
1014 
1015 template<sp_d T> std::mutex& Factory<T>::get_trial_constraint_resistance_mutex() { return global_mutex.at(4); }
1016 
1017 template<sp_d T> std::mutex& Factory<T>::get_trial_load_mutex() { return global_mutex.at(5); }
1018 
1019 template<sp_d T> std::mutex& Factory<T>::get_trial_settlement_mutex() { return global_mutex.at(6); }
1020 
1021 template<sp_d T> std::mutex& Factory<T>::get_reference_load_mutex() { return global_mutex.at(7); }
1022 
1023 template<sp_d T> std::mutex& Factory<T>::get_mass_mutex() { return global_mutex.at(8); }
1024 
1025 template<sp_d T> std::mutex& Factory<T>::get_damping_mutex() { return global_mutex.at(9); }
1026 
1027 template<sp_d T> std::mutex& Factory<T>::get_nonviscous_mutex() { return global_mutex.at(10); }
1028 
1029 template<sp_d T> std::mutex& Factory<T>::get_stiffness_mutex() { return global_mutex.at(11); }
1030 
1031 template<sp_d T> std::mutex& Factory<T>::get_geometry_mutex() { return global_mutex.at(12); }
1032 
1033 template<sp_d T> const Col<T>& Factory<T>::get_eigenvalue() const { return eigenvalue; }
1034 
1035 template<sp_d T> const Mat<T>& Factory<T>::get_eigenvector() const { return eigenvector; }
1036 
1037 template<sp_d T> void Factory<T>::commit_energy() {
1038  auto se = std::async([&] { if(!trial_resistance.empty() && !incre_displacement.empty()) strain_energy += .5 * dot(trial_resistance + current_resistance, incre_displacement); });
1039  auto ke = std::async([&] { if(!trial_inertial_force.empty() && !trial_velocity.empty()) kinetic_energy = .5 * dot(global_mass * trial_velocity, trial_velocity); });
1040  auto ve = std::async([&] { if(!trial_damping_force.empty() && !incre_displacement.empty()) viscous_energy += .5 * dot(trial_damping_force + current_damping_force, incre_displacement); });
1041  auto ne = std::async([&] { if(!trial_nonviscous_force.empty() && !incre_displacement.empty()) nonviscous_energy += .5 * dot(trial_nonviscous_force + current_nonviscous_force, incre_displacement); });
1042  auto ce = std::async([&] { if(!trial_displacement.empty() && !incre_resistance.empty()) complementary_energy += .5 * dot(trial_displacement + current_displacement, incre_resistance); });
1043  auto mm = std::async([&] { if(!trial_inertial_force.empty() && !trial_velocity.empty()) momentum = global_mass * trial_velocity; });
1044 
1045  se.get();
1046  ke.get();
1047  ve.get();
1048  ne.get();
1049  ce.get();
1050  mm.get();
1051 }
1052 
1053 template<sp_d T> void Factory<T>::clear_energy() {
1054  strain_energy = T(0);
1055  kinetic_energy = T(0);
1056  viscous_energy = T(0);
1057  nonviscous_energy = T(0);
1058  complementary_energy = T(0);
1059  momentum.zeros();
1060 }
1061 
1062 template<sp_d T> void Factory<T>::commit_status() {
1063  ninja.zeros();
1064 
1065  commit_energy();
1066 
1067  commit_time();
1068  commit_load_factor();
1069  commit_load();
1070  commit_settlement();
1071  commit_resistance();
1072  commit_damping_force();
1073  commit_nonviscous_force();
1074  commit_inertial_force();
1075  commit_displacement();
1076  commit_velocity();
1077  commit_acceleration();
1078  commit_temperature();
1079  commit_auxiliary_resistance();
1080 }
1081 
1082 template<sp_d T> void Factory<T>::commit_time() {
1083  current_time = trial_time;
1084  incre_time = T(0);
1085 }
1086 
1087 template<sp_d T> void Factory<T>::commit_load_factor() {
1088  if(trial_load_factor.is_empty()) return;
1089  current_load_factor = trial_load_factor;
1090  incre_load_factor.zeros();
1091 }
1092 
1093 template<sp_d T> void Factory<T>::commit_load() {
1094  if(trial_load.is_empty()) return;
1095  current_load = trial_load;
1096  incre_load.zeros();
1097 }
1098 
1099 template<sp_d T> void Factory<T>::commit_settlement() {
1100  if(trial_settlement.is_empty()) return;
1101  current_settlement = trial_settlement;
1102  incre_settlement.zeros();
1103 }
1104 
1105 template<sp_d T> void Factory<T>::commit_resistance() {
1106  if(trial_resistance.is_empty()) return;
1107  current_resistance = trial_resistance;
1108  incre_resistance.zeros();
1109 }
1110 
1111 template<sp_d T> void Factory<T>::commit_damping_force() {
1112  if(trial_damping_force.is_empty()) return;
1113  current_damping_force = trial_damping_force;
1114  incre_damping_force.zeros();
1115 }
1116 
1117 template<sp_d T> void Factory<T>::commit_nonviscous_force() {
1118  if(trial_nonviscous_force.is_empty()) return;
1119  current_nonviscous_force = trial_nonviscous_force;
1120  incre_nonviscous_force.zeros();
1121 }
1122 
1123 template<sp_d T> void Factory<T>::commit_inertial_force() {
1124  if(trial_inertial_force.is_empty()) return;
1125  current_inertial_force = trial_inertial_force;
1126  incre_inertial_force.zeros();
1127 }
1128 
1129 template<sp_d T> void Factory<T>::commit_displacement() {
1130  if(trial_displacement.is_empty()) return;
1131  current_displacement = trial_displacement;
1132  incre_displacement.zeros();
1133 }
1134 
1135 template<sp_d T> void Factory<T>::commit_velocity() {
1136  if(trial_velocity.is_empty()) return;
1137  current_velocity = trial_velocity;
1138  incre_velocity.zeros();
1139 }
1140 
1141 template<sp_d T> void Factory<T>::commit_acceleration() {
1142  if(trial_acceleration.is_empty()) return;
1143  current_acceleration = trial_acceleration;
1144  incre_acceleration.zeros();
1145 }
1146 
1147 template<sp_d T> void Factory<T>::commit_temperature() {
1148  if(trial_temperature.is_empty()) return;
1149  current_temperature = trial_temperature;
1150  incre_temperature.zeros();
1151 }
1152 
1154  if(trial_constraint_resistance.is_empty()) return;
1155  current_constraint_resistance = trial_constraint_resistance;
1156 }
1157 
1158 template<sp_d T> void Factory<T>::commit_pre_status() {
1159  commit_pre_time();
1160  commit_pre_load_factor();
1161  commit_pre_load();
1162  commit_pre_settlement();
1163  commit_pre_resistance();
1164  commit_pre_damping_force();
1165  commit_pre_nonviscous_force();
1166  commit_pre_inertial_force();
1167  commit_pre_displacement();
1168  commit_pre_velocity();
1169  commit_pre_acceleration();
1170  commit_pre_temperature();
1171 }
1172 
1173 template<sp_d T> void Factory<T>::commit_pre_time() { pre_time = current_time; }
1174 
1175 template<sp_d T> void Factory<T>::commit_pre_load_factor() { if(!current_load_factor.is_empty()) pre_load_factor = current_load_factor; }
1176 
1177 template<sp_d T> void Factory<T>::commit_pre_load() { if(!current_load.is_empty()) pre_load = current_load; }
1178 
1179 template<sp_d T> void Factory<T>::commit_pre_settlement() { if(!current_settlement.is_empty()) pre_settlement = current_settlement; }
1180 
1181 template<sp_d T> void Factory<T>::commit_pre_resistance() { if(!current_resistance.is_empty()) pre_resistance = current_resistance; }
1182 
1183 template<sp_d T> void Factory<T>::commit_pre_damping_force() { if(!current_damping_force.is_empty()) pre_damping_force = current_damping_force; }
1184 
1185 template<sp_d T> void Factory<T>::commit_pre_nonviscous_force() { if(!current_nonviscous_force.is_empty()) pre_nonviscous_force = current_nonviscous_force; }
1186 
1187 template<sp_d T> void Factory<T>::commit_pre_inertial_force() { if(!current_inertial_force.is_empty()) pre_inertial_force = current_inertial_force; }
1188 
1189 template<sp_d T> void Factory<T>::commit_pre_displacement() { if(!current_displacement.is_empty()) pre_displacement = current_displacement; }
1190 
1191 template<sp_d T> void Factory<T>::commit_pre_velocity() { if(!current_velocity.is_empty()) pre_velocity = current_velocity; }
1192 
1193 template<sp_d T> void Factory<T>::commit_pre_acceleration() { if(!current_acceleration.is_empty()) pre_acceleration = current_acceleration; }
1194 
1195 template<sp_d T> void Factory<T>::commit_pre_temperature() { if(!current_temperature.is_empty()) pre_temperature = current_temperature; }
1196 
1197 template<sp_d T> void Factory<T>::clear_status() {
1198  access::rw(initialized) = false;
1199 
1200  ninja.zeros();
1201 
1202  clear_energy();
1203 
1204  clear_time();
1205  clear_load_factor();
1206  clear_load();
1207  clear_settlement();
1208  clear_resistance();
1209  clear_damping_force();
1210  clear_nonviscous_force();
1211  clear_inertial_force();
1212  clear_displacement();
1213  clear_velocity();
1214  clear_acceleration();
1215  clear_temperature();
1216  clear_auxiliary_resistance();
1217 }
1218 
1219 template<sp_d T> void Factory<T>::clear_time() { trial_time = incre_time = current_time = 0.; }
1220 
1221 template<sp_d T> void Factory<T>::clear_load_factor() {
1222  if(!pre_load_factor.is_empty()) pre_load_factor.zeros();
1223  if(!trial_load_factor.is_empty()) trial_load_factor.zeros();
1224  if(!incre_load_factor.is_empty()) incre_load_factor.zeros();
1225  if(!current_load_factor.is_empty()) current_load_factor.zeros();
1226 }
1227 
1228 template<sp_d T> void Factory<T>::clear_load() {
1229  if(!pre_load.is_empty()) pre_load.zeros();
1230  if(!trial_load.is_empty()) trial_load.zeros();
1231  if(!incre_load.is_empty()) incre_load.zeros();
1232  if(!current_load.is_empty()) current_load.zeros();
1233 }
1234 
1235 template<sp_d T> void Factory<T>::clear_settlement() {
1236  if(!pre_settlement.is_empty()) pre_settlement.zeros();
1237  if(!trial_settlement.is_empty()) trial_settlement.zeros();
1238  if(!incre_settlement.is_empty()) incre_settlement.zeros();
1239  if(!current_settlement.is_empty()) current_settlement.zeros();
1240 }
1241 
1242 template<sp_d T> void Factory<T>::clear_resistance() {
1243  if(!pre_resistance.is_empty()) pre_resistance.zeros();
1244  if(!trial_resistance.is_empty()) trial_resistance.zeros();
1245  if(!incre_resistance.is_empty()) incre_resistance.zeros();
1246  if(!current_resistance.is_empty()) current_resistance.zeros();
1247 }
1248 
1249 template<sp_d T> void Factory<T>::clear_damping_force() {
1250  if(!pre_damping_force.is_empty()) pre_damping_force.zeros();
1251  if(!trial_damping_force.is_empty()) trial_damping_force.zeros();
1252  if(!incre_damping_force.is_empty()) incre_damping_force.zeros();
1253  if(!current_damping_force.is_empty()) current_damping_force.zeros();
1254 }
1255 
1256 template<sp_d T> void Factory<T>::clear_nonviscous_force() {
1257  if(!pre_nonviscous_force.is_empty()) pre_nonviscous_force.zeros();
1258  if(!trial_nonviscous_force.is_empty()) trial_nonviscous_force.zeros();
1259  if(!incre_nonviscous_force.is_empty()) incre_nonviscous_force.zeros();
1260  if(!current_nonviscous_force.is_empty()) current_nonviscous_force.zeros();
1261 }
1262 
1263 template<sp_d T> void Factory<T>::clear_inertial_force() {
1264  if(!pre_inertial_force.is_empty()) pre_inertial_force.zeros();
1265  if(!trial_inertial_force.is_empty()) trial_inertial_force.zeros();
1266  if(!incre_inertial_force.is_empty()) incre_inertial_force.zeros();
1267  if(!current_inertial_force.is_empty()) current_inertial_force.zeros();
1268 }
1269 
1270 template<sp_d T> void Factory<T>::clear_displacement() {
1271  if(!pre_displacement.is_empty()) pre_displacement.zeros();
1272  if(!trial_displacement.is_empty()) trial_displacement.zeros();
1273  if(!incre_displacement.is_empty()) incre_displacement.zeros();
1274  if(!current_displacement.is_empty()) current_displacement.zeros();
1275 }
1276 
1277 template<sp_d T> void Factory<T>::clear_velocity() {
1278  if(!pre_velocity.is_empty()) pre_velocity.zeros();
1279  if(!trial_velocity.is_empty()) trial_velocity.zeros();
1280  if(!incre_velocity.is_empty()) incre_velocity.zeros();
1281  if(!current_velocity.is_empty()) current_velocity.zeros();
1282 }
1283 
1284 template<sp_d T> void Factory<T>::clear_acceleration() {
1285  if(!pre_acceleration.is_empty()) pre_acceleration.zeros();
1286  if(!trial_acceleration.is_empty()) trial_acceleration.zeros();
1287  if(!incre_acceleration.is_empty()) incre_acceleration.zeros();
1288  if(!current_acceleration.is_empty()) current_acceleration.zeros();
1289 }
1290 
1291 template<sp_d T> void Factory<T>::clear_temperature() {
1292  if(!pre_temperature.is_empty()) pre_temperature.zeros();
1293  if(!trial_temperature.is_empty()) trial_temperature.zeros();
1294  if(!incre_temperature.is_empty()) incre_temperature.zeros();
1295  if(!current_temperature.is_empty()) current_temperature.zeros();
1296 }
1297 
1299  if(!trial_constraint_resistance.is_empty()) trial_constraint_resistance.zeros();
1300  if(!current_constraint_resistance.is_empty()) current_constraint_resistance.zeros();
1301 }
1302 
1303 template<sp_d T> void Factory<T>::reset_status() {
1304  ninja.zeros();
1305 
1306  reset_time();
1307  reset_load_factor();
1308  reset_load();
1309  reset_settlement();
1310  reset_resistance();
1311  reset_damping_force();
1312  reset_nonviscous_force();
1313  reset_inertial_force();
1314  reset_displacement();
1315  reset_velocity();
1316  reset_acceleration();
1317  reset_temperature();
1318  reset_auxiliary_resistance();
1319 }
1320 
1321 template<sp_d T> void Factory<T>::reset_time() {
1322  trial_time = current_time;
1323  incre_time = T(0);
1324 }
1325 
1326 template<sp_d T> void Factory<T>::reset_load_factor() {
1327  if(trial_load_factor.is_empty()) return;
1328  trial_load_factor = current_load_factor;
1329  incre_load_factor.zeros();
1330 }
1331 
1332 template<sp_d T> void Factory<T>::reset_load() {
1333  if(trial_load.is_empty()) return;
1334  trial_load = current_load;
1335  incre_load.zeros();
1336 }
1337 
1338 template<sp_d T> void Factory<T>::reset_settlement() {
1339  if(trial_settlement.is_empty()) return;
1340  trial_settlement = current_settlement;
1341  incre_settlement.zeros();
1342 }
1343 
1344 template<sp_d T> void Factory<T>::reset_resistance() {
1345  if(trial_resistance.is_empty()) return;
1346  trial_resistance = current_resistance;
1347  incre_resistance.zeros();
1348 }
1349 
1350 template<sp_d T> void Factory<T>::reset_damping_force() {
1351  if(trial_damping_force.is_empty()) return;
1352  trial_damping_force = current_damping_force;
1353  incre_damping_force.zeros();
1354 }
1355 
1356 template<sp_d T> void Factory<T>::reset_nonviscous_force() {
1357  if(trial_nonviscous_force.is_empty()) return;
1358  trial_nonviscous_force = current_nonviscous_force;
1359  incre_nonviscous_force.zeros();
1360 }
1361 
1362 template<sp_d T> void Factory<T>::reset_inertial_force() {
1363  if(trial_inertial_force.is_empty()) return;
1364  trial_inertial_force = current_inertial_force;
1365  incre_inertial_force.zeros();
1366 }
1367 
1368 template<sp_d T> void Factory<T>::reset_displacement() {
1369  if(trial_displacement.is_empty()) return;
1370  trial_displacement = current_displacement;
1371  incre_displacement.zeros();
1372 }
1373 
1374 template<sp_d T> void Factory<T>::reset_velocity() {
1375  if(trial_velocity.is_empty()) return;
1376  trial_velocity = current_velocity;
1377  incre_velocity.zeros();
1378 }
1379 
1380 template<sp_d T> void Factory<T>::reset_acceleration() {
1381  if(trial_acceleration.is_empty()) return;
1382  trial_acceleration = current_acceleration;
1383  incre_acceleration.zeros();
1384 }
1385 
1386 template<sp_d T> void Factory<T>::reset_temperature() {
1387  if(trial_temperature.is_empty()) return;
1388  trial_temperature = current_temperature;
1389  incre_temperature.zeros();
1390 }
1391 
1393  if(trial_constraint_resistance.is_empty()) return;
1394  trial_constraint_resistance = current_constraint_resistance;
1395 }
1396 
1397 template<sp_d T> void Factory<T>::clear_eigen() {
1398  if(!eigenvalue.is_empty()) eigenvalue.zeros();
1399  if(!eigenvector.is_empty()) eigenvector.zeros();
1400 }
1401 
1402 template<sp_d T> void Factory<T>::clear_mass() { if(global_mass != nullptr) global_mass->zeros(); }
1403 
1404 template<sp_d T> void Factory<T>::clear_damping() { if(global_damping != nullptr) global_damping->zeros(); }
1405 
1406 template<sp_d T> void Factory<T>::clear_nonviscous() { if(global_nonviscous != nullptr) global_nonviscous->zeros(); }
1407 
1408 template<sp_d T> void Factory<T>::clear_stiffness() { if(global_stiffness != nullptr) global_stiffness->zeros(); }
1409 
1410 template<sp_d T> void Factory<T>::clear_geometry() { if(global_geometry != nullptr) global_geometry->zeros(); }
1411 
1412 template<sp_d T> void Factory<T>::clear_auxiliary() {
1413  n_mpc = 0;
1414  auxiliary_load.reset();
1415  auxiliary_stiffness.set_size(n_size, 0);
1416  auxiliary_resistance.reset();
1417  auxiliary_encoding.reset();
1418 }
1419 
1420 template<sp_d T> void Factory<T>::reset() {
1421  global_mass = nullptr;
1422  global_damping = nullptr;
1423  global_nonviscous = nullptr;
1424  global_stiffness = nullptr;
1425  global_geometry = nullptr;
1426 }
1427 
1428 template<sp_d T> void Factory<T>::assemble_resistance(const Mat<T>& ER, const uvec& EI) {
1429  if(ER.is_empty()) return;
1430  for(auto I = 0llu; I < EI.n_elem; ++I) trial_resistance(EI(I)) += ER(I);
1431 }
1432 
1433 template<sp_d T> void Factory<T>::assemble_damping_force(const Mat<T>& ER, const uvec& EI) {
1434  if(ER.is_empty()) return;
1435  for(auto I = 0llu; I < EI.n_elem; ++I) trial_damping_force(EI(I)) += ER(I);
1436 }
1437 
1438 template<sp_d T> void Factory<T>::assemble_nonviscous_force(const Mat<T>& ER, const uvec& EI) {
1439  if(ER.is_empty()) return;
1440  for(auto I = 0llu; I < EI.n_elem; ++I) trial_nonviscous_force(EI(I)) += ER(I);
1441 }
1442 
1443 template<sp_d T> void Factory<T>::assemble_inertial_force(const Mat<T>& ER, const uvec& EI) {
1444  if(ER.is_empty()) return;
1445  for(auto I = 0llu; I < EI.n_elem; ++I) trial_inertial_force(EI(I)) += ER(I);
1446 }
1447 
1448 template<sp_d T> void Factory<T>::assemble_matrix_helper(shared_ptr<MetaMat<T>>& GM, const Mat<T>& EM, const uvec& EI, const std::vector<MappingDOF>& MAP) {
1449  if(EM.is_empty()) return;
1450 
1451  if(StorageScheme::BANDSYMM == storage_type || StorageScheme::SYMMPACK == storage_type) for(const auto [g_row, g_col, l_row, l_col] : MAP) GM->unsafe_at(g_row, g_col) += EM(l_row, l_col);
1452  else for(auto I = 0llu; I < EI.n_elem; ++I) for(auto J = 0llu; J < EI.n_elem; ++J) GM->unsafe_at(EI(J), EI(I)) += EM(J, I);
1453 }
1454 
1455 template<sp_d T> void Factory<T>::assemble_mass(const Mat<T>& EM, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_mass, EM, EI, MAP); }
1456 
1457 template<sp_d T> void Factory<T>::assemble_damping(const Mat<T>& EC, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_damping, EC, EI, MAP); }
1458 
1459 template<sp_d T> void Factory<T>::assemble_nonviscous(const Mat<T>& EC, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_nonviscous, EC, EI, MAP); }
1460 
1461 template<sp_d T> void Factory<T>::assemble_stiffness(const Mat<T>& EK, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_stiffness, EK, EI, MAP); }
1462 
1463 template<sp_d T> void Factory<T>::assemble_geometry(const Mat<T>& EG, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_geometry, EG, EI, MAP); }
1464 
1465 template<sp_d T> void Factory<T>::assemble_stiffness(const SpMat<T>& EK, const uvec& EI) {
1466  if(EK.is_empty()) return;
1467  for(auto I = EK.begin(); I != EK.end(); ++I) global_stiffness->at(EI(I.row()), EI(I.col())) += *I;
1468 }
1469 
1470 template<sp_d T> void Factory<T>::print() const {
1471  suanpan_info("A Factory object with size of {}.\n", n_size);
1472 }
1473 
1474 template<sp_d T> unique_ptr<MetaMat<T>> Factory<T>::get_basic_container() {
1475  switch(storage_type) {
1476  case StorageScheme::FULL:
1477 #ifdef SUANPAN_CUDA
1478  if(SolverType::CUDA == solver) return std::make_unique<FullMatCUDA<T>>(n_size, n_size);
1479 #endif
1480  return std::make_unique<FullMat<T>>(n_size, n_size);
1481  case StorageScheme::BAND:
1482  if(SolverType::SPIKE == solver) return std::make_unique<BandMatSpike<T>>(n_size, n_lobw, n_upbw);
1483  return std::make_unique<BandMat<T>>(n_size, n_lobw, n_upbw);
1485  return std::make_unique<BandSymmMat<T>>(n_size, n_lobw);
1487  return std::make_unique<SymmPackMat<T>>(n_size);
1488  case StorageScheme::SPARSE:
1489  if(SolverType::MUMPS == solver) return std::make_unique<SparseMatMUMPS<T>>(n_size, n_size, n_elem);
1490  if(SolverType::LIS == solver) return std::make_unique<SparseMatLis<T>>(n_size, n_size, n_elem);
1491  if(SolverType::SUPERLU == solver) return std::make_unique<SparseMatSuperLU<T>>(n_size, n_size, n_elem);
1492 #ifdef SUANPAN_MKL
1493  if(SolverType::PARDISO == solver) return std::make_unique<SparseMatPARDISO<T>>(n_size, n_size, n_elem);
1494  if(SolverType::FGMRES == solver) return std::make_unique<SparseMatFGMRES<T>>(n_size, n_size, n_elem);
1495 #endif
1496 #ifdef SUANPAN_CUDA
1497  if(SolverType::CUDA == solver) return std::make_unique<SparseMatCUDA<T>>(n_size, n_size, n_elem);
1498 #ifdef SUANPAN_MAGMA
1499  if(SolverType::MAGMA == solver) return std::make_unique<SparseMatMAGMA<T>>(n_size, n_size, magma_setting);
1500 #endif
1501 #endif
1502  return std::make_unique<SparseMatSuperLU<T>>(n_size, n_size, n_elem);
1504 #ifdef SUANPAN_MKL
1505  if(SolverType::FGMRES == solver) return std::make_unique<SparseSymmMatFGMRES<T>>(n_size, n_size, n_elem);
1506 #endif
1507  return std::make_unique<SparseSymmMatMUMPS<T>>(n_size, n_size, n_elem);
1508  default:
1509  throw invalid_argument("need a proper storage scheme");
1510  }
1511 }
1512 
1513 template<sp_d T> unique_ptr<MetaMat<T>> Factory<T>::get_matrix_container() {
1514  auto global_mat = get_basic_container();
1515 
1516  global_mat->set_solver_setting(setting);
1517 
1518  return global_mat;
1519 }
1520 
1521 template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_mass() { return global_mass; }
1522 
1523 template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_damping() { return global_damping; }
1524 
1525 template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_nonviscous() { return global_nonviscous; }
1526 
1527 template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_stiffness() { return global_stiffness; }
1528 
1529 template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_geometry() { return global_geometry; }
1530 
1531 template<sp_d T> Col<T>& Factory<T>::modify_ninja() { return ninja; }
1532 
1533 template<sp_d T> Col<T>& Factory<T>::modify_sushi() { return sushi; }
1534 
1535 template<sp_d T> suanpan::set<uword>& Factory<T>::modify_reference_dof() { return reference_dof; }
1536 
1537 template<sp_d T> SpMat<T>& Factory<T>::modify_reference_load() { return reference_load; }
1538 
1539 template<sp_d T> uvec& Factory<T>::modify_auxiliary_encoding() { return auxiliary_encoding; }
1540 
1541 template<sp_d T> Col<T>& Factory<T>::modify_auxiliary_lambda() { return auxiliary_lambda; }
1542 
1543 template<sp_d T> Col<T>& Factory<T>::modify_auxiliary_resistance() { return auxiliary_resistance; }
1544 
1545 template<sp_d T> Col<T>& Factory<T>::modify_auxiliary_load() { return auxiliary_load; }
1546 
1547 template<sp_d T> SpMat<T>& Factory<T>::modify_auxiliary_stiffness() { return auxiliary_stiffness; }
1548 
1549 template<sp_d T> SpCol<T>& Factory<T>::modify_trial_constraint_resistance() { return trial_constraint_resistance; }
1550 
1551 template<sp_d T> SpCol<T>& Factory<T>::modify_current_constraint_resistance() { return current_constraint_resistance; }
1552 
1553 template<sp_d T> Col<T>& Factory<T>::modify_eigenvalue() { return eigenvalue; }
1554 
1555 template<sp_d T> Mat<T>& Factory<T>::modify_eigenvector() { return eigenvector; }
1556 
1557 template<sp_d T> void Factory<T>::set_trial_time(const T M) { trial_time = M; }
1558 
1559 template<sp_d T> void Factory<T>::set_trial_load_factor(const Col<T>& L) { trial_load_factor = L; }
1560 
1561 template<sp_d T> void Factory<T>::set_trial_load(const Col<T>& L) { trial_load = L; }
1562 
1563 template<sp_d T> void Factory<T>::set_trial_settlement(const Col<T>& S) { trial_settlement = S; }
1564 
1565 template<sp_d T> void Factory<T>::set_trial_resistance(const Col<T>& R) { trial_resistance = R; }
1566 
1567 template<sp_d T> void Factory<T>::set_trial_damping_force(const Col<T>& R) { trial_damping_force = R; }
1568 
1569 template<sp_d T> void Factory<T>::set_trial_nonviscous_force(const Col<T>& R) { trial_nonviscous_force = R; }
1570 
1571 template<sp_d T> void Factory<T>::set_trial_inertial_force(const Col<T>& R) { trial_inertial_force = R; }
1572 
1573 template<sp_d T> void Factory<T>::set_trial_displacement(const Col<T>& D) { trial_displacement = D; }
1574 
1575 template<sp_d T> void Factory<T>::set_trial_velocity(const Col<T>& V) { trial_velocity = V; }
1576 
1577 template<sp_d T> void Factory<T>::set_trial_acceleration(const Col<T>& A) { trial_acceleration = A; }
1578 
1579 template<sp_d T> void Factory<T>::set_trial_temperature(const Col<T>& M) { trial_temperature = M; }
1580 
1581 template<sp_d T> void Factory<T>::set_incre_time(const T M) { incre_time = M; }
1582 
1583 template<sp_d T> void Factory<T>::set_incre_load_factor(const Col<T>& L) { incre_load_factor = L; }
1584 
1585 template<sp_d T> void Factory<T>::set_incre_load(const Col<T>& L) { incre_load = L; }
1586 
1587 template<sp_d T> void Factory<T>::set_incre_settlement(const Col<T>& S) { incre_settlement = S; }
1588 
1589 template<sp_d T> void Factory<T>::set_incre_resistance(const Col<T>& R) { incre_resistance = R; }
1590 
1591 template<sp_d T> void Factory<T>::set_incre_damping_force(const Col<T>& R) { incre_damping_force = R; }
1592 
1593 template<sp_d T> void Factory<T>::set_incre_nonviscous_force(const Col<T>& R) { incre_nonviscous_force = R; }
1594 
1595 template<sp_d T> void Factory<T>::set_incre_inertial_force(const Col<T>& R) { incre_inertial_force = R; }
1596 
1597 template<sp_d T> void Factory<T>::set_incre_displacement(const Col<T>& D) { incre_displacement = D; }
1598 
1599 template<sp_d T> void Factory<T>::set_incre_velocity(const Col<T>& V) { incre_velocity = V; }
1600 
1601 template<sp_d T> void Factory<T>::set_incre_acceleration(const Col<T>& A) { incre_acceleration = A; }
1602 
1603 template<sp_d T> void Factory<T>::set_incre_temperature(const Col<T>& M) { incre_temperature = M; }
1604 
1605 template<sp_d T> void Factory<T>::set_current_time(const T M) { current_time = M; }
1606 
1607 template<sp_d T> void Factory<T>::set_current_load_factor(const Col<T>& L) { current_load_factor = L; }
1608 
1609 template<sp_d T> void Factory<T>::set_current_load(const Col<T>& L) { current_load = L; }
1610 
1611 template<sp_d T> void Factory<T>::set_current_settlement(const Col<T>& S) { current_settlement = S; }
1612 
1613 template<sp_d T> void Factory<T>::set_current_resistance(const Col<T>& R) { current_resistance = R; }
1614 
1615 template<sp_d T> void Factory<T>::set_current_damping_force(const Col<T>& R) { current_damping_force = R; }
1616 
1617 template<sp_d T> void Factory<T>::set_current_nonviscous_force(const Col<T>& R) { current_nonviscous_force = R; }
1618 
1619 template<sp_d T> void Factory<T>::set_current_inertial_force(const Col<T>& R) { current_inertial_force = R; }
1620 
1621 template<sp_d T> void Factory<T>::set_current_displacement(const Col<T>& D) { current_displacement = D; }
1622 
1623 template<sp_d T> void Factory<T>::set_current_velocity(const Col<T>& V) { current_velocity = V; }
1624 
1625 template<sp_d T> void Factory<T>::set_current_acceleration(const Col<T>& A) { current_acceleration = A; }
1626 
1627 template<sp_d T> void Factory<T>::set_current_temperature(const Col<T>& M) { current_temperature = M; }
1628 
1629 template<sp_d T> void Factory<T>::set_pre_time(const T M) { pre_time = M; }
1630 
1631 template<sp_d T> void Factory<T>::set_pre_load_factor(const Col<T>& L) { pre_load_factor = L; }
1632 
1633 template<sp_d T> void Factory<T>::set_pre_load(const Col<T>& L) { pre_load = L; }
1634 
1635 template<sp_d T> void Factory<T>::set_pre_settlement(const Col<T>& S) { pre_settlement = S; }
1636 
1637 template<sp_d T> void Factory<T>::set_pre_resistance(const Col<T>& R) { pre_resistance = R; }
1638 
1639 template<sp_d T> void Factory<T>::set_pre_damping_force(const Col<T>& R) { pre_damping_force = R; }
1640 
1641 template<sp_d T> void Factory<T>::set_pre_nonviscous_force(const Col<T>& R) { pre_nonviscous_force = R; }
1642 
1643 template<sp_d T> void Factory<T>::set_pre_inertial_force(const Col<T>& R) { pre_inertial_force = R; }
1644 
1645 template<sp_d T> void Factory<T>::set_pre_displacement(const Col<T>& D) { pre_displacement = D; }
1646 
1647 template<sp_d T> void Factory<T>::set_pre_velocity(const Col<T>& V) { pre_velocity = V; }
1648 
1649 template<sp_d T> void Factory<T>::set_pre_acceleration(const Col<T>& A) { pre_acceleration = A; }
1650 
1651 template<sp_d T> void Factory<T>::set_pre_temperature(const Col<T>& M) { pre_temperature = M; }
1652 
1653 template<sp_d T> T Factory<T>::get_trial_time() const { return trial_time; }
1654 
1655 template<sp_d T> const Col<T>& Factory<T>::get_trial_load_factor() const { return trial_load_factor; }
1656 
1657 template<sp_d T> const Col<T>& Factory<T>::get_trial_load() const { return trial_load; }
1658 
1659 template<sp_d T> const Col<T>& Factory<T>::get_trial_settlement() const { return trial_settlement; }
1660 
1661 template<sp_d T> const Col<T>& Factory<T>::get_trial_resistance() const { return trial_resistance; }
1662 
1663 template<sp_d T> const Col<T>& Factory<T>::get_trial_damping_force() const { return trial_damping_force; }
1664 
1665 template<sp_d T> const Col<T>& Factory<T>::get_trial_nonviscous_force() const { return trial_nonviscous_force; }
1666 
1667 template<sp_d T> const Col<T>& Factory<T>::get_trial_inertial_force() const { return trial_inertial_force; }
1668 
1669 template<sp_d T> const Col<T>& Factory<T>::get_trial_displacement() const { return trial_displacement; }
1670 
1671 template<sp_d T> const Col<T>& Factory<T>::get_trial_velocity() const { return trial_velocity; }
1672 
1673 template<sp_d T> const Col<T>& Factory<T>::get_trial_acceleration() const { return trial_acceleration; }
1674 
1675 template<sp_d T> const Col<T>& Factory<T>::get_trial_temperature() const { return trial_temperature; }
1676 
1677 template<sp_d T> T Factory<T>::get_incre_time() const { return incre_time; }
1678 
1679 template<sp_d T> const Col<T>& Factory<T>::get_incre_load_factor() const { return incre_load_factor; }
1680 
1681 template<sp_d T> const Col<T>& Factory<T>::get_incre_load() const { return incre_load; }
1682 
1683 template<sp_d T> const Col<T>& Factory<T>::get_incre_settlement() const { return incre_settlement; }
1684 
1685 template<sp_d T> const Col<T>& Factory<T>::get_incre_resistance() const { return incre_resistance; }
1686 
1687 template<sp_d T> const Col<T>& Factory<T>::get_incre_damping_force() const { return incre_damping_force; }
1688 
1689 template<sp_d T> const Col<T>& Factory<T>::get_incre_nonviscous_force() const { return incre_nonviscous_force; }
1690 
1691 template<sp_d T> const Col<T>& Factory<T>::get_incre_inertial_force() const { return incre_inertial_force; }
1692 
1693 template<sp_d T> const Col<T>& Factory<T>::get_incre_displacement() const { return incre_displacement; }
1694 
1695 template<sp_d T> const Col<T>& Factory<T>::get_incre_velocity() const { return incre_velocity; }
1696 
1697 template<sp_d T> const Col<T>& Factory<T>::get_incre_acceleration() const { return incre_acceleration; }
1698 
1699 template<sp_d T> const Col<T>& Factory<T>::get_incre_temperature() const { return incre_temperature; }
1700 
1701 template<sp_d T> T Factory<T>::get_current_time() const { return current_time; }
1702 
1703 template<sp_d T> const Col<T>& Factory<T>::get_current_load_factor() const { return current_load_factor; }
1704 
1705 template<sp_d T> const Col<T>& Factory<T>::get_current_load() const { return current_load; }
1706 
1707 template<sp_d T> const Col<T>& Factory<T>::get_current_settlement() const { return current_settlement; }
1708 
1709 template<sp_d T> const Col<T>& Factory<T>::get_current_resistance() const { return current_resistance; }
1710 
1711 template<sp_d T> const Col<T>& Factory<T>::get_current_damping_force() const { return current_damping_force; }
1712 
1713 template<sp_d T> const Col<T>& Factory<T>::get_current_nonviscous_force() const { return current_nonviscous_force; }
1714 
1715 template<sp_d T> const Col<T>& Factory<T>::get_current_inertial_force() const { return current_inertial_force; }
1716 
1717 template<sp_d T> const Col<T>& Factory<T>::get_current_displacement() const { return current_displacement; }
1718 
1719 template<sp_d T> const Col<T>& Factory<T>::get_current_velocity() const { return current_velocity; }
1720 
1721 template<sp_d T> const Col<T>& Factory<T>::get_current_acceleration() const { return current_acceleration; }
1722 
1723 template<sp_d T> const Col<T>& Factory<T>::get_current_temperature() const { return current_temperature; }
1724 
1725 template<sp_d T> T Factory<T>::get_pre_time() const { return pre_time; }
1726 
1727 template<sp_d T> const Col<T>& Factory<T>::get_pre_load_factor() const { return pre_load_factor; }
1728 
1729 template<sp_d T> const Col<T>& Factory<T>::get_pre_load() const { return pre_load; }
1730 
1731 template<sp_d T> const Col<T>& Factory<T>::get_pre_settlement() const { return pre_settlement; }
1732 
1733 template<sp_d T> const Col<T>& Factory<T>::get_pre_resistance() const { return pre_resistance; }
1734 
1735 template<sp_d T> const Col<T>& Factory<T>::get_pre_damping_force() const { return pre_damping_force; }
1736 
1737 template<sp_d T> const Col<T>& Factory<T>::get_pre_nonviscous_force() const { return pre_nonviscous_force; }
1738 
1739 template<sp_d T> const Col<T>& Factory<T>::get_pre_inertial_force() const { return pre_inertial_force; }
1740 
1741 template<sp_d T> const Col<T>& Factory<T>::get_pre_displacement() const { return pre_displacement; }
1742 
1743 template<sp_d T> const Col<T>& Factory<T>::get_pre_velocity() const { return pre_velocity; }
1744 
1745 template<sp_d T> const Col<T>& Factory<T>::get_pre_acceleration() const { return pre_acceleration; }
1746 
1747 template<sp_d T> const Col<T>& Factory<T>::get_pre_temperature() const { return pre_temperature; }
1748 
1749 template<sp_d T> T& Factory<T>::modify_trial_time() { return trial_time; }
1750 
1751 template<sp_d T> Col<T>& Factory<T>::modify_trial_load_factor() { return trial_load_factor; }
1752 
1753 template<sp_d T> Col<T>& Factory<T>::modify_trial_load() { return trial_load; }
1754 
1755 template<sp_d T> Col<T>& Factory<T>::modify_trial_settlement() { return trial_settlement; }
1756 
1757 template<sp_d T> Col<T>& Factory<T>::modify_trial_resistance() { return trial_resistance; }
1758 
1759 template<sp_d T> Col<T>& Factory<T>::modify_trial_damping_force() { return trial_damping_force; }
1760 
1761 template<sp_d T> Col<T>& Factory<T>::modify_trial_nonviscous_force() { return trial_nonviscous_force; }
1762 
1763 template<sp_d T> Col<T>& Factory<T>::modify_trial_inertial_force() { return trial_inertial_force; }
1764 
1765 template<sp_d T> Col<T>& Factory<T>::modify_trial_displacement() { return trial_displacement; }
1766 
1767 template<sp_d T> Col<T>& Factory<T>::modify_trial_velocity() { return trial_velocity; }
1768 
1769 template<sp_d T> Col<T>& Factory<T>::modify_trial_acceleration() { return trial_acceleration; }
1770 
1771 template<sp_d T> Col<T>& Factory<T>::modify_trial_temperature() { return trial_temperature; }
1772 
1773 template<sp_d T> T& Factory<T>::modify_incre_time() { return incre_time; }
1774 
1775 template<sp_d T> Col<T>& Factory<T>::modify_incre_load_factor() { return incre_load_factor; }
1776 
1777 template<sp_d T> Col<T>& Factory<T>::modify_incre_load() { return incre_load; }
1778 
1779 template<sp_d T> Col<T>& Factory<T>::modify_incre_settlement() { return incre_settlement; }
1780 
1781 template<sp_d T> Col<T>& Factory<T>::modify_incre_resistance() { return incre_resistance; }
1782 
1783 template<sp_d T> Col<T>& Factory<T>::modify_incre_damping_force() { return incre_damping_force; }
1784 
1785 template<sp_d T> Col<T>& Factory<T>::modify_incre_nonviscous_force() { return incre_nonviscous_force; }
1786 
1787 template<sp_d T> Col<T>& Factory<T>::modify_incre_inertial_force() { return incre_inertial_force; }
1788 
1789 template<sp_d T> Col<T>& Factory<T>::modify_incre_displacement() { return incre_displacement; }
1790 
1791 template<sp_d T> Col<T>& Factory<T>::modify_incre_velocity() { return incre_velocity; }
1792 
1793 template<sp_d T> Col<T>& Factory<T>::modify_incre_acceleration() { return incre_acceleration; }
1794 
1795 template<sp_d T> Col<T>& Factory<T>::modify_incre_temperature() { return incre_temperature; }
1796 
1797 template<sp_d T> T& Factory<T>::modify_current_time() { return current_time; }
1798 
1799 template<sp_d T> Col<T>& Factory<T>::modify_current_load_factor() { return current_load_factor; }
1800 
1801 template<sp_d T> Col<T>& Factory<T>::modify_current_load() { return current_load; }
1802 
1803 template<sp_d T> Col<T>& Factory<T>::modify_current_settlement() { return current_settlement; }
1804 
1805 template<sp_d T> Col<T>& Factory<T>::modify_current_resistance() { return current_resistance; }
1806 
1807 template<sp_d T> Col<T>& Factory<T>::modify_current_damping_force() { return current_damping_force; }
1808 
1809 template<sp_d T> Col<T>& Factory<T>::modify_current_nonviscous_force() { return current_nonviscous_force; }
1810 
1811 template<sp_d T> Col<T>& Factory<T>::modify_current_inertial_force() { return current_inertial_force; }
1812 
1813 template<sp_d T> Col<T>& Factory<T>::modify_current_displacement() { return current_displacement; }
1814 
1815 template<sp_d T> Col<T>& Factory<T>::modify_current_velocity() { return current_velocity; }
1816 
1817 template<sp_d T> Col<T>& Factory<T>::modify_current_acceleration() { return current_acceleration; }
1818 
1819 template<sp_d T> Col<T>& Factory<T>::modify_current_temperature() { return current_temperature; }
1820 
1821 template<sp_d T> T& Factory<T>::modify_pre_time() { return pre_time; }
1822 
1823 template<sp_d T> Col<T>& Factory<T>::modify_pre_load_factor() { return pre_load_factor; }
1824 
1825 template<sp_d T> Col<T>& Factory<T>::modify_pre_load() { return pre_load; }
1826 
1827 template<sp_d T> Col<T>& Factory<T>::modify_pre_settlement() { return pre_settlement; }
1828 
1829 template<sp_d T> Col<T>& Factory<T>::modify_pre_resistance() { return pre_resistance; }
1830 
1831 template<sp_d T> Col<T>& Factory<T>::modify_pre_damping_force() { return pre_damping_force; }
1832 
1833 template<sp_d T> Col<T>& Factory<T>::modify_pre_nonviscous_force() { return pre_nonviscous_force; }
1834 
1835 template<sp_d T> Col<T>& Factory<T>::modify_pre_inertial_force() { return pre_inertial_force; }
1836 
1837 template<sp_d T> Col<T>& Factory<T>::modify_pre_displacement() { return pre_displacement; }
1838 
1839 template<sp_d T> Col<T>& Factory<T>::modify_pre_velocity() { return pre_velocity; }
1840 
1841 template<sp_d T> Col<T>& Factory<T>::modify_pre_acceleration() { return pre_acceleration; }
1842 
1843 template<sp_d T> Col<T>& Factory<T>::modify_pre_temperature() { return pre_temperature; }
1844 
1845 template<sp_d T> void Factory<T>::update_trial_time(const T M) {
1846  trial_time = M;
1847  incre_time = trial_time - current_time;
1848 }
1849 
1850 template<sp_d T> void Factory<T>::update_trial_load_factor(const Col<T>& L) {
1851  trial_load_factor = L;
1852  incre_load_factor = trial_load_factor - current_load_factor;
1853 }
1854 
1855 template<sp_d T> void Factory<T>::update_trial_load(const Col<T>& L) {
1856  trial_load = L;
1857  incre_load = trial_load - current_load;
1858 }
1859 
1860 template<sp_d T> void Factory<T>::update_trial_settlement(const Col<T>& S) {
1861  trial_settlement = S;
1862  incre_settlement = trial_settlement - current_settlement;
1863 }
1864 
1865 template<sp_d T> void Factory<T>::update_trial_resistance(const Col<T>& R) {
1866  trial_resistance = R;
1867  incre_resistance = trial_resistance - current_resistance;
1868 }
1869 
1870 template<sp_d T> void Factory<T>::update_trial_damping_force(const Col<T>& R) {
1871  trial_damping_force = R;
1872  incre_damping_force = trial_damping_force - current_damping_force;
1873 }
1874 
1875 template<sp_d T> void Factory<T>::update_trial_nonviscous_force(const Col<T>& R) {
1876  trial_nonviscous_force = R;
1877  incre_nonviscous_force = trial_nonviscous_force - current_nonviscous_force;
1878 }
1879 
1880 template<sp_d T> void Factory<T>::update_trial_inertial_force(const Col<T>& R) {
1881  trial_inertial_force = R;
1882  incre_inertial_force = trial_inertial_force - current_inertial_force;
1883 }
1884 
1885 template<sp_d T> void Factory<T>::update_trial_displacement(const Col<T>& D) {
1886  trial_displacement = D;
1887  incre_displacement = trial_displacement - current_displacement;
1888 }
1889 
1890 template<sp_d T> void Factory<T>::update_trial_velocity(const Col<T>& V) {
1891  trial_velocity = V;
1892  incre_velocity = trial_velocity - current_velocity;
1893 }
1894 
1895 template<sp_d T> void Factory<T>::update_trial_acceleration(const Col<T>& A) {
1896  trial_acceleration = A;
1897  incre_acceleration = trial_acceleration - current_acceleration;
1898 }
1899 
1900 template<sp_d T> void Factory<T>::update_trial_temperature(const Col<T>& M) {
1901  trial_temperature = M;
1902  incre_temperature = trial_temperature - current_temperature;
1903 }
1904 
1905 template<sp_d T> void Factory<T>::update_incre_time(const T M) {
1906  incre_time = M;
1907  trial_time = current_time + incre_time;
1908 }
1909 
1910 template<sp_d T> void Factory<T>::update_incre_load_factor(const Col<T>& L) {
1911  incre_load_factor = L;
1912  trial_load_factor = current_load_factor + incre_load_factor;
1913 }
1914 
1915 template<sp_d T> void Factory<T>::update_incre_load(const Col<T>& L) {
1916  incre_load = L;
1917  trial_load = current_load + incre_load;
1918 }
1919 
1920 template<sp_d T> void Factory<T>::update_incre_settlement(const Col<T>& S) {
1921  incre_settlement = S;
1922  trial_settlement = current_settlement + incre_settlement;
1923 }
1924 
1925 template<sp_d T> void Factory<T>::update_incre_resistance(const Col<T>& R) {
1926  incre_resistance = R;
1927  trial_resistance = current_resistance + incre_resistance;
1928 }
1929 
1930 template<sp_d T> void Factory<T>::update_incre_damping_force(const Col<T>& R) {
1931  incre_damping_force = R;
1932  trial_damping_force = current_damping_force + incre_damping_force;
1933 }
1934 
1935 template<sp_d T> void Factory<T>::update_incre_nonviscous_force(const Col<T>& R) {
1936  incre_nonviscous_force = R;
1937  trial_nonviscous_force = current_nonviscous_force + incre_nonviscous_force;
1938 }
1939 
1940 template<sp_d T> void Factory<T>::update_incre_inertial_force(const Col<T>& R) {
1941  incre_inertial_force = R;
1942  trial_inertial_force = current_inertial_force + incre_inertial_force;
1943 }
1944 
1945 template<sp_d T> void Factory<T>::update_incre_displacement(const Col<T>& D) {
1946  incre_displacement = D;
1947  trial_displacement = current_displacement + incre_displacement;
1948 }
1949 
1950 template<sp_d T> void Factory<T>::update_incre_velocity(const Col<T>& V) {
1951  incre_velocity = V;
1952  trial_velocity = current_velocity + incre_velocity;
1953 }
1954 
1955 template<sp_d T> void Factory<T>::update_incre_acceleration(const Col<T>& A) {
1956  incre_acceleration = A;
1957  trial_acceleration = current_acceleration + incre_acceleration;
1958 }
1959 
1960 template<sp_d T> void Factory<T>::update_incre_temperature(const Col<T>& M) {
1961  incre_temperature = M;
1962  trial_temperature = current_temperature + incre_temperature;
1963 }
1964 
1965 template<sp_d T> void Factory<T>::update_current_time(const T M) {
1966  trial_time = current_time = M;
1967  incre_time = T(0);
1968 }
1969 
1970 template<sp_d T> void Factory<T>::update_current_load_factor(const Col<T>& L) {
1971  trial_load_factor = current_load_factor = L;
1972  incre_load_factor.zeros();
1973 }
1974 
1975 template<sp_d T> void Factory<T>::update_current_load(const Col<T>& L) {
1976  trial_load = current_load = L;
1977  incre_load.zeros();
1978 }
1979 
1980 template<sp_d T> void Factory<T>::update_current_settlement(const Col<T>& S) {
1981  trial_settlement = current_settlement = S;
1982  incre_settlement.zeros();
1983 }
1984 
1985 template<sp_d T> void Factory<T>::update_current_resistance(const Col<T>& R) {
1986  trial_resistance = current_resistance = R;
1987  incre_resistance.zeros();
1988 }
1989 
1990 template<sp_d T> void Factory<T>::update_current_damping_force(const Col<T>& R) {
1991  trial_damping_force = current_damping_force = R;
1992  incre_damping_force.zeros();
1993 }
1994 
1995 template<sp_d T> void Factory<T>::update_current_nonviscous_force(const Col<T>& R) {
1996  trial_nonviscous_force = current_nonviscous_force = R;
1997  incre_nonviscous_force.zeros();
1998 }
1999 
2000 template<sp_d T> void Factory<T>::update_current_inertial_force(const Col<T>& R) {
2001  trial_inertial_force = current_inertial_force = R;
2002  incre_inertial_force.zeros();
2003 }
2004 
2005 template<sp_d T> void Factory<T>::update_current_displacement(const Col<T>& D) {
2006  trial_displacement = current_displacement = D;
2007  incre_displacement.zeros();
2008 }
2009 
2010 template<sp_d T> void Factory<T>::update_current_velocity(const Col<T>& V) {
2011  trial_velocity = current_velocity = V;
2012  incre_velocity.zeros();
2013 }
2014 
2015 template<sp_d T> void Factory<T>::update_current_acceleration(const Col<T>& A) {
2016  trial_acceleration = current_acceleration = A;
2017  incre_acceleration.zeros();
2018 }
2019 
2020 template<sp_d T> void Factory<T>::update_current_temperature(const Col<T>& M) {
2021  trial_temperature = current_temperature = M;
2022  incre_temperature.zeros();
2023 }
2024 
2025 template<sp_d T> void Factory<T>::update_trial_time_by(const T M) {
2026  trial_time += M;
2027  incre_time = trial_time - current_time;
2028 }
2029 
2030 template<sp_d T> void Factory<T>::update_trial_load_factor_by(const Col<T>& L) {
2031  trial_load_factor += L;
2032  incre_load_factor = trial_load_factor - current_load_factor;
2033 }
2034 
2035 template<sp_d T> void Factory<T>::update_trial_load_by(const Col<T>& L) {
2036  trial_load += L;
2037  incre_load = trial_load - current_load;
2038 }
2039 
2040 template<sp_d T> void Factory<T>::update_trial_settlement_by(const Col<T>& S) {
2041  trial_settlement += S;
2042  incre_settlement = trial_settlement - current_settlement;
2043 }
2044 
2045 template<sp_d T> void Factory<T>::update_trial_resistance_by(const Col<T>& R) {
2046  trial_resistance += R;
2047  incre_resistance = trial_resistance - current_resistance;
2048 }
2049 
2050 template<sp_d T> void Factory<T>::update_trial_damping_force_by(const Col<T>& R) {
2051  trial_damping_force += R;
2052  incre_damping_force = trial_damping_force - current_damping_force;
2053 }
2054 
2055 template<sp_d T> void Factory<T>::update_trial_nonviscous_force_by(const Col<T>& R) {
2056  trial_nonviscous_force += R;
2057  incre_nonviscous_force = trial_nonviscous_force - current_nonviscous_force;
2058 }
2059 
2060 template<sp_d T> void Factory<T>::update_trial_inertial_force_by(const Col<T>& R) {
2061  trial_inertial_force += R;
2062  incre_inertial_force = trial_inertial_force - current_inertial_force;
2063 }
2064 
2065 template<sp_d T> void Factory<T>::update_trial_displacement_by(const Col<T>& D) {
2066  trial_displacement += D;
2067  incre_displacement = trial_displacement - current_displacement;
2068 }
2069 
2070 template<sp_d T> void Factory<T>::update_trial_velocity_by(const Col<T>& V) {
2071  trial_velocity += V;
2072  incre_velocity = trial_velocity - current_velocity;
2073 }
2074 
2075 template<sp_d T> void Factory<T>::update_trial_acceleration_by(const Col<T>& A) {
2076  trial_acceleration += A;
2077  incre_acceleration = trial_acceleration - current_acceleration;
2078 }
2079 
2080 template<sp_d T> void Factory<T>::update_trial_temperature_by(const Col<T>& M) {
2081  trial_temperature += M;
2082  incre_temperature = trial_temperature - current_temperature;
2083 }
2084 
2085 template<sp_d T> void Factory<T>::update_incre_time_by(const T M) {
2086  incre_time += M;
2087  trial_time = current_time + incre_time;
2088 }
2089 
2090 template<sp_d T> void Factory<T>::update_incre_load_factor_by(const Col<T>& L) {
2091  incre_load_factor += L;
2092  trial_load_factor = current_load_factor + incre_load_factor;
2093 }
2094 
2095 template<sp_d T> void Factory<T>::update_incre_load_by(const Col<T>& L) {
2096  incre_load += L;
2097  trial_load = current_load + incre_load;
2098 }
2099 
2100 template<sp_d T> void Factory<T>::update_incre_settlement_by(const Col<T>& S) {
2101  incre_settlement += S;
2102  trial_settlement = current_settlement + incre_settlement;
2103 }
2104 
2105 template<sp_d T> void Factory<T>::update_incre_resistance_by(const Col<T>& R) {
2106  incre_resistance += R;
2107  trial_resistance = current_resistance + incre_resistance;
2108 }
2109 
2110 template<sp_d T> void Factory<T>::update_incre_damping_force_by(const Col<T>& R) {
2111  incre_damping_force += R;
2112  trial_damping_force = current_damping_force + incre_damping_force;
2113 }
2114 
2115 template<sp_d T> void Factory<T>::update_incre_nonviscous_force_by(const Col<T>& R) {
2116  incre_nonviscous_force += R;
2117  trial_nonviscous_force = current_nonviscous_force + incre_nonviscous_force;
2118 }
2119 
2120 template<sp_d T> void Factory<T>::update_incre_inertial_force_by(const Col<T>& R) {
2121  incre_inertial_force += R;
2122  trial_inertial_force = current_inertial_force + incre_inertial_force;
2123 }
2124 
2125 template<sp_d T> void Factory<T>::update_incre_displacement_by(const Col<T>& D) {
2126  incre_displacement += D;
2127  trial_displacement = current_displacement + incre_displacement;
2128 }
2129 
2130 template<sp_d T> void Factory<T>::update_incre_velocity_by(const Col<T>& V) {
2131  incre_velocity += V;
2132  trial_velocity = current_velocity + incre_velocity;
2133 }
2134 
2135 template<sp_d T> void Factory<T>::update_incre_acceleration_by(const Col<T>& A) {
2136  incre_acceleration += A;
2137  trial_acceleration = current_acceleration + incre_acceleration;
2138 }
2139 
2140 template<sp_d T> void Factory<T>::update_incre_temperature_by(const Col<T>& M) {
2141  incre_temperature += M;
2142  trial_temperature = current_temperature + incre_temperature;
2143 }
2144 
2145 template<sp_d T> void Factory<T>::update_current_time_by(const T M) {
2146  trial_time = current_time += M;
2147  incre_time = 0.;
2148 }
2149 
2150 template<sp_d T> void Factory<T>::update_current_load_factor_by(const Col<T>& L) {
2151  trial_load_factor = current_load_factor += L;
2152  incre_load_factor.zeros();
2153 }
2154 
2155 template<sp_d T> void Factory<T>::update_current_load_by(const Col<T>& L) {
2156  trial_load = current_load += L;
2157  incre_load.zeros();
2158 }
2159 
2160 template<sp_d T> void Factory<T>::update_current_settlement_by(const Col<T>& S) {
2161  trial_settlement = current_settlement += S;
2162  incre_settlement.zeros();
2163 }
2164 
2165 template<sp_d T> void Factory<T>::update_current_resistance_by(const Col<T>& R) {
2166  trial_resistance = current_resistance += R;
2167  incre_resistance.zeros();
2168 }
2169 
2170 template<sp_d T> void Factory<T>::update_current_damping_force_by(const Col<T>& R) {
2171  trial_damping_force = current_damping_force += R;
2172  incre_damping_force.zeros();
2173 }
2174 
2175 template<sp_d T> void Factory<T>::update_current_nonviscous_force_by(const Col<T>& R) {
2176  trial_nonviscous_force = current_nonviscous_force += R;
2177  incre_nonviscous_force.zeros();
2178 }
2179 
2180 template<sp_d T> void Factory<T>::update_current_inertial_force_by(const Col<T>& R) {
2181  trial_inertial_force = current_inertial_force += R;
2182  incre_inertial_force.zeros();
2183 }
2184 
2185 template<sp_d T> void Factory<T>::update_current_displacement_by(const Col<T>& D) {
2186  trial_displacement = current_displacement += D;
2187  incre_displacement.zeros();
2188 }
2189 
2190 template<sp_d T> void Factory<T>::update_current_velocity_by(const Col<T>& V) {
2191  trial_velocity = current_velocity += V;
2192  incre_velocity.zeros();
2193 }
2194 
2195 template<sp_d T> void Factory<T>::update_current_acceleration_by(const Col<T>& A) {
2196  trial_acceleration = current_acceleration += A;
2197  incre_acceleration.zeros();
2198 }
2199 
2200 template<sp_d T> void Factory<T>::update_current_temperature_by(const Col<T>& M) {
2201  trial_temperature = current_temperature += M;
2202  incre_temperature.zeros();
2203 }
2204 
2205 #endif // FACTORY_HPP
2206 
void reset(ExternalMaterialData *data, int *info)
Definition: ElasticExternal.cpp:74
A Factory class.
Definition: Factory.hpp:73
const bool initialized
Definition: Factory.hpp:190
A MetaMat class that holds matrices.
Definition: MetaMat.hpp:72
void set_pre_inertial_force(const Col< T > &)
Definition: Factory.hpp:1643
T & modify_trial_time()
Definition: Factory.hpp:1749
void set_eigenvalue(const Col< T > &)
Definition: Factory.hpp:959
void clear_acceleration()
Definition: Factory.hpp:1284
void initialize_load()
Definition: Factory.hpp:845
void set_mpc(unsigned)
Definition: Factory.hpp:939
shared_ptr< MetaMat< T > > & modify_stiffness()
Definition: Factory.hpp:1527
Col< T > & modify_trial_inertial_force()
Definition: Factory.hpp:1763
Col< T > & modify_trial_nonviscous_force()
Definition: Factory.hpp:1761
void commit_temperature()
Definition: Factory.hpp:1147
bool is_sparse() const
Definition: Factory.hpp:751
void update_current_nonviscous_force(const Col< T > &)
Definition: Factory.hpp:1995
void set_incre_velocity(const Col< T > &)
Definition: Factory.hpp:1599
uword get_entry() const
Definition: Factory.hpp:709
void set_incre_temperature(const Col< T > &)
Definition: Factory.hpp:1603
void assemble_stiffness(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition: Factory.hpp:1461
void print() const
Definition: Factory.hpp:1470
void clear_resistance()
Definition: Factory.hpp:1242
void initialize_auxiliary_resistance()
Definition: Factory.hpp:905
void clear_load_factor()
Definition: Factory.hpp:1221
void set_incre_acceleration(const Col< T > &)
Definition: Factory.hpp:1601
int initialize()
Definition: Factory.hpp:782
Col< T > & modify_incre_load()
Definition: Factory.hpp:1777
const Col< T > & get_current_acceleration() const
Definition: Factory.hpp:1721
void update_incre_damping_force_by(const Col< T > &)
Definition: Factory.hpp:2110
const Col< T > & get_pre_damping_force() const
Definition: Factory.hpp:1735
void reset_displacement()
Definition: Factory.hpp:1368
void update_incre_temperature_by(const Col< T > &)
Definition: Factory.hpp:2140
void set_ninja(const Col< T > &)
Definition: Factory.hpp:933
Col< T > & modify_trial_damping_force()
Definition: Factory.hpp:1759
void update_trial_temperature_by(const Col< T > &)
Definition: Factory.hpp:2080
bool is_nlgeom() const
Definition: Factory.hpp:717
void update_incre_load_factor_by(const Col< T > &)
Definition: Factory.hpp:2090
const Mat< T > & get_eigenvector() const
Definition: Factory.hpp:1035
const Col< T > & get_trial_velocity() const
Definition: Factory.hpp:1671
const Col< T > & get_trial_resistance() const
Definition: Factory.hpp:1661
void update_trial_load_factor_by(const Col< T > &)
Definition: Factory.hpp:2030
const Col< T > & get_ninja() const
Definition: Factory.hpp:963
Col< T > & modify_trial_temperature()
Definition: Factory.hpp:1771
void set_pre_velocity(const Col< T > &)
Definition: Factory.hpp:1647
const Col< T > & get_current_load() const
Definition: Factory.hpp:1705
void set_trial_settlement(const Col< T > &)
Definition: Factory.hpp:1563
void set_incre_load_factor(const Col< T > &)
Definition: Factory.hpp:1583
void set_trial_time(T)
Definition: Factory.hpp:1557
Col< T > & modify_incre_resistance()
Definition: Factory.hpp:1781
void update_incre_load_by(const Col< T > &)
Definition: Factory.hpp:2095
void set_solver_type(SolverType)
Definition: Factory.hpp:727
const Col< T > & get_pre_load() const
Definition: Factory.hpp:1729
void commit_acceleration()
Definition: Factory.hpp:1141
AnalysisType get_analysis_type() const
Definition: Factory.hpp:741
void clear_nonviscous()
Definition: Factory.hpp:1406
void set_current_resistance(const Col< T > &)
Definition: Factory.hpp:1613
void update_incre_displacement(const Col< T > &)
Definition: Factory.hpp:1945
void set_incre_displacement(const Col< T > &)
Definition: Factory.hpp:1597
void update_current_acceleration_by(const Col< T > &)
Definition: Factory.hpp:2195
const uvec & get_auxiliary_encoding() const
Definition: Factory.hpp:971
void reset_auxiliary_resistance()
Definition: Factory.hpp:1392
Col< T > & modify_pre_displacement()
Definition: Factory.hpp:1837
void commit_resistance()
Definition: Factory.hpp:1105
const Col< T > & get_incre_inertial_force() const
Definition: Factory.hpp:1691
void set_storage_scheme(StorageScheme)
Definition: Factory.hpp:743
void update_trial_nonviscous_force_by(const Col< T > &)
Definition: Factory.hpp:2055
void initialize_damping()
Definition: Factory.hpp:912
void update_incre_velocity(const Col< T > &)
Definition: Factory.hpp:1950
Col< T > & modify_incre_acceleration()
Definition: Factory.hpp:1793
void initialize_geometry()
Definition: Factory.hpp:922
void update_incre_displacement_by(const Col< T > &)
Definition: Factory.hpp:2125
void update_current_damping_force_by(const Col< T > &)
Definition: Factory.hpp:2170
void assemble_mass(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition: Factory.hpp:1455
void set_trial_velocity(const Col< T > &)
Definition: Factory.hpp:1575
void update_current_resistance(const Col< T > &)
Definition: Factory.hpp:1985
void clear_inertial_force()
Definition: Factory.hpp:1263
Col< T > & modify_ninja()
Definition: Factory.hpp:1531
void set_pre_nonviscous_force(const Col< T > &)
Definition: Factory.hpp:1641
void update_trial_damping_force_by(const Col< T > &)
Definition: Factory.hpp:2050
unsigned get_mpc() const
Definition: Factory.hpp:967
Col< T > & modify_eigenvalue()
Definition: Factory.hpp:1553
const Col< T > & get_current_settlement() const
Definition: Factory.hpp:1707
const SolverSetting< double > & get_solver_setting() const
Definition: Factory.hpp:733
void commit_load_factor()
Definition: Factory.hpp:1087
std::mutex & get_damping_mutex()
Definition: Factory.hpp:1025
void update_trial_velocity_by(const Col< T > &)
Definition: Factory.hpp:2070
void update_current_load_by(const Col< T > &)
Definition: Factory.hpp:2155
SpCol< T > & modify_trial_constraint_resistance()
Definition: Factory.hpp:1549
const shared_ptr< MetaMat< T > > & get_damping() const
Definition: Factory.hpp:999
void set_pre_settlement(const Col< T > &)
Definition: Factory.hpp:1635
void set_nlgeom(bool)
Definition: Factory.hpp:711
StorageScheme get_storage_scheme() const
Definition: Factory.hpp:749
void commit_inertial_force()
Definition: Factory.hpp:1123
void update_current_time(T)
Definition: Factory.hpp:1965
void initialize_resistance()
Definition: Factory.hpp:857
const Col< T > & get_current_damping_force() const
Definition: Factory.hpp:1711
void set_solver_setting(const SolverSetting< double > &)
Definition: Factory.hpp:731
void update_trial_displacement_by(const Col< T > &)
Definition: Factory.hpp:2065
void reset_velocity()
Definition: Factory.hpp:1374
const shared_ptr< MetaMat< T > > & get_nonviscous() const
Definition: Factory.hpp:1001
const Col< T > & get_trial_load_factor() const
Definition: Factory.hpp:1655
void set_current_load(const Col< T > &)
Definition: Factory.hpp:1609
std::mutex & get_trial_load_mutex()
Definition: Factory.hpp:1017
void initialize_stiffness()
Definition: Factory.hpp:920
void initialize_damping_force()
Definition: Factory.hpp:863
T get_current_time() const
Definition: Factory.hpp:1701
const Col< T > & get_pre_settlement() const
Definition: Factory.hpp:1731
void update_trial_temperature(const Col< T > &)
Definition: Factory.hpp:1900
void commit_velocity()
Definition: Factory.hpp:1135
void set_pre_time(T)
Definition: Factory.hpp:1629
T get_complementary_energy()
Definition: Factory.hpp:993
void reset()
Definition: Factory.hpp:1420
const Col< T > & get_incre_displacement() const
Definition: Factory.hpp:1693
void set_current_damping_force(const Col< T > &)
Definition: Factory.hpp:1615
Col< T > & modify_current_nonviscous_force()
Definition: Factory.hpp:1809
void update_trial_acceleration(const Col< T > &)
Definition: Factory.hpp:1895
void set_current_time(T)
Definition: Factory.hpp:1605
const Col< T > & get_auxiliary_load() const
Definition: Factory.hpp:977
void commit_pre_settlement()
Definition: Factory.hpp:1179
Col< T > & modify_pre_settlement()
Definition: Factory.hpp:1827
const shared_ptr< MetaMat< T > > & get_mass() const
Definition: Factory.hpp:997
void commit_status()
Definition: Factory.hpp:1062
void commit_pre_load_factor()
Definition: Factory.hpp:1175
const suanpan::set< uword > & get_reference_dof() const
Definition: Factory.hpp:776
std::pair< unsigned, unsigned > get_bandwidth() const
Definition: Factory.hpp:761
void set_incre_time(T)
Definition: Factory.hpp:1581
shared_ptr< MetaMat< T > > & modify_damping()
Definition: Factory.hpp:1523
Col< T > & modify_pre_velocity()
Definition: Factory.hpp:1839
std::mutex & get_reference_load_mutex()
Definition: Factory.hpp:1021
void assemble_damping(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition: Factory.hpp:1457
void update_incre_nonviscous_force(const Col< T > &)
Definition: Factory.hpp:1935
void update_trial_inertial_force(const Col< T > &)
Definition: Factory.hpp:1880
void set_current_velocity(const Col< T > &)
Definition: Factory.hpp:1623
void commit_displacement()
Definition: Factory.hpp:1129
void update_current_velocity_by(const Col< T > &)
Definition: Factory.hpp:2190
const Col< T > & get_incre_damping_force() const
Definition: Factory.hpp:1687
void update_trial_load_factor(const Col< T > &)
Definition: Factory.hpp:1850
void reset_load_factor()
Definition: Factory.hpp:1326
void clear_velocity()
Definition: Factory.hpp:1277
void update_current_load(const Col< T > &)
Definition: Factory.hpp:1975
void update_incre_time(T)
Definition: Factory.hpp:1905
void update_incre_settlement_by(const Col< T > &)
Definition: Factory.hpp:2100
unsigned get_reference_size() const
Definition: Factory.hpp:770
void update_current_acceleration(const Col< T > &)
Definition: Factory.hpp:2015
StorageScheme
Definition: Factory.hpp:52
void commit_pre_temperature()
Definition: Factory.hpp:1195
void commit_pre_status()
Definition: Factory.hpp:1158
void update_current_time_by(T)
Definition: Factory.hpp:2145
suanpan::set< uword > & modify_reference_dof()
Definition: Factory.hpp:1535
void set_trial_acceleration(const Col< T > &)
Definition: Factory.hpp:1577
void update_trial_damping_force(const Col< T > &)
Definition: Factory.hpp:1870
const Col< T > & get_incre_settlement() const
Definition: Factory.hpp:1683
void reset_inertial_force()
Definition: Factory.hpp:1362
void set_reference_dof(const suanpan::set< uword > &)
Definition: Factory.hpp:774
const Col< T > & get_momentum()
Definition: Factory.hpp:995
void update_current_displacement(const Col< T > &)
Definition: Factory.hpp:2005
void clear_settlement()
Definition: Factory.hpp:1235
void set_pre_load(const Col< T > &)
Definition: Factory.hpp:1633
std::mutex & get_auxiliary_stiffness_mutex()
Definition: Factory.hpp:1013
T get_incre_time() const
Definition: Factory.hpp:1677
void commit_pre_damping_force()
Definition: Factory.hpp:1183
Col< T > & modify_current_inertial_force()
Definition: Factory.hpp:1811
const Col< T > & get_current_nonviscous_force() const
Definition: Factory.hpp:1713
T get_strain_energy()
Definition: Factory.hpp:985
void commit_pre_inertial_force()
Definition: Factory.hpp:1187
void update_trial_settlement(const Col< T > &)
Definition: Factory.hpp:1860
void set_trial_inertial_force(const Col< T > &)
Definition: Factory.hpp:1571
const Col< T > & get_eigenvalue() const
Definition: Factory.hpp:1033
void commit_pre_displacement()
Definition: Factory.hpp:1189
std::mutex & get_geometry_mutex()
Definition: Factory.hpp:1031
void set_incre_nonviscous_force(const Col< T > &)
Definition: Factory.hpp:1593
void commit_pre_acceleration()
Definition: Factory.hpp:1193
const Col< T > & get_incre_velocity() const
Definition: Factory.hpp:1695
Col< T > & modify_current_settlement()
Definition: Factory.hpp:1803
T get_trial_time() const
Definition: Factory.hpp:1653
void update_current_temperature_by(const Col< T > &)
Definition: Factory.hpp:2200
void set_incre_load(const Col< T > &)
Definition: Factory.hpp:1585
void update_reference_size()
Definition: Factory.hpp:763
void set_pre_acceleration(const Col< T > &)
Definition: Factory.hpp:1649
void set_incre_inertial_force(const Col< T > &)
Definition: Factory.hpp:1595
const Col< T > & get_trial_nonviscous_force() const
Definition: Factory.hpp:1665
T get_kinetic_energy()
Definition: Factory.hpp:987
void set_current_settlement(const Col< T > &)
Definition: Factory.hpp:1611
void update_trial_settlement_by(const Col< T > &)
Definition: Factory.hpp:2040
void initialize_inertial_force()
Definition: Factory.hpp:875
void set_pre_resistance(const Col< T > &)
Definition: Factory.hpp:1637
void commit_load()
Definition: Factory.hpp:1093
void update_current_nonviscous_force_by(const Col< T > &)
Definition: Factory.hpp:2175
void update_trial_acceleration_by(const Col< T > &)
Definition: Factory.hpp:2075
const Col< T > & get_pre_acceleration() const
Definition: Factory.hpp:1745
void update_incre_acceleration_by(const Col< T > &)
Definition: Factory.hpp:2135
void commit_pre_nonviscous_force()
Definition: Factory.hpp:1185
void commit_pre_load()
Definition: Factory.hpp:1177
void commit_pre_velocity()
Definition: Factory.hpp:1191
void update_trial_resistance_by(const Col< T > &)
Definition: Factory.hpp:2045
Col< T > & modify_current_resistance()
Definition: Factory.hpp:1805
AnalysisType
Definition: Factory.hpp:43
bool is_nonviscous() const
Definition: Factory.hpp:725
void commit_pre_time()
Definition: Factory.hpp:1173
shared_ptr< MetaMat< T > > & modify_geometry()
Definition: Factory.hpp:1529
void commit_auxiliary_resistance()
Definition: Factory.hpp:1153
void update_current_settlement(const Col< T > &)
Definition: Factory.hpp:1980
void set_stiffness(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:955
void update_trial_load(const Col< T > &)
Definition: Factory.hpp:1855
Col< T > & modify_pre_nonviscous_force()
Definition: Factory.hpp:1833
void update_incre_load(const Col< T > &)
Definition: Factory.hpp:1915
void update_trial_time(T)
Definition: Factory.hpp:1845
void clear_nonviscous_force()
Definition: Factory.hpp:1256
const shared_ptr< MetaMat< T > > & get_stiffness() const
Definition: Factory.hpp:1003
void commit_settlement()
Definition: Factory.hpp:1099
void set_trial_nonviscous_force(const Col< T > &)
Definition: Factory.hpp:1569
const Col< T > & get_auxiliary_lambda() const
Definition: Factory.hpp:973
void reset_load()
Definition: Factory.hpp:1332
void reset_damping_force()
Definition: Factory.hpp:1350
void update_current_inertial_force_by(const Col< T > &)
Definition: Factory.hpp:2180
void set_trial_temperature(const Col< T > &)
Definition: Factory.hpp:1579
std::mutex & get_auxiliary_encoding_mutex()
Definition: Factory.hpp:1007
Col< T > & modify_trial_load()
Definition: Factory.hpp:1753
SolverType
Definition: Factory.hpp:61
void assemble_damping_force(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1433
void set_current_nonviscous_force(const Col< T > &)
Definition: Factory.hpp:1617
const Col< T > & get_pre_displacement() const
Definition: Factory.hpp:1741
const Col< T > & get_current_inertial_force() const
Definition: Factory.hpp:1715
std::mutex & get_mass_mutex()
Definition: Factory.hpp:1023
void set_nonviscous(bool)
Definition: Factory.hpp:719
const Col< T > & get_pre_load_factor() const
Definition: Factory.hpp:1727
Col< T > & modify_current_acceleration()
Definition: Factory.hpp:1817
const Col< T > & get_trial_displacement() const
Definition: Factory.hpp:1669
void update_trial_velocity(const Col< T > &)
Definition: Factory.hpp:1890
void update_current_temperature(const Col< T > &)
Definition: Factory.hpp:2020
void reset_temperature()
Definition: Factory.hpp:1386
std::mutex & get_auxiliary_resistance_mutex()
Definition: Factory.hpp:1009
void set_current_acceleration(const Col< T > &)
Definition: Factory.hpp:1625
T get_nonviscous_energy()
Definition: Factory.hpp:991
void update_current_settlement_by(const Col< T > &)
Definition: Factory.hpp:2160
void initialize_settlement()
Definition: Factory.hpp:851
void clear_damping()
Definition: Factory.hpp:1404
const Col< T > & get_incre_nonviscous_force() const
Definition: Factory.hpp:1689
const Col< T > & get_incre_load_factor() const
Definition: Factory.hpp:1679
void clear_energy()
Definition: Factory.hpp:1053
const SpMat< T > & get_reference_load() const
Definition: Factory.hpp:969
void initialize_displacement()
Definition: Factory.hpp:881
Col< T > & modify_auxiliary_load()
Definition: Factory.hpp:1545
void reset_time()
Definition: Factory.hpp:1321
const Col< T > & get_current_displacement() const
Definition: Factory.hpp:1717
Col< T > & modify_current_velocity()
Definition: Factory.hpp:1815
unsigned get_size() const
Definition: Factory.hpp:702
T & modify_current_time()
Definition: Factory.hpp:1797
void reset_settlement()
Definition: Factory.hpp:1338
void clear_damping_force()
Definition: Factory.hpp:1249
void set_analysis_type(AnalysisType)
Definition: Factory.hpp:735
const Col< T > & get_incre_load() const
Definition: Factory.hpp:1681
void initialize_nonviscous_force()
Definition: Factory.hpp:869
void set_entry(uword)
Definition: Factory.hpp:704
Col< T > & modify_pre_damping_force()
Definition: Factory.hpp:1831
const Col< T > & get_incre_acceleration() const
Definition: Factory.hpp:1697
T get_pre_time() const
Definition: Factory.hpp:1725
void clear_geometry()
Definition: Factory.hpp:1410
T get_viscous_energy()
Definition: Factory.hpp:989
std::mutex & get_trial_constraint_resistance_mutex()
Definition: Factory.hpp:1015
const Col< T > & get_current_temperature() const
Definition: Factory.hpp:1723
void update_incre_time_by(T)
Definition: Factory.hpp:2085
void assemble_nonviscous(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition: Factory.hpp:1459
void set_reference_size(unsigned)
Definition: Factory.hpp:765
Col< T > & modify_trial_settlement()
Definition: Factory.hpp:1755
void initialize_temperature()
Definition: Factory.hpp:899
void set_current_displacement(const Col< T > &)
Definition: Factory.hpp:1621
void update_reference_dof(const uvec &)
Definition: Factory.hpp:772
const Col< T > & get_trial_inertial_force() const
Definition: Factory.hpp:1667
void set_current_inertial_force(const Col< T > &)
Definition: Factory.hpp:1619
SpCol< T > & modify_current_constraint_resistance()
Definition: Factory.hpp:1551
Col< T > & modify_incre_load_factor()
Definition: Factory.hpp:1775
void set_trial_load(const Col< T > &)
Definition: Factory.hpp:1561
Col< T > & modify_auxiliary_resistance()
Definition: Factory.hpp:1543
Col< T > & modify_pre_resistance()
Definition: Factory.hpp:1829
Col< T > & modify_incre_temperature()
Definition: Factory.hpp:1795
T & modify_pre_time()
Definition: Factory.hpp:1821
void update_incre_temperature(const Col< T > &)
Definition: Factory.hpp:1960
shared_ptr< MetaMat< T > > & modify_nonviscous()
Definition: Factory.hpp:1525
void update_current_velocity(const Col< T > &)
Definition: Factory.hpp:2010
void set_error(T)
Definition: Factory.hpp:778
Col< T > & modify_incre_nonviscous_force()
Definition: Factory.hpp:1785
void assemble_geometry(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition: Factory.hpp:1463
const Col< T > & get_trial_settlement() const
Definition: Factory.hpp:1659
void set_eigenvector(const Mat< T > &)
Definition: Factory.hpp:961
void initialize_load_factor()
Definition: Factory.hpp:835
Col< T > & modify_incre_inertial_force()
Definition: Factory.hpp:1787
const Col< T > & get_trial_load() const
Definition: Factory.hpp:1657
uvec & modify_auxiliary_encoding()
Definition: Factory.hpp:1539
void update_current_displacement_by(const Col< T > &)
Definition: Factory.hpp:2185
void set_trial_resistance(const Col< T > &)
Definition: Factory.hpp:1565
void commit_nonviscous_force()
Definition: Factory.hpp:1117
void set_pre_load_factor(const Col< T > &)
Definition: Factory.hpp:1631
void update_current_damping_force(const Col< T > &)
Definition: Factory.hpp:1990
Col< T > & modify_current_damping_force()
Definition: Factory.hpp:1807
void update_current_load_factor(const Col< T > &)
Definition: Factory.hpp:1970
void commit_energy()
Definition: Factory.hpp:1037
void commit_damping_force()
Definition: Factory.hpp:1111
T get_error() const
Definition: Factory.hpp:780
void update_trial_time_by(T)
Definition: Factory.hpp:2025
void update_incre_acceleration(const Col< T > &)
Definition: Factory.hpp:1955
const SpMat< T > & get_auxiliary_stiffness() const
Definition: Factory.hpp:979
void assemble_inertial_force(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1443
void clear_status()
Definition: Factory.hpp:1197
void set_bandwidth(unsigned, unsigned)
Definition: Factory.hpp:753
void update_incre_velocity_by(const Col< T > &)
Definition: Factory.hpp:2130
const Col< T > & get_current_load_factor() const
Definition: Factory.hpp:1703
void commit_time()
Definition: Factory.hpp:1082
void clear_mass()
Definition: Factory.hpp:1402
const Col< T > & get_pre_inertial_force() const
Definition: Factory.hpp:1739
void commit_pre_resistance()
Definition: Factory.hpp:1181
void update_trial_nonviscous_force(const Col< T > &)
Definition: Factory.hpp:1875
void reset_resistance()
Definition: Factory.hpp:1344
SolverType get_solver_type() const
Definition: Factory.hpp:729
Col< T > & modify_pre_temperature()
Definition: Factory.hpp:1843
std::mutex & get_auxiliary_load_mutex()
Definition: Factory.hpp:1011
const Col< T > & get_trial_acceleration() const
Definition: Factory.hpp:1673
void assemble_resistance(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1428
void set_incre_damping_force(const Col< T > &)
Definition: Factory.hpp:1591
Col< T > & modify_incre_displacement()
Definition: Factory.hpp:1789
Col< T > & modify_pre_load()
Definition: Factory.hpp:1825
const SpCol< T > & get_trial_constraint_resistance() const
Definition: Factory.hpp:981
void update_incre_damping_force(const Col< T > &)
Definition: Factory.hpp:1930
const Col< T > & get_pre_velocity() const
Definition: Factory.hpp:1743
const Col< T > & get_current_resistance() const
Definition: Factory.hpp:1709
void assemble_nonviscous_force(const Mat< T > &, const uvec &)
Definition: Factory.hpp:1438
void initialize_velocity()
Definition: Factory.hpp:887
std::mutex & get_stiffness_mutex()
Definition: Factory.hpp:1029
Col< T > & modify_pre_acceleration()
Definition: Factory.hpp:1841
void update_incre_inertial_force_by(const Col< T > &)
Definition: Factory.hpp:2120
const Col< T > & get_current_velocity() const
Definition: Factory.hpp:1719
void update_incre_load_factor(const Col< T > &)
Definition: Factory.hpp:1910
SpMat< T > & modify_auxiliary_stiffness()
Definition: Factory.hpp:1547
void set_incre_resistance(const Col< T > &)
Definition: Factory.hpp:1589
void initialize_nonviscous()
Definition: Factory.hpp:914
Col< T > & modify_trial_velocity()
Definition: Factory.hpp:1767
Col< T > & modify_trial_displacement()
Definition: Factory.hpp:1765
const Col< T > & get_pre_temperature() const
Definition: Factory.hpp:1747
Col< T > & modify_trial_acceleration()
Definition: Factory.hpp:1769
void clear_load()
Definition: Factory.hpp:1228
const Col< T > & get_incre_temperature() const
Definition: Factory.hpp:1699
Col< T > & modify_pre_load_factor()
Definition: Factory.hpp:1823
Factory(unsigned=0, AnalysisType=AnalysisType::NONE, StorageScheme=StorageScheme::FULL)
Definition: Factory.hpp:691
Col< T > & modify_trial_load_factor()
Definition: Factory.hpp:1751
void set_size(unsigned)
Definition: Factory.hpp:696
void update_trial_displacement(const Col< T > &)
Definition: Factory.hpp:1885
const Col< T > & get_incre_resistance() const
Definition: Factory.hpp:1685
Col< T > & modify_trial_resistance()
Definition: Factory.hpp:1757
const Col< T > & get_sushi() const
Definition: Factory.hpp:965
const Col< T > & get_auxiliary_resistance() const
Definition: Factory.hpp:975
const Col< T > & get_trial_damping_force() const
Definition: Factory.hpp:1663
Mat< T > & modify_eigenvector()
Definition: Factory.hpp:1555
void update_current_resistance_by(const Col< T > &)
Definition: Factory.hpp:2165
void set_damping(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:951
void set_pre_temperature(const Col< T > &)
Definition: Factory.hpp:1651
void update_trial_load_by(const Col< T > &)
Definition: Factory.hpp:2035
Col< T > & modify_auxiliary_lambda()
Definition: Factory.hpp:1541
const SpCol< T > & get_current_constraint_resistance() const
Definition: Factory.hpp:983
void clear_auxiliary_resistance()
Definition: Factory.hpp:1298
const Col< T > & get_pre_nonviscous_force() const
Definition: Factory.hpp:1737
void update_incre_resistance_by(const Col< T > &)
Definition: Factory.hpp:2105
void update_sushi_by(const Col< T > &)
Definition: Factory.hpp:937
void set_incre_settlement(const Col< T > &)
Definition: Factory.hpp:1587
Col< T > & modify_pre_inertial_force()
Definition: Factory.hpp:1835
void set_sushi(const Col< T > &)
Definition: Factory.hpp:935
shared_ptr< MetaMat< T > > & modify_mass()
Definition: Factory.hpp:1521
void update_incre_settlement(const Col< T > &)
Definition: Factory.hpp:1920
void set_pre_displacement(const Col< T > &)
Definition: Factory.hpp:1645
std::mutex & get_nonviscous_mutex()
Definition: Factory.hpp:1027
void set_current_load_factor(const Col< T > &)
Definition: Factory.hpp:1607
void set_trial_displacement(const Col< T > &)
Definition: Factory.hpp:1573
void update_trial_resistance(const Col< T > &)
Definition: Factory.hpp:1865
Col< T > & modify_current_load_factor()
Definition: Factory.hpp:1799
void clear_eigen()
Definition: Factory.hpp:1397
void clear_auxiliary()
Definition: Factory.hpp:1412
void set_geometry(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:957
void set_reference_load(const SpMat< T > &)
Definition: Factory.hpp:947
T & modify_incre_time()
Definition: Factory.hpp:1773
void update_incre_resistance(const Col< T > &)
Definition: Factory.hpp:1925
SpMat< T > & modify_reference_load()
Definition: Factory.hpp:1537
Col< T > & modify_incre_damping_force()
Definition: Factory.hpp:1783
Col< T > & modify_incre_settlement()
Definition: Factory.hpp:1779
void update_current_load_factor_by(const Col< T > &)
Definition: Factory.hpp:2150
void set_trial_load_factor(const Col< T > &)
Definition: Factory.hpp:1559
void initialize_eigen()
Definition: Factory.hpp:928
void update_incre_inertial_force(const Col< T > &)
Definition: Factory.hpp:1940
void clear_time()
Definition: Factory.hpp:1219
void clear_displacement()
Definition: Factory.hpp:1270
void set_mass(const shared_ptr< MetaMat< T >> &)
Definition: Factory.hpp:949
const shared_ptr< MetaMat< T > > & get_geometry() const
Definition: Factory.hpp:1005
Col< T > & modify_current_load()
Definition: Factory.hpp:1801
void reset_status()
Definition: Factory.hpp:1303
Col< T > & modify_current_temperature()
Definition: Factory.hpp:1819
void update_current_inertial_force(const Col< T > &)
Definition: Factory.hpp:2000
const Col< T > & get_trial_temperature() const
Definition: Factory.hpp:1675
void set_trial_damping_force(const Col< T > &)
Definition: Factory.hpp:1567
std::mutex & get_trial_settlement_mutex()
Definition: Factory.hpp:1019
void initialize_mass()
Definition: Factory.hpp:910
Col< T > & modify_sushi()
Definition: Factory.hpp:1533
void set_pre_damping_force(const Col< T > &)
Definition: Factory.hpp:1639
void clear_temperature()
Definition: Factory.hpp:1291
Col< T > & modify_current_displacement()
Definition: Factory.hpp:1813
void initialize_acceleration()
Definition: Factory.hpp:893
void update_trial_inertial_force_by(const Col< T > &)
Definition: Factory.hpp:2060
void set_current_temperature(const Col< T > &)
Definition: Factory.hpp:1627
Col< T > & modify_incre_velocity()
Definition: Factory.hpp:1791
void reset_acceleration()
Definition: Factory.hpp:1380
void reset_nonviscous_force()
Definition: Factory.hpp:1356
const Col< T > & get_pre_resistance() const
Definition: Factory.hpp:1733
void clear_stiffness()
Definition: Factory.hpp:1408
void update_incre_nonviscous_force_by(const Col< T > &)
Definition: Factory.hpp:2115
void error(const std::string_view file_name, const int line, const std::string_view format_str, const T &... args)
Definition: suanPan.h:234
std::set< T > set
Definition: container.h:54
#define suanpan_info
Definition: suanPan.h:305