File: Thermostat.h

package info (click to toggle)
lammps 20220106.git7586adbb6a%2Bds1-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 348,064 kB
  • sloc: cpp: 831,421; python: 24,896; xml: 14,949; f90: 10,845; ansic: 7,967; sh: 4,226; perl: 4,064; fortran: 2,424; makefile: 1,501; objc: 238; lisp: 163; csh: 16; awk: 14; tcl: 6
file content (1061 lines) | stat: -rw-r--r-- 30,675 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
#ifndef THERMOSTAT_H
#define THERMOSTAT_H

#include "AtomicRegulator.h"
#include "PerAtomQuantityLibrary.h"
#include <map>
#include <set>
#include <string>

namespace ATC {

  static const int myLambdaMaxIterations = 50;

  // forward declarations
  class ThermalTimeIntegrator;
  class AtfShapeFunctionRestriction;
  class FundamentalAtomQuantity;
  class PrescribedDataManager;

  /**
   *  @class  Thermostat
   *  @brief  Manager class for atom-continuum control of thermal energy
   */
  class Thermostat : public AtomicRegulator {

  public:

    // constructor
    Thermostat(ATC_Coupling * atc,
               const std::string & regulatorPrefix = "");

    // destructor
    virtual ~Thermostat(){};

    /** parser/modifier */
    virtual bool modify(int narg, char **arg);

    /** instantiate up the desired method(s) */
    virtual void construct_methods();

    // data access, intended for method objects
    /** reset the nodal power to a prescribed value */
    virtual void reset_lambda_contribution(const DENS_MAT & target);

    /** return value for the correction maximum number of iterations */
    int lambda_max_iterations() {return lambdaMaxIterations_;};

  protected:

    // data regarding fixed nodes and applied fluxes
    /** set of all fixed nodes */
    std::set<int> fixedNodes_;
    /** set of all nodes which have a flux applied */
    std::set<int> fluxNodes_;

    /** maximum number of iterations used in iterative solve for lambda */
    int lambdaMaxIterations_;

  private:

    // DO NOT define this
    Thermostat();

  };

  /**
   *  @class  ThermostatShapeFunction
   *  @brief  Class for thermostat algorithms using the shape function matrices
   *          (thermostats have general for of N^T w N lambda = rhs)
   */

  class ThermostatShapeFunction : public RegulatorShapeFunction {


  public:

    ThermostatShapeFunction(AtomicRegulator * thermostat,
                            const std::string & regulatorPrefix = "");

    virtual ~ThermostatShapeFunction() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

  protected:

    // methods
    /** set weighting factor for in matrix Nhat^T * weights * Nhat */
    virtual void set_weights();

    // member data

    /** MD mass matrix */
    DIAG_MAN & mdMassMatrix_;

    /** pointer to atom velocities */
    FundamentalAtomQuantity * atomVelocities_;

    /** workspace variables */
    DENS_VEC _weightVector_, _maskedWeightVector_;

  private:

    // DO NOT define this
    ThermostatShapeFunction();

  };

  /**
   *  @class  ThermostatRescale
   *  @brief  Enforces constraint on atomic kinetic energy based on FE temperature
   */

  class ThermostatRescale : public ThermostatShapeFunction {

  public:

    friend class KinetoThermostatRescale; // since this is sometimes used as a set of member functions for friend

    ThermostatRescale(AtomicRegulator * thermostat);

    virtual ~ThermostatRescale() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** applies thermostat to atoms in the post-corrector phase */
    virtual void apply_post_corrector(double dt);

    /** compute boundary flux, requires thermostat input since it is part of the coupling scheme */
    virtual void compute_boundary_flux(FIELDS & /* fields */)
      {boundaryFlux_[TEMPERATURE] = 0.;};

    /** get data for output */
    virtual void output(OUTPUT_LIST & outputData);

  protected:

    /** set weighting factor for in matrix Nhat^T * weights * Nhat */
    virtual void set_weights();

    /** sets up and solves thermostat equations */
    virtual void compute_thermostat(double dt);

    /** apply solution to atomic quantities */
    void apply_to_atoms(PerAtomQuantity<double> * atomVelocities);

    /** construct the RHS vector */
    virtual void set_rhs(DENS_MAT & rhs);

    /** FE temperature field */
    DENS_MAN & nodalTemperature_;

