Source code for Mordicus.Modules.Safran.Containers.Loadings.Dirichlet

# -*- coding: utf-8 -*-
import os
from mpi4py import MPI
if MPI.COMM_WORLD.Get_size() > 1: # pragma: no cover
    os.environ["OMP_NUM_THREADS"] = "1"
    os.environ["OPENBLAS_NUM_THREADS"] = "1"
    os.environ["MKL_NUM_THREADS"] = "1"
    os.environ["VECLIB_MAXIMUM_THREADS"] = "1"
    os.environ["NUMEXPR_NUM_THREADS"] = "1"
import numpy as np

from Mordicus.Core.Containers.Loadings.LoadingBase import LoadingBase
#import collections


[docs]class Dirichlet(LoadingBase): """ Class containing a Loading of type Dirichlet boundary condition. Attributes ---------- function : function returning the value being imposed for the loading, with respect to the spatial position assembledBC : numpy.ndarray of size (numberOfModes,) containing the precomputed vector "red" for deriving the Dirichlet reduced external forces contribution is the form: tgv*red tgv : float very large value, used for weak enforcement of the loading using penalisation """ def __init__(self, solutionName, set): assert isinstance(set, str) assert isinstance(solutionName, str) super(Dirichlet, self).__init__(solutionName, set, "dirichlet") self.function = "" self.assembledBC = None self.tgv = 1.e3
[docs] def SetFunction(self, function): """ Sets function Parameters ---------- function : function returning the value being imposed for the loading, with respect to the spatial position """ self.function = function
[docs] def GetAssembledReducedFieldAtTime(self): """ Returns assembledBC Returns ------- numpy.ndarray of size (numberOfModes,) containing the precomputed vector "red" for deriving the Dirichlet reduced external forces contribution is the form: tgv*red """ return self.assembledBC
[docs] def ReduceLoading(self, mesh, problemData, reducedOrderBases, operatorCompressionData): """ Computes and sets the reduced representation of the loading Parameters ---------- mesh : BasicTools.Containers.UnstructuredMesh mesh of the high-fidelity model problemData : ProblemData problemData containing the loading reducedOrderBases : dict(str: np.ndarray) dictionary with solutionNames (str) as keys and reducedOrderBases (np.ndarray of size (numberOfModes, numberOfDOFs)) as values operatorCompressionData : dict(str: custom_data_structure) not used in this loading dictionary with solutionNames (str) as keys and data structure generated by the operator compression step as values """ from Mordicus.Modules.Safran.FE import FETools as FT integrationWeightsSet, phiAtIntegPointSet = FT.ComputePhiAtIntegPoint(mesh, [self.set], -1) numberOfNodes = mesh.GetNumberOfNodes() positionIntegPointsOnSet = phiAtIntegPointSet.dot(mesh.GetNodes()) #print("mesh.GetNodes() =", mesh.GetNodes()) numberOfComponents = reducedOrderBases[self.solutionName].shape[1]//numberOfNodes numberOfModes = reducedOrderBases[self.solutionName].shape[0] numberOfIntegrationPointsSet = phiAtIntegPointSet.shape[0] assert numberOfComponents == self.function(positionIntegPointsOnSet[0]).shape[0], 'dirichlet condition do not have the same number of components as the provided reducedOrderBasis' valueAtIntegPointSet = np.empty((numberOfIntegrationPointsSet,numberOfComponents)) for i in range(numberOfIntegrationPointsSet): valueAtIntegPointSet[i,:] = self.function(positionIntegPointsOnSet[i]) componentReducedOrderBasis = [] for i in range(numberOfComponents): componentReducedOrderBasis.append(reducedOrderBases[self.solutionName][:,i*numberOfNodes:(i+1)*numberOfNodes].T) reducedPhiAtIntegPointSet = np.empty((numberOfComponents,numberOfIntegrationPointsSet,numberOfModes)) for i in range(numberOfComponents): reducedPhiAtIntegPointSet[i] = phiAtIntegPointSet.dot(componentReducedOrderBasis[i]) #print(np.einsum("ti,ti,t", valueAtIntegPointSet, valueAtIntegPointSet, integrationWeightsSet, optimize = True)) #print(np.einsum("ti,itj,t->j", valueAtIntegPointSet, reducedPhiAtIntegPointSet, integrationWeightsSet, optimize = True)) self.assembledBC = np.einsum("ti,itj,t->j", valueAtIntegPointSet, reducedPhiAtIntegPointSet, integrationWeightsSet, optimize = True)
#def HyperReduceLoading(self, mesh, problemData, reducedOrderBases, operatorCompressionData): # return
[docs] def ComputeContributionToReducedExternalForces(self, time): """ Computes and returns the reduced external forces contribution of the loading Parameters ---------- time : float Returns ------- np.ndarray of size (numberOfModes,) """ # assert type of time assert isinstance(time, (float, np.float64)) return self.tgv*self.GetAssembledReducedFieldAtTime()
def __getstate__(self): state = {} state["solutionName"] = self.solutionName state["set"] = self.set state["type"] = self.type state["function"] = self.function return state def __str__(self): res = "Dirichlet Loading with set "+self.GetSet() return res
if __name__ == "__main__":# pragma: no cover from Mordicus import RunTestFile RunTestFile(__file__)