예제 #1
0
def test_weights(comm, pnp, nx, ny, basenpoints):
    """compares the MPI-Implemtation of the halfspace with the serial one"""
    nx += basenpoints
    ny += basenpoints
    sx = 30.0
    sy = 1.0
    # equivalent Young's modulus
    E_s = 1.0

    substrate = PeriodicFFTElasticHalfSpace((nx, ny),
                                            E_s, (sx, sy),
                                            fft='mpi',
                                            communicator=comm)
    reference = PeriodicFFTElasticHalfSpace((nx, ny),
                                            E_s, (sx, sy),
                                            fft="fftw",
                                            communicator=MPI.COMM_SELF)
    np.testing.assert_allclose(
        reference.greens_function[substrate.fourier_slices],
        substrate.greens_function,
        rtol=0,
        atol=1e-16,
        err_msg="weights are different")
    np.testing.assert_allclose(
        reference.surface_stiffness[substrate.fourier_slices],
        substrate.surface_stiffness,
        rtol=0,
        atol=1e-16,
        err_msg="iweights are different")
    def test_unit_neutrality(self):
        tol = 1e-7
        # runs the same problem in two unit sets and checks whether results are
        # changed

        # Conversion factors
        length_c = 1. + np.random.rand()
        force_c = 1. + np.random.rand()
        pressure_c = force_c / length_c**2
        # energy_c = force_c * length_c

        # length_rc = (1., 1. / length_c)
        force_rc = (1., 1. / force_c)
        # pressure_rc = (1., 1. / pressure_c)
        # energy_rc = (1., 1. / energy_c)
        nb_grid_pts = (32, 32)
        young = (self.young, pressure_c * self.young)
        size = self.physical_sizes[0], 2 * self.physical_sizes[1]
        size = (size, tuple((length_c * s for s in size)))
        # print('SELF.SIZE = {}'.format(self.physical_sizes))

        disp = np.random.random(nb_grid_pts)
        disp -= disp.mean()
        disp = (disp, disp * length_c)

        forces = list()
        for i in range(2):
            sub = PeriodicFFTElasticHalfSpace(nb_grid_pts, young[i], size[i])
            force = sub.evaluate_force(disp[i])
            forces.append(force * force_rc[i])
        error = Tools.mean_err(forces[0], forces[1])
        self.assertTrue(error < tol,
                        "error = {} ≥ tol = {}".format(error, tol))