    /** construction for prolongation of lambda to atoms */
    AtomicVelocityRescaleFactor * atomVelocityRescalings_;

    /** workspace variables */
    DENS_MAT _rhs_;

  private:

    // DO NOT define this
    ThermostatRescale();

  };

  /**
   *  @class  ThermostatRescaleMixedKePe
   *  @brief  Enforces constraint on atomic kinetic energy based on FE temperature
   *          when the temperature is a mix of the KE and PE
   */

  class ThermostatRescaleMixedKePe : public ThermostatRescale {

  public:

    ThermostatRescaleMixedKePe(AtomicRegulator * thermostat);

    virtual ~ThermostatRescaleMixedKePe() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** pre-run initialization of method data */
    virtual void initialize();

  protected:

    /** set weighting factor for in matrix Nhat^T * weights * Nhat */
    virtual void set_weights();

    /** construct the RHS vector */
    virtual void set_rhs(DENS_MAT & rhs);

    /** nodal fluctuating potential energy */
    DENS_MAN * nodalAtomicFluctuatingPotentialEnergy_;

    /** fraction of temperature from KE */
    double keMultiplier_;

    /** fraction of temperature from PE */
    double peMultiplier_;

  private:

    // DO NOT define this
    ThermostatRescaleMixedKePe();

  };

    /**
   *  @class  ThermostatFsSolver
   *  @brief  Class for solving the linear system for lambda
   *          (thermostats have general for of N^T w N lambda = rhs)
   */

  class ThermostatFsSolver : public RegulatorShapeFunction {


  public:

    ThermostatFsSolver(AtomicRegulator * thermostat,
                       int lambdaMaxIterations,
                       const std::string & regulatorPrefix = "");

    virtual ~ThermostatFsSolver() {};

    /** pre-run initialization of method data */
    virtual void initialize();

    /* sets up and solves the linear system for lambda */
    virtual void compute_lambda(const DENS_MAT & rhs,
                                bool iterateSolution = true);

    /* scales lambda */
    virtual void scale_lambda(double factor) {*lambda_ *= factor;};

    /** change the time step factor */
    virtual void set_timestep_factor(double dtFactor) {dtFactor_ = dtFactor;};

  protected:

    // methods
    /** solves the non-linear equation for lambda iteratively */
    void iterate_lambda(const MATRIX & rhs);

    /** set weighting factor for in matrix Nhat^T * weights * Nhat */
    virtual void set_weights();

    // data
    /** mapping from all to regulated nodes */
    DENS_MAT rhsMap_;

    /** maximum number of iterations used in iterative solve for lambda */
    int lambdaMaxIterations_;

    /** pointer to the values of lambda interpolated to atoms */
    DENS_MAN * rhsLambdaSquared_;

    /** fraction of timestep over which constraint is exactly enforced */
    double dtFactor_;

    // workspace
    DENS_MAT _lambdaOld_; // lambda from previous iteration
    DENS_MAT _rhsOverlap_; // normal RHS vector mapped to overlap nodes
    DENS_VEC _rhsTotal_; // normal + 2nd order RHS for the iteration loop
    DENS_VEC _weightVector_, _maskedWeightVector_;

  private:

    // DO NOT define this
    ThermostatFsSolver();

  };

  /**
   *  @class  ThermostatGlcFs
   *  @brief  Class for thermostat algorithms which perform the time-integration component of the fractional step method
   */

  class ThermostatGlcFs : public RegulatorMethod {

  public:

    ThermostatGlcFs(AtomicRegulator * thermostat,
                    int lambdaMaxIterations,
                    const std::string & regulatorPrefix = "");

    virtual ~ThermostatGlcFs() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** pre-run initialization of method data */
    virtual void initialize();

    /** applies thermostat to atoms in the predictor phase */
    virtual void apply_pre_predictor(double dt);

    /** applies thermostat to atoms in the pre-corrector phase */
    virtual void apply_pre_corrector(double dt);

    /** applies thermostat to atoms in the post-corrector phase */
    virtual void apply_post_corrector(double dt);

    /** compute boundary flux, requires regulator input since it is part of the coupling scheme */
    virtual void compute_boundary_flux(FIELDS & fields);

    /** get data for output */
    virtual void output(OUTPUT_LIST & outputData);

