File: Stress.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 (212 lines) | stat: -rw-r--r-- 7,627 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
#ifndef STRESS_H
#define STRESS_H

#include <map>
#include <string>
#include <vector>
#include <fstream>
#include "MatrixLibrary.h"
#include "ATC_TypeDefs.h"
#include "ATC_TypeDefs.h"
#include "NonLinearSolver.h"

namespace ATC {
  enum ElasticityTensorType {FIRST_ELASTICITY_TENSOR=0, SECOND_ELASTICITY_TENSOR};
  /**
   * @class Stress
   * @brief Base class that defines interface for a constitutive law
   * @brief that computes stress given all field and gradient information.
   */
  class Stress
  {
    public:
      Stress()  {};
      virtual ~Stress() {};
      virtual void initialize(void){};
      //* Returns parameter values, (Nothing uses this).
      virtual void parameters(std::map<std::string,double> & /* parameters */) {}
      //* Computes stress given a displacement gradient.
      //* Units: mvv/L^3 (i.e. for units Real: g/(mol ps^2 A^2) )
      virtual void stress(const FIELD_MATS &fields,
                          const GRAD_FIELD_MATS &gradFields,
                          DENS_MAT_VEC &stress)=0;
      //* Computes free (T>0)/potential(T=0) energy density
      //* Units: mvv/L^3 (i.e. for units Real: g/(mol ps^2 A^2) )
      virtual void elastic_energy(const FIELD_MATS &fields,
                                  const GRAD_FIELD_MATS &gradFields,
                                  DENS_MAT &energy) const;
      //* Returns the material tangent at a given deformation gradient.
      virtual void tangent(const MATRIX & /* F */, MATRIX & /* C */) const
        {throw ATC_Error("Stress::tangent: unimplemented function");}
  };


  /**
   *  @class  StressCubicElastic
   *  @brief  Class for computing stress for a cubic elastic material
   */

  class StressCubicElastic : public Stress
  {
    public:
      StressCubicElastic():c11_(0),c12_(0),c44_(0){};
      StressCubicElastic(std::fstream &matfile);
      StressCubicElastic(double c11, double c12, double c44)
        : c11_(c11), c12_(c12), c44_(c44) { }
      void stress(const FIELD_MATS &fields,
                  const GRAD_FIELD_MATS &gradFields,
                  DENS_MAT_VEC &flux);
      virtual void elastic_energy(const FIELD_MATS &fields,
                          const GRAD_FIELD_MATS &gradFields,
                          DENS_MAT &energy) const;
      virtual void tangent(const MATRIX & /* F */, MATRIX &C) const {C=C_;}
    protected:
      double c11_, c12_, c44_;
      DENS_MAT C_;
      void set_tangent();
  };

  /**
   *  @class  StressCubicElasticDamped
   *  @brief  Class for computing stress for a cubic elastic material w/ damping
   */

  class StressCubicElasticDamped : public StressCubicElastic
  {
    public:
      StressCubicElasticDamped(std::fstream &matfile);
      StressCubicElasticDamped(double c11, double c12, double c44, double gamma)
        : StressCubicElastic(c11,c12,c44), gamma_(gamma) { }
      void stress(const FIELD_MATS &fields,
                  const GRAD_FIELD_MATS &gradFields,
                  DENS_MAT_VEC &flux);
    protected:
      double gamma_;
  };

  /**
   *  @class  StressLinearElastic
   *  @brief  Class for computing stress for a linear elastic material
   */

  class StressLinearElastic : public StressCubicElastic
  {
    public:
      StressLinearElastic(std::fstream &matfile);
      void stress(const FIELD_MATS &fields,
                  const GRAD_FIELD_MATS &gradFields,
                  DENS_MAT_VEC &flux);
    protected:
      double E_, nu_;
      double mu_, lambda_;
  };

  // forward declarations needed by StressCauchyBorn
  class CbPotential;
  class CBLattice;

