We discretize the problem with conforming Lagrangian FE spaces. For this formulation, the nonlinear weak form reads: find such that for all . As in previous tutorials, the space is the set of functions in that fulfill the Dirichlet boundary conditions, whereas is composed by functions in that vanish at the Dirichlet boundary. The weak residual evaluated at a function is the linear form defined as
In order to solve this nonlinear weak equation, we consider a Newton-Raphson method, which is associated with a linearization of the problem in an arbitrary direction , namely . In previous formula, is the Jacobian evaluated at , which is the bilinear form
Note that the solution of this nonlinear PDE with a Newton-Raphson method, will require to discretize both the residual and the Jacobian . In Gridap, this is done by following an approach similar to the one already shown in previous tutorials for discretizing the bilinear and linear forms associated with a linear FE problem. The specific details are discussed now.
As in previous tutorials, the first step to solve the PDE is to load a discretization of the computational domain. In this case, we load the model from the same file as in the first tutorial
using Gridap
using DrWatson
model = DiscreteModelFromFile("meshes/poisson.json")
As stated before, we want to impose Dirichlet boundary conditions on and , but none of these boundaries is identified in the model. E.g., you can easily see by writing the model in vtk format
writevtk(model,"model")
and by opening the file "model_0"
in ParaView that the boundary identified as "sides"
only includes the vertices in the interior of , but here we want to impose Dirichlet boundary conditions in the closure of , i.e., also on the vertices on the contour of . Fortunately, the objects on the contour of are identified with the tag "sides_c"
(see next figure). Thus, the Dirichlet boundary can be built as the union of the objects identified as "sides"
and "sides_c"
.
Gridap provides a convenient way to create new object identifiers (referred to as "tags") from existing ones. First, we need to extract from the model, the object that holds the information about the boundary identifiers (referred to as FaceLabeling
):
labels = get_face_labeling(model)
Then, we can add new identifiers (aka "tags") to it. In the next line, we create a new tag called "diri0"
as the union of the objects identified as "sides"
and "sides_c"
, which is precisely what we need to represent the closure of the Dirichlet boundary .
add_tag_from_tags!(labels,"diri0",["sides", "sides_c"])
We follow the same approach to build a new identifier for the closure of the Dirichlet boundary . In this case, the boundary is expressed as the union of the objects identified with the tags "circle"
, "circle_c"
, "triangle"
, "triangle_c"
, "square"
, "square_c"
. Thus, we create a new tag for , called "dirig"
simply as follows:
add_tag_from_tags!(labels,"dirig",
["circle","circle_c", "triangle", "triangle_c", "square", "square_c"])
Now, we can build the FE space by using the newly defined boundary tags.
reffe = ReferenceFE(lagrangian,Float64,1)
V0 = TestFESpace(model,reffe,conformity=:H1,labels=labels,dirichlet_tags=["diri0", "dirig"])
The construction of this space is essentially the same as in the first tutorial (we build a continuous scalar-valued Lagrangian interpolation of first order). However, we also pass here the labels
object (that contains the newly created boundary tags). From this FE space, we define the trial FE spaces
g = 1
Ug = TrialFESpace(V0,[0,g])
At this point, we are ready to build the nonlinear FE problem. To this end, we need to define the weak residual and also its corresponding Jacobian. This is done following a similar procedure to the one considered in previous tutorials to define the bilinear and linear forms associated with linear FE problems. We first need to define the usual objects for numerical integration:
degree=2
Ω = Triangulation(model)
dΩ = Measure(Ω,degree)
On the one hand, the weak residual is built as follows
const p = 3
flux(∇u) = norm(∇u)^(p-2) * ∇u
f(x) = 1
res(u,v) = ∫( ∇(v)⊙(flux∘∇(u)) - v*f )*dΩ
Function res
is the one representing the integrand of the weak residual . The first argument of function res
stands for the function , where the residual is evaluated, and the second argument stands for a generic test function .
On the other hand, we (optionally) implement a function jac
representing the Jacobian.
dflux(∇du,∇u) = (p-2)*norm(∇u)^(p-4)*(∇u⊙∇du)*∇u+norm(∇u)^(p-2)*∇du
jac(u,du,v) = ∫( ∇(v)⊙(dflux∘(∇(du),∇(u))) )*dΩ
The first argument of function jac
stands for function , where the Jacobian is evaluated. The second argument is a test function , and the third argument represents an arbitrary direction .
We finally construct the nonlinear FE problem
op = FEOperator(res,jac,Ug,V0)
Here, we have constructed an instance of FEOperator
, which is the type that represents a general nonlinear FE problem in Gridap. The constructor takes the functions representing the weak residual and Jacobian, and the test and trial spaces. If only the function for the residual is provided, the Jacobian is computed internally with automatic differentiation:
op_AD = FEOperator(res,Ug,V0)
We have already built the nonlinear FE problem. Now, the remaining step is to solve it. In Gridap, nonlinear (and also linear) FE problems can be solved with instances of the type FESolver
.
We construct an instance of FESolver
as follows:
using LineSearches: BackTracking
nls = NLSolver(
show_trace=true, method=:newton, linesearch=BackTracking())
solver = FESolver(nls)
Note that the NLSolver
function used above internally calls the nlsolve
function of the NLsolve package with the provided key-word arguments. Thus, one can use any of the nonlinear methods available via the function nlsolve
to solve the nonlinear FE problem. Here, we have selected a Newton-Raphson method with a back-tracking line-search from the LineSearches package.
We are finally in place to solve the nonlinear FE problem. The initial guess is a FEFunction
, which we build from a vector of random (free) nodal values:
import Random
Random.seed!(1234)
x = rand(Float64,num_free_dofs(Ug))
uh0 = FEFunction(Ug,x)
uh, = solve!(uh0,solver,op)
We finish this tutorial by writing the computed solution for visualization (see next figure).
writevtk(Ω,datadir("p_laplacian"),cellfields=["uh"=>uh])