    /* flag for performing the full lambda prediction calculation */
    bool full_prediction();

    /** set up atom to material identification */
    virtual void reset_atom_materials(const Array<int> & elementToMaterialMap,
                                      const MatrixDependencyManager<DenseMatrix, int> * atomElement);

  protected:

    // methods
    /** sets up appropriate rhs for thermostat equations */
    virtual void set_thermostat_rhs(DENS_MAT & rhs,
                                    double dt) = 0;

    /** apply forces to atoms */
    virtual void apply_to_atoms(PerAtomQuantity<double> * atomicVelocity,
                                const DENS_MAN * nodalAtomicEnergy,
                                const DENS_MAT & lambdaForce,
                                DENS_MAT & nodalAtomicLambdaPower,
                                double dt);

    /** add contributions from thermostat to FE energy */
    virtual void add_to_energy(const DENS_MAT & nodalLambdaPower,
                               DENS_MAT & deltaEnergy,
                               double dt) = 0;

    /* sets up and solves the linear system for lambda */
    virtual void compute_lambda(double dt,
                                bool iterateSolution = true);

    // member data
    /** solver for lambda linear system */
    ThermostatFsSolver * lambdaSolver_;


    /** MD mass matrix */
    DIAG_MAN & mdMassMatrix_;

    /** pointer to atom velocities */
    FundamentalAtomQuantity * atomVelocities_;

    /** reference to AtC FE temperature */
    DENS_MAN & temperature_;

    /** pointer to a time filtering object */
    TimeFilter * timeFilter_;

    /** power induced by lambda */
    DENS_MAN * nodalAtomicLambdaPower_;

    /** filtered lambda power */
    DENS_MAN * lambdaPowerFiltered_;

    /** lambda at atomic locations */
    PerAtomQuantity<double> * atomLambdas_;

    /** atomic force induced by lambda */
    AtomicThermostatForce * atomThermostatForces_;

    /** pointer to atom masses */
    FundamentalAtomQuantity * atomMasses_;

    /** hack to determine if first timestep has been passed */
    bool isFirstTimestep_;

    /** nodal atomic energy */
    DENS_MAN * nodalAtomicEnergy_;

    /** local version of velocity used as predicted final veloctiy */
    PerAtomQuantity<double> * atomPredictedVelocities_;

    /** predicted nodal atomic energy */
    AtfShapeFunctionRestriction * nodalAtomicPredictedEnergy_;

    /** pointer for force applied in first time step */
    DENS_MAN * firstHalfAtomForces_;

    /** FE temperature change from thermostat during predictor phase in second half of timestep */
    DENS_MAT deltaEnergy1_;

    /** FE temperature change from thermostat during corrector phase in second half of timestep */
    DENS_MAT deltaEnergy2_;

    /** right-hand side data for thermostat equation */
    DENS_MAT rhs_;

    // workspace
    DENS_MAT _lambdaPowerOutput_; // power applied by lambda in output format
    DENS_MAT _velocityDelta_; // change in velocity when lambda force is applied
    DENS_VEC _lambdaOverlap_; // lambda in MD overlapping FE nodes

  private:

    // DO NOT define this
    ThermostatGlcFs();

  };

  /**
   *  @class  ThermostatSolverFlux
   *  @brief  Class enforces GLC on atomic forces based on FE power when using fractional step time integration
   */

  class ThermostatSolverFlux : public ThermostatFsSolver {

  public:

    ThermostatSolverFlux(AtomicRegulator * thermostat,
                         int lambdaMaxIterations,
                         const std::string & regulatorPrefix = "");

    virtual ~ThermostatSolverFlux() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

  protected:

    // methods
    /** sets up the transfer which is the set of nodes being regulated */
    virtual void construct_regulated_nodes();

  private:

    // DO NOT define this
    ThermostatSolverFlux();

  };

  /**
   *  @class  ThermostatIntegratorFlux
   *  @brief  Class integrates GLC on atomic forces based on FE power when using fractional step time integration
   */

  class ThermostatIntegratorFlux : public ThermostatGlcFs {

  public:

    ThermostatIntegratorFlux(AtomicRegulator * thermostat,
                             int lambdaMaxIterations,
                             const std::string & regulatorPrefix = "");

    virtual ~ThermostatIntegratorFlux() {};

    /** pre-run initialization of method data */
    virtual void initialize();