예제 #3
0
def test_sineWave_disp_rotation_invariance(comm, pnp, nx, ny, basenpoints):
    """
    for a sinusoidal displacement, test if the energy depends on if the wave is
    oriented in x or y direction

    Parameters
    ----------
    comm
    pnp
    fftengine_class
    nx
    ny
    basenpoints

    Returns
    -------

    """
    nx += basenpoints
    ny += basenpoints
    sx = 3.  # 30.0
    sy = 3.

    # equivalent Young's modulus
    E_s = 1.0

    computedenergies = []
    computedenergies_kspace = []
    for k in [(min(nx, ny) // 2, 0), (0, min(nx, ny) // 2)]:
        qx = k[0] * np.pi * 2 / sx
        qy = k[1] * np.pi * 2 / sy

        Y, X = np.meshgrid(
            np.linspace(0, sy, ny + 1)[:-1],
            np.linspace(0, sx, nx + 1)[:-1])
        # At the Nyquist frequency for even number of points, the energy
        # computation can only be exact for this point
        disp = np.cos(qx * X + qy * Y) + np.sin(qx * X + qy * Y)

        substrate = PeriodicFFTElasticHalfSpace((nx, ny),
                                                E_s, (sx, sy),
                                                fft='mpi',
                                                communicator=comm)

        computedenergies_kspace += [
            substrate.evaluate(disp[substrate.subdomain_slices],
                               pot=True,
                               forces=False)[0]
        ]
        computedenergies += [
            substrate.evaluate(disp[substrate.subdomain_slices],
                               pot=True,
                               forces=True)[0]
        ]

    # np.testing.assert_allclose(computedpressures[0],computedpressures[1].T)
    np.testing.assert_allclose(*computedenergies, rtol=1e-10)
    np.testing.assert_allclose(*computedenergies_kspace, rtol=1e-10)
 def test_uniform_displacement(self):
     """ test whether uniform displacement returns stiffness_q0"""
     sq0 = 1.43
     hs = PeriodicFFTElasticHalfSpace(self.res,
                                      self.young,
                                      self.physical_sizes,
                                      stiffness_q0=sq0)
     force = hs.evaluate_force(-np.ones(self.res))
     self.assertAlmostEqual(force.sum() / np.prod(self.physical_sizes), sq0)
 def test_same_energy(self):
     """
     Asserts that the energies computed in the real space and in the fourier
     space are the same
     """
     for res in [(16, 16), (16, 15), (15, 16), (15, 9)]:
         disp = np.random.normal(size=res)
         hs = PeriodicFFTElasticHalfSpace(res, self.young,
                                          self.physical_sizes)
         np.testing.assert_allclose(
             hs.evaluate(disp, pot=True, forces=True)[0],
             hs.evaluate(disp, pot=True, forces=False)[0])
 def test_uniform_displacement_finite_height(self):
     """ test whether uniform displacement returns stiffness_q0"""
     h0 = 3.45
     hs = PeriodicFFTElasticHalfSpace(self.res,
                                      self.young,
                                      self.physical_sizes,
                                      thickness=h0,
                                      poisson=self.poisson)
     force = hs.evaluate_force(-np.ones(self.res))
     M = (1 - self.poisson) / ((1 - 2 * self.poisson) *
                               (1 + self.poisson)) * self.young
     self.assertAlmostEqual(force.sum() / np.prod(self.physical_sizes),
                            M / h0)
 def test_limit_of_large_thickness(self):
     hs = PeriodicFFTElasticHalfSpace(self.res,
                                      self.young,
                                      self.physical_sizes,
                                      poisson=self.poisson)
     hsf = PeriodicFFTElasticHalfSpace(self.res,
                                       self.young,
                                       self.physical_sizes,
                                       thickness=20,
                                       poisson=self.poisson)
     # diff = hs.weights - hsf.weights
     np.testing.assert_allclose(hs.greens_function.ravel()[1:],
                                hsf.greens_function.ravel()[1:],
                                atol=1e-6)
 def test_parabolic_shape_disp(self):
     """ test whether the Elastic energy is a quadratic function of the
         applied displacement"""
     hs = PeriodicFFTElasticHalfSpace(self.res, self.young,
                                      self.physical_sizes)
     disp = random(self.res)
     disp -= disp.mean()
     nb_tests = 4
     El = np.zeros(nb_tests)
     for i in range(nb_tests):
         force = hs.evaluate_force(i * disp)
         El[i] = hs.evaluate_elastic_energy(i * force, disp)
     tol = 1e-10
     error = norm(El / El[1] - np.arange(nb_tests)**2)
     self.assertTrue(error < tol)
예제 #9
0
def test_sineWave_force(comm, pnp, nx, ny, basenpoints):
    """
    for  a given sinusoidal force, compares displacement with a reference
    solution

    Parameters
    ----------
    comm
    pnp
    fftengine_class
    nx
    ny
    basenpoints

    Returns
    -------

    """
    nx += basenpoints
    ny += basenpoints
    sx = 2  # 30.0
    sy = 1.0

    # equivalent Young's modulus
    E_s = 1.0

    Y, X = np.meshgrid(
        np.linspace(0, sy, ny + 1)[:-1],
        np.linspace(0, sx, nx + 1)[:-1])

    qx = 1 * np.pi * 2 / sx
    qy = 4 * np.pi * 2 / sy

    q = np.sqrt(qx**2 + qy**2)
    p = np.cos(qx * X + qy * Y)

    refdisp = -p / E_s * 2 / q

    substrate = PeriodicFFTElasticHalfSpace((nx, ny),
                                            E_s, (sx, sy),
                                            fft='mpi',
                                            communicator=comm)
    computeddisp = substrate.evaluate_disp(p[substrate.subdomain_slices] *
                                           substrate.area_per_pt)
    np.testing.assert_allclose(computeddisp,
                               refdisp[substrate.subdomain_slices],
                               atol=1e-7,
                               rtol=1e-10)
    def test_unit_neutrality1D(self):
        tol = 1e-7
        # runs the same problem in two unit sets and checks whether results are
        # changed

        # Conversion factors
        length_c = 1. + np.random.rand()
        force_c = 2. + np.random.rand()
        pressure_c = force_c / length_c**2
        # energy_c = force_c * length_c
        force_per_length_c = force_c / length_c

        # length_rc = (1., 1. / length_c)
        # force_rc = (1., 1. / force_c)
        # pressure_rc = (1., 1. / pressure_c)
        # energy_rc = (1., 1. / energy_c)
        force_per_length_rc = (1., 1. / force_per_length_c)

        nb_grid_pts = (32, )
        young = (self.young, pressure_c * self.young)
        size = (self.physical_sizes[0], length_c * self.physical_sizes[0])

        disp = np.random.random(nb_grid_pts)
        disp -= disp.mean()
        disp = (disp, disp * length_c)

        subs = tuple((PeriodicFFTElasticHalfSpace(nb_grid_pts, y, s)
                      for y, s in zip(young, size)))
        forces = tuple(
            (s.evaluate_force(d) * f_p_l
             for s, d, f_p_l in zip(subs, disp, force_per_length_rc)))
        error = Tools.mean_err(forces[0], forces[1])
        self.assertTrue(error < tol,
                        "error = {} ≥ tol = {}".format(error, tol))
 def test_no_nans(self):
     hs = PeriodicFFTElasticHalfSpace(self.res,
                                      self.young,
                                      self.physical_sizes,
                                      thickness=100,
                                      poisson=self.poisson)
     self.assertTrue(np.count_nonzero(np.isnan(hs.greens_function)) == 0)
    def test_consistency(self):
        pressure = list()
        base_res = 128
        tol = 1e-4
        for i in (1, 2):
            s_res = base_res * i
            test_res = (s_res, s_res)
            hs = PeriodicFFTElasticHalfSpace(test_res, self.young,
                                             self.physical_sizes)
            forces = np.zeros(test_res)
            forces[:s_res // 2, :s_res // 2] = 1.

            pressure.append(
                hs.evaluate_disp(forces)[::i, ::i] * hs.area_per_pt)
        error = ((pressure[0] - pressure[1])**2).sum().sum() / base_res**2
        self.assertTrue(error < tol, "error = {}".format(error))
    def test_gradient(self):
        res = size = (2, 2)
        disp = random(res)
        disp -= disp.mean()
        hs = PeriodicFFTElasticHalfSpace(res, self.young, size)
        hs.compute(disp, forces=True)
        f = hs.energy
        g = -hs.force
        approx_g = Tools.evaluate_gradient(
            lambda x: hs.evaluate(x, forces=True)[0], disp, 1e-5)

        tol = 1e-8
        error = Tools.mean_err(g, approx_g)
        msg = []
        msg.append("f = {}".format(f))
        msg.append("g = {}".format(g))
        msg.append('approx = {}'.format(approx_g))
        msg.append("error = {}".format(error))
        msg.append("tol = {}".format(tol))
        self.assertTrue(error < tol, ", ".join(msg))
    def test_force_disp_reversibility(self):
        # since only the zero-frequency is rejected, any force/disp field with
        # zero mean should be fully reversible
        tol = 1e-10
        for res in ((self.res[0], ), self.res, (self.res[0] + 1, self.res[1]),
                    (self.res[0], self.res[1] + 1)):
            hs = PeriodicFFTElasticHalfSpace(res, self.young,
                                             self.physical_sizes)
            disp = random(res)
            disp -= disp.mean()

            error = Tools.mean_err(disp,
                                   hs.evaluate_disp(hs.evaluate_force(disp)))
            self.assertTrue(
                error < tol,
                "for nb_grid_pts = {}, error = {} > tol = {}".format(
                    res, error, tol))

            force = random(res)
            force -= force.mean()

            error = Tools.mean_err(force,
                                   hs.evaluate_force(hs.evaluate_disp(force)))
            self.assertTrue(
                error < tol,
                "for nb_grid_pts = {}, error = {} > tol = {}".format(
                    res, error, tol))
예제 #15
0
def test_constrained_conjugate_gradients():
    nb_grid_pts = (512, 512)
    physical_sizes = (1., 1.)
    hurst = 0.8
    rms_slope = 0.1
    modulus = 1

    np.random.seed(999)
    topography = fourier_synthesis(nb_grid_pts,
                                   physical_sizes,
                                   hurst,
                                   rms_slope=rms_slope)

    substrate = PeriodicFFTElasticHalfSpace(nb_grid_pts, modulus,
                                            physical_sizes)
    system = make_system(substrate, topography)
    system.minimize_proxy(offset=0.1)
def test_constrained_conjugate_gradients():
    # system physical_sizes
    sx = 30.0
    sy = 1.0
    # equivalent Young's modulus
    E_s = 3.56
    for nx, ny in [(256, 16)]:  # , (256, 15), (255, 16)]:
        for disp0, normal_force in [(-0.9, None),
                                    (-0.1, None)]:  # (0.1, None),
            substrate = PeriodicFFTElasticHalfSpace((nx, ny), E_s,
                                                    (sx, sy))
            profile = np.resize(np.cos(2 * np.pi * np.arange(nx) / nx),
                                (ny, nx))
            surface = Topography(profile.T, (sx, sy))
            system = make_system(substrate, surface)

            result = system.minimize_proxy(offset=disp0,
                                           external_force=normal_force,
                                           pentol=1e-9)
            # offset = result.offset
            forces = result.jac
            # displ = result.x[:forces.shape[0], :forces.shape[1]]
            converged = result.success
            assert converged

            x = np.arange(nx) * sx / nx
            mean_pressure = np.mean(forces) / substrate.area_per_pt
            pth = mean_pressure * _pressure(
                x / sx,
                mean_pressure=sx * mean_pressure / E_s)

            # import matplotlib.pyplot as plt
            # plt.figure()
            # plt.plot(np.arange(nx)*self.sx/nx, profile)
            # plt.plot(x, displ[:, 0], 'r-')
            # plt.plot(x, surface[:, 0]+offset, 'k-')
            # plt.figure()
            # plt.plot(x, forces[:, 0]/substrate.area_per_pt, 'k-')
            # plt.plot(x, pth, 'r-')
            # plt.show()
            assert (np.allclose(
                forces[:nx // 2, 0] /
                substrate.area_per_pt, pth[:nx // 2], atol=1e-2))
    def test_energy(self):
        tol = 1e-10
        L = 2 + rand()  # domain length
        a = 3 + rand()  # amplitude of force
        E = 4 + rand()  # Young's Mod
        for res in [4, 8, 16]:
            area_per_pt = L / res
            x = np.arange(res) * area_per_pt
            force = a * np.cos(2 * np.pi / L * x)

            # theoretical FFT of force
            Fforce = np.zeros_like(x)
            Fforce[1] = Fforce[-1] = res / 2. * a

            # theoretical FFT of disp
            Fdisp = np.zeros_like(x)
            Fdisp[1] = Fdisp[-1] = res / 2. * a / E * L / (np.pi)

            # verify consistency
            hs = PeriodicFFTElasticHalfSpace(res, E, L)
            fforce = rfftn(force.T).T
            fdisp = hs.greens_function * fforce
            self.assertTrue(
                Tools.mean_err(fforce, Fforce, rfft=True) < tol,
                "fforce = \n{},\nFforce = \n{}".format(fforce.real, Fforce))
            self.assertTrue(
                Tools.mean_err(fdisp, Fdisp, rfft=True) < tol,
                "fdisp = \n{},\nFdisp = \n{}".format(fdisp.real, Fdisp))

            # Fourier energy
            E = .5 * np.dot(Fforce / area_per_pt, Fdisp) / res

            disp = hs.evaluate_disp(force)
            e = hs.evaluate_elastic_energy(force, disp)
            kdisp = hs.evaluate_k_disp(force)
            self.assertTrue(
                abs(disp - irfftn(kdisp.T).T).sum() < tol,
                ("disp   = {}\n"
                 "ikdisp = {}").format(disp,
                                       irfftn(kdisp.T).T))
            ee = hs.evaluate_elastic_energy_k_space(fforce, kdisp)
            self.assertTrue(
                abs(e - ee) < tol,
                "violate Parseval: e = {}, ee = {}, ee/e = {}".format(
                    e, ee, ee / e))

            self.assertTrue(
                abs(E - e) < tol,
                "theoretical E = {}, computed e = {}, diff(tol) = {}({})".
                format(E, e, E - e, tol))
def test_hard_wall_bearing_area(comm):
    # Test that at very low hardness we converge to (almost) the bearing
    # area geometry
    pnp = Reduction(comm)
    fullsurface = open_topography(os.path.join(FIXTURE_DIR, 'surface1.out'))
    fullsurface = fullsurface.topography(
        physical_sizes=fullsurface.channels[0].nb_grid_pts)
    nb_domain_grid_pts = fullsurface.nb_grid_pts
    substrate = PeriodicFFTElasticHalfSpace(nb_domain_grid_pts,
                                            1.0,
                                            fft='mpi',
                                            communicator=comm)
    surface = Topography(
        fullsurface.heights(),
        physical_sizes=nb_domain_grid_pts,
        decomposition='domain',
        subdomain_locations=substrate.topography_subdomain_locations,
        nb_subdomain_grid_pts=substrate.topography_nb_subdomain_grid_pts,
        communicator=substrate.communicator)

    plastic_surface = PlasticTopography(surface, 1e-12)
    system = PlasticNonSmoothContactSystem(substrate, plastic_surface)
    offset = -0.002
    if comm.rank == 0:

        def cb(it, p_r, d):
            print("{0}: area = {1}".format(it, d["area"]))
    else:

        def cb(it, p_r, d):
            pass

    result = system.minimize_proxy(offset=offset, callback=cb)
    assert result.success
    c = result.jac > 0.0
    ncontact = pnp.sum(c)
    assert plastic_surface.plastic_area == ncontact * surface.area_per_pt
    bearing_area = bisect(
        lambda x: pnp.sum((surface.heights() > x)) - ncontact, -0.03, 0.03)
    cba = surface.heights() > bearing_area
    # print(comm.Get_rank())
    assert pnp.sum(np.logical_not(c == cba)) < 25
예제 #19
0
print('RMS heights of surfaces = {} {}'.format(
    surface1.rms_height_from_area(), surface2.rms_height_from_area()))

# This is the grid nb_grid_pts of the two surfaces.
nx, ny = surface1.nb_grid_pts

# TranslatedSurface knows how to translate a surface into some direction.
translated_surface1 = TranslatedTopography(surface1)

# This is the compound of the two surfaces, effectively creating a surface
# that is the difference between the two profiles.
compound_surface = CompoundTopography(translated_surface1, surface2)

# Periodic substrate and hard-wall interactions.
substrate = PeriodicFFTElasticHalfSpace((nx, ny), E_s, (sx, sx))

# This creates a "System" object that knows about substrate, interaction and
# surface.
system = make_system(substrate, compound_surface)

# Initial displacement field.
disp = np.zeros(surface1.shape)

# Dump some information to this NetCDF file. Inspect the NetCDF with the
# 'ncdump' command.
container = NetCDFContainer('traj.nc', mode='w', double=True)
# NetCDF needs to know the nb_grid_pts/shape
container.set_shape(surface2)
# This creates a field called 'surface2' inside the NetCDF file.
container.surface2 = surface2.heights()
예제 #20
0
def test_multipleSineWaves_evaluate(comm, pnp, nx, ny, basenpoints):
    """
    displacements: superposition of sinwaves, compares forces and energes with
    analytical solution

    Parameters
    ----------
    comm
    pnp
    fftengine_class
    nx
    ny
    basenpoints

    Returns
    -------

    """
    nx += basenpoints
    ny += basenpoints
    sx = 2  # 30.0
    sy = 1.0
    # equivalent Young's modulus
    E_s = 1.0

    Y, X = np.meshgrid(
        np.linspace(0, sy, ny + 1)[:-1],
        np.linspace(0, sx, nx + 1)[:-1])

    disp = np.zeros((nx, ny))
    refForce = np.zeros((nx, ny))

    refEnergy = 0
    for qx, qy in zip((1, 0, 5, nx // 2 - 1), (4, 4, 0, ny // 2 - 2)):
        qx = qx * np.pi * 2 / sx
        qy = qy * np.pi * 2 / sy

        q = np.sqrt(qx**2 + qy**2)
        h = 1  # q**(-0.8)
        disp += h * (np.cos(qx * X + qy * Y) + np.sin(qx * X + qy * Y))
        refForce += h * (np.cos(qx * X + qy * Y) +
                         np.sin(qx * X + qy * Y)) * E_s / 2 * q
        refEnergy += E_s / 8 * q * 2 * h**2
        # * 2 because the amplitude of cos(x) + sin(x) is sqrt(2)

    # max possible Wavelengths at the edge

    for qx, qy in zip((nx // 2, nx // 2, 0), (ny // 2, 0, ny // 2)):
        qx = qx * np.pi * 2 / sx
        qy = qy * np.pi * 2 / sy

        q = np.sqrt(qx**2 + qy**2)
        h = 1  # q**(-0.8)
        disp += h * (np.cos(qx * X + qy * Y) + np.sin(qx * X + qy * Y))
        refForce += h * (np.cos(qx * X + qy * Y) +
                         np.sin(qx * X + qy * Y)) * E_s / 2 * q

        refEnergy += E_s / 8 * q * h**2 * 2
        # * 2 because the amplitude of cos(x) + sin(x) is sqrt(2)

    refEnergy *= sx * sy
    refForce *= -sx * sy / (nx * ny)

    substrate = PeriodicFFTElasticHalfSpace((nx, ny),
                                            E_s, (sx, sy),
                                            fft='mpi',
                                            communicator=comm)
    computed_E_k_space = substrate.evaluate(disp[substrate.subdomain_slices],
                                            pot=True,
                                            forces=False)[0]
    # If force is not queried this computes the energy using kspace
    computed_E_realspace, computed_force = substrate.evaluate(
        disp[substrate.subdomain_slices], pot=True, forces=True)

    # print("{}: Local: E_kspace: {}, E_realspace: {}"
    # .format(substrate.fftengine.comm.Get_rank(),computed_E_k_space,computed_E_realspace))
    # print(computed_E_k_space)
    # print(refEnergy)

    # if substrate.fftengine.comm.Get_rank() == 0 :
    #    print(computed_E_k_space)
    #    print(computed_E_realspace)

    # print("{}: Global: E_kspace: {}, E_realspace: {}"
    # .format(substrate.fftengine.comm.Get_rank(),
    # computed_E_k_space, computed_E_realspace))

    # Make an MPI-Reduce of the Energies !
    # print(substrate.evaluate_elastic_energy(refForce, disp))
    # print(0.5*np.vdot(refForce,disp))
    # print(substrate.evaluate_elastic_energy(substrate.evaluate_force(disp),disp))
    # print(computed_E_k_space)
    # print(computed_E_realspace)
    # print(refEnergy)

    np.testing.assert_almost_equal(computed_E_k_space, refEnergy)
    np.testing.assert_almost_equal(computed_E_realspace, refEnergy)
    np.testing.assert_allclose(computed_force,
                               refForce[substrate.subdomain_slices],
                               atol=1e-7,
                               rtol=1e-10)
예제 #21
0
def test_sineWave_disp(comm, pnp, nx, ny, basenpoints):
    """
    for given sinusoidal displacements, compares the pressures and the energies
    to the analytical solutions

    Special cases at the edges of the fourier domain are done

    Parameters
    ----------
    comm
    pnp
    fftengine_class
    nx
    ny
    basenpoints

    Returns
    -------

    """
    nx += basenpoints
    ny += basenpoints
    sx = 2.45  # 30.0
    sy = 1.0

    # equivalent Young's modulus
    E_s = 1.0

    for k in [(1, 0), (0, 1), (1, 2), (nx // 2, 0), (1, ny // 2), (0, 2),
              (nx // 2, ny // 2), (0, ny // 2)]:
        # print("testing wavevector ({}* np.pi * 2 / sx,
        # {}* np.pi * 2 / sy) ".format(*k))
        qx = k[0] * np.pi * 2 / sx
        qy = k[1] * np.pi * 2 / sy
        q = np.sqrt(qx**2 + qy**2)

        Y, X = np.meshgrid(
            np.linspace(0, sy, ny + 1)[:-1],
            np.linspace(0, sx, nx + 1)[:-1])
        disp = np.cos(qx * X + qy * Y) + np.sin(qx * X + qy * Y)

        refpressure = -disp * E_s / 2 * q

        substrate = PeriodicFFTElasticHalfSpace((nx, ny),
                                                E_s, (sx, sy),
                                                fft='mpi',
                                                communicator=comm)
        fftengine = FFT((nx, ny), fft='mpi', communicator=comm)
        fftengine.create_plan(1)

        kpressure = substrate.evaluate_k_force(
            disp[substrate.subdomain_slices]) / substrate.area_per_pt / (nx *
                                                                         ny)
        expected_k_disp = np.zeros((nx // 2 + 1, ny), dtype=complex)
        expected_k_disp[k[0], k[1]] += .5 - .5j

        # add the symetrics
        if k[0] == 0:
            expected_k_disp[0, -k[1]] += .5 + .5j
        if k[0] == nx // 2 and nx % 2 == 0:
            expected_k_disp[k[0], -k[1]] += .5 + .5j

        fft_disp = np.zeros(substrate.nb_fourier_grid_pts,
                            order='f',
                            dtype=complex)
        fftengine.fft(disp[substrate.subdomain_slices], fft_disp)
        np.testing.assert_allclose(fft_disp / (nx * ny),
                                   expected_k_disp[substrate.fourier_slices],
                                   rtol=1e-7,
                                   atol=1e-10)

        expected_k_pressure = -E_s / 2 * q * expected_k_disp
        np.testing.assert_allclose(
            kpressure,
            expected_k_pressure[substrate.fourier_slices],
            rtol=1e-7,
            atol=1e-10)

        computedpressure = substrate.evaluate_force(
            disp[substrate.subdomain_slices]) / substrate.area_per_pt
        np.testing.assert_allclose(computedpressure,
                                   refpressure[substrate.subdomain_slices],
                                   atol=1e-10,
                                   rtol=1e-7)

        computedenergy_kspace = \
            substrate.evaluate(disp[substrate.subdomain_slices], pot=True,
                               forces=False)[0]
        computedenergy = \
            substrate.evaluate(disp[substrate.subdomain_slices], pot=True,
                               forces=True)[0]
        refenergy = E_s / 8 * 2 * q * sx * sy

        # print(substrate.nb_domain_grid_pts[-1] % 2)
        # print(substrate.nb_fourier_grid_pts)
        # print(substrate.fourier_locations[-1] +
        # substrate.nb_fourier_grid_pts[-1] - 1)
        # print(substrate.nb_domain_grid_pts[-1] // 2 )
        # print(computedenergy)
        # print(computedenergy_kspace)
        # print(refenergy)
        np.testing.assert_allclose(
            computedenergy,
            refenergy,
            rtol=1e-10,
            err_msg="wavevektor {} for nb_domain_grid_pts {}, "
            "subdomain nb_grid_pts {}, nb_fourier_grid_pts {}".format(
                k, substrate.nb_domain_grid_pts,
                substrate.nb_subdomain_grid_pts,
                substrate.nb_fourier_grid_pts))
        np.testing.assert_allclose(
            computedenergy_kspace,
            refenergy,
            rtol=1e-10,
            err_msg="wavevektor {} for nb_domain_grid_pts {}, "
            "subdomain nb_grid_pts {}, nb_fourier_grid_pts {}".format(
                k, substrate.nb_domain_grid_pts,
                substrate.nb_subdomain_grid_pts,
                substrate.nb_fourier_grid_pts))
def test_lbfgsb_1D(dx, n):
    # test the 1D periodic objective is working

    s = n * dx

    Es = 1.

    substrate = PeriodicFFTElasticHalfSpace((n,), young=Es,
                                            physical_sizes=(s,), fft='serial')

    surface = UniformLineScan(np.cos(2 * np.pi * np.arange(n) / n), (s,))
    system = make_system(substrate, surface)

    offset = 0.005
    lbounds = np.zeros(n)
    bnds = system._reshape_bounds(lbounds, )
    init_gap = np.zeros(n, )  # .flatten()
    disp = init_gap + surface.heights() + offset

    res = optim.minimize(system.primal_objective(offset, gradient=True),
                         disp,
                         method='L-BFGS-B', jac=True,
                         bounds=bnds,
                         options=dict(gtol=1e-5 * Es * surface.rms_slope_from_profile()
                                      * surface.area_per_pt,
                                      ftol=0))

    forces = res.jac
    gap = res.x
    converged = res.success
    assert converged
    if hasattr(res.message, "decode"):
        decoded_message = res.message.decode()
    else:
        decoded_message = res.message

    assert decoded_message == \
        'CONVERGENCE: NORM_OF_PROJECTED_GRADIENT_<=_PGTOL'

    x = np.arange(n) * s / n
    mean_pressure = np.mean(forces) / substrate.area_per_pt

    assert mean_pressure > 0

    pth = mean_pressure * _pressure(
        x / s,
        mean_pressure=s * mean_pressure / Es)

    if False:
        import matplotlib.pyplot as plt
        plt.figure()
        # plt.plot(np.arange(n)*s/n, surface.heights())
        plt.plot(x, gap + (surface.heights()[:] + offset), 'r-')
        plt.plot(x, (surface.heights()[:] + offset), 'k-')
        plt.figure()
        plt.plot(x, forces / substrate.area_per_pt, 'k-')
        plt.plot(x, pth, 'r-')
        plt.show()

    assert (np.allclose(
        forces /
        substrate.area_per_pt, pth, atol=1e-2))
예제 #23
0
logger.pr('SurfaceTopography has dimension of {} and physical_sizes of {} {}.'
          .format(surface.nb_grid_pts, surface.physical_sizes, surface.info['unit']))
logger.pr('RMS height = {}, RMS gradient = {}'.format(surface.rms_height_from_area(),
                                                   surface.rms_gradient()))
if arguments.detrend is not None:
    surface = surface.detrend(detrend_mode=arguments.detrend)
    logger.pr('After detrending: RMS height = {}, RMS gradient = {}' \
              .format(surface.rms_height_from_area(), surface.rms_gradient()))

if arguments.hardness is not None:
    surface = PlasticTopography(surface, arguments.hardness)

# Initialize elastic half-space.
if arguments.boundary == 'periodic':
    substrate = PeriodicFFTElasticHalfSpace(surface.nb_grid_pts, arguments.modulus,
                                            surface.physical_sizes,
                                            thickness=arguments.thickness,
                                            poisson=arguments.poisson)
elif arguments.boundary == 'nonperiodic':
    if arguments.thickness is not None:
        raise ValueError('"thickness" arguments cannot be used with '
                         'nonperiodic boundaries.')
    substrate = FreeFFTElasticHalfSpace(
        surface.nb_grid_pts,
        arguments.modulus/(1-arguments.poisson**2),
        surface.physical_sizes
        )
else:
    raise ValueError('Unknown boundary conditions: '
                     '{}'.format(arguments.boundary))

# Piece the full system together. In particular the PyCo.System.SystemBase
def test_constrained_conjugate_gradients(comm,):
    sx = 30.0
    sy = 1.0
    # equivalent Young's modulus
    E_s = 3.56
    pnp = Reduction(comm)

    for nx, ny in [(16384, 8 * comm.Get_size())]:
        # if ny is too small, one processor will end
        # with one empty fourier subdomain, what is not supported
        # 16384 points are indeed needed to reach the relative error of 1e-2

        for disp0, normal_force in [(-0.9, None),
                                    (-0.1, None)]:  # (0.1, None),
            subs = PeriodicFFTElasticHalfSpace(
                (nx, ny), E_s, (sx, sy),
                fft="serial" if comm.Get_size() == 1 else "mpi",
                communicator=comm)
            profile = np.resize(np.cos(2 * np.pi * np.arange(nx) / nx),
                                (ny, nx))
            surface = Topography(
                profile.T, physical_sizes=(sx, sy),
                # nb_grid_pts=substrate.nb_grid_pts,
                decomposition='domain',
                subdomain_locations=subs.topography_subdomain_locations,
                nb_subdomain_grid_pts=subs.topography_nb_subdomain_grid_pts,
                communicator=subs.communicator)
            system = make_system(subs, surface)

            result = system.minimize_proxy(offset=disp0,
                                           external_force=normal_force,
                                           pentol=1e-12)

            forces = result.jac

            converged = result.success
            assert converged

            # print(forces)
            # print(displ)

            x = np.arange(nx) * sx / nx
            mean_pressure = pnp.sum(forces) / np.prod(subs.physical_sizes)
            analytical_pressures = mean_pressure * \
                _pressure(x / sx, mean_pressure=sx * mean_pressure / E_s)

            # symetrize the Profile
            analytical_pressures[1:] = analytical_pressures[1:] \
                + analytical_pressures[:0:-1]

            if False:
                offset = result.offset
                displ = result.x[:forces.shape[0], :forces.shape[1]]
                import matplotlib.pyplot as plt
                plt.figure()
                plt.plot(np.arange(nx)*sx/nx, profile[0, :], "--k")
                plt.plot(x, displ[:, 0], 'r-', label="displacements")
                plt.plot(x, surface[:, 0] + offset, 'k-', label="indenter")
                plt.legend()
                plt.figure()
                plt.plot(x, forces[:, 0] / system.substrate.area_per_pt, 'k-',
                         label="numerical")
                plt.plot(x, analytical_pressures, 'r--', label="analytical")
                plt.legend()
                plt.show()

            # boolean array indicating where tolerance is not reached
            error_mask = np.abs(
                (forces[:, 0] / subs.area_per_pt - analytical_pressures[
                    subs.subdomain_slices[0]]) >= 1e-12 + 1e-2 * np.abs(
                    analytical_pressures[subs.subdomain_slices[0]]))

            # parallelise the check, this leads to cleaner failures
            assert pnp.sum(np.count_nonzero(
                error_mask)) == 0, "max relative diff at index {} with " \
                                   "ref = {}, computed= {}".format(
                np.arange(subs.nb_subdomain_grid_pts[0])[error_mask],
                analytical_pressures[subs.subdomain_slices[0]][error_mask],
                forces[:, 0][error_mask] / subs.area_per_pt)
예제 #25
0
    return result


###

# Read the topography from file.
topography = read_topography('surface1.out', physical_sizes=(sx, sy))

print('RMS height of topography = {}'.format(
    topography.rms_height_from_area()))
print('RMS slope of topography = {}'.format(topography.rms_gradient()))

# This is the grid nb_grid_pts of the topography.
nx, ny = topography.nb_grid_pts

# Periodic substrate, i.e. the elastic half-space.
substrate = PeriodicFFTElasticHalfSpace((nx, ny), E_s,
                                        (sx, sx))  # physical physical_sizes

# Contact pressure: 0.05 h_rms' E_s / kappa with kappa = 2, which means ~ 5%
# contact area
res = constrained_conjugate_gradients(
    substrate,
    topography.heights(),
    external_force=0.05 * topography.rms_slope() * E_s * sx * sy / 2,
    logger=screen)

# Show contact area
plt.pcolormesh(res.jac > 0)
plt.show()