Skip to content

Operations

General Operations¤

Operations are combination of residual operators for given PDEs.

A operation can be designed with inheriting the FieldOperations class which provide all the available operators:

ConvDO.operations.FieldOperations ¤

A class that performs various operations on fields.

Parameters:

Name Type Description Default
order int

The order of the operations.

required
device str

The device to perform the operations on. Defaults to "cpu".

'cpu'
dtype dtype

The data type of the operations. Defaults to torch.float32.

float32

Attributes:

Name Type Description
nabla ConvNabla

The gradient operator.

nabla2 ConvLaplacian

The Laplacian operator.

grad_x ConvGrad

The gradient operator in the x direction.

grad_y ConvGrad

The gradient operator in the y direction.

Source code in ConvDO/operations.py
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
class FieldOperations():
    r"""
    A class that performs various operations on fields.

    Args:
        order (int): The order of the operations.
        device (str, optional): The device to perform the operations on. Defaults to "cpu".
        dtype (torch.dtype, optional): The data type of the operations. Defaults to torch.float32.

    Attributes:
        nabla (ConvNabla): The gradient operator.
        nabla2 (ConvLaplacian): The Laplacian operator.
        grad_x (ConvGrad): The gradient operator in the x direction.
        grad_y (ConvGrad): The gradient operator in the y direction.
    """

    def __init__(self, order:int, device="cpu", dtype=torch.float32) -> None:
        self.nabla = ConvNabla(order, device=device, dtype=dtype)
        self.nabla2 = ConvLaplacian(order, device=device, dtype=dtype)
        self.grad_x = ConvGrad(order, direction='x', device=device, dtype=dtype)
        self.grad_y = ConvGrad(order, direction='y', device=device, dtype=dtype)
        self.grad2_x = ConvGrad2(order, direction='x', device=device, dtype=dtype)
        self.grad2_y = ConvGrad2(order, direction='y', device=device, dtype=dtype)
nabla instance-attribute ¤
nabla = ConvNabla(order, device=device, dtype=dtype)
nabla2 instance-attribute ¤
nabla2 = ConvLaplacian(order, device=device, dtype=dtype)
grad_x instance-attribute ¤
grad_x = ConvGrad(
    order, direction="x", device=device, dtype=dtype
)
grad_y instance-attribute ¤
grad_y = ConvGrad(
    order, direction="y", device=device, dtype=dtype
)
grad2_x instance-attribute ¤
grad2_x = ConvGrad2(
    order, direction="x", device=device, dtype=dtype
)
grad2_y instance-attribute ¤
grad2_y = ConvGrad2(
    order, direction="y", device=device, dtype=dtype
)
__init__ ¤
__init__(
    order: int, device="cpu", dtype=torch.float32
) -> None
Source code in ConvDO/operations.py
22
23
24
25
26
27
28
def __init__(self, order:int, device="cpu", dtype=torch.float32) -> None:
    self.nabla = ConvNabla(order, device=device, dtype=dtype)
    self.nabla2 = ConvLaplacian(order, device=device, dtype=dtype)
    self.grad_x = ConvGrad(order, direction='x', device=device, dtype=dtype)
    self.grad_y = ConvGrad(order, direction='y', device=device, dtype=dtype)
    self.grad2_x = ConvGrad2(order, direction='x', device=device, dtype=dtype)
    self.grad2_y = ConvGrad2(order, direction='y', device=device, dtype=dtype)

Pre-defined Operations¤

We also give some pre-defined operations which mainly focus on the flow dynamics:

ConvDO.operations.TransientNS ¤

Bases: FieldOperations

Class representing the transient Navier-Stokes equations. Returns a three channel tensor containing the residual of momentum equation in x and y directions and the divergence of the velocity field. Defining \(\mathbf{u}_{int}=(\mathbf{u}_0+\mathbf{u}_1)/2=((u_{x,0}+u_{x,1})/2,(u_{y,0}+u_{y,1})/2)\), the residual is given by: \(\begin{matrix}\frac{\mathbf{u}_1-\mathbf{u}_0}{\Delta t}+\mathbf{u}_{int}\cdot\nabla\mathbf{u}_{int}+\nabla\left(\frac{p_0+p_1}{2}\right)-\nu\nabla^2\mathbf{u}_{int}\\(\nabla\cdot\mathbf{u}_{int}+\nabla\cdot\mathbf{u}_1)/2\end{matrix}\)

Parameters:

Name Type Description Default
domain_u Domain

The domain for the x-velocity component.

required
domain_v Domain

The domain for the y-velocity component.

required
domain_p Domain

The domain for the pressure component.

required
viscosity float

The viscosity coefficient.

required
dt float

The time step size.

required
order int

The order of accuracy for the finite difference scheme.