  protected:

    /** sets up appropriate rhs for thermostat equations */
    virtual void set_thermostat_rhs(DENS_MAT & rhs,
                                    double dt);

    /** add contributions from thermostat to FE energy */
    virtual void add_to_energy(const DENS_MAT & nodalLambdaPower,
                               DENS_MAT & deltaEnergy,
                               double dt);

    // data
    /** reference to ATC sources coming from prescribed data, AtC coupling, and extrinsic coupling */
    DENS_MAN & heatSource_;

  private:

    // DO NOT define this
    ThermostatIntegratorFlux();

  };

  /**
   *  @class  ThermostatSolverFixed
   *  @brief  Class enforces GLC on atomic forces based on FE power when using fractional step time integration
   */

  class ThermostatSolverFixed : public ThermostatFsSolver {

  public:

    ThermostatSolverFixed(AtomicRegulator * thermostat,
                          int lambdaMaxIterations,
                          const std::string & regulatorPrefix = "");

    virtual ~ThermostatSolverFixed() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

  protected:

    // methods
    /** sets up the transfer which is the set of nodes being regulated */
    virtual void construct_regulated_nodes();

   private:

    // DO NOT define this
    ThermostatSolverFixed();

  };

  /**
   *  @class  ThermostatIntegratorFixed
   *  @brief  Class integrates GLC on atomic forces based on FE power when using fractional step time integration
   */

  class ThermostatIntegratorFixed : public ThermostatGlcFs {

  public:

    ThermostatIntegratorFixed(AtomicRegulator * thermostat,
                              int lambdaMaxIterations,
                              const std::string & regulatorPrefix = "");

    virtual ~ThermostatIntegratorFixed() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** pre-run initialization of method data */
    virtual void initialize();

    /** applies thermostat to atoms in the predictor phase */
    virtual void apply_pre_predictor(double dt);

    /** applies thermostat to atoms in the pre-corrector phase */
    virtual void apply_pre_corrector(double dt);

    /** applies thermostat to atoms in the post-corrector phase */
    virtual void apply_post_corrector(double dt);

    /** compute boundary flux, requires thermostat input since it is part of the coupling scheme */
    virtual void compute_boundary_flux(FIELDS & /* fields */)
      {boundaryFlux_[TEMPERATURE] = 0.;};

    /** determine if local shape function matrices are needed */
    virtual bool use_local_shape_functions() const {return atomicRegulator_->use_localized_lambda();};

  protected:

    // methods
    /** initialize data for tracking the change in nodal atomic temperature */
    virtual void initialize_delta_nodal_atomic_energy(double dt);

    /** compute the change in nodal atomic temperature */
    virtual void compute_delta_nodal_atomic_energy(double dt);

    /** sets up appropriate rhs for thermostat equations */
    virtual void set_thermostat_rhs(DENS_MAT & rhs,
                                    double dt);

    /** add contributions from thermostat to FE energy */
    virtual void add_to_energy(const DENS_MAT & nodalLambdaPower,
                               DENS_MAT & deltaEnergy,
                               double dt);

    /* sets up and solves the linear system for lambda */
    virtual void compute_lambda(double dt,
                                bool iterateSolution = true);

    /** flag for halving the applied force to mitigate numerical errors */
    bool halve_force();

    // data
    /** change in FE energy over a timestep */
    DENS_MAT deltaFeEnergy_;

    /** initial FE energy used to compute change */
    DENS_MAT initialFeEnergy_;

    /** change in restricted atomic FE energy over a timestep */
    DENS_MAT deltaNodalAtomicEnergy_;

    /** initial restricted atomic FE energy used to compute change */
    DENS_MAT initialNodalAtomicEnergy_;

    /** filtered nodal atomic energy */
    DENS_MAN nodalAtomicEnergyFiltered_;

    /** forces depending on predicted velocities for correct updating with fixed nodes */
    AtomicThermostatForce * atomThermostatForcesPredVel_;

    /** coefficient to account for effect of time filtering on rhs terms */
    double filterCoefficient_;

    /** kinetic energy multiplier in total energy (used for temperature expression) */
    double keMultiplier_;

    // workspace
    DENS_MAT _tempNodalAtomicEnergyFiltered_; // stores filtered energy change in atoms for persistence during predictor

  private:

