warp.fem#
warp.fem
模块旨在促进求解描述为微分方程的物理系统。例如,它可以求解基于有限元 (FEM) Galerkin 方法的扩散、对流、流体流动和弹性问题的 PDE,并允许用户快速试验各种 FEM 公式和离散化方案。
积分算子#
FEM 工具包的核心功能是能够在各种域上并使用任意插值基函数积分常数、线性和双线性形式。
主要机制是 integrand()
装饰器,例如
@integrand
def linear_form(
s: Sample,
domain: Domain,
v: Field,
):
x = domain(s)
return v(s) * wp.max(0.0, 1.0 - wp.length(x))
@integrand
def diffusion_form(s: Sample, u: Field, v: Field, nu: float):
return nu * wp.dot(
grad(u, s),
grad(v, s),
)
积分算子是普通的 Warp 内核,这意味着它们可以包含任意 Warp 函数。但是,它们接受一些特殊的参数
Sample
包含有关当前积分采样点的信息,例如单元索引和单元坐标。
Field
指定一个抽象字段,它将在调用时被实际的字段类型替换,例如离散字段、field.TestField
或field.TrialField
,定义在某个FunctionSpace
上,一个包装任意函数的ImplicitField
,或者任何其他可用的 字段。然后可以使用通常的调用运算符u(s)
在给定的样本 s 处评估字段 u。还有其他几个运算符可用,例如梯度grad()
;请参阅 运算符 部分。
Domain
指定一个抽象积分域。在样本 s 处评估域,如domain(s)
产生相应的世界位置,并且还为域提供了一些运算符,例如在给定样本处评估法线@integrand def boundary_form( s: Sample, domain: Domain, u: Field, ): nor = normal(domain, s) return wp.dot(u(s), nor)
积分算子不能直接与 warp.launch()
一起使用,而必须通过 integrate()
或 interpolate()
调用。Sample
和 Domain
根积分算子的参数(传递给 integrate()
或 interpolate()
调用的 integrand 参数)将自动填充。Field
参数必须作为字典在启动器函数的 fields 参数中传递,所有其他标准 Warp 类型参数必须作为字典在启动器函数的 values 参数中传递,例如
integrate(diffusion_form, fields={"u": trial, "v": test}, values={"nu": viscosity})
基本工作流程#
使用 warp.fem
求解线性化 PDE 的典型步骤如下
定义一个
Geometry
(网格、网状结构等)。目前,支持 2D 和 3D 规则网格、NanoVDB 体以及三角形、四边形、四面体和六面体非结构化网格。定义一个或多个
FunctionSpace
,通过为几何元素配备形函数。请参阅make_polynomial_space()
。目前,支持连续/非连续拉格朗日 (\(P_{k[d]}, Q_{k[d]}\)) 和 Serendipity (\(S_k\)) 阶数为 \(k \leq 3\) 的形函数,以及线性 Nédélec(第一类)和 Raviart-Thomas 向量值形函数。定义一个积分域,例如几何体的单元格 (
Cells
) 或边界侧 (BoundarySides
)。积分线性形式以构建系统的右侧。使用
make_test()
在函数空间上定义一个测试函数,使用Quadrature
公式(或让模块根据函数空间度数选择一个),并使用线性形式积分算子调用integrate()
。结果是一个warp.array
,其中包含每个函数空间自由度的积分结果。积分双线性形式以构建系统的左侧。使用
make_trial()
在函数空间上定义一个试验函数,然后使用双线性形式积分算子调用integrate()
。结果是一个warp.sparse.BsrMatrix
,其中包含每对测试和试验函数空间自由度的积分结果。请注意,试验函数和测试函数不必定义在同一函数空间上,因此支持混合 FEM。使用您选择的求解器求解生成的线性系统,例如内置的 迭代线性求解器 之一。
以下摘录自介绍性示例 warp/examples/fem/example_diffusion.py
概述了此过程
# Grid geometry
geo = Grid2D(n=50, cell_size=2)
# Domain and function spaces
domain = Cells(geometry=geo)
scalar_space = make_polynomial_space(geo, degree=3)
# Right-hand-side (forcing term)
test = make_test(space=scalar_space, domain=domain)
rhs = integrate(linear_form, fields={"v": test})
# Weakly-imposed boundary conditions on Y sides
boundary = BoundarySides(geo)
bd_test = make_test(space=scalar_space, domain=boundary)
bd_trial = make_trial(space=scalar_space, domain=boundary)
bd_matrix = integrate(y_mass_form, fields={"u": bd_trial, "v": bd_test})
# Diffusion form
trial = make_trial(space=scalar_space, domain=domain)
matrix = integrate(diffusion_form, fields={"u": trial, "v": test}, values={"nu": viscosity})
# Assemble linear system (add diffusion and boundary condition matrices)
matrix += bd_matrix * boundary_strength
# Solve linear system using Conjugate Gradient
x = wp.zeros_like(rhs)
bsr_cg(matrix, b=rhs, x=x)
注意
integrate()
函数不检查传递的积分算子是否真的是线性或双线性形式;由用户来确保它们是。为了求解非线性 PDE,可以使用迭代过程,并将研究函数的当前值 DiscreteField
参数传递给积分算子,其中允许任意操作。但是,形式的结果必须在测试和试验字段中保持线性。该策略在 example_mixed_elasticity.py
示例中得到了演示。
入门示例#
warp.fem
附带了 warp/examples/fem
目录中的示例列表,演示了如何解决经典模型问题。
example_diffusion.py
:具有齐次 Neumann 和 Dirichlet 边界条件的 2D 扩散
example_diffusion_3d.py
:扩散问题的 3D 变体
example_convection_diffusion.py
:使用半拉格朗日平流的 2D 对流扩散
example_convection_diffusion_dg.py
:使用非连续伽辽金与逆风传输和对称内部惩罚的 2D 对流扩散
example_burgers.py
:使用非连续伽辽金与逆风传输和斜率限制器的 2D 无粘性 Burgers 方程
example_stokes.py
:使用混合 \(P_k/P_{k-1}\) 或 \(Q_k/P_{(k-1)d}\) 元素的 2D 不可压缩 Stokes 流
example_navier_stokes.py
:使用混合 \(P_k/P_{k-1}\) 元素的 2D Navier-Stokes 流
example_mixed_elasticity.py
:使用混合连续/非连续 \(S_k/P_{(k-1)d}\) 元素的 2D 非线性弹性
example_distortion_energy.py
: 3D 表面参数化,使其 2D 非线性畸变能量最小化
example_magnetostatics.py
: 使用旋度-旋度公式的 2D 静磁学
高级用法#
高阶(弯曲)几何体#
可以将任何 Geometry
(网格和显式网格)转换为弯曲的高阶变体,方法是使用 make_deformed_geometry()
方法,用任意阶次的位移场对其进行变形。该过程如下所示
# Define a base geometry
base_geo = fem.Grid3D(res=resolution)
# Define a displacement field on the base geometry
deformation_space = fem.make_polynomial_space(base_geo, degree=deformation_degree, dtype=wp.vec3)
deformation_field = deformation_space.make_field()
# Populate the field value by interpolating an expression
fem.interpolate(deformation_field_expr, dest=deformation_field)
# Construct the deformed geometry from the displacement field
deform_geo = deformation_field.make_deformed_geometry()
# Define new function spaces on the deformed geometry
scalar_space = fem.make_polynomial_space(deformed_geo, degree=scalar_space_degree)
有关完整示例,请参见 example_deformed_geometry.py
。也可以从 ImplicitField
定义变形场,如 example_magnetostatics.py
中所示。
基于粒子的正交积分和位置查找#
全局 lookup()
算子允许从任意位置生成一个 Sample
; example_streamlines.py
示例演示了通过追踪速度场来生成 3D 流线。
此算子也被 PicQuadrature
利用,以提供一种从一组任意粒子定义 Particle-In-Cell 正交积分的方法,从而可以实现 MPM 类型的方法。在初始化正交积分时,粒子会自动分组到几何单元中。example_stokes_transfer.py
和 example_apic_fluid.py
示例说明了这一点。
注意
当前,全局 lookup()
算子不支持 Quadmesh2D
、Hexmesh
和变形的几何体。
非一致性场#
在给定的 Geometry
上定义的场不能直接用于对不同的几何体进行积分;但是,可以将它们包装在 NonconformingField
中以达到此目的。example_nonconforming_contact.py
利用这一点来模拟单独离散化的接触体。
注意
当前,NonconformingField
不支持包装试验场,因此尚无法在不同的几何体上定义双线性形式。
注意
不同几何体之间的映射是基于位置的,因此 NonconformingField
无法准确捕获不连续函数空间。此外,积分域必须支持 lookup()
算子。
分区#
FEM 工具包可以对域元素的子集执行积分,并可能重新索引自由度,以便线性系统仅包含局部自由度。这对于分布式计算(请参见 warp/examples/fem/example_diffusion_mgpu.py
)或仅将仿真域限制为活动单元的子集(请参见 warp/examples/fem/example_stokes_transfer.py
)非常有用。
可以使用 GeometryPartition
的子类(例如 LinearGeometryPartition
或 ExplicitGeometryPartition
)来定义仿真几何体的分区。
然后可以使用 make_space_partition()
根据几何分区对函数空间进行分区。生成的 SpacePartition
对象允许在空间范围和分区范围的节点索引之间进行转换,并区分内部节点、边界节点和外部节点。
可以使用 Subdomain
类在元素的子集上进行积分,同时保留完整的自由度集,即,无需重新索引; example_streamlines.py
示例说明了这一点,用于定义流入和流出边界。
自适应性#
虽然非结构化网格细化当前不在范围内,但 warp.fem
提供了稀疏网格几何体的自适应版本 AdaptiveNanogrid
,具有二的幂次体素比例。还提供了用于从网格层次结构或细化预言构建此类几何体的助手,请参见 adaptive_nanogrid_from_field()
和 adaptive_nanogrid_from_hierarchy()
。在 warp/examples/fem/example_adaptive_grid.py
中提供了一个示例。
注意
分辨率边界处“T 形结”的存在意味着通常的三项式形状函数将不再全局连续。可以使用不连续伽辽金或类似技术来考虑多分辨率面上的“跳跃”。
内存管理#
一些 warp.fem
函数需要分配临时缓冲区才能执行其计算。如果在紧密循环中多次调用此类函数,则许多分配和取消分配可能会降低性能,尽管当使用 流排序内存池分配器 时,这种情况的严重程度会大大降低。为了解决这个问题,可以创建一个 cache.TemporaryStore
对象来持久化和重用调用之间的临时分配,方法是全局使用 set_default_temporary_store()
,或在每个函数的粒度上使用相应的参数。
可视化#
大多数函数空间都定义了一个 FunctionSpace.cells_to_vtk()
方法,该方法返回 VTK 兼容的单元类型和节点索引的列表。这可以用于在支持 VTK 的查看器(例如 pyvista
)中可视化离散场,例如
import numpy as np
import pyvista
import warp as wp
import warp.fem as fem
@fem.integrand
def ackley(s: fem.Sample, domain: fem.Domain):
x = domain(s)
return (
-20.0 * wp.exp(-0.2 * wp.sqrt(0.5 * wp.length_sq(x)))
- wp.exp(0.5 * (wp.cos(2.0 * wp.pi * x[0]) + wp.cos(2.0 * wp.pi * x[1])))
+ wp.e
+ 20.0
)
# Define field
geo = fem.Grid2D(res=wp.vec2i(64, 64), bounds_lo=wp.vec2(-4.0, -4.0), bounds_hi=wp.vec2(4.0, 4.0))
space = fem.make_polynomial_space(geo, degree=3)
field = space.make_field()
fem.interpolate(ackley, dest=field)
# Extract cells, nodes and values
cells, types = field.space.vtk_cells()
nodes = field.space.node_positions().numpy()
values = field.dof_values.numpy()
positions = np.hstack((nodes, values[:, np.newaxis]))
# Visualize with pyvista
grid = pyvista.UnstructuredGrid(cells, types, positions)
grid.point_data["scalars"] = values
plotter = pyvista.Plotter()
plotter.add_mesh(grid)
plotter.show()
算子#
域算子#
- warp.fem.normal(domain: Domain, s: Sample)#
计算样本点 s 处的单元法线。如果该单元是边,或者几何体嵌入到更高维度的空间中(例如
Trimesh3D
),则非零
- warp.fem.lookup(domain: Domain, x)#
查找与世界位置 x 对应的样本点,并投影到域上最接近的点。
- 参数:
x – 要在几何体中查找的点的世界位置
guess – (可选)
Sample
初始猜测,可能有助于执行查询
注意
目前,此运算符不支持
Hexmesh
、Quadmesh2D
、Quadmesh3D
和变形几何体。
场算子#
积分#
- warp.fem.integrate(
- 被积函数,
- domain=None,
- quadrature=None,
- nodal=False,
- fields=None,
- values=None,
- accumulate_dtype=wp.float64,
- output_dtype=None,
- output=None,
- device=None,
- temporary_store=None,
- kernel_options=None,
- assembly=None,
- add=False,
- bsr_options=None,
对常数、线性或双线性形式进行积分,并分别返回标量、数组或稀疏矩阵。
- 参数:
integrand (Integrand) – 要积分的形式,必须具有
integrand()
装饰器domain (GeometryDomain | None) – 积分域。如果为 None,则从场推导
quadrature (Quadrature | None) – 正交公式。如果为 None,则从域和场阶数推导。
nodal (bool) – 已弃用。请改用等效的 assembly=”nodal”。
fields (Dict[str, FieldLike] | None) – 要传递给被积函数的离散场、测试场和试探场。字典中的键必须与被积函数参数名称匹配。
values (Dict[str, Any] | None) – 要传递给被积函数的其他变量值,可以是 warp 内核启动接受的任何类型。字典中的键必须与被积函数参数名称匹配。
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
accumulate_dtype (type) – 用于累积积分样本的标量类型
output_dtype (type | None) – 如果未提供 output,则返回结果的标量类型。如果为 None,则默认为 accumulate_dtype
device – 执行积分的设备
kernel_options (Dict[str, Any] | None) – 要传递给内核构建器的重载选项(例如,
{"enable_backward": True}
)assembly (str | None) – 指定组装积分向量或矩阵的策略: - “nodal”: 对于线性或双线性形式,使用测试函数节点作为正交点。假定使用拉格朗日插值函数,并且不在测试或试探函数上评估微分或 DG 算子。 - “generic”: 单遍积分和形状函数评估。不假设被积函数的内容,但可能导致许多冗余计算。 - “dispatch”: 对于线性或双线性形式,首先在正交点评估该形式,然后在第二遍中分派到节点。对于评估成本高的被积函数更有效。与测试或试探函数上的 at_node 算子不兼容。 - None(默认):自动选择合适的组装策略(“generic” 或 “dispatch”)
add (bool) – 如果为 True 且提供了 output,则将积分结果添加到 output 而不是替换其内容
bsr_options (Dict[str, Any] | None) – 要传递给稀疏矩阵构造算法的其他选项。请参阅
warp.sparse.bsr_set_from_triplets()
- warp.fem.interpolate(
- 被积函数,
- dest=None,
- quadrature=None,
- dim=0,
- domain=None,
- fields=None,
- values=None,
- device=None,
- kernel_options=None,
- temporary_store=None,
- bsr_options=None,
在有限的采样点集上插值一个函数,并可选择将结果分配给离散场或原始 warp 数组。
- 参数:
integrand (Integrand | FieldLike) – 要插值的函数:具有
warp.fem.integrand()
装饰器的函数或字段dest (DiscreteField | FieldRestriction | array | None) –
存储插值结果的位置。可以是
一个
DiscreteField
,或者离散场到域的限制(来自make_restriction()
)。在这种情况下,将在每个节点上执行插值。一个普通的 warp
array
,或者None
。在这种情况下,插值样本将由 quadrature 或 dim 参数决定,按此顺序。
quadrature (Quadrature | None) – 如果 dest 不是离散场或场限制,则定义插值样本的正交公式。
dim (int) – 如果 dest 不是离散场或限制并且 quadrature 是
None
,则插值样本的数量。在这种情况下,传递给 integrand 的Sample
将无效,但样本点索引s.qp_index
可用于定义自定义插值逻辑。domain (Domain | None) – 插值域,仅当 dest 不是场限制且 quadrature 为
None
时使用fields (Dict[str, FieldLike] | None) – 要传递给被积函数的离散场。字典中的键必须与被积函数参数名称匹配。
values (Dict[str, Any] | None) – 要传递给被积函数的其他变量值,可以是 warp 内核启动接受的任何类型。字典中的键必须与被积函数参数名称匹配。
device – 执行插值的设备
kernel_options (Dict[str, Any] | None) – 要传递给内核构建器的重载选项(例如,
{"enable_backward": True}
)temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
bsr_options (Dict[str, Any] | None) – 要传递给稀疏矩阵构造算法的其他选项。请参阅
warp.sparse.bsr_set_from_triplets()
- class warp.fem.Sample#
用于在被积函数中评估字段和相关算子的每个样本点上下文。
几何体#
- class warp.fem.Grid2D(res, bounds_lo=None, bounds_hi=None)[source]#
基类:
Geometry
二维规则网格几何体
- 参数:
res (vec2i)
bounds_lo (vec2f | None)
bounds_hi (vec2f | None)
- class warp.fem.Trimesh2D(
- tri_vertex_indices,
- positions,
- build_bvh=False,
- temporary_store=None,
基类:
Trimesh
2D 三角网格几何体
- 参数:
tri_vertex_indices (array)
positions (array)
build_bvh (bool)
temporary_store (TemporaryStore | None)
- __init__(
- tri_vertex_indices,
- positions,
- build_bvh=False,
- temporary_store=None,
构造一个 D 维三角形网格。
- 参数:
tri_vertex_indices (array) – 形状为 (num_tris, 3) 的 warp 数组,包含每个三角形的顶点索引
positions (array) – 形状为 (num_vertices, D) 的 warp 数组,包含每个顶点的位置
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
build_bvh (bool) – 是否同时构建三角形 BVH,这是全局 fem.lookup 算子在没有初始猜测的情况下运行所必需的
- class warp.fem.Trimesh3D(
- tri_vertex_indices,
- positions,
- build_bvh=False,
- temporary_store=None,
基类:
Trimesh
3D 三角网格几何体
- 参数:
tri_vertex_indices (array)
positions (array)
build_bvh (bool)
temporary_store (TemporaryStore | None)
- __init__(
- tri_vertex_indices,
- positions,
- build_bvh=False,
- temporary_store=None,
构造一个 D 维三角形网格。
- 参数:
tri_vertex_indices (array) – 形状为 (num_tris, 3) 的 warp 数组,包含每个三角形的顶点索引
positions (array) – 形状为 (num_vertices, D) 的 warp 数组,包含每个顶点的位置
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
build_bvh (bool) – 是否同时构建三角形 BVH,这是全局 fem.lookup 算子在没有初始猜测的情况下运行所必需的
- class warp.fem.Quadmesh2D(quad_vertex_indices, positions, temporary_store=None)[source]#
基类:
Quadmesh
二维四边形网格
- 参数:
quad_vertex_indices (array)
positions (array)
temporary_store (TemporaryStore | None)
- __init__(
- quad_vertex_indices,
- positions,
- temporary_store=None,
构造一个 D 维四边形网格。
- 参数:
quad_vertex_indices (array) – 形状为 (num_tris, 4) 的 warp 数组,包含每个四边形的顶点索引,按逆时针顺序排列
positions (array) – 形状为 (num_vertices, D) 的 warp 数组,包含每个顶点的位置
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
- class warp.fem.Quadmesh3D(quad_vertex_indices, positions, temporary_store=None)[source]#
基类:
Quadmesh
三维四边形网格
- 参数:
quad_vertex_indices (array)
positions (array)
temporary_store (TemporaryStore | None)
- __init__(
- quad_vertex_indices,
- positions,
- temporary_store=None,
构造一个 D 维四边形网格。
- 参数:
quad_vertex_indices (array) – 形状为 (num_tris, 4) 的 warp 数组,包含每个四边形的顶点索引,按逆时针顺序排列
positions (array) – 形状为 (num_vertices, D) 的 warp 数组,包含每个顶点的位置
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
- class warp.fem.Grid3D(res, bounds_lo=None, bounds_hi=None)[source]#
基类:
Geometry
三维规则网格几何体
- 参数:
res (vec3i)
bounds_lo (vec3f | None)
bounds_hi (vec3f | None)
- class warp.fem.Tetmesh(
- tet_vertex_indices,
- positions,
- build_bvh=False,
- temporary_store=None,
基类:
Geometry
四面体网格几何体
- 参数:
tet_vertex_indices (array)
positions (array)
build_bvh (bool)
temporary_store (TemporaryStore | None)
- __init__(
- tet_vertex_indices,
- positions,
- build_bvh=False,
- temporary_store=None,
构造一个四面体网格。
- 参数:
tet_vertex_indices (array) – 形状为 (num_tets, 4) 的 warp 数组,包含每个四面体的顶点索引
positions (array) – 形状为 (num_vertices, 3) 的 warp 数组,包含每个顶点的 3d 位置
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
build_bvh (bool) – 是否同时构建四面体 BVH,这是全局 fem.lookup 算子在没有初始猜测的情况下运行所必需的
- class warp.fem.Hexmesh(hex_vertex_indices, positions, temporary_store=None)[source]#
基类:
Geometry
六面体网格几何体
- 参数:
hex_vertex_indices (array)
positions (array)
temporary_store (TemporaryStore | None)
- __init__(hex_vertex_indices, positions, temporary_store=None)[source]#
构造一个四面体网格。
- 参数:
hex_vertex_indices (array) – 形状为 (num_hexes, 8) 的 warp 数组,包含每个六面体的顶点索引,遵循标准顺序(底面顶点按逆时针顺序排列,然后是顶面顶点)
positions (array) – 形状为 (num_vertices, 3) 的 warp 数组,包含每个顶点的 3d 位置
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
- class warp.fem.Nanogrid(grid, temporary_store=None)[source]#
基类:
Geometry
稀疏网格几何体
- 参数:
grid (Volume)
temporary_store (TemporaryStore | None)
- __init__(grid, temporary_store=None)[source]#
从内存中的 NanoVDB 体积构造一个稀疏网格几何体。
- 参数:
grid (Volume) – NanoVDB 体积。接受任何类型,但为了索引效率,建议使用索引网格。如果 grid 是一个 “on” 索引网格,则仅为活动体素创建单元格,否则将为所有叶子体素创建单元格。
temporary_store (TemporaryStore | None) – 用于分配临时数组的共享池
- class warp.fem.AdaptiveNanogrid(cell_grid, cell_level, level_count, temporary_store)[source]#
基类:
Geometry
自适应稀疏网格
- 参数:
cell_grid (Volume)
cell_level (array)
level_count (int)
temporary_store (TemporaryStore)
- __init__(
- cell_grid,
- cell_level,
- level_count,
- temporary_store,
从内存中的 NanoVDB 体积和级别列表构造一个自适应稀疏网格几何体。
不建议直接使用此构造函数;请参阅辅助函数
warp.fem.adaptive_nanogrid_from_field()
和warp.fem.adaptive_nanogrid_from_hierarchy()
- 参数:
cell_grid (Volume) – 一个 warp 体积(理想情况下由索引网格支持),其体素坐标对应于每个单元格的最低精细分辨率体素。单元格的范围由 cell_level 数组给出。例如,坐标
ijk
处的体素和级别0
对应于相同坐标处的精细单元格,坐标2*ijk
处的体素和级别1
对应于一个单元格,该单元格跨越从2*ijk
到2*ijk + (1,1,1)
的2^3
体素,依此类推。cell_level (array) – 体积中每个体素的细化级别。级别 0 是最精细的,级别
level_count-1
是最粗糙的。level_count (int) – 网格中的级别数
temporary_store (TemporaryStore)
- class warp.fem.LinearGeometryPartition(
- geometry,
- partition_rank,
- partition_count,
- device=None,
- temporary_store=None,
- 参数:
geometry (Geometry)
partition_rank (int)
partition_count (int)
temporary_store (TemporaryStore)
- class warp.fem.ExplicitGeometryPartition(geometry, cell_mask, temporary_store=None)[source]#
- 参数:
geometry (Geometry)
cell_mask (wp.array(dtype=int))
temporary_store (TemporaryStore)
- __init__(
- geometry,
- cell_mask,
- temporary_store=None,
通过均匀分割单元索引创建几何分区
- 参数:
geometry (Geometry) – 要分区的几何体
cell_mask (array(ndim=1, dtype=int32)) – 长度为
geometry.cell_count()
的 warp 数组,指示选择了哪些单元格。数组值必须为1
(已选择)或0
(未选择)。temporary_store (TemporaryStore)
- class warp.fem.Cells(geometry)[source]#
基类:
GeometryDomain
一个包含几何体或几何分区的所有单元格的域
- 参数:
geometry (Geometry | GeometryPartition)
- __init__(geometry)[source]#
- 参数:
geometry (Geometry | GeometryPartition)
- class warp.fem.Sides(geometry)[source]#
基类:
GeometryDomain
一个包含几何体或几何分区的所有(内部和边界)面的域
- 参数:
geometry (Geometry | GeometryPartition)
- __init__(geometry)[source]#
- 参数:
geometry (Geometry | GeometryPartition)
- class warp.fem.BoundarySides(geometry)[source]#
基类:
Sides
一个包含几何体或几何分区的边界面的域
- 参数:
geometry (Geometry | GeometryPartition)
- __init__(geometry)[source]#
- 参数:
geometry (Geometry | GeometryPartition)
- class warp.fem.FrontierSides(geometry)[source]#
基类:
Sides
一个包含几何分区的前沿面的域(与至少另一个分区共享的面)
- 参数:
geometry (Geometry | GeometryPartition)
- __init__(geometry)[source]#
- 参数:
geometry (Geometry | GeometryPartition)
- class warp.fem.Subdomain(
- domain,
- element_mask=None,
- element_indices=None,
- temporary_store=None,
基类:
GeometryDomain
子域 – 将域限制为其元素的子集
- 参数:
domain (GeometryDomain)
element_mask (array | None)
element_indices (array | None)
temporary_store (TemporaryStore | None)
- __init__(
- domain,
- element_mask=None,
- element_indices=None,
- temporary_store=None,
从元素的子集创建子域。
应提供 element_mask 和 element_indices 中的一个。
- 参数:
domain (GeometryDomain) – 包含域
element_mask (array | None) – 长度为
domain.element_count()
的数组,指示应包括哪些元素。数组值必须为1
(已选择)或0
(未选择)。element_indices (array | None) – 要包括的元素索引的显式数组
temporary_store (TemporaryStore | None)
- class warp.fem.Polynomial(*values)[source]#
定义区间上插值节点的Polynomial族
- GAUSS_LEGENDRE = 'GL'#
Gauss–Legendre 1D polynomial族 (不包括端点)
- LOBATTO_GAUSS_LEGENDRE = 'LGL'#
Lobatto–Gauss–Legendre 1D polynomial族 (包括端点)
- EQUISPACED_CLOSED = 'closed'#
闭合 1D polynomial 族,具有均匀分布的节点 (包括端点)
- EQUISPACED_OPEN = 'open'#
开放 1D polynomial 族,具有均匀分布的节点 (不包括端点)
- class warp.fem.RegularQuadrature(domain, order, family=None)[source]#
基类:
_QuadratureWithRegularEvaluationPoints
规则正交公式,每个单元使用一组固定的正交点
- 参数:
domain (GeometryDomain)
order (int)
family (Polynomial)
- __init__(domain, order, family=None)[source]#
- 参数:
domain (GeometryDomain)
order (int)
family (Polynomial)
- class warp.fem.NodalQuadrature(domain, space)[source]#
基类:
Quadrature
使用空间节点作为正交点的正交方法
请注意,与 nodal=True 标志用于
integrate()
不同,使用此正交方法并不意味着对形状函数的正交性有任何假设,因此可以安全地用于任意被积函数。- 参数:
domain (GeometryDomain | None)
space (FunctionSpace)
- __init__(domain, space)[source]#
- 参数:
domain (GeometryDomain | None)
space (FunctionSpace)
- class warp.fem.ExplicitQuadrature(domain, points, weights)[source]#
基类:
_QuadratureWithRegularEvaluationPoints
使用显式逐单元点和权重进行正交。
假设每个单元的正交点数量是恒定的,并从点和权重数组的形状推断出来。可以为整个几何体或仅为域的元素提供正交点。
- 参数:
domain (GeometryDomain) – 正交公式的定义域
points (wp.array2d(dtype=Coords)) – 形状为
(domain.element_count(), points_per_cell)
或(domain.geometry_element_count(), points_per_cell)
的二维数组,包含每个正交点的坐标。weights (wp.array2d(dtype=float)) – 形状为
(domain.element_count(), points_per_cell)
或(domain.geometry_element_count(), points_per_cell)
的二维数组,包含每个正交点的权重。
另请参阅:
PicQuadrature
- __init__(domain, points, weights)[source]#
- 参数:
domain (GeometryDomain)
points (array(ndim=2, dtype=vec3f))
weights (array(ndim=2, dtype=float32))
- class warp.fem.PicQuadrature(
- domain,
- positions,
- measures=None,
- requires_grad=False,
- temporary_store=None,
基类:
Quadrature
基于粒子的正交公式,使用全局的一组点不均匀地分布在几何元素上。
对于 Particle-In-Cell 和派生方法很有用。
- 参数:
domain (GeometryDomain) – 正交的底层域
positions (wp.array(dtype=wp.vecXd) | Tuple[wp.array(dtype=ElementIndex), wp.array(dtype=Coords)]) – 包含所有粒子世界位置的数组,或包含每个粒子的单元索引和坐标的数组的元组。请注意,前者要求底层几何体定义一个全局
Geometry.cell_lookup()
方法;目前仅适用于Grid2D
和Grid3D
。measures (wp.array(dtype=float) | None) – 包含每个粒子的度量(面积/体积)的数组,用于定义积分权重。如果
None
,则默认为单元格度量除以单元格中的粒子数。requires_grad (bool) – 是否应为计算量分配梯度
temporary_store (TemporaryStore) – 用于分配临时数组的共享池
- __init__(
- domain,
- positions,
- measures=None,
- requires_grad=False,
- temporary_store=None,
- 参数:
domain (GeometryDomain)
positions (wp.array(dtype=wp.vecXd) | Tuple[wp.array(dtype=ElementIndex), wp.array(dtype=Coords)])
measures (wp.array(dtype=float) | None)
requires_grad (bool)
temporary_store (TemporaryStore)
函数空间#
- warp.fem.make_polynomial_space(
- geo,
- dtype=float,
- dof_mapper=None,
- degree=1,
- element_basis=None,
- discontinuous=False,
- family=None,
为几何体配备一个并置的多项式函数空间。相当于依次调用
make_polynomial_basis_space()
然后调用 make_collocated_function_space, make_covariant_function_space 或 make_contravariant_function_space。- 参数:
geo (Geometry) – 构建空间的几何体
dtype (type) – 函数空间的值类型。如果提供了
dof_mapper
,则将使用来自 DofMapper 的值类型。dof_mapper ( DofMapper | None ) – 从节点自由度到函数值的映射,默认为 Identity。适用于降阶坐标,例如
SymmetricTensorMapper
将 2x2(或 3x3)对称张量映射到 3(或 6)个自由度。degree ( int ) – 每个元素的形函数的多项式阶数
discontinuous ( bool ) – 如果为 True,则使用不连续伽辽金形函数。如果阶数为 0,即分段常数形函数,则隐含不连续。
element_basis ( ElementBasis | None ) – 单个元素的基函数类型
family ( Polynomial | None ) – 用于生成形函数基的多项式族。如果未提供,将选择合理的基。
- 返回:
构造的函数空间
- 返回类型:
CollocatedFunctionSpace
- warp.fem.make_polynomial_basis_space(
- geo,
- degree=1,
- element_basis=None,
- discontinuous=False,
- family=None,
为几何体配备多项式基。
- 参数:
geo (Geometry) – 构建空间的几何体
degree ( int ) – 每个元素的形函数的多项式阶数
discontinuous ( bool ) – 如果为 True,则使用不连续伽辽金形函数。如果阶数为 0,即分段常数形函数,则隐含不连续。
element_basis ( ElementBasis | None ) – 单个元素的基函数类型
family ( Polynomial | None ) – 用于生成形函数基的多项式族。如果未提供,将选择合理的基。
- 返回:
构造的基空间
- 返回类型:
- warp.fem.make_collocated_function_space(
- basis_space,
- dtype=float,
- dof_mapper=None,
从标量值基空间和值类型构造函数空间,使得值类型的所有自由度都存储在每个基节点上。
- 参数:
geo – 构建空间的几何体
dtype (type) – 函数空间的值类型。如果提供了
dof_mapper
,则将使用来自 DofMapper 的值类型。dof_mapper ( DofMapper | None ) – 从节点自由度到函数值的映射,默认为 Identity。适用于降阶坐标,例如
SymmetricTensorMapper
将 2x2(或 3x3)对称张量映射到 3(或 6)个自由度。basis_space ( BasisSpace )
- 返回:
构造的函数空间
- 返回类型:
CollocatedFunctionSpace
- warp.fem.make_covariant_function_space(basis_space)[source]#
从向量值基空间构造协变函数空间
- 参数:
basis_space ( BasisSpace )
- 返回类型:
CovariantFunctionSpace
- warp.fem.make_contravariant_function_space(basis_space)[source]#
从向量值基空间构造逆变函数空间
- 参数:
basis_space ( BasisSpace )
- 返回类型:
ContravariantFunctionSpace
- warp.fem.make_space_partition(
- space=None,
- geometry_partition=None,
- space_topology=None,
- with_halo=True,
- device=None,
- temporary_store=None,
计算函数空间拓扑中与几何体分区接触的节点子集
必须提供 space_topology 或 space 之一(并将按该顺序考虑)。
- 参数:
space ( FunctionSpace | None ) – (已弃用)如果 space_topology 为
None
,则定义拓扑的函数空间。geometry_partition ( GeometryPartition | None ) – 空间几何体的子集。如果未提供,则使用整个几何体。
space_topology ( SpaceTopology | None ) – 要考虑的函数空间拓扑。如果
None
,则从 space 推导。with_halo ( bool ) – 如果为 True,则包含 halo 节点(从外部边界单元到分区的节点)
device – 执行和存储计算的 Warp 设备
temporary_store (TemporaryStore)
- 返回:
生成的空间分区
- 返回类型:
- warp.fem.make_space_restriction(
- space=None,
- space_partition=None,
- domain=None,
- space_topology=None,
- device=None,
- temporary_store=None,
将函数空间分区限制到域,即其元素的子集。
必须提供 space_partition、space_topology 或 space 之一(并将按该顺序考虑)。
- 参数:
space ( FunctionSpace | None ) – (已弃用)如果未提供 space_partition 和 space_topology,则定义要限制的拓扑的空间
space_partition ( SpacePartition | None ) – 要考虑的空间拓扑中的节点子集
domain ( GeometryDomain | None ) – 将空间限制到的域,默认为空间几何体或分区的所有单元格。
space_topology ( SpaceTopology | None ) – 要限制的空间拓扑,如果 space_partition 为
None
。device – 执行和存储计算的设备
temporary_store ( Optional[ warp.fem.cache.TemporaryStore ] ) – 用于分配临时数组的共享池
- 返回类型:
- class warp.fem.ElementBasis(*values)[source]#
为单个元素配备的基函数的选择
- LAGRANGE = 'P'#
单纯形的拉格朗日基函数 \(P_k\),正方形和立方体的张量积 \(Q_k\)
- SERENDIPITY = 'S'#
Serendipity 元素 \(S_k\),对应于去除内部点的拉格朗日节点(对于阶数 <= 3)
- NONCONFORMING_POLYNOMIAL = 'dP'#
嵌入到非一致参考元素(例如正方形或立方体)中的单纯形拉格朗日基函数 \(P_{kd}\)。仅不连续。
- NEDELEC_FIRST_KIND = 'N1'#
Nédélec(第一类)H(curl) 形函数。应与协变函数空间一起使用。
- RAVIART_THOMAS = 'RT'#
Raviart-Thomas H(div) 形函数。应与逆变函数空间一起使用。
- class warp.fem.SymmetricTensorMapper(dtype, mapping=Mapping.VOIGT)[source]#
基类:
DofMapper
从 R^{n (n+1)} 到 nxn 对称张量的正交同构,对向量使用通常的 L2 范数,对张量使用半 Frobenius 范数,(tau : tau)/2。
- 参数:
dtype ( type )
mapping ( Mapping )
- class warp.fem.SkewSymmetricTensorMapper(dtype)[source]#
基类:
DofMapper
从 R^{n (n-1)} 到 nxn 反对称张量的正交同构,向量使用通常的 L2 范数,张量使用一半 Frobenius 范数,(tau : tau)/2。
- 参数:
dtype ( type )
- class warp.fem.PointBasisSpace(quadrature)[source]#
基类:
BasisSpace
一个非结构化的
BasisSpace
,仅在有限的点集处非零。节点位置和节点积分权重由
Quadrature
公式定义。- 参数:
quadrature (Quadrature)
- __init__(quadrature)[source]#
- 参数:
quadrature (Quadrature)
字段#
- warp.fem.make_test(
- space,
- space_restriction=None,
- space_partition=None,
- domain=None,
- device=None,
构造一个函数空间或其限制上的测试场
- 参数:
space (FunctionSpace) – 函数空间
space_restriction (SpaceRestriction | None) – 空间拓扑对域的限制
space_partition (SpacePartition | None) – 如果 space_restriction 为
None
,则要考虑的节点索引的可选子集domain (GeometryDomain | None) – 如果 space_restriction 为
None
,则要考虑的元素的可选子集device – 执行和存储计算的 Warp 设备
- 返回:
测试场
- 返回类型:
- warp.fem.make_trial(
- space,
- space_restriction=None,
- space_partition=None,
- domain=None,
构造函数空间或分区上的试验场
- 参数:
space (FunctionSpace) – 函数空间或函数空间限制
space_restriction (SpaceRestriction | None) – 空间拓扑对域的限制
space_partition (SpacePartition | None) – 如果 space_restriction 为
None
,则要考虑的节点索引的可选子集domain (GeometryDomain | None) – 如果 space_restriction 为
None
,则要考虑的元素的可选子集device – 执行和存储计算的 Warp 设备
- 返回:
试验场
- 返回类型:
- warp.fem.make_discrete_field(space, space_partition=None)[source]#
构造函数空间或分区上的零初始化离散场
另请参阅:
warp.fem.FunctionSpace.make_field()
- 参数:
space (FunctionSpace)
space_partition (SpacePartition | None)
- 返回类型:
- class warp.fem.ImplicitField(
- domain,
- func,
- values=None,
- grad_func=None,
- div_func=None,
- degree=0,
基类:
GeometryField
从域上的任意函数定义的场。尚不支持自动微分,因此如果需要梯度/散度评估,则必须提供相应的函数。
- 参数:
domain (GeometryDomain) – 定义场的域
func (Function) – 计算给定位置处场的 Warp 函数。必须接受至少一个参数,第一个参数是评估位置 (
wp.vec2
或wp.vec3
)。grad_func (Function | None) – 可选的梯度评估函数;必须采用与 func 相同的参数
div_func (Function | None) – 可选的散度评估函数;必须采用与 func 相同的参数
degree – 用于自动确定积分阶数的可选提示,用于积分该场
- __init__(
- domain,
- func,
- values=None,
- grad_func=None,
- div_func=None,
- degree=0,
- 参数:
domain (GeometryDomain)
func (Function)
grad_func (Function | None)
div_func (Function | None)
- class warp.fem.UniformField(domain, value)[source]#
基类:
GeometryField
定义为域上常量值的场。
- 参数:
domain (GeometryDomain) – 定义场的域
value (Any) – 域上的统一值
- __init__(domain, value)[source]#
- 参数:
domain (GeometryDomain)
value (Any)
- class warp.fem.NonconformingField(domain, field, background=0.0)[source]#
基类:
GeometryField
定义为离散场在非协调几何上的映射的场。
- 参数:
domain (GeometryDomain) – 将评估非协调场的新域
field (DiscreteField) – 非协调离散场
background (Any) – 确定 field 定义几何体之外的值的统一值或域协调场
- __init__(domain, field, background=0.0)[source]#
- 参数:
domain (GeometryDomain)
field (DiscreteField)
background (Any)
- warp.fem.make_restriction(
- field,
- space_restriction=None,
- domain=None,
- device=None,
将离散场约束到元素的子集。
- 参数:
field (DiscreteField) – 要约束的离散场
space_restriction (SpaceRestriction | None) – 定义要考虑的元素子集的功能空间约束
domain (GeometryDomain | None) – 如果未提供
space_restriction
,则定义要考虑的元素子集的Domain
device – 执行和存储计算的 Warp 设备
- 返回:
场约束
- 返回类型:
边界条件#
- warp.fem.normalize_dirichlet_projector(projector_matrix, fixed_value=None)[source]#
缩放投影算子,使其变为幂等的,如果提供了 fixed_value,则应用相同的缩放
自适应性#
- warp.fem.adaptive_nanogrid_from_hierarchy(
- grids,
- grading=None,
- temporary_store=None,
从非重叠的网格层次结构构造一个
warp.fem.AdaptiveNanogrid
。警告:如果级别之间存在部分重叠,则结果是未定义的,也就是说,如果级别 l 的单元格仅被级别 l-1 或更低的单元格部分覆盖。
- 参数:
grading (str | None) – 补充分级条件,可以是
None
、“face”或“vertex”;请参阅enforce_nanogrid_grading()
temporary_store (TemporaryStore | None) – 临时分配的存储
- 返回类型:
- warp.fem.adaptive_nanogrid_from_field(
- coarse_grid,
- level_count,
- refinement_field,
- samples_per_voxel=64,
- grading=None,
- temporary_store=None,
从粗糙网格和细化场构造一个
warp.fem.AdaptiveNanogrid
。- 参数:
coarse_grid (Volume) – 从中开始细化的基础网格。不会在基础网格之外添加任何体素。
level_count (int) – 最大细化级别数
refinement_field (GeometryField) – 用作细化 oracle 的标量场。如果返回的值为负,则会雕刻出相应的体素。正值表示所需的细化程度,其中 0.0 对应于最精细的级别,而 1.0 对应于最粗糙的级别。
samples_per_voxel (int) – 用于评估每个体素内的细化场的样本数
grading (str | None) – 补充分级条件,可以是
None
、“face”或“vertex”;请参阅enforce_nanogrid_grading()
temporary_store (TemporaryStore | None) – 临时分配的存储
- 返回类型:
内存管理#
- warp.fem.set_default_temporary_store(temporary_store)[source]#
全局设置默认的
TemporaryStore
实例,用于warp.fem
函数中的临时分配。如果默认临时存储设置为
None
,则除非在每个函数的粒度上提供TemporaryStore
,否则临时分配不会持久保存。- 参数:
temporary_store (TemporaryStore | None)
接口#
接口类不应直接构造,但可以派生自接口类以扩展内置功能。
- class warp.fem.GeometryDomain(geometry)[源代码]#
域的接口类,即 Geometry 中元素的(部分)视图
- 参数:
geometry (Geometry | GeometryPartition)
- __init__(geometry)[源代码]#
- 参数:
geometry (Geometry | GeometryPartition)
- class warp.fem.Quadrature(domain)[源代码]#
积分规则的接口类
- 参数:
domain (GeometryDomain)
- __init__(domain)[源代码]#
- 参数:
domain (GeometryDomain)
- class warp.fem.FunctionSpace(topology)[源代码]#
函数空间的接口类,即几何体 + 插值基
- 函数 f 在位置 x 的值通常计算为
f(x) = L(x)[sum_i f_i N_i(x)]
- 其中
f_i
第 i 个节点的自由度 (dof) 的值N_i(x)
与 x 处节点相关的权重L(x)
从节点空间到世界空间的局部线性变换
- 参数:
topology (SpaceTopology)
- __init__(topology)[源代码]#
- 参数:
topology (SpaceTopology)
- property topology: SpaceTopology[源代码]#
底层几何体
- make_field(space_partition=None)[source]#
在函数空间上创建一个零初始化的离散场,该场保存空间分区中节点的所有自由度的值
- 参数:
space_partition – 如果提供,则为要考虑的节点子集
另请参见:
make_space_partition()
- class warp.fem.SpaceTopology(geometry, max_nodes_per_element)[source]#
用于定义函数空间拓扑的接口类。
拓扑仅考虑每个元素中节点的索引,因此,考虑函数空间的连接模式。 它不指定节点在元素中的实际位置或估值函数。
- class warp.fem.BasisSpace(topology)[source]#
用于定义几何体上的形函数空间的接口类。
基空间使定义共享相同基(因此共享节点)但具有不同估值函数的多个函数空间变得容易; 但是,它不是函数空间所必需的组件。
另请参见:
make_polynomial_basis_space()
,make_collocated_function_space()
- 参数:
topology (SpaceTopology)
- __init__(topology)[source]#
- 参数:
topology (SpaceTopology)
- property topology: SpaceTopology[source]#
基空间的底层拓扑
- class warp.fem.SpacePartition(space_topology, geo_partition)[source]#
- 参数:
space_topology (SpaceTopology)
geo_partition (GeometryPartition)
- __init__(space_topology, geo_partition)[source]#
- 参数:
space_topology (SpaceTopology)
geo_partition (GeometryPartition)
- class warp.fem.SpaceRestriction(
- space_partition,
- domain,
- device=None,
- temporary_store=None,
将空间分区限制到给定的 GeometryDomain
- 参数:
space_partition (SpacePartition)
domain (GeometryDomain)
temporary_store (TemporaryStore)
- __init__(
- space_partition,
- domain,
- device=None,
- temporary_store=None,
- 参数:
space_partition (SpacePartition)
domain (GeometryDomain)
temporary_store (TemporaryStore)
- class warp.fem.DiscreteField(space, space_partition)[source]#
基类:
SpaceField
在离散函数空间的分区上定义的显式值场
- 参数:
space (FunctionSpace)
space_partition (SpacePartition)
- __init__(space, space_partition)[source]#
- 参数:
space (FunctionSpace)
space_partition (SpacePartition)
- class warp.fem.field.FieldRestriction(space_restriction, field)[source]#
离散场到给定GeometryDomain的限制
- 参数:
space_restriction (SpaceRestriction)
field (DiscreteField)
- __init__(space_restriction, field)[source]#
- 参数:
space_restriction (SpaceRestriction)
field (DiscreteField)
- class warp.fem.field.GeometryField[source]#
基类:
FieldLike
在几何体上定义的场的基类
- make_deformed_geometry(relative=True)[source]#
返回底层几何体的变形版本,其位置根据此场的值进行位移。
- 参数:
relative – 如果
True
,则该场被解释为相对于原始几何体的相对位移。 如果False
,则场值被解释为绝对位置。- 返回类型:
- __init__()#
- class warp.fem.field.SpaceField(space, space_partition)[source]#
基类:
GeometryField
在函数空间上定义的场的基类
- 参数:
space (FunctionSpace)
space_partition (SpacePartition)
- __init__(space, space_partition)[source]#
- 参数:
space (FunctionSpace)
space_partition (SpacePartition)
- class warp.fem.field.TestField(space_restriction, space)[source]#
基类:
AdjointField
在空间限制上定义的场,可以用作测试函数。
为了重用计算,可以使用为与测试函数值类型不同的值类型定义的SpaceRestriction来定义测试字段,只要节点拓扑相似即可。
- 参数:
space_restriction (SpaceRestriction)
space (FunctionSpace)
- __init__(space_restriction, space)[source]#
- 参数:
space_restriction (SpaceRestriction)
space (FunctionSpace)
- class warp.fem.field.TrialField(space, space_partition, domain)[source]#
基类:
AdjointField
在可以用作试验函数的域上定义的场
- 参数:
space (FunctionSpace)
space_partition (SpacePartition)
domain (GeometryDomain)
- __init__(space, space_partition, domain)[source]#
- 参数:
space (FunctionSpace)
space_partition (SpacePartition)
domain (GeometryDomain)
- class warp.fem.TemporaryStore[source]#
临时数组的共享池,这些临时数组将在
warp.fem
函数的调用之间持久存在并重用。TemporaryStore
实例可以显式传递给接受此类参数的warp.fem
函数,例如integrate.integrate()
,也可以使用set_default_temporary_store()
全局设置为默认存储。默认情况下,没有默认的临时存储,因此临时分配不会持久存在。
- class warp.fem.cache.Temporary(*args, **kwargs)[source]#
来自
TemporaryStore
的临时数组的句柄。销毁此对象后,该数组将自动返回到临时池以供重用,除非使用
detach()
将临时数组显式地从池中分离。 也可以在使用release()
销毁之前显式地将其返回到池中。