Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #pragma once
11 :
12 : #include "RadialReturnCreepStressUpdateBase.h"
13 : #include "nlohmann/json.h"
14 :
15 : template <bool is_ad>
16 : class LAROMANCEStressUpdateBaseTempl : public RadialReturnCreepStressUpdateBaseTempl<is_ad>
17 : {
18 :
19 : public:
20 : static InputParameters validParams();
21 :
22 : LAROMANCEStressUpdateBaseTempl(const InputParameters & parameters);
23 :
24 : virtual void resetIncrementalMaterialProperties() override;
25 : virtual void
26 : storeIncrementalMaterialProperties(const unsigned int total_number_substeps) override;
27 :
28 : protected:
29 : virtual void exportJSON();
30 :
31 : virtual bool substeppingCapabilityEnabled() override;
32 :
33 : enum class ROMInputTransform
34 : {
35 : LINEAR,
36 : LOG,
37 : EXP
38 : };
39 :
40 : virtual void initialSetup() override;
41 :
42 : // Setup unit conversion factors. The required units in the ROM are:
43 : // Cell dislocation density: m^-2
44 : // Wall dislocation density: m^-2
45 : // MX phase fracture: nondim.
46 : // stress: MPa
47 : // strain: nondim.
48 : // temperature: K
49 : virtual void setupUnitConversionFactors(const InputParameters & parameters);
50 :
51 : virtual void initQpStatefulProperties() override;
52 : virtual void
53 : computeStressInitialize(const GenericReal<is_ad> & effective_trial_stress,
54 : const GenericRankFourTensor<is_ad> & elasticity_tensor) override;
55 : virtual GenericReal<is_ad> computeResidual(const GenericReal<is_ad> & effective_trial_stress,
56 : const GenericReal<is_ad> & scalar) override;
57 :
58 : virtual GenericReal<is_ad>
59 341396 : computeDerivative(const GenericReal<is_ad> & /*effective_trial_stress*/,
60 : const GenericReal<is_ad> & /*scalar*/) override
61 : {
62 341396 : return _derivative;
63 : }
64 :
65 : virtual void
66 : computeStressFinalize(const GenericRankTwoTensor<is_ad> & plastic_strain_increment) override;
67 : virtual GenericReal<is_ad>
68 : maximumPermissibleValue(const GenericReal<is_ad> & effective_trial_stress) const override;
69 : virtual Real computeTimeStepLimit() override;
70 :
71 : void outputIterationSummary(std::stringstream * iter_output,
72 : const unsigned int total_it) override;
73 :
74 : virtual void outputIterationStep(std::stringstream * iter_output,
75 : const GenericReal<is_ad> & effective_trial_stress,
76 : const GenericReal<is_ad> & scalar,
77 : const Real reference_residual) override;
78 :
79 : /// Enum to error, warn, ignore, or extrapolate if input is outside of window of applicability
80 : enum class WindowFailure
81 : {
82 : ERROR,
83 : EXCEPTION,
84 : WARN,
85 : IGNORE,
86 : DONOTHING,
87 : USELIMIT,
88 : EXTRAPOLATE
89 : };
90 :
91 : /**
92 : * Precompute the ROM strain rate information for all inputs except for strain. Strain will be
93 : * computed in the radial return algorithm several times, while the remainder of the inputs remain
94 : * constant.
95 : * @param out_index Output index
96 : */
97 : void precomputeROM(const unsigned out_index);
98 :
99 : /**
100 : * Computes the ROM calculated increment for a given output and tile.
101 : * @param tile Tile index
102 : * @param partition Partition index
103 : * @param out_index Output index
104 : * @param derivative Optional flag to return derivative of ROM increment with respect to stress.
105 : * @return ROM computed increment
106 : */
107 : GenericReal<is_ad> computeROM(const unsigned int tile,
108 : const unsigned int partition,
109 : const unsigned out_index,
110 : const bool derivative = false);
111 :
112 : /**
113 : * Checks if the input combination is in a specific tile
114 : * @param p Partition index
115 : * @param t Tile index
116 : * @return bool if in tile
117 : */
118 : bool checkInTile(const unsigned int p, const unsigned int t) const;
119 :
120 : /**
121 : * Checks if two tile domains are equal
122 : * @param p Partition index
123 : * @param t Tile 1 index
124 : * @param tt Tile 2 index
125 : * @param in_index input index
126 : * @return integer 0 (false) or 1 (true)
127 : */
128 : bool areTilesNotIdentical(const unsigned int p,
129 : const unsigned int t,
130 : const unsigned int tt,
131 : const unsigned int in_index);
132 :
133 : /**
134 : * Convert the input variables into the form expected by the ROM Legendre polynomials to have a
135 : * normalized space from [-1, 1] so that every variable has equal weight
136 : * @param input Input value
137 : * @param transform ROMInputTransform enum indicating how the input is to be transformed
138 : * @param transform_coef Transform coefficient for the given input
139 : * @param transformed_limits Transformed limits for the given input
140 : * @param derivative Optional flag to return derivative of converted input with respect to stress.
141 : * @return Converted input
142 : */
143 : GenericReal<is_ad> normalizeInput(const GenericReal<is_ad> & input,
144 : const ROMInputTransform transform,
145 : const Real transform_coef,
146 : const std::vector<Real> & transformed_limits,
147 : const bool derivative = false);
148 :
149 : /**
150 : * Assemble the array of Legendre polynomials to be multiplied by the ROM coefficients
151 : * @param p Partition index
152 : * @param rom_input ROM input
153 : * @param polynomial_inputs Vector of transformed Legendre polynomials
154 : * @param drom_input Optional derivative of ROM input with respect to stress
155 : * @param derivative Optional flag to return derivative of converted input with respect to stress.
156 : */
157 : void buildPolynomials(const unsigned int p,
158 : const GenericReal<is_ad> & rom_input,
159 : std::vector<GenericReal<is_ad>> & polynomial_inputs,
160 2182784 : const GenericReal<is_ad> & drom_input = 0,
161 : const bool derivative = false);
162 :
163 : /**
164 : * Arranges the calculated Legendre polynomials into the proper oder and multiplies the Legendre
165 : * polynomials by the ROM coefficients to compute the predicted output values. This method works
166 : * with all inputs besides stress, while stress is handled by computeValues
167 : * @param p Partition index
168 : * @param coefs Legendre polynomial coefficients
169 : * @param polynomial_inputs Vector of transformed Legendre polynomial
170 : * @param precomputed Vector that holds the precomputed ROM values
171 : */
172 : void precomputeValues(const unsigned int p,
173 : const std::vector<Real> & coefs,
174 : const std::vector<std::vector<GenericReal<is_ad>>> & polynomial_inputs,
175 : std::vector<GenericReal<is_ad>> & precomputed);
176 :
177 : /**
178 : * Arranges the calculated Legendre polynomials into the proper oder and multiplies the Legendre
179 : * polynomials by the ROM coefficients to compute the predicted output values. This method only
180 : * manipulates the stress input, with all others handled in precomputeValues
181 : * @param p Partition index
182 : * @param precomputed Precomputed multiplication of polynomials
183 : * @param polynomial_inputs Vector of Legendre polynomial transformation
184 : * @param dpolynomial_inputs Vector of derivative of Legendre polynomial transformation with
185 : * respect to stress
186 : * @param derivative Optional flag to return derivative of converted computed values with respect
187 : * to stress.
188 : * @return ROM output
189 : */
190 : GenericReal<is_ad>
191 : computeValues(const unsigned int p,
192 : const std::vector<GenericReal<is_ad>> & precomputed,
193 : const std::vector<std::vector<GenericReal<is_ad>>> & polynomial_inputs,
194 : const std::vector<GenericReal<is_ad>> & dpolynomial_inputs = {},
195 : const bool derivative = false);
196 :
197 : /**
198 : * Computes the output variable increments from the ROM predictions by bringing out of the
199 : * normalized map to the actual physical values
200 : * @param p Partition index
201 : * @param old_input_values Previous timestep values of ROM inputs
202 : * @param rom_output Outputs from ROM
203 : * @param out_index Output index
204 : * @param drom_output Derivative of output with respect to stress
205 : * @param derivative Optional flag to return derivative of output with respect to stress.
206 : * @return Converted ROM output
207 : */
208 : virtual GenericReal<is_ad> convertOutput(const unsigned int p,
209 : const std::vector<Real> & old_input_values,
210 : const GenericReal<is_ad> & rom_output,
211 : const unsigned out_index,
212 872563 : const GenericReal<is_ad> & drom_output = 0.0,
213 : const bool derivative = false);
214 :
215 : /**
216 : * Calculate the value or derivative of Legendre polynomial up to 3rd order
217 : * @param value Input to Legendre polynomial
218 : * @param degree Degree of Legendre polynomial
219 : * @param derivative Optional flag to return derivative of Legendre polynomial Legendre
220 : * @return Computed value from Legendre polynomial
221 : */
222 : GenericReal<is_ad> computePolynomial(const GenericReal<is_ad> & value,
223 : const unsigned int degree,
224 : const bool derivative = false);
225 :
226 : /**
227 : * Calculate the sigmoid function weighting for the input based on the limits
228 : * @param lower Lower limit
229 : * @param upper Upper limit
230 : * @param val Input value
231 : * @param derivative Optional flag to return derivative of the sigmoid w.r.t. the input
232 : * @return weight
233 : */
234 : GenericReal<is_ad> sigmoid(const Real lower,
235 : const Real upper,
236 : const GenericReal<is_ad> & val,
237 : const bool derivative = false);
238 :
239 : /**
240 : * Compute the contribution (weight) of each tile in each partition,
241 : * based on the input and tile boundaries (in terms of input domain).
242 : * @param weights Weights for each tile
243 : * @param input Input value
244 : * @param in_index Input index
245 : * @param derivative Optional flag to return derivative of the sigmoid w.r.t. the input
246 : */
247 : void computeTileWeight(std::vector<std::vector<GenericReal<is_ad>>> & weights,
248 : GenericReal<is_ad> & input,
249 : const unsigned int in_index,
250 : const bool derivative = false);
251 :
252 : /**
253 : * Compute the weight of the different partitions
254 : * @param weights Weights for each partition
255 : * @param derivative Optional flag to return derivative of the sigmoid w.r.t. the input
256 : */
257 286664 : virtual void computePartitionWeights(std::vector<GenericReal<is_ad>> & weights,
258 : std::vector<GenericReal<is_ad>> & dweights_dstress)
259 : {
260 286664 : if (_num_partitions != 1)
261 0 : mooseError("Internal error: If number of partitions is not one, then computePartitionWeights "
262 : "must be defined");
263 286664 : weights[0] = 1.0;
264 286664 : dweights_dstress[0] = 0.0;
265 286664 : }
266 :
267 : /**
268 : * Convert input based on the transform type
269 : * @param x Input value
270 : * @param transform Enum declaring transform to be performed
271 : * @param coef Coefficient applied during transformation
272 : * @param derivative Optional flag to return derivative of the sigmoid w.r.t. the input
273 : */
274 : template <typename T>
275 4263193 : void convertValue(T & x,
276 : const ROMInputTransform transform,
277 : const Real coef,
278 : const bool derivative = false)
279 : {
280 : using std::exp, std::log;
281 :
282 4263193 : if (transform == ROMInputTransform::EXP)
283 : {
284 : mooseAssert(coef != 0, "Coefficient must not be zero.");
285 178476 : if (derivative)
286 0 : x = exp(x / coef) / coef;
287 : else
288 178476 : x = exp(x / coef);
289 : }
290 4084717 : else if (transform == ROMInputTransform::LOG)
291 : {
292 : mooseAssert(x + coef > 0, "Sum must be greater than 0.");
293 2478312 : if (derivative)
294 343936 : x = 1.0 / (x + coef);
295 : else
296 2199256 : x = log(x + coef);
297 : }
298 1606405 : else if (transform == ROMInputTransform::LINEAR)
299 : {
300 : mooseAssert(coef == 0.0, "Coefficient cannot be supplied with linear transformation");
301 1606405 : if (derivative)
302 29475 : x = 1.0;
303 : }
304 4263193 : }
305 :
306 : /*
307 : * Calculates and returns vector utilized in assign values
308 : * @param p Partition index
309 : * @return Vector that preallocates indexing calculations for polynomial calculation
310 : */
311 : std::vector<unsigned int> getMakeFrameHelper(const unsigned int p) const;
312 :
313 : /*
314 : * Calculates and returns the transformed limits for the ROM calculations
315 : * Indexes are [partition][tile][output][input].
316 : * Inputs ordering is
317 : * input[0]: cell_old
318 : * input[1]: wall_old
319 : * input[2]: trial stress,
320 : * input[3]: effective strain old,
321 : * input[4]: temperature
322 : * input[5]: environmental factor (optional)
323 : * output ordering is:
324 : * output[0]: cell dislocations increment
325 : * output[1]: wall dislocations increment
326 : * output[2]: strain increment
327 : * @param p Partition index
328 : * @param limits Human readable limits
329 : * @return Multi-dimentional vector of transformed limits
330 : */
331 : std::vector<std::vector<std::vector<std::vector<Real>>>>
332 : getTransformedLimits(const unsigned int p,
333 : const std::vector<std::vector<std::vector<Real>>> limits);
334 :
335 : /*
336 : * Returns vector of the functions to use for the conversion of input variables.
337 : * Indexes are [partition][tile][output][input].
338 : * Inputs ordering is
339 : * input[0]: cell_old
340 : * input[1]: wall_old
341 : * input[2]: trial stress,
342 : * input[3]: effective strain old,
343 : * input[4]: temperature
344 : * input[5]: environmental factor (optional)
345 : * output ordering is:
346 : * output[0]: cell dislocations increment
347 : * output[1]: wall dislocations increment
348 : * output[2]: strain increment
349 : * @return vector of the functions to use for the conversion of input variables.
350 : */
351 63 : virtual std::vector<std::vector<std::vector<std::vector<ROMInputTransform>>>> getTransform()
352 : {
353 63 : checkJSONKey("transform");
354 : return _json["transform"]
355 126 : .template get<std::vector<std::vector<std::vector<std::vector<ROMInputTransform>>>>>();
356 : }
357 :
358 : /*
359 : * Returns factors for the functions for the conversion functions given in getTransform
360 : * Indexes are [partition][tile][output][input].
361 : * Inputs ordering is
362 : * input[0]: cell_old
363 : * input[1]: wall_old
364 : * input[2]: trial stress,
365 : * input[3]: effective strain old,
366 : * input[4]: temperature
367 : * input[5]: environmental factor (optional)
368 : * output ordering is:
369 : * output[0]: cell dislocations increment
370 : * output[1]: wall dislocations increment
371 : * output[2]: strain increment
372 : * @return factors for the functions for the conversion functions given in getTransform
373 : */
374 63 : virtual std::vector<std::vector<std::vector<std::vector<Real>>>> getTransformCoefs()
375 : {
376 63 : checkJSONKey("transform_coefs");
377 : return _json["transform_coefs"]
378 126 : .template get<std::vector<std::vector<std::vector<std::vector<Real>>>>>();
379 : }
380 :
381 : /* Optional method that returns human-readable limits used for normalization. Default is to just
382 : * use the input limits.
383 : * Indexes are [partition][tile][input][upper/lower].
384 : * Inputs ordering is
385 : * input[0]: cell_old
386 : * input[1]: wall_old
387 : * input[2]: trial stress,
388 : * input[3]: effective strain old,
389 : * input[4]: temperature
390 : * input[5]: environmental factor (optional)
391 : * @return human-readable limits for the normalization limits
392 : */
393 432 : virtual std::vector<std::vector<std::vector<std::vector<Real>>>> getNormalizationLimits()
394 : {
395 432 : if (_json.contains("normalization_limits"))
396 : return _json["normalization_limits"]
397 63 : .template get<std::vector<std::vector<std::vector<std::vector<Real>>>>>();
398 :
399 369 : return getInputLimits();
400 : }
401 :
402 : /* Returns human-readable limits for the inputs.
403 : * Indexes are [partition][tile][input][upper/lower].
404 : * Inputs ordering is
405 : * input[0]: cell_old
406 : * input[1]: wall_old
407 : * input[2]: trial stress,
408 : * input[3]: effective strain old,
409 : * input[4]: temperature
410 : * input[5]: environmental factor (optional)
411 : * @return human-readable limits for the input limits
412 : */
413 63 : virtual std::vector<std::vector<std::vector<std::vector<Real>>>> getInputLimits()
414 : {
415 63 : checkJSONKey("input_limits");
416 : return _json["input_limits"]
417 126 : .template get<std::vector<std::vector<std::vector<std::vector<Real>>>>>();
418 : }
419 :
420 : /*
421 : * Material specific coefficients multiplied by the Legendre polynomials for each of the input
422 : * variables
423 : * @return Legendre polynomial coefficients
424 : */
425 63 : virtual std::vector<std::vector<std::vector<std::vector<Real>>>> getCoefs()
426 : {
427 63 : checkJSONKey("coefs");
428 126 : return _json["coefs"].template get<std::vector<std::vector<std::vector<std::vector<Real>>>>>();
429 : }
430 :
431 : /*
432 : * Material specific orientations of tiling
433 : * variables. Indexing is partition, then input
434 : * @return Vector of a vector declaring tiling orientation
435 : */
436 432 : virtual std::vector<std::vector<unsigned int>> getTilings()
437 : {
438 432 : if (_json.contains("tiling"))
439 63 : return _json["tiling"].template get<std::vector<std::vector<unsigned int>>>();
440 :
441 369 : if (_environmental)
442 0 : return {{1, 1, 1, 1, 1, 1}};
443 738 : return {{1, 1, 1, 1, 1}};
444 738 : };
445 :
446 : /*
447 : * Minimum strain value allowed by the ROM. This is material specific, and needs to be overwritten
448 : * by individual roms and each partition
449 : * @return Vector of material specific ROM low strain value for each partition
450 : */
451 63 : virtual std::vector<Real> getStrainCutoff()
452 : {
453 63 : checkJSONKey("cutoff");
454 126 : return _json["cutoff"].template get<std::vector<Real>>();
455 : }
456 :
457 : /// Coupled temperature variable
458 : const GenericVariableValue<is_ad> & _temperature;
459 :
460 : /// Optionally coupled environmental factor
461 : const GenericMaterialProperty<Real, is_ad> * _environmental;
462 :
463 : /*
464 : * Vector of vectors WindowFailure enum that informs how to handle input that is outside of the
465 : * limits. Shape is number of inputs by 2 (lower and upper window enum)
466 : */
467 : std::vector<std::pair<WindowFailure, WindowFailure>> _window_failure;
468 :
469 : /// Flag to output verbose infromation
470 : const bool _verbose;
471 :
472 : ///@{Material properties for cell (glissile) dislocation densities (1/m^2)
473 : GenericMaterialProperty<Real, is_ad> & _cell_dislocations;
474 : const MaterialProperty<Real> & _cell_dislocations_old;
475 : ///@}
476 :
477 : /// Maximum cell dislocation increment
478 : const Real _max_cell_increment;
479 :
480 : /// Optional cell dislocation forcing function
481 : const Function * const _cell_function;
482 :
483 : /// Container for cell dislocation increment
484 : GenericReal<is_ad> _cell_dislocation_increment;
485 :
486 : ///@{Material properties for wall (locked) dislocation densities (1/m^2)
487 : GenericMaterialProperty<Real, is_ad> & _wall_dislocations;
488 : const MaterialProperty<Real> & _wall_dislocations_old;
489 : ///@}
490 :
491 : /// Maximum wall dislocation increment
492 : const Real _max_wall_increment;
493 :
494 : /// Optional wall dislocation forcing function
495 : const Function * const _wall_function;
496 :
497 : /// Optiontal effective stress forcing function
498 : const Function * const _stress_function;
499 :
500 : /// Container for wall dislocation increment
501 : GenericReal<is_ad> _wall_dislocation_increment;
502 :
503 : /// Index corresponding to the position for the dislocations with in the cell in the input vector
504 : const unsigned int _cell_input_index;
505 :
506 : /// Index corresponding to the position for the dislocations within the cell wall in the input vector
507 : const unsigned int _wall_input_index;
508 :
509 : /// Index corresponding to the position for the stress in the input vector
510 : const unsigned int _stress_input_index;
511 :
512 : /// Index corresponding to the position for the old strain in the input vector
513 : const unsigned int _old_strain_input_index;
514 :
515 : /// Index corresponding to the position for the tempeature in the input vector
516 : const unsigned int _temperature_input_index;
517 :
518 : /// Index corresponding to the position for the environmental factor in the input vector
519 : const unsigned int _environmental_input_index;
520 :
521 : /// Index corresponding to the position for cell dislocations increment in the output vector
522 : const unsigned int _cell_output_index;
523 :
524 : /// Index corresponding to the position for wall dislocations increment in the output vector
525 : const unsigned int _wall_output_index;
526 :
527 : /// Index corresponding to the position for strain increment in the output vector
528 : const unsigned int _strain_output_index;
529 :
530 : /// Optional old creep strain forcing function
531 : const Function * const _creep_strain_old_forcing_function;
532 :
533 : /// Number of partitions
534 : unsigned int _num_partitions;
535 :
536 : /// Number of ROM tiles per partition
537 : std::vector<unsigned int> _num_tiles;
538 :
539 : /// Number of inputs for the ROM data set
540 : unsigned int _num_inputs;
541 :
542 : /// Number of inputs to the ROM data set
543 : unsigned int _num_outputs;
544 :
545 : /// Legendre polynomial degree for the ROM data set for each partition
546 : std::vector<unsigned int> _degree;
547 :
548 : /// Total number of Legendre polynomial coefficients for the ROM data set in each parition
549 : std::vector<unsigned int> _num_coefs;
550 :
551 : /// Transform rules defined by the ROM data set for each partition
552 : std::vector<std::vector<std::vector<std::vector<ROMInputTransform>>>> _transform;
553 :
554 : /// Transform coefficients defined by the ROM data set for each partition
555 : std::vector<std::vector<std::vector<std::vector<Real>>>> _transform_coefs;
556 :
557 : /// Input limits defined by the ROM data set for each partition
558 : std::vector<std::vector<std::vector<std::vector<Real>>>> _input_limits;
559 :
560 : /// Normalization limits defined by the ROM data set for each partition
561 : std::vector<std::vector<std::vector<std::vector<Real>>>> _normalization_limits;
562 :
563 : /// Coefficients used with Legendre polynomials defined by the ROM data set for each partition
564 : std::vector<std::vector<std::vector<std::vector<Real>>>> _coefs;
565 :
566 : /// Limits transformed from readabile input to ROM readable limits for normalization
567 : std::vector<std::vector<std::vector<std::vector<std::vector<Real>>>>>
568 : _transformed_normalization_limits;
569 :
570 : /// Helper container defined by the ROM data set
571 : std::vector<std::vector<unsigned int>> _makeframe_helper;
572 :
573 : /// Creep rate material property
574 : GenericMaterialProperty<Real, is_ad> & _creep_rate;
575 :
576 : /// Cell dislocations rate of change
577 : GenericMaterialProperty<Real, is_ad> & _cell_rate;
578 :
579 : /// Wall dislocations rate of change
580 : GenericMaterialProperty<Real, is_ad> & _wall_rate;
581 :
582 : /// Material property to hold smootherstep applied in order to extrapolate.
583 : MaterialProperty<Real> & _extrapolation;
584 :
585 : /// Material property to store partition weight.
586 : GenericMaterialProperty<Real, is_ad> & _second_partition_weight;
587 :
588 : /// Container for derivative of creep increment with respect to strain
589 : GenericReal<is_ad> _derivative;
590 :
591 : /// Container for input values
592 : std::vector<GenericReal<is_ad>> _input_values;
593 :
594 : /// Container for old input values
595 : std::vector<Real> _old_input_values;
596 :
597 : /// Container for converted rom_inputs
598 : std::vector<std::vector<std::vector<GenericReal<is_ad>>>> _rom_inputs;
599 :
600 : /// Container for ROM polynomial inputs
601 : std::vector<std::vector<std::vector<std::vector<GenericReal<is_ad>>>>> _polynomial_inputs;
602 :
603 : /// Container for ROM precomputed values
604 : std::vector<std::vector<std::vector<GenericReal<is_ad>>>> _precomputed_vals;
605 :
606 : /// Container for global limits
607 : std::vector<std::pair<Real, Real>> _global_limits;
608 :
609 : /// Container for weights for each tile as computed for all input values beside stress
610 : std::vector<std::vector<GenericReal<is_ad>>> _non_stress_weights;
611 :
612 : /// Container for weights for each tile as computed for all input values beside stress
613 : std::vector<std::vector<GenericReal<is_ad>>> _weights;
614 :
615 : /// Container for weights for each tile as computed for all input values beside stress
616 : std::vector<GenericReal<is_ad>> _partition_weights;
617 :
618 : /// Container for d_parition_weights d_stress
619 : std::vector<GenericReal<is_ad>> _dpartition_weight_dstress;
620 :
621 : /// Container for tiling orientations
622 : std::vector<std::vector<unsigned int>> _tiling;
623 :
624 : /// Container for strain cutoff
625 : std::vector<Real> _cutoff;
626 :
627 : /// Unit conversion factors required to convert from the specified unit to MPa
628 : Real _stress_ucf;
629 :
630 : ///@{Material properties accumulated at substeps
631 : GenericMaterialProperty<Real, is_ad> & _wall_dislocations_step;
632 : GenericMaterialProperty<Real, is_ad> & _cell_dislocations_step;
633 : ///@}
634 :
635 : /// Total plastic strain increment in step (summing substep contributions)
636 : RankTwoTensor _plastic_strain_increment;
637 :
638 : /// Material property capturing number of substeps for output purposes (defaults to one if substepping isn't used)
639 : MaterialProperty<Real> & _number_of_substeps;
640 :
641 : /// index names for error output
642 : std::vector<std::string> _index_name;
643 :
644 : /// check if a JSON file was loaded and if the specified key exists
645 : void checkJSONKey(const std::string & key);
646 :
647 : /// JSON object constructed from the datafile
648 : nlohmann::json _json;
649 :
650 : usingTransientInterfaceMembers;
651 : using Material::_name;
652 : using Material::_q_point;
653 : using Material::_qp;
654 : using Material::coupledGenericValue;
655 : using RadialReturnCreepStressUpdateBaseTempl<is_ad>::computeResidual;
656 : using RadialReturnCreepStressUpdateBaseTempl<is_ad>::computeDerivative;
657 : using RadialReturnCreepStressUpdateBaseTempl<is_ad>::_apply_strain;
658 : using RadialReturnCreepStressUpdateBaseTempl<is_ad>::initQpStatefulProperties;
659 : using RadialReturnCreepStressUpdateBaseTempl<is_ad>::outputIterationStep;
660 : using RadialReturnCreepStressUpdateBaseTempl<is_ad>::outputIterationSummary;
661 : };
662 :
663 : typedef LAROMANCEStressUpdateBaseTempl<false> LAROMANCEStressUpdateBase;
664 : typedef LAROMANCEStressUpdateBaseTempl<true> ADLAROMANCEStressUpdateBase;
|