required
device str

The device to use for computation (default: "cpu").

'cpu'
dtype dtype

The data type to use for computation (default: torch.float32).

float32
Source code in ConvDO/operations.py
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
class TransientNS(FieldOperations):
    r"""
    Class representing the transient Navier-Stokes equations.
        Returns a three channel tensor containing the residual of momentum equation in x and y directions and the divergence of the velocity field.
        Defining $\mathbf{u}_{int}=(\mathbf{u}_0+\mathbf{u}_1)/2=((u_{x,0}+u_{x,1})/2,(u_{y,0}+u_{y,1})/2)$, the residual is given by:
        $\begin{matrix}\frac{\mathbf{u}_1-\mathbf{u}_0}{\Delta t}+\mathbf{u}_{int}\cdot\nabla\mathbf{u}_{int}+\nabla\left(\frac{p_0+p_1}{2}\right)-\nu\nabla^2\mathbf{u}_{int}\\(\nabla\cdot\mathbf{u}_{int}+\nabla\cdot\mathbf{u}_1)/2\end{matrix}$

    Args:
        domain_u (Domain): The domain for the x-velocity component.
        domain_v (Domain): The domain for the y-velocity component.
        domain_p (Domain): The domain for the pressure component.
        viscosity (float): The viscosity coefficient.
        dt (float): The time step size.
        order (int): The order of accuracy for the finite difference scheme.
        device (str, optional): The device to use for computation (default: "cpu").
        dtype (torch.dtype, optional): The data type to use for computation (default: torch.float32).
    """

    def __init__(self, 
                 domain_u: Domain,
                 domain_v: Domain, 
                 domain_p: Domain, 
                 viscosity: float, 
                 dt: float,
                 order:int,
                 device="cpu", 
                 dtype=torch.float32,
                 ) -> None:
        super().__init__(order, device=device, dtype=dtype)
        self.p_0 = ScalarField(domain=domain_p)
        self.p_1 = ScalarField(domain=domain_p)
        self.velocity_0 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
        self.velocity_1 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
        self.viscosity = viscosity
        self.dt = dt


    def __call__(self, u_0, v_0, p_0, u_1, v_1, p_1):
        """
        Compute the solution of the transient Navier-Stokes equations with external force.

        Args:
            u_0 (torch.Tensor): The x-velocity component at time step t.
            v_0 (torch.Tensor): The y-velocity component at time step t.
            p_0 (torch.Tensor): The pressure component at time step t.
            u_1 (torch.Tensor): The x-velocity component at time step t+1.
            v_1 (torch.Tensor): The y-velocity component at time step t+1.
            p_1 (torch.Tensor): The pressure component at time step t+1.

        Returns:
            residual (torch.Tensor): The concatenated tensor of the x-velocity, y-velocity, and divergence components.
        """
        self.velocity_0.ux.register_value(u_0)
        self.velocity_0.uy.register_value(v_0)
        self.p_0.register_value(p_0)
        self.velocity_1.ux.register_value(u_1)
        self.velocity_1.uy.register_value(v_1)
        self.p_1.register_value(p_1)
        u_inter = (self.velocity_0 + self.velocity_1) * 0.5
        transient = (self.velocity_1 - self.velocity_0) / self.dt
        advection = u_inter @ (self.nabla * u_inter)
        pressure = self.nabla * ((self.p_0 + self.p_1) * 0.5)
        vis = -1 * self.viscosity * (self.nabla2 * u_inter)
        ns_res = transient + advection + pressure + vis
        divergence = ((self.nabla @ u_inter) + (self.nabla @ self.velocity_1)) * 0.5
        return torch.cat([ns_res.ux.value, ns_res.uy.value, divergence.value], dim=-3)
