Page MenuHomec4science

plasticity.py
No OneTemporary

File Metadata

Created
Sat, Jul 6, 22:58

plasticity.py

import numpy as np
import tamaas as tm
from tamaas.dumpers import UVWDumper as Dumper
from tamaas.nonlinear_solvers import DFSANESolver as Solver
tm.initialize(2)
# Definition of modeled domain
model_type = tm.model_type.volume_2d
discretization = [32, 51, 51]
flat_domain = [1, 1]
system_size = [0.5] + flat_domain
# Creation of model
model = tm.ModelFactory.createModel(model_type,
system_size,
discretization)
model.E = 1.
model.nu = 0.3
# Setup for plasticity
residual = tm.ModelFactory.createResidual(model, 0.1 * model.E, 0.01 * model.E)
epsolver = Solver(residual, model)
# Setup for contact
x = np.linspace(0, system_size[1], discretization[1], endpoint=False)
y = np.linspace(0, system_size[2], discretization[2], endpoint=False)
xx, yy = np.meshgrid(x, y, indexing='ij')
R = 0.2
surface = -((xx - flat_domain[0] / 2)**2 +
(yy - flat_domain[1] / 2)**2) / (2 * R)
csolver = tm.PolonskyKeerRey(model, surface, 1e-12,
tm.PolonskyKeerRey.pressure,
tm.PolonskyKeerRey.pressure)
# EPIC setup
epic = tm.EPICSolver(csolver, epsolver, 1e-9)
# Dumper
dumper_helper = Dumper(residual, 'hertz')
model.setDumper(dumper_helper)
loads = np.linspace(0.001, 0.005, 3)
for i, load in enumerate(loads):
epic.solve(load)
print("---> Solved load step {}/{}".format(i+1, len(loads)))

Event Timeline