Coef
Descriptions of the contents of the Coef class, which cotntains matrices to perform operation on SEM data.
- class pysemtools.datatypes.coef.Coef(msh, comm, get_area=False, apply_1d_operators=True, bckend='numpy')
Class that contains arrays like mass matrix, jacobian, jacobian inverse, etc.
This class can be used when mathematical operations such as derivation and integration is needed on the sem mesh.
- Parameters:
- mshMesh
Mesh object.
- commComm
MPI comminicator object.
- get_areabool, optional
If True, the area integration weight and normal vectors will be calculated. (Default value = True).
- apply_1d_operatorsbool, optional
If True, the 1D operators will be applied instead of building 3D operators. (Default value = True).
- Attributes:
- drdxndarray
component [0,0] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- drdyndarray
component [0,1] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- drdzndarray
component [0,2] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- dsdxndarray
component [1,0] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- dsdyndarray
component [1,1] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- dsdzndarray
component [1,2] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- dtdxndarray
component [2,0] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- dtdyndarray
component [2,1] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- dtdzndarray
component [2,2] of the jacobian inverse tensor for each point. shape is (nelv, lz, ly, lx).
- Bndarray
Mass matrix for each point. shape is (nelv, lz, ly, lx).
- areandarray
Area integration weight for each point in the facets. shape is (nelv, 6, ly, lx).
- nxndarray
x component of the normal vector for each point in the facets. shape is (nelv, 6, ly, lx).
- nyndarray
y component of the normal vector for each point in the facets. shape is (nelv, 6, ly, lx).
- nzndarray
z component of the normal vector for each point in the facets. shape is (nelv, 6, ly, lx).
Methods
dssum
(field, msh)Peform average of given field over shared points in each rank.
dudrst
(field[, direction])Perform derivative with respect to reference coordinate r/s/t.
dudrst_1d_operator
(field, dr, ds[, dt])Perform derivative applying the 1d operators provided as inputs.
dudrst_1d_operator_torch
(field, dr, ds[, dt])Perform derivative applying the 1D operators provided as inputs.
dudrst_3d_operator
(field, dr)Perform derivative with respect to reference coordinate r.
dudxyz
(field, drdx, dsdx[, dtdx])Perform derivative with respect to physical coordinate x,y,z.
glsum
(a, comm[, dtype])Peform global summatin of given qunaitity a using MPI.
- Returns:
Examples
Assuming you have a mesh object and MPI communicator object, you can initialize the Coef object as follows:
>>> from pysemtools import Coef >>> coef = Coef(msh, comm)
- dssum(field, msh)
Peform average of given field over shared points in each rank.
This method averages the field over shared points in the same rank. It uses the connectivity data in the mesh object. dssum might be a missleading name.
- Parameters:
- fieldndarray
Field to average over shared points.
- mshMesh
pySEMTools Mesh object.
- Returns:
- ndarray
Input field with shared points averaged with shared points in the SAME rank.
Examples
Assuming you have a Coef object and are working on a 3d field:
>>> dudx = coef.dssum(dudx, msh)
- dudrst(field, direction='r')
Perform derivative with respect to reference coordinate r/s/t.
Used to perform the derivative in the reference coordinates
- Parameters:
- fieldndarray
Field to take derivative of. Shape should be (nelv, lz, ly, lx).
- directionstr
Direction to take the derivative. Can be ‘r’, ‘s’, or ‘t’. (Default value = ‘r’).
- Returns:
- ndarray
Derivative of the field with respect to r/s/t. Shape is the same as the input field.
- dudrst_1d_operator(field, dr, ds, dt=None)
Perform derivative applying the 1d operators provided as inputs.
This method uses the 1D operators to apply the derivative. To apply derivative in r. you mush provide the 1d differenciation matrix in that direction and identity in the others.
- Parameters:
- fieldndarray
Field to take derivative of. Shape should be (nelv, lz, ly, lx).
- drndarray
Derivative matrix in the r direction to apply to each element. Shape should be (lx, lx).
- dsndarray
Derivative matrix in the s direction to apply to each element. Shape should be (ly, ly).
- dtndarray
Derivative matrix in the t direction to apply to each element. Shape should be (lz, lz). This is optional. If none is passed, it is assumed that the field is 2D.
- Returns:
- ndarray
Derivative of the field with respect to r/s/t. Shape is the same as the input field.
Examples
Assuming you have a Coef object
>>> dxdr = coef.dudrst(x, coef.dr, np.eye(ly, dtype=coef.dtype), np.eye(lz, dtype=coef.dtype))
- dudrst_1d_operator_torch(field, dr, ds, dt=None)
Perform derivative applying the 1D operators provided as inputs.
- Parameters:
- fieldtorch.Tensor
Field to take derivative of. Shape should be (nelv, lz, ly, lx).
- drtorch.Tensor
Derivative matrix in the r direction to apply to each element. Shape should be (lx, lx).
- dstorch.Tensor
Derivative matrix in the s direction to apply to each element. Shape should be (ly, ly).
- dttorch.Tensor (optional)
Derivative matrix in the t direction. Shape should be (lz, lz).
- Returns:
- torch.Tensor
Derivative of the field with respect to r/s/t. Shape is the same as the input field.
- dudrst_3d_operator(field, dr)
Perform derivative with respect to reference coordinate r.
This method uses derivation matrices from the lagrange polynomials at the GLL points.
- Parameters:
- fieldndarray
Field to take derivative of. Shape should be (nelv, lz, ly, lx).
- drndarray
Derivative matrix in the r/s/t direction to apply to each element. Shape should be (lx*ly*lz, lx*ly*lz).
- Returns:
- ndarray
Derivative of the field with respect to r/s/t. Shape is the same as the input field.
Examples
Assuming you have a Coef object
>>> dxdr = coef.dudrst(x, coef.dr)
- dudxyz(field, drdx, dsdx, dtdx=None)
Perform derivative with respect to physical coordinate x,y,z.
This method uses the chain rule, first evaluating derivatives with respect to rst, then multiplying by the inverse of the jacobian to map to xyz.
- Parameters:
- fieldndarray
Field to take derivative of. Shape should be (nelv, lz, ly, lx).
- drdxndarray
Derivative of the reference coordinates with respect to x, i.e., first entry in the appropiate row of the jacobian inverse. Shape should be the same as the field.
- dsdxndarray
Derivative of the reference coordinates with respect to y, i.e., second entry in the appropiate row of the jacobian inverse. Shape should be the same as the field.
- dtdxndarray
Derivative of the reference coordinates with respect to z, i.e., third entry in the appropiate row of the jacobian inverse. Shape should be the same as the field. (Default value = None) Only valid for 3D fields.
- Returns:
- ndarray
Derivative of the field with respect to x,y,z. Shape is the same as the input field.
Examples
Assuming you have a Coef object and are working on a 3d field:
>>> dudx = coef.dudxyz(u, coef.drdx, coef.dsdx, coef.dtdx)
- glsum(a, comm, dtype=<class 'numpy.float64'>)
Peform global summatin of given qunaitity a using MPI.
This method uses MPI to sum over all MPI ranks. It works with any numpy array shape and returns one value.
- Parameters:
- andarray
Quantity to sum over all mpiranks.
- commComm
MPI communicator object.
- datypenumpy.dtype
(Default value = np.double).
- Returns:
- float
Sum of the quantity a over all MPI ranks.
Examples
Assuming you have a Coef object and are working on a 3d field:
>>> volume = coef.glsum(coef.B, comm)