    // DO NOT define this
    ThermostatIntegratorFixed();

  };

  /**
   *  @class  ThermostatIntegratorFluxFiltered
   *  @brief  Class integrates GLC on atomic forces based on FE power when using fractional step time integration
   *          in conjunction with time filtering
   */

  class ThermostatIntegratorFluxFiltered : public ThermostatIntegratorFlux {

  public:

    ThermostatIntegratorFluxFiltered(AtomicRegulator * thermostat,
                                     int lambdaMaxIterations,
                                     const std::string & regulatorPrefix = "");

    virtual ~ThermostatIntegratorFluxFiltered() {};

    /** pre-run initialization of method data */
    virtual void initialize();

    /** applies thermostat to atoms in the post-corrector phase */
    virtual void apply_post_corrector(double dt);

    /** get data for output */
    virtual void output(OUTPUT_LIST & outputData);

  protected:

    /** sets up appropriate rhs for thermostat equations */
    virtual void set_thermostat_rhs(DENS_MAT & rhs,
                                    double dt);

    /** add contributions from thermostat to FE energy */
    virtual void add_to_energy(const DENS_MAT & nodalLambdaPower,
                               DENS_MAT & deltaEnergy,
                               double dt);

    // data
    /** heat source time history required to get instantaneous heat sources */
    DENS_MAT heatSourceOld_;
    DENS_MAT instantHeatSource_;
    DENS_MAT timeStepSource_;

  private:

    // DO NOT define this
    ThermostatIntegratorFluxFiltered();

  };

  /**
   *  @class  ThermostatIntegratorFixedFiltered
   *  @brief  Class for thermostatting using the temperature matching constraint and is compatible with
 the fractional step time-integration with time filtering
   */

  class ThermostatIntegratorFixedFiltered : public ThermostatIntegratorFixed {

  public:

    ThermostatIntegratorFixedFiltered(AtomicRegulator * thermostat,
                                      int lambdaMaxIterations,
                                      const std::string & regulatorPrefix = "");

    virtual ~ThermostatIntegratorFixedFiltered() {};

    /** get data for output */
    virtual void output(OUTPUT_LIST & outputData);

  protected:

    // methods
    /** initialize data for tracking the change in nodal atomic temperature */
    virtual void initialize_delta_nodal_atomic_energy(double dt);

    /** compute the change in nodal atomic temperature */
    virtual void compute_delta_nodal_atomic_energy(double dt);

    /** sets up appropriate rhs for thermostat equations */
    virtual void set_thermostat_rhs(DENS_MAT & rhs,
                                    double dt);

    /** add contributions from thermostat to temperature for uncoupled nodes */
    virtual void add_to_energy(const DENS_MAT & nodalLambdaPower,
                               DENS_MAT & deltaEnergy,
                               double dt);

  private:

    // DO NOT define this
    ThermostatIntegratorFixedFiltered();

  };

  /**
   *  @class  ThermostatFluxFixed
   *  @brief  Class for thermostatting using the temperature matching constraint one one set of nodes and the flux matching constraint on another
   */

  class ThermostatFluxFixed : public RegulatorMethod {

  public:

    ThermostatFluxFixed(AtomicRegulator * thermostat,
                        int lambdaMaxIterations,
                        bool constructThermostats = true);

    virtual ~ThermostatFluxFixed();

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** pre-run initialization of method data */
    virtual void initialize();

    /** applies thermostat to atoms in the predictor phase */
    virtual void apply_pre_predictor(double dt);

    /** applies thermostat to atoms in the pre-corrector phase */
    virtual void apply_pre_corrector(double dt);

    /** applies thermostat to atoms in the post-corrector phase */
    virtual void apply_post_corrector(double dt);

    /** get data for output */
    virtual void output(OUTPUT_LIST & outputData);

    /** compute boundary flux, requires thermostat input since it is part of the coupling scheme */
    virtual void compute_boundary_flux(FIELDS & fields)
      {thermostatBcs_->compute_boundary_flux(fields);};

  protected:

    // data
    /** thermostat for imposing the fluxes */
    ThermostatIntegratorFlux * thermostatFlux_;

    /** thermostat for imposing fixed nodes */
    ThermostatIntegratorFixed * thermostatFixed_;

    /** pointer to whichever thermostat should compute the flux, based on coupling method */
    ThermostatGlcFs * thermostatBcs_;