nabla instance-attribute ¤
nabla = ConvNabla(order, device=device, dtype=dtype)
nabla2 instance-attribute ¤
nabla2 = ConvLaplacian(order, device=device, dtype=dtype)
grad_x instance-attribute ¤
grad_x = ConvGrad(
    order, direction="x", device=device, dtype=dtype
)
grad_y instance-attribute ¤
grad_y = ConvGrad(
    order, direction="y", device=device, dtype=dtype
)
grad2_x instance-attribute ¤
grad2_x = ConvGrad2(
    order, direction="x", device=device, dtype=dtype
)
grad2_y instance-attribute ¤
grad2_y = ConvGrad2(
    order, direction="y", device=device, dtype=dtype
)
p_0 instance-attribute ¤
p_0 = ScalarField(domain=domain_p)
p_1 instance-attribute ¤
p_1 = ScalarField(domain=domain_p)
velocity_0 instance-attribute ¤
velocity_0 = VectorValue(
    ScalarField(domain=domain_u),
    ScalarField(domain=domain_v),
)
velocity_1 instance-attribute ¤
velocity_1 = VectorValue(
    ScalarField(domain=domain_u),
    ScalarField(domain=domain_v),
)
viscosity instance-attribute ¤
viscosity = viscosity
dt instance-attribute ¤
dt = dt
__init__ ¤
__init__(
    domain_u: Domain,
    domain_v: Domain,
    domain_p: Domain,
    viscosity: float,
    dt: float,
    order: int,
    device="cpu",
    dtype=torch.float32,
) -> None
Source code in ConvDO/operations.py
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
def __init__(self, 
             domain_u: Domain,
             domain_v: Domain, 
             domain_p: Domain, 
             viscosity: float, 
             dt: float,
             order:int,
             device="cpu", 
             dtype=torch.float32,
             ) -> None:
    super().__init__(order, device=device, dtype=dtype)
    self.p_0 = ScalarField(domain=domain_p)
    self.p_1 = ScalarField(domain=domain_p)
    self.velocity_0 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
    self.velocity_1 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
    self.viscosity = viscosity
    self.dt = dt
__call__ ¤
__call__(u_0, v_0, p_0, u_1, v_1, p_1)

Compute the solution of the transient Navier-Stokes equations with external force.

Parameters:

Name Type Description Default
u_0 Tensor

The x-velocity component at time step t.

required
v_0 Tensor

The y-velocity component at time step t.

required
p_0 Tensor

The pressure component at time step t.

required
u_1 Tensor

The x-velocity component at time step t+1.

required
v_1 Tensor

The y-velocity component at time step t+1.

required
p_1 Tensor

The pressure component at time step t+1.

required

Returns:

Name Type Description
residual Tensor

The concatenated tensor of the x-velocity, y-velocity, and divergence components.

Source code in ConvDO/operations.py
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
def __call__(self, u_0, v_0, p_0, u_1, v_1, p_1):
    """
    Compute the solution of the transient Navier-Stokes equations with external force.

    Args:
        u_0 (torch.Tensor): The x-velocity component at time step t.
        v_0 (torch.Tensor): The y-velocity component at time step t.
        p_0 (torch.Tensor): The pressure component at time step t.
        u_1 (torch.Tensor): The x-velocity component at time step t+1.
        v_1 (torch.Tensor): The y-velocity component at time step t+1.
        p_1 (torch.Tensor): The pressure component at time step t+1.

    Returns:
        residual (torch.Tensor): The concatenated tensor of the x-velocity, y-velocity, and divergence components.
    """
    self.velocity_0.ux.register_value(u_0)
    self.velocity_0.uy.register_value(v_0)
    self.p_0.register_value(p_0)
    self.velocity_1.ux.register_value(u_1)
    self.velocity_1.uy.register_value(v_1)
    self.p_1.register_value(p_1)
    u_inter = (self.velocity_0 + self.velocity_1) * 0.5
    transient = (self.velocity_1 - self.velocity_0) / self.dt
    advection = u_inter @ (self.nabla * u_inter)
    pressure = self.nabla * ((self.p_0 + self.p_1) * 0.5)
    vis = -1 * self.viscosity * (self.nabla2 * u_inter)
    ns_res = transient + advection + pressure + vis
    divergence = ((self.nabla @ u_inter) + (self.nabla @ self.velocity_1)) * 0.5
    return torch.cat([ns_res.ux.value, ns_res.uy.value, divergence.value], dim=-3)

ConvDO.operations.TransientNSWithForce ¤

Bases: FieldOperations

Class representing the transient Navier-Stokes equations with external force. Returns a three channel tensor containing the residual of momentum equation in x and y directions and the divergence of the velocity field.
Defining \(\mathbf{u}_{int}=(\mathbf{u}_0+\mathbf{u}_1)/2=((u_{x,0}+u_{x,1})/2,(u_{y,0}+u_{y,1})/2)\), the residual is given by: \(\begin{matrix}\frac{\mathbf{u}_1-\mathbf{u}_0}{\Delta t}+\mathbf{u}_{int}\cdot\nabla\mathbf{u}_{int}+\nabla\left(\frac{p_0+p_1}{2}\right)-\nu\nabla^2\mathbf{u}_{int}-\mathbf{f}\\(\nabla\cdot\mathbf{u}_{int}+\nabla\cdot\mathbf{u}_1)/2\end{matrix}\)

Parameters:

Name Type Description Default
domain_u Domain

The domain for the x-velocity component.

required
domain_v Domain

The domain for the y-velocity component.

required
domain_p Domain

The domain for the pressure component.

required
force_x Tensor

The x-component of the external force.

required
force_y Tensor

The y-component of the external force.

