def test_rfftn(self): force = np.zeros([2 * r for r in self.res]) force[:self.res[0], :self.res[1]] = np.random.random(self.res) from muFFT import FFT ref = rfftn(force.T).T fftengine = FFT([2 * r for r in self.res], fft="serial") fftengine.create_plan(1) tested = np.zeros(fftengine.nb_fourier_grid_pts, order='f', dtype=complex) fftengine.fft(force, tested) np.testing.assert_allclose(ref.real, tested.real) np.testing.assert_allclose(ref.imag, tested.imag)
class PeriodicFFTElasticHalfSpace(ElasticSubstrate): """ Uses the FFT to solve the displacements and stresses in an elastic Halfspace due to a given array of point forces. This halfspace implementation cheats somewhat: since a net pressure would result in infinite displacement, the first term of the FFT is systematically dropped. The implementation follows the description in Stanley & Kato J. Tribol. 119(3), 481-485 (Jul 01, 1997) """ name = "periodic_fft_elastic_halfspace" _periodic = True def __init__(self, nb_grid_pts, young, physical_sizes=2 * np.pi, stiffness_q0=None, thickness=None, poisson=0.0, superclass=True, fft="serial", communicator=None): """ Parameters ---------- nb_grid_pts : int tuple containing number of points in spatial directions. The length of the tuple determines the spatial dimension of the problem. young : float Young's modulus, if poisson is not specified it is the contact modulus as defined in Johnson, Contact Mechanics physical_sizes : float or float tuple (default 2π) domain size. For multidimensional problems, a tuple can be provided to specify the lengths per dimension. If the tuple has less entries than dimensions, the last value in repeated. stiffness_q0 : float, optional Substrate stiffness at the Gamma-point (wavevector q=0). If None, this is taken equal to the lowest nonvanishing stiffness. Cannot be used in combination with thickness. thickness : float, optional Thickness of the elastic half-space. If None, this models an infinitely deep half-space. Cannot be used in combination with stiffness_q0. poisson : float Default 0 Poisson number. Need only be specified for substrates of finite thickness. If left unspecified for substrates of infinite thickness, then young is the contact modulus. superclass : bool (default True) client software never uses this. Only inheriting subclasses use this. fft: string Default: 'serial' FFT engine to use. Options are 'fftw', 'fftwmpi', 'pfft' and 'p3dfft'. 'serial' and 'mpi' can also be specified, where the choice of the appropriate fft is made by muFFT communicator : mpi4py communicator or NuMPI stub communicator MPI communicator object. """ super().__init__() if not hasattr(nb_grid_pts, "__iter__"): nb_grid_pts = (nb_grid_pts, ) if not hasattr(physical_sizes, "__iter__"): physical_sizes = (physical_sizes, ) self.__dim = len(nb_grid_pts) if self.dim not in (1, 2): raise self.Error( ("Dimension of this problem is {}. Only 1 and 2-dimensional " "problems are supported").format(self.dim)) if stiffness_q0 is not None and thickness is not None: raise self.Error("Please specify either stiffness_q0 or thickness " "or neither.") self._nb_grid_pts = nb_grid_pts tmpsize = list() for i in range(self.dim): tmpsize.append(physical_sizes[min(i, len(physical_sizes) - 1)]) self._physical_sizes = tuple(tmpsize) try: self._steps = tuple( float(size) / res for size, res in zip(self.physical_sizes, self.nb_grid_pts)) except ZeroDivisionError as err: raise ZeroDivisionError( ("{}, when trying to handle " " self._steps = tuple(" " float(physical_sizes)/res for physical_sizes, res in" " zip(self.physical_sizes, self.nb_grid_pts))" "Parameters: self.physical_sizes = {}, self.nb_grid_pts = {}" "").format(err, self.physical_sizes, self.nb_grid_pts)) self.young = young self.poisson = poisson self.contact_modulus = young / (1 - poisson**2) self.stiffness_q0 = stiffness_q0 self.thickness = thickness self.fftengine = FFT(self.nb_domain_grid_pts, fft=fft, communicator=communicator, allow_temporary_buffer=False, allow_destroy_input=True) # Allocate buffers and create plan for one degree of freedom self.real_buffer = self.fftengine.register_real_space_field( "real-space", 1) self.fourier_buffer = self.fftengine.register_fourier_space_field( "fourier-space", 1) self.greens_function = None self.surface_stiffness = None self._communicator = communicator self.pnp = Reduction(communicator) if superclass: self.greens_function = self._compute_greens_function() self.surface_stiffness = self._compute_surface_stiffness() @property def dim(self, ): "return the substrate's physical dimension" return self.__dim @property def nb_grid_pts(self): return self._nb_grid_pts @property def area_per_pt(self): return np.prod(self.physical_sizes) / np.prod(self.nb_grid_pts) @property def physical_sizes(self): return self._physical_sizes @property def nb_domain_grid_pts(self, ): """ usually, the nb_grid_pts of the system is equal to the geometric nb_grid_pts (of the surface). For example free boundary conditions, require the computational nb_grid_pts to differ from the geometric one, see FreeFFTElasticHalfSpace. """ return self.nb_grid_pts @property def nb_subdomain_grid_pts(self): """ When working in Parallel one processor holds only Part of the Data :return: """ return self.fftengine.nb_subdomain_grid_pts @property def topography_nb_subdomain_grid_pts(self): return self.nb_subdomain_grid_pts @property def subdomain_locations(self): """ When working in Parallel one processor holds only Part of the Data :return: """ return self.fftengine.subdomain_locations @property def topography_subdomain_locations(self): return self.subdomain_locations @property def subdomain_slices(self): """ When working in Parallel one processor holds only Part of the Data :return: """ return self.fftengine.subdomain_slices @property def topography_subdomain_slices(self): return tuple([ slice(s, s + n) for s, n in zip(self.topography_subdomain_locations, self.topography_nb_subdomain_grid_pts) ]) @property def local_topography_subdomain_slices(self): """ slice representing the local subdomain without the padding area """ return tuple( [slice(0, n) for n in self.topography_nb_subdomain_grid_pts]) @property def nb_fourier_grid_pts(self): """ When working in Parallel one processor holds only Part of the Data :return: """ return self.fftengine.nb_fourier_grid_pts @property def fourier_locations(self): """ When working in Parallel one processor holds only Part of the Data :return: """ return self.fftengine.fourier_locations @property def fourier_slices(self): """ When working in Parallel one processor holds only Part of the Data :return: """ return self.fftengine.fourier_slices @property def communicator(self): """Return the MPI communicator""" return self._communicator def __repr__(self): dims = 'x', 'y', 'z' size_str = ', '.join('{}: {}({})'.format(dim, size, nb_grid_pts) for dim, size, nb_grid_pts in zip( dims, self.physical_sizes, self.nb_grid_pts)) return "{0.dim}-dimensional halfspace '{0.name}', " \ "physical_sizes(nb_grid_pts) in {1}, E' = {0.young}" \ .format(self, size_str) def _compute_greens_function(self): r""" Compute the weights w relating fft(displacement) to fft(pressure): fft(u) = w*fft(p), see (6) Stanley & Kato J. Tribol. 119(3), 481-485 (Jul 01, 1997). For the infinite halfspace, .. math :: w = q E^* / 2 q is the wavevector (:math:`2 \pi / wavelength`) WARNING: the paper is dimensionally *incorrect*. see for the correct 1D formulation: Section 13.2 in K. L. Johnson. (1985). Contact Mechanics. [Online]. Cambridge: Cambridge University Press. Available from: Cambridge Books Online <http://dx.doi.org/10.1017/CBO9781139171731> [Accessed 16 February 2015] for correct 2D formulation: Appendix 1, eq A.2 in Johnson, Greenwood and Higginson, "The Contact of Elastic Regular Wavy surfaces", Int. J. Mech. Sci. Vol. 27 No. 6, pp. 383-396, 1985 <http://dx.doi.org/10.1016/0020-7403(85)90029-3> [Accessed 18 March 2015] """ if self.dim == 1: nx, = self.nb_grid_pts sx, = self.physical_sizes # Note: q-values from 0 to 1, not from 0 to 2*pi qx = np.arange(self.fourier_locations[0], self.fourier_locations[0] + self.nb_fourier_grid_pts[0], dtype=np.float64) qx = np.where(qx <= nx // 2, qx / sx, (nx - qx) / sx) surface_stiffness = np.pi * self.contact_modulus * qx if self.stiffness_q0 is None: surface_stiffness[0] = surface_stiffness[1].real elif self.stiffness_q0 == 0.0: surface_stiffness[0] = 1.0 else: surface_stiffness[0] = self.stiffness_q0 greens_function = 1 / surface_stiffness if self.fourier_locations == (0, ): if self.stiffness_q0 == 0.0: greens_function[0, 0] = 0.0 elif self.dim == 2: if np.prod(self.nb_fourier_grid_pts) == 0: greens_function = np.zeros(self.nb_fourier_grid_pts, order='f', dtype=complex) else: nx, ny = self.nb_grid_pts sx, sy = self.physical_sizes # Note: q-values from 0 to 1, not from 0 to 2*pi qx = np.arange(self.fourier_locations[0], self.fourier_locations[0] + self.nb_fourier_grid_pts[0], dtype=np.float64) qx = np.where(qx <= nx // 2, qx / sx, (nx - qx) / sx) qy = np.arange(self.fourier_locations[1], self.fourier_locations[1] + self.nb_fourier_grid_pts[1], dtype=np.float64) qy = np.where(qy <= ny // 2, qy / sy, (ny - qy) / sy) q = np.sqrt((qx * qx).reshape(-1, 1) + (qy * qy).reshape(1, -1)) if self.fourier_locations == (0, 0): q[0, 0] = np.NaN # q[0,0] has no Impact on the end result, # but q[0,0] = 0 produces runtime Warnings # (because corr[0,0]=inf) surface_stiffness = np.pi * self.contact_modulus * q # E* / 2 (2 \pi / \lambda) # (q is 1 / lambda, here) if self.thickness is not None: # Compute correction for finite thickness q *= 2 * np.pi * self.thickness fac = 3 - 4 * self.poisson off = 4 * self.poisson * (2 * self.poisson - 3) + 5 with np.errstate(over="ignore", invalid="ignore", divide="ignore"): corr = (fac * np.cosh(2 * q) + 2 * q ** 2 + off) / \ (fac * np.sinh(2 * q) - 2 * q) # The expression easily overflows numerically. These are # then q-values that are converged to the infinite system # expression. corr[np.isnan(corr)] = 1.0 surface_stiffness *= corr if self.fourier_locations == (0, 0): surface_stiffness[0, 0] = \ self.young / self.thickness * \ (1 - self.poisson) / ((1 - 2 * self.poisson) * (1 + self.poisson)) else: if self.fourier_locations == (0, 0): if self.stiffness_q0 is None: surface_stiffness[0, 0] = \ (surface_stiffness[1, 0].real + surface_stiffness[0, 1].real) / 2 elif self.stiffness_q0 == 0.0: surface_stiffness[0, 0] = 1.0 else: surface_stiffness[0, 0] = self.stiffness_q0 greens_function = 1 / surface_stiffness if self.fourier_locations == (0, 0): if self.stiffness_q0 == 0.0: greens_function[0, 0] = 0.0 return greens_function def _compute_surface_stiffness(self): """ Invert the weights w relating fft(displacement) to fft(pressure): """ surface_stiffness = np.zeros(self.nb_fourier_grid_pts, order='f', dtype=complex) surface_stiffness[self.greens_function != 0] = \ 1. / self.greens_function[self.greens_function != 0] return surface_stiffness def evaluate_disp(self, forces): """ Computes the displacement due to a given force array Keyword Arguments: forces -- a numpy array containing point forces (*not* pressures) """ if forces.shape != self.nb_subdomain_grid_pts: raise self.Error( ("force array has a different shape ({0}) than this " "halfspace's nb_grid_pts ({1})").format( forces.shape, self.nb_subdomain_grid_pts)) self.real_buffer.array()[...] = -forces self.fftengine.fft(self.real_buffer, self.fourier_buffer) self.fourier_buffer.array()[...] *= self.greens_function self.fftengine.ifft(self.fourier_buffer, self.real_buffer) return self.real_buffer.array().real / \ self.area_per_pt * self.fftengine.normalisation def evaluate_force(self, disp): """ Computes the force (*not* pressures) due to a given displacement array. Keyword Arguments: disp -- a numpy array containing point displacements """ if disp.shape != self.nb_subdomain_grid_pts: raise self.Error( ("displacements array has a different shape ({0}) than " "this halfspace's nb_grid_pts ({1})").format( disp.shape, self.nb_subdomain_grid_pts)) self.real_buffer.array()[...] = disp self.fftengine.fft(self.real_buffer, self.fourier_buffer) self.fourier_buffer.array()[...] *= self.surface_stiffness self.fftengine.ifft(self.fourier_buffer, self.real_buffer) return -self.real_buffer.array().real * \ self.area_per_pt * self.fftengine.normalisation def evaluate_k_disp(self, forces): """ Computes the K-space displacement due to a given force array Keyword Arguments: forces -- a numpy array containing point forces (*not* pressures) """ if forces.shape != self.nb_subdomain_grid_pts: raise self.Error( ("force array has a different shape ({0}) than this halfspace'" "s nb_grid_pts ({1})").format( forces.shape, self.nb_subdomain_grid_pts)) # nopep8 self.real_buffer.array()[...] = -forces self.fftengine.fft(self.real_buffer, self.fourier_buffer) return self.greens_function * \ self.fourier_buffer.array() / self.area_per_pt def evaluate_k_force(self, disp): """ Computes the K-space forces (*not* pressures) due to a given displacement array. Keyword Arguments: disp -- a numpy array containing point displacements """ if disp.shape != self.nb_subdomain_grid_pts: raise self.Error( ("displacements array has a different shape ({0}) than this " "halfspace's nb_grid_pts ({1})").format( disp.shape, self.nb_subdomain_grid_pts)) # nopep8 self.real_buffer.array()[...] = disp self.fftengine.fft(self.real_buffer, self.fourier_buffer) return -self.surface_stiffness * \ self.fourier_buffer.array() * self.area_per_pt def evaluate_k_force_k(self, disp_k): """ Computes the K-space forces (*not* pressures) due to a given displacement array. Parameters: ----------- disp_k: complex nd_array a numpy array containing the rfft of point displacements """ return -self.surface_stiffness * disp_k * self.area_per_pt def evaluate_elastic_energy(self, forces, disp): """ computes and returns the elastic energy due to forces and displacements Arguments: forces -- array of forces disp -- array of displacements """ # pylint: disable=no-self-use return .5 * self.pnp.dot(np.ravel(disp), np.ravel(-forces)) def evaluate_scalar_product_k_space(self, ka, kb): r""" Computes the scalar product, i.e. the power, between the `a` and `b`, given their fourier representation. `Power theorem <https://ccrma.stanford.edu/~jos/mdft/Power_Theorem.html>`_: .. math :: P = \sum_{ij} a_{ij} b_{ij} = \frac{1}{n_x n_y}\sum_{ij} \tilde a_{ij} \overline{\tilde b_{ij}} Note that for `a`, `b` real, .. math :: P = \sum_{kl} Re(\tilde a_{kl}) Re(\tilde b_{kl}) + Im(\tilde a_{kl}) Im(\tilde b_{kl}) Parameters ---------- ka, kb: arrays of complex type and of size substrate.nb_fourier_grid_pts Fourier representation (output of a 2D rfftn) `a` (resp. `b`) (`nx, ny` real array) Returns ------- P The scalar product of a and b """ # ka and kb are the output of the 2D rfftn, that means the a # part of the transform is omitted because of the symetry along the # last dimension # # That's why the components whose symetrics have been omitted are # weighted with a factor of 2. # # The first column (indexes [...,0], wavevector 0 along the last # dimension) has no symetric # # When the number of points in the last dimension is even, the last # column (Nyquist Frequency) has also no symetric. # # The serial code implementation would look like this # if (self.nb_domain_grid_pts[-1] % 2 == 0) # return .5*(np.vdot(ka, kb).real + # # adding the data that has been omitted by rfftn # np.vdot(ka[..., 1:-1], kb[..., 1:-1]).real # # because of symetry # )/self.nb_pts # else : # return .5 * (np.vdot(ka, kb).real + # # adding the data that has been omitted by rfftn # # np.vdot(ka[..., 1:], kb[..., 1:]).real # # # because of symetry # # )/self.nb_pts # # Parallelized Version # The inner part of the fourier data should always be symetrized (i.e. # multiplied by 2). When the fourier subdomain contains boundary values # (wavevector 0 (even and odd) and ny//2 (only for odd)) these values # should only be added once if ka.size > 0: if self.fourier_locations[0] == 0: # First row of this fourier data is first of global data fact0 = 1 elif self.nb_fourier_grid_pts[0] > 1: # local first row is not the first in the global data fact0 = 2 else: fact0 = 0 if self.fourier_locations[0] == 0 and \ self.nb_fourier_grid_pts[0] == 1: factend = 0 elif (self.nb_domain_grid_pts[0] % 2 == 1): # odd number of points, last row have always to be symmetrized factend = 2 elif self.fourier_locations[0] + \ self.nb_fourier_grid_pts[0] - 1 == \ self.nb_domain_grid_pts[0] // 2: # last row of the global rfftn already contains it's symmetric factend = 1 # print("last Element of the even data has to be accounted # only once") else: factend = 2 # print("last element of this local slice is not last element # of the total global data") # print("fact0={}".format(fact0)) # print("factend={}".format(factend)) if self.nb_fourier_grid_pts[0] > 2: factmiddle = 2 else: factmiddle = 0 # vdot(a, b) = conj(a) . b locsum = (factmiddle * np.vdot(ka[1:-1, ...], kb[1:-1, ...]).real + fact0 * np.vdot(ka[0, ...], kb[0, ...]).real + factend * np.vdot(ka[-1, ...], kb[-1, ...]).real) / np.prod( self.nb_domain_grid_pts) # nopep8 # We divide by the total number of points to get the appropriate # normalisation of the Fourier transform (in numpy the division by # happens only at the inverse transform) else: # This handles the case where the processor holds an empty # subdomain locsum = np.array([], dtype=ka.real.dtype) # print(locsum) return self.pnp.sum(locsum) def evaluate_elastic_energy_k_space(self, kforces, kdisp): r""" Computes the Energy due to forces and displacements using their Fourier representation. .. math :: E_{el} &= - \frac{1}{2} \sum_{ij} u_{ij} f_{ij} &= - \frac{1}{2} \frac{1}{n_x n_y} \sum_{kl} \tilde u{kl} \overline{\tilde f_{kl}} (:math:`\tilde f_{ij} = - \tilde K_{ijkl} u`) In a parallelized code kforces and kdisp contain only the slice attributed to this processor Parameters ---------- kforces: array of complex type and of size substrate.nb_fourier_grid_pts Fourier representation (output of a 2D rfftn) of the forces acting on the grid points kdisp: array of complex type and of physical_sizes substrate.nb_fourier_grid_pts Fourier representation (output of a 2D rfftn) of the displacements of the grid points Returns ------- E The elastic energy due to the forces and displacements """ # noqa: E501, W291, W293 return -0.5 * self.evaluate_scalar_product_k_space(kdisp, kforces) def evaluate(self, disp, pot=True, forces=False): """Evaluates the elastic energy and the point forces Keyword Arguments: disp -- array of distances pot -- (default True) if true, returns potential energy forces -- (default False) if true, returns forces """ force = potential = None if forces: force = self.evaluate_force(disp) if pot: potential = self.evaluate_elastic_energy(force, disp) elif pot: kforce = self.evaluate_k_force(disp) # TODO: OPTIMISATION: here kdisp is computed twice, because it's # needed in kforce self.real_buffer.array()[...] = disp self.fftengine.fft(self.real_buffer, self.fourier_buffer) potential = self.evaluate_elastic_energy_k_space( kforce, self.fourier_buffer.array()) return potential, force def evaluate_k(self, disp_k, pot=True, forces=False): """Evaluates the elastic energy and the point forces Keyword Arguments: disp -- array of distances pot -- (default True) if true, returns potential energy forces -- (default False) if true, returns forces """ potential = None if forces: force_k = self.evaluate_k_force_k(disp_k) if pot: potential = self.evaluate_elastic_energy_k_space( force_k, disp_k) elif pot: force_k = self.evaluate_k_force_k(disp_k) potential = self.evaluate_elastic_energy_k_space(force_k, disp_k) return potential, force_k
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))