  private:

    // DO NOT define this
    ThermostatFluxFixed();
  };

  /**
   *  @class  ThermostatFluxFixedFiltered
   *  @brief  Class for thermostatting using the temperature matching constraint one one set of nodes and the flux matching constraint on another with time filtering
   */

  class ThermostatFluxFixedFiltered : public ThermostatFluxFixed {

  public:

    ThermostatFluxFixedFiltered(AtomicRegulator * thermostat,
                                int lambdaMaxIterations);

    virtual ~ThermostatFluxFixedFiltered(){};

  private:

    // DO NOT define this
    ThermostatFluxFixedFiltered();

  };

  /**
   *  @class  ThermostatGlc
   *  @brief  Class for thermostat algorithms based on Gaussian least constraints (GLC)
   */

  class ThermostatGlc : public ThermostatShapeFunction {

  public:

    ThermostatGlc(AtomicRegulator * thermostat);

    virtual ~ThermostatGlc() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

  protected:

    // methods
    /** apply forces to atoms */
    virtual void apply_to_atoms(PerAtomQuantity<double> * atomicVelocity,
                                const DENS_MAT & lambdaForce,
                                double dt);

    // member data
    /** pointer to a time filtering object */
    TimeFilter * timeFilter_;

    /** filtered lambda power */
    DENS_MAN * lambdaPowerFiltered_;

    /** atomic force induced by lambda */
    PerAtomQuantity<double> * atomThermostatForces_;

    /** pointer to access prescribed data for fixed nodes */
    PrescribedDataManager * prescribedDataMgr_;

    /** pointer to atom masses */
    FundamentalAtomQuantity * atomMasses_;

    /** workspace variables */
    DENS_MAT _velocityDelta_;

  private:

    // DO NOT define this
    ThermostatGlc();

  };

  /**
   *  @class  ThermostatPowerVerlet
   *  @brief  Class for thermostatting using the heat flux matching constraint and is compatible with
 the Gear time-integration
   */

  class ThermostatPowerVerlet : public ThermostatGlc {

  public:

    ThermostatPowerVerlet(AtomicRegulator * thermostat);

    virtual ~ThermostatPowerVerlet() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** pre-run initialization of method data */
    virtual void initialize();

    /** applies thermostat to atoms in the predictor phase */
    virtual void apply_pre_predictor(double dt);

    /** applies thermostat to atoms in the pre-corrector phase */
    virtual void apply_pre_corrector(double dt);

    /** get data for output */
    virtual void output(OUTPUT_LIST & outputData);

    /** final tasks of a run */
    virtual void finish();

    /** determine if local shape function matrices are needed */
    virtual bool use_local_shape_functions() const {return (!(atomicRegulator_->use_lumped_lambda_solve()) && atomicRegulator_->use_localized_lambda());};

  protected:

    /** nodal temperature rate of change */
    DENS_MAN & nodalTemperatureRoc_;

    /** reference to ATC sources coming from prescribed data, AtC coupling, and extrinsic coupling */
    DENS_MAN & heatSource_;

    /** pointer to nodal atomic power */
    DENS_MAN * nodalAtomicPower_;

    /** power applied to each atom by lambda force */
    AtfShapeFunctionRestriction * nodalAtomicLambdaPower_;

    /** workspace variables */
    DENS_MAT _rhs_;

    /** sets up and solves thermostat equations */
    virtual void compute_thermostat(double dt);

    /** sets up appropriate rhs for thermostat equations */
    virtual void set_thermostat_rhs(DENS_MAT & rhs_nodes);

    /** add contributions (if any) to the finite element right-hand side */
    virtual void add_to_rhs(FIELDS & rhs);

    // workspace
    DENS_MAT _nodalAtomicLambdaPowerOut_; // power induced by lambda in output format

  private:

    // DO NOT define this
    ThermostatPowerVerlet();

  };

  /**
   *  @class  ThermostatHooverVerlet
   *  @brief  Classfor thermostatting using the temperature matching constraint and is compatible with
 Gear time-integration
   */

  class ThermostatHooverVerlet : public ThermostatPowerVerlet {

  public:

    ThermostatHooverVerlet(AtomicRegulator * thermostat);

    virtual ~ThermostatHooverVerlet() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** final tasks of a run */
    virtual void finish() {};