required
domain_force_x Domain

The domain for the x-component of the external force.

required
domain_force_y Domain

The domain for the y-component of the external force.

required
viscosity float

The viscosity coefficient.

required
dt float

The time step size.

required
order int

The order of accuracy for the finite difference scheme.

required
device str

The device to use for computation (default: "cpu").

'cpu'
dtype dtype

The data type to use for computation (default: torch.float32).

float32
Source code in ConvDO/operations.py
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
class TransientNSWithForce(FieldOperations):
    r"""
    Class representing the transient Navier-Stokes equations with external force.
        Returns a three channel tensor containing the residual of momentum equation in x and y directions and the divergence of the velocity field.   
        Defining $\mathbf{u}_{int}=(\mathbf{u}_0+\mathbf{u}_1)/2=((u_{x,0}+u_{x,1})/2,(u_{y,0}+u_{y,1})/2)$, the residual is given by:
        $\begin{matrix}\frac{\mathbf{u}_1-\mathbf{u}_0}{\Delta t}+\mathbf{u}_{int}\cdot\nabla\mathbf{u}_{int}+\nabla\left(\frac{p_0+p_1}{2}\right)-\nu\nabla^2\mathbf{u}_{int}-\mathbf{f}\\(\nabla\cdot\mathbf{u}_{int}+\nabla\cdot\mathbf{u}_1)/2\end{matrix}$


    Args:
        domain_u (Domain): The domain for the x-velocity component.
        domain_v (Domain): The domain for the y-velocity component.
        domain_p (Domain): The domain for the pressure component.
        force_x (torch.Tensor): The x-component of the external force.
        force_y (torch.Tensor): The y-component of the external force.
        domain_force_x (Domain): The domain for the x-component of the external force.
        domain_force_y (Domain): The domain for the y-component of the external force.
        viscosity (float): The viscosity coefficient.
        dt (float): The time step size.
        order (int): The order of accuracy for the finite difference scheme.
        device (str, optional): The device to use for computation (default: "cpu").
        dtype (torch.dtype, optional): The data type to use for computation (default: torch.float32).
    """

    def __init__(self, 
                 domain_u: Domain,
                 domain_v: Domain, 
                 domain_p: Domain, 
                 force_x: torch.Tensor,
                 force_y: torch.Tensor,
                 domain_force_x: Domain, 
                 domain_force_y: Domain, 
                 viscosity: float, 
                 dt: float,
                 order:int,
                 device="cpu", 
                 dtype=torch.float32,
                 ) -> None:
        super().__init__(order, device=device, dtype=dtype)
        self.p_0 = ScalarField(domain=domain_p)
        self.p_1 = ScalarField(domain=domain_p)
        self.velocity_0 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
        self.velocity_1 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
        self.force = VectorValue(ScalarField(force_x, domain=domain_force_x),
                                 ScalarField(force_y, domain=domain_force_y))
        self.viscosity = viscosity
        self.dt = dt


    def __call__(self, u_0, v_0, p_0, u_1, v_1, p_1):
        """
        Compute the solution of the transient Navier-Stokes equations with external force.

        Args:
            u_0 (torch.Tensor): The x-velocity component at time step t.
            v_0 (torch.Tensor): The y-velocity component at time step t.
            p_0 (torch.Tensor): The pressure component at time step t.
            u_1 (torch.Tensor): The x-velocity component at time step t+1.
            v_1 (torch.Tensor): The y-velocity component at time step t+1.
            p_1 (torch.Tensor): The pressure component at time step t+1.

        Returns:
            residual (torch.Tensor): The concatenated tensor of the x-velocity, y-velocity, and divergence components.
        """
        self.velocity_0.ux.register_value(u_0)
        self.velocity_0.uy.register_value(v_0)
        self.p_0.register_value(p_0)
        self.velocity_1.ux.register_value(u_1)
        self.velocity_1.uy.register_value(v_1)
        self.p_1.register_value(p_1)
        u_inter = (self.velocity_0 + self.velocity_1) * 0.5
        transient = (self.velocity_1 - self.velocity_0) / self.dt
        advection = u_inter @ (self.nabla * u_inter)
        pressure = self.nabla * ((self.p_0 + self.p_1) * 0.5)
        vis = -1 * self.viscosity * (self.nabla2 * u_inter)
        ns_res = transient + advection + pressure + vis - self.force
        divergence = ((self.nabla @ u_inter) + (self.nabla @ self.velocity_1)) * 0.5
        return torch.cat([ns_res.ux.value, ns_res.uy.value, divergence.value], dim=-3)