  /**
   * Structure of lattice properties needed by StressCauchyBorn.
   */
  struct CbData {
    double e2mvv;           //*> Energy conversion factor (1/mvv2e).
    double boltzmann;       //*> Boltzmann constant (in LAMMPS units)
    double hbar;            //*> Planck's constant (in LAMMPS units)
    double inv_atom_volume; //*> Volume of atom.
    double atom_mass;       //*> Mass of an atom.
    DENS_MAT cell_vectors;  //*> Unit vectors for lattice cells.
    DENS_MAT basis_vectors; //*> Positions of atoms within a lattice cell.
  };

  /**
   * @class StressCauchyBorn
   * @brief Class for computing the stress and elastic constants for a
   * @brief Cauchy-Born material.
   */


  class StressCauchyBorn : public Stress
  {
    public:
      StressCauchyBorn(std::fstream &matfile, CbData &cb);
      virtual ~StressCauchyBorn();
      virtual void initialize(void);
      //* Returns the stress computed from a 0K Cauchy-Born approxmation.
      virtual void stress(const FIELD_MATS &fields, const GRAD_FIELD_MATS &gradFields,
                  DENS_MAT_VEC &flux);
      //* Computes free (T>0)/potential(T=0) energy density
      //* Units: mvv/L^3 (i.e. for units Real: g/(mol ps^2 A^2) )
      virtual void elastic_energy(const FIELD_MATS &fields,
                          const GRAD_FIELD_MATS &gradFields,
                          DENS_MAT &energy) const;
      //* Computes entropic energy density
      void entropic_energy(const FIELD_MATS &fields, const GRAD_FIELD_MATS &gradFields,
                          DENS_MAT &energy) const;
      //* Returns the material tangent at a given deformation gradient.
      virtual void tangent(const MATRIX &F, MATRIX &C) const;
      double stiffness() const;
      //* Creates a linearization for a deformation gradient.
      DENS_VEC elasticity_tensor(const VECTOR &Fv, MATRIX &C, const ElasticityTensorType type=FIRST_ELASTICITY_TENSOR) const;
      DENS_VEC elasticity_tensor(const MATRIX &F, MATRIX &C, const ElasticityTensorType type=FIRST_ELASTICITY_TENSOR) const;
    protected:
      void linearize(MATRIX *F=nullptr);
      CBLattice   *cblattice_;       //*> CbLattice -> makes atom clusters.
      CbPotential *potential_;       //*> CbPotential -> interatomic forces.
      bool makeLinear_;
      StressCubicElastic *cubicMat_; //*> Stores optional linear elastic law.
      bool initialized_;
      double fixed_temperature_;     //*> Specifies a uniform temperature.
      CbData cbdata_;                //*> Lattice & atom volume/mass.
  };


    // adaptor to NonLinearSolver
    class CBElasticTangentOperator : public TangentOperator {
      public:
        CBElasticTangentOperator (StressCauchyBorn * cauchyBornStress,
                                  DENS_VEC & targetP) :
          TangentOperator(),
          cauchyBornStress_(cauchyBornStress),
          targetP_(targetP) {};
        void function(const VECTOR & F, DENS_VEC & R)
        {
          DENS_MAT B;
          tangent(F,R,B);
        }
        void tangent(const VECTOR & F, DENS_VEC & R, MATRIX & B)
        {
          cbP_ = cauchyBornStress_->elasticity_tensor(F, B);

          R = cbP_ - targetP_;
        }
      private:
        StressCauchyBorn * cauchyBornStress_;
        DENS_VEC targetP_, cbP_;
    };

    // adaptor to NonLinearSolver
    class CB2ndElasticTangentOperator : public TangentOperator {
      public:
        CB2ndElasticTangentOperator (StressCauchyBorn * cauchyBornStress,
                                  DENS_VEC & targetS) :
          TangentOperator(),
          cauchyBornStress_(cauchyBornStress),
          targetS_(targetS) {};
        void function(const VECTOR & U, DENS_VEC & r)
        {
          DENS_MAT C;
          tangent(U,r,C);
        }
        void tangent(const VECTOR & U, DENS_VEC & r, MATRIX & C)
        {
          cbS_ = cauchyBornStress_->elasticity_tensor(U, C, SECOND_ELASTICITY_TENSOR);

          r = cbS_ - targetS_;
        }
      private:
        StressCauchyBorn * cauchyBornStress_;
        DENS_VEC targetS_, cbS_;
    };

}
#endif