    /** compute boundary flux, requires thermostat input since it is part of the coupling scheme */
    virtual void compute_boundary_flux(FIELDS & /* fields */)
      {boundaryFlux_[TEMPERATURE] = 0.;};

  protected:

    /** lambda coupling parameter for hoover thermostat */
    DENS_MAN * lambdaHoover_;

    /** workspace variables */
    DENS_MAT _myNodalLambdaPower_;

    /** sets up and solves thermostat equations */
    virtual void compute_thermostat(double dt);

    /** sets up Hoover component of the thermostat */
    void set_hoover_rhs(DENS_MAT & rhs);

    /** add Hoover contributions to lambda power */
    void add_to_lambda_power(const DENS_MAT & myLambdaForce,
                             double dt);

    /** power applied to each atom by hoover lambda force */
    AtfShapeFunctionRestriction * nodalAtomicHooverLambdaPower_;

    /** add contributions (if any) to the finite element right-hand side */
    virtual void add_to_rhs(FIELDS & rhs);

  private:

    // DO NOT implement this
    ThermostatHooverVerlet();

  };

  /**
   *  @class  ThermostatPowerVerletFiltered
   *  @brief  Class for thermostatting using the heat flux matching constraint and is compatible with
 Gear time-integration with time filtering
   */

  class ThermostatPowerVerletFiltered : public ThermostatPowerVerlet {

  public:

    ThermostatPowerVerletFiltered(AtomicRegulator * thermostat);

    virtual ~ThermostatPowerVerletFiltered(){};

    /** get data for output */
    virtual void output(OUTPUT_LIST & outputData);

    /** compute boundary flux, requires thermostat input since it is part of the coupling scheme */
    virtual void compute_boundary_flux(FIELDS & fields);

  protected:

    /** sets up appropriate rhs for thermostat equations */
    virtual void set_thermostat_rhs(DENS_MAT & rhs_nodes);

    /** add contributions (if any) to the finite element right-hand side */
    virtual void add_to_rhs(FIELDS & rhs);

    /** nodal temperature 2nd rate of change (i.e. second time derivative) */
    DENS_MAN & nodalTemperature2Roc_;

    /** reference to ATC rate of change sources coming from prescribed data, AtC coupling, and extrinsic coupling */
    DENS_MAN heatSourceRoc_;

    /** references to ATC field rates of changing for inverting the filtered heat sources */
    FIELDS & fieldsRoc_;

    /** flux rate of changes for inverting filtered fluxes */
    FIELDS fluxRoc_;

    /** time scale for the time filter */
    double filterScale_;

  private:

    // DO NOT define this
    ThermostatPowerVerletFiltered();

  };

  /**
   *  @class  ThermostatHooverVerletFiltered
   *  @brief  Class for thermostatting using the temperature matching constraint and is compatible with
 Gear time-integration with time filtering
   */

  class ThermostatHooverVerletFiltered : public ThermostatPowerVerletFiltered {

  public:

    ThermostatHooverVerletFiltered(AtomicRegulator * thermostat);

    virtual ~ThermostatHooverVerletFiltered() {};

    /** instantiate all needed data */
    virtual void construct_transfers();

    /** final tasks of a run */
    virtual void finish() {};

    /** compute boundary flux, requires thermostat input since it is part of the coupling scheme */
    virtual void compute_boundary_flux(FIELDS & /* fields */)
      {boundaryFlux_[TEMPERATURE] = 0.;};

  protected:

    /** lambda coupling parameter for hoover thermostat */
    DENS_MAN * lambdaHoover_;

    /** workspace variables */
    DENS_MAT _myNodalLambdaPower_;

    /** sets up and solves thermostat equations */
    virtual void compute_thermostat(double dt);

    /** sets up Hoover component of the thermostat */
    void set_hoover_rhs(DENS_MAT & rhs);

    /** add Hoover contributions to lambda power */
    void add_to_lambda_power(const DENS_MAT & myLambdaForce,
                             double dt);

    /** power applied to each atom by hoover lambda force */
    DENS_MAN * nodalAtomicHooverLambdaPower_;

    /** add contributions (if any) to the finite element right-hand side */
    virtual void add_to_rhs(FIELDS & rhs);

  private:

    // DO NOT implement this
    ThermostatHooverVerletFiltered();

  };

};

#endif