nabla instance-attribute ¤
nabla = ConvNabla(order, device=device, dtype=dtype)
nabla2 instance-attribute ¤
nabla2 = ConvLaplacian(order, device=device, dtype=dtype)
grad_x instance-attribute ¤
grad_x = ConvGrad(
    order, direction="x", device=device, dtype=dtype
)
grad_y instance-attribute ¤
grad_y = ConvGrad(
    order, direction="y", device=device, dtype=dtype
)
grad2_x instance-attribute ¤
grad2_x = ConvGrad2(
    order, direction="x", device=device, dtype=dtype
)
grad2_y instance-attribute ¤
grad2_y = ConvGrad2(
    order, direction="y", device=device, dtype=dtype
)
p_0 instance-attribute ¤
p_0 = ScalarField(domain=domain_p)
p_1 instance-attribute ¤
p_1 = ScalarField(domain=domain_p)
velocity_0 instance-attribute ¤
velocity_0 = VectorValue(
    ScalarField(domain=domain_u),
    ScalarField(domain=domain_v),
)
velocity_1 instance-attribute ¤
velocity_1 = VectorValue(
    ScalarField(domain=domain_u),
    ScalarField(domain=domain_v),
)
force instance-attribute ¤
force = VectorValue(
    ScalarField(force_x, domain=domain_force_x),
    ScalarField(force_y, domain=domain_force_y),
)
viscosity instance-attribute ¤
viscosity = viscosity
dt instance-attribute ¤
dt = dt
__init__ ¤
__init__(
    domain_u: Domain,
    domain_v: Domain,
    domain_p: Domain,
    force_x: torch.Tensor,
    force_y: torch.Tensor,
    domain_force_x: Domain,
    domain_force_y: Domain,
    viscosity: float,
    dt: float,
    order: int,
    device="cpu",
    dtype=torch.float32,
) -> None
Source code in ConvDO/operations.py
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
def __init__(self, 
             domain_u: Domain,
             domain_v: Domain, 
             domain_p: Domain, 
             force_x: torch.Tensor,
             force_y: torch.Tensor,
             domain_force_x: Domain, 
             domain_force_y: Domain, 
             viscosity: float, 
             dt: float,
             order:int,
             device="cpu", 
             dtype=torch.float32,
             ) -> None:
    super().__init__(order, device=device, dtype=dtype)
    self.p_0 = ScalarField(domain=domain_p)
    self.p_1 = ScalarField(domain=domain_p)
    self.velocity_0 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
    self.velocity_1 = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
    self.force = VectorValue(ScalarField(force_x, domain=domain_force_x),
                             ScalarField(force_y, domain=domain_force_y))
    self.viscosity = viscosity
    self.dt = dt
__call__ ¤
__call__(u_0, v_0, p_0, u_1, v_1, p_1)

Compute the solution of the transient Navier-Stokes equations with external force.

Parameters:

Name Type Description Default
u_0 Tensor

The x-velocity component at time step t.

required
v_0 Tensor

The y-velocity component at time step t.

required
p_0 Tensor

The pressure component at time step t.

required
u_1 Tensor

The x-velocity component at time step t+1.

required
v_1 Tensor

The y-velocity component at time step t+1.

required
p_1 Tensor

The pressure component at time step t+1.

required

Returns:

Name Type Description
residual Tensor

The concatenated tensor of the x-velocity, y-velocity, and divergence components.

Source code in ConvDO/operations.py
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
def __call__(self, u_0, v_0, p_0, u_1, v_1, p_1):
    """
    Compute the solution of the transient Navier-Stokes equations with external force.

    Args:
        u_0 (torch.Tensor): The x-velocity component at time step t.
        v_0 (torch.Tensor): The y-velocity component at time step t.
        p_0 (torch.Tensor): The pressure component at time step t.
        u_1 (torch.Tensor): The x-velocity component at time step t+1.
        v_1 (torch.Tensor): The y-velocity component at time step t+1.
        p_1 (torch.Tensor): The pressure component at time step t+1.

    Returns:
        residual (torch.Tensor): The concatenated tensor of the x-velocity, y-velocity, and divergence components.
    """
    self.velocity_0.ux.register_value(u_0)
    self.velocity_0.uy.register_value(v_0)
    self.p_0.register_value(p_0)
    self.velocity_1.ux.register_value(u_1)
    self.velocity_1.uy.register_value(v_1)
    self.p_1.register_value(p_1)
    u_inter = (self.velocity_0 + self.velocity_1) * 0.5
    transient = (self.velocity_1 - self.velocity_0) / self.dt
    advection = u_inter @ (self.nabla * u_inter)
    pressure = self.nabla * ((self.p_0 + self.p_1) * 0.5)
    vis = -1 * self.viscosity * (self.nabla2 * u_inter)
    ns_res = transient + advection + pressure + vis - self.force
    divergence = ((self.nabla @ u_inter) + (self.nabla @ self.velocity_1)) * 0.5
    return torch.cat([ns_res.ux.value, ns_res.uy.value, divergence.value], dim=-3)

ConvDO.operations.PoissonDivergence ¤

Bases: FieldOperations

Class representing the Poisson equation for pressure and the divergence of the velocity field. Returns a two channel tensor containing the residual of the Poisson equation and the divergence of the velocity field. Residual of Poisson equation:\(\left( {{\partial^2 p} \over {\partial x^2}} + {{\partial^2 p} \over {\partial y^2}} \right) = \left( {{\partial u} \over {\partial x}} \right)^2 + 2 {{\partial u} \over {\partial y}} {{\partial v} \over {\partial x}}+ \left( {{\partial v} \over {\partial y}} \right)^2\) Residual of continuity:\(\nabla\cdot\mathbf{u}\)

Parameters:

Name Type Description Default
domain_u Domain

The domain for the x-velocity component.

required
domain_v Domain

The domain for the y-velocity component.

required
domain_p Domain

The domain for the pressure component.

required
order int

The order of accuracy for the finite difference scheme.

required
device str

The device to use for computation (default: "cpu").

'cpu'
dtype dtype

The data type to use for computation (default: torch.float32).

float32
Source code in ConvDO/operations.py
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
class PoissonDivergence(FieldOperations):
    r"""
    Class representing the Poisson equation for pressure and the divergence of the velocity field.
        Returns a two channel tensor containing the residual of the Poisson equation and the divergence of the velocity field.
        Residual of Poisson equation:$\left( {{\partial^2 p} \over {\partial x^2}} + {{\partial^2 p} \over {\partial y^2}} \right) = \left( {{\partial u} \over {\partial x}} \right)^2 + 2 {{\partial u} \over {\partial y}} {{\partial v} \over {\partial x}}+ \left( {{\partial v} \over {\partial y}} \right)^2$
        Residual of continuity:$\nabla\cdot\mathbf{u}$

    Args:
        domain_u (Domain): The domain for the x-velocity component.
        domain_v (Domain): The domain for the y-velocity component.
        domain_p (Domain): The domain for the pressure component.
        order (int): The order of accuracy for the finite difference scheme.
        device (str, optional): The device to use for computation (default: "cpu").
        dtype (torch.dtype, optional): The data type to use for computation (default: torch.float32).
    """

    def __init__(self, 
                 domain_u: Domain,
                 domain_v: Domain, 
                 domain_p: Domain, 
                 order: int,
                 device="cpu", 
                 dtype=torch.float32) -> None:
        super().__init__(order, device=device, dtype=dtype)
        self.pressure = ScalarField(domain=domain_p)
        self.velocity = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))

    def __call__(self, u, v, p):
        """
        Compute the solution of the Poisson equation for pressure and the divergence of the velocity field.

        Args:
            u (torch.Tensor): The x-velocity component.
            v (torch.Tensor): The y-velocity component.
            p (torch.Tensor): The pressure component.

        Returns:
            residual (torch.Tensor): The concatenated tensor of the Poisson equation and the divergence components.
        """
        self.velocity.ux.register_value(u)
        self.velocity.uy.register_value(v)
        self.pressure.register_value(p)
        poisson = (self.nabla2*self.pressure)+(self.grad_x*self.velocity.ux)**2+2*(self.grad_y*self.velocity.ux)*(self.grad_x*self.velocity.uy)+(self.grad_y*self.velocity.uy)**2
        divergence = self.nabla@self.velocity
        return torch.cat([poisson.value, divergence.value], dim=-3)
nabla instance-attribute ¤
nabla = ConvNabla(order, device=device, dtype=dtype)
nabla2 instance-attribute ¤
nabla2 = ConvLaplacian(order, device=device, dtype=dtype)
grad_x instance-attribute ¤
grad_x = ConvGrad(
    order, direction="x", device=device, dtype=dtype
)
grad_y instance-attribute ¤
grad_y = ConvGrad(
    order, direction="y", device=device, dtype=dtype
)
grad2_x instance-attribute ¤
grad2_x = ConvGrad2(
    order, direction="x", device=device, dtype=dtype
)
grad2_y instance-attribute ¤
grad2_y = ConvGrad2(
    order, direction="y", device=device, dtype=dtype
)
pressure instance-attribute ¤
pressure = ScalarField(domain=domain_p)
velocity instance-attribute ¤
velocity = VectorValue(
    ScalarField(domain=domain_u),
    ScalarField(domain=domain_v),
)
__init__ ¤
__init__(
    domain_u: Domain,
    domain_v: Domain,
    domain_p: Domain,
    order: int,
    device="cpu",
    dtype=torch.float32,
) -> None
Source code in ConvDO/operations.py
251
252
253
254
255
256
257
258
259
260
def __init__(self, 
             domain_u: Domain,
             domain_v: Domain, 
             domain_p: Domain, 
             order: int,
             device="cpu", 
             dtype=torch.float32) -> None:
    super().__init__(order, device=device, dtype=dtype)
    self.pressure = ScalarField(domain=domain_p)
    self.velocity = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
__call__ ¤
__call__(u, v, p)

Compute the solution of the Poisson equation for pressure and the divergence of the velocity field.

Parameters:

Name Type Description Default
u Tensor

The x-velocity component.

required
v Tensor

The y-velocity component.

required
p Tensor

The pressure component.

required

Returns:

Name Type Description
residual Tensor

The concatenated tensor of the Poisson equation and the divergence components.

Source code in ConvDO/operations.py
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
def __call__(self, u, v, p):
    """
    Compute the solution of the Poisson equation for pressure and the divergence of the velocity field.

    Args:
        u (torch.Tensor): The x-velocity component.
        v (torch.Tensor): The y-velocity component.
        p (torch.Tensor): The pressure component.

    Returns:
        residual (torch.Tensor): The concatenated tensor of the Poisson equation and the divergence components.
    """
    self.velocity.ux.register_value(u)
    self.velocity.uy.register_value(v)
    self.pressure.register_value(p)
    poisson = (self.nabla2*self.pressure)+(self.grad_x*self.velocity.ux)**2+2*(self.grad_y*self.velocity.ux)*(self.grad_x*self.velocity.uy)+(self.grad_y*self.velocity.uy)**2
    divergence = self.nabla@self.velocity
    return torch.cat([poisson.value, divergence.value], dim=-3)

ConvDO.operations.PoissonDivergenceWithForce ¤

Bases: FieldOperations

Class representing the Poisson equation for pressure and the divergence of the velocity field. Returns a two channel tensor containing the residual of the Poisson equation and the divergence of the velocity field. Residual of Poisson equation: \(\left( {{\partial^2 p} \over {\partial x^2}} + {{\partial^2 p} \over {\partial y^2}} \right) = \left( {{\partial u} \over {\partial x}} \right)^2 + 2 {{\partial u} \over {\partial y}} {{\partial v} \over {\partial x}}+ \left( {{\partial v} \over {\partial y}} \right)^2\) Residual of continuity:\(\nabla\cdot\mathbf{u}\)

Parameters:

Name Type Description Default
domain_u Domain

The domain for the x-velocity component.

required
domain_v Domain

The domain for the y-velocity component.

required
domain_p Domain

The domain for the pressure component.

required
force_x Tensor

The x-component of the external force.

required
force_y Tensor

The y-component of the external force.

required
domain_force_x Domain

The domain for the x-component of the external force.

required
domain_force_y Domain

The domain for the y-component of the external force.

required
order int

The order of accuracy for the finite difference scheme.

required
device str

The device to use for computation (default: "cpu").

'cpu'
dtype dtype

The data type to use for computation (default: torch.float32).

float32
Source code in ConvDO/operations.py
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
class PoissonDivergenceWithForce(FieldOperations):
    r"""
    Class representing the Poisson equation for pressure and the divergence of the velocity field.
        Returns a two channel tensor containing the residual of the Poisson equation and the divergence of the velocity field.
        Residual of Poisson equation:
        $\left( {{\partial^2 p} \over {\partial x^2}} + {{\partial^2 p} \over {\partial y^2}} \right) = \left( {{\partial u} \over {\partial x}} \right)^2 + 2 {{\partial u} \over {\partial y}} {{\partial v} \over {\partial x}}+ \left( {{\partial v} \over {\partial y}} \right)^2$
        Residual of continuity:$\nabla\cdot\mathbf{u}$

    Args:
        domain_u (Domain): The domain for the x-velocity component.
        domain_v (Domain): The domain for the y-velocity component.
        domain_p (Domain): The domain for the pressure component.
        force_x (torch.Tensor): The x-component of the external force.
        force_y (torch.Tensor): The y-component of the external force.
        domain_force_x (Domain): The domain for the x-component of the external force.
        domain_force_y (Domain): The domain for the y-component of the external force.
        order (int): The order of accuracy for the finite difference scheme.
        device (str, optional): The device to use for computation (default: "cpu").
        dtype (torch.dtype, optional): The data type to use for computation (default: torch.float32).
    """

    def __init__(self, 
                 domain_u: Domain,
                 domain_v: Domain, 
                 domain_p: Domain, 
                 force_x: torch.Tensor,
                 force_y: torch.Tensor,
                 domain_force_x: Domain, 
                 domain_force_y: Domain,
                 order: int,
                 device="cpu", 
                 dtype=torch.float32) -> None:
        super().__init__(order, device=device, dtype=dtype)
        self.pressure = ScalarField(domain=domain_p)
        self.velocity = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
        self.force = VectorValue(ScalarField(force_x, domain=domain_force_x),ScalarField(force_y, domain=domain_force_y))

    def __call__(self, u, v, p):
        """
        Compute the solution of the Poisson equation for pressure and the divergence of the velocity field.

        Args:
            u (torch.Tensor): The x-velocity component.
            v (torch.Tensor): The y-velocity component.
            p (torch.Tensor): The pressure component.

        Returns:
            residual (torch.Tensor): The concatenated tensor of the Poisson equation and the divergence components.
        """
        self.velocity.ux.register_value(u)
        self.velocity.uy.register_value(v)
        self.pressure.register_value(p)
        poisson = (self.nabla2*self.pressure)+(self.grad_x*self.velocity.ux)**2+2*(self.grad_y*self.velocity.ux)*(self.grad_x*self.velocity.uy)+(self.grad_y*self.velocity.uy)**2-self.nabla@self.force
        divergence = self.nabla@self.velocity
        return torch.cat([poisson.value, divergence.value], dim=-3)
nabla instance-attribute ¤
nabla = ConvNabla(order, device=device, dtype=dtype)
nabla2 instance-attribute ¤
nabla2 = ConvLaplacian(order, device=device, dtype=dtype)
grad_x instance-attribute ¤
grad_x = ConvGrad(
    order, direction="x", device=device, dtype=dtype
)
grad_y instance-attribute ¤
grad_y = ConvGrad(
    order, direction="y", device=device, dtype=dtype
)
grad2_x instance-attribute ¤
grad2_x = ConvGrad2(
    order, direction="x", device=device, dtype=dtype
)
grad2_y instance-attribute ¤
grad2_y = ConvGrad2(
    order, direction="y", device=device, dtype=dtype
)
pressure instance-attribute ¤
pressure = ScalarField(domain=domain_p)
velocity instance-attribute ¤
velocity = VectorValue(
    ScalarField(domain=domain_u),
    ScalarField(domain=domain_v),
)
force instance-attribute ¤
force = VectorValue(
    ScalarField(force_x, domain=domain_force_x),
    ScalarField(force_y, domain=domain_force_y),
)
__init__ ¤
__init__(
    domain_u: Domain,
    domain_v: Domain,
    domain_p: Domain,
    force_x: torch.Tensor,
    force_y: torch.Tensor,
    domain_force_x: Domain,
    domain_force_y: Domain,
    order: int,
    device="cpu",
    dtype=torch.float32,
) -> None
Source code in ConvDO/operations.py
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
def __init__(self, 
             domain_u: Domain,
             domain_v: Domain, 
             domain_p: Domain, 
             force_x: torch.Tensor,
             force_y: torch.Tensor,
             domain_force_x: Domain, 
             domain_force_y: Domain,
             order: int,
             device="cpu", 
             dtype=torch.float32) -> None:
    super().__init__(order, device=device, dtype=dtype)
    self.pressure = ScalarField(domain=domain_p)
    self.velocity = VectorValue(ScalarField(domain=domain_u), ScalarField(domain=domain_v))
    self.force = VectorValue(ScalarField(force_x, domain=domain_force_x),ScalarField(force_y, domain=domain_force_y))
__call__ ¤
__call__(u, v, p)

Compute the solution of the Poisson equation for pressure and the divergence of the velocity field.

Parameters:

Name Type Description Default
u Tensor

The x-velocity component.

required
v Tensor

The y-velocity component.

required
p Tensor

The pressure component.

required

Returns:

Name Type Description
residual Tensor

The concatenated tensor of the Poisson equation and the divergence components.

Source code in ConvDO/operations.py
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
def __call__(self, u, v, p):
    """
    Compute the solution of the Poisson equation for pressure and the divergence of the velocity field.

    Args:
        u (torch.Tensor): The x-velocity component.
        v (torch.Tensor): The y-velocity component.
        p (torch.Tensor): The pressure component.

    Returns:
        residual (torch.Tensor): The concatenated tensor of the Poisson equation and the divergence components.
    """
    self.velocity.ux.register_value(u)
    self.velocity.uy.register_value(v)
    self.pressure.register_value(p)
    poisson = (self.nabla2*self.pressure)+(self.grad_x*self.velocity.ux)**2+2*(self.grad_y*self.velocity.ux)*(self.grad_x*self.velocity.uy)+(self.grad_y*self.velocity.uy)**2-self.nabla@self.force
    divergence = self.nabla@self.velocity
    return torch.cat([poisson.value, divergence.value], dim=-3)