xref: /libCEED/examples/solids/src/matops.c (revision ccaff0309dc399f656ea11018b919b30feb8b0fa)
1*ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2*ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3*ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4*ccaff030SJeremy L Thompson //
5*ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6*ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7*ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8*ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9*ccaff030SJeremy L Thompson //
10*ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11*ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12*ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13*ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14*ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15*ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16*ccaff030SJeremy L Thompson 
17*ccaff030SJeremy L Thompson /// @file
18*ccaff030SJeremy L Thompson /// Matrix shell operations for solid mechanics example using PETSc
19*ccaff030SJeremy L Thompson 
20*ccaff030SJeremy L Thompson #include "../elasticity.h"
21*ccaff030SJeremy L Thompson 
22*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
23*ccaff030SJeremy L Thompson // libCEED Operators for MatShell
24*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
25*ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator
26*ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) {
27*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
28*ccaff030SJeremy L Thompson   PetscScalar *x, *y;
29*ccaff030SJeremy L Thompson 
30*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
31*ccaff030SJeremy L Thompson 
32*ccaff030SJeremy L Thompson   // Global-to-local
33*ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
34*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
35*ccaff030SJeremy L Thompson 
36*ccaff030SJeremy L Thompson   // Setup CEED vectors
37*ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
38*ccaff030SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
39*ccaff030SJeremy L Thompson   CeedVectorSetArray(user->Xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
40*ccaff030SJeremy L Thompson   CeedVectorSetArray(user->Yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
41*ccaff030SJeremy L Thompson 
42*ccaff030SJeremy L Thompson   // Apply CEED operator
43*ccaff030SJeremy L Thompson   CeedOperatorApply(user->op, user->Xceed, user->Yceed, CEED_REQUEST_IMMEDIATE);
44*ccaff030SJeremy L Thompson   CeedVectorSyncArray(user->Yceed, CEED_MEM_HOST);
45*ccaff030SJeremy L Thompson 
46*ccaff030SJeremy L Thompson   // Restore PETSc vectors
47*ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
48*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
49*ccaff030SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
50*ccaff030SJeremy L Thompson 
51*ccaff030SJeremy L Thompson   // Local-to-global
52*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
53*ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, user->Yloc, ADD_VALUES, Y); CHKERRQ(ierr);
54*ccaff030SJeremy L Thompson 
55*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
56*ccaff030SJeremy L Thompson };
57*ccaff030SJeremy L Thompson 
58*ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual
59*ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
60*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
61*ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
62*ccaff030SJeremy L Thompson 
63*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
64*ccaff030SJeremy L Thompson 
65*ccaff030SJeremy L Thompson   // Use computed BCs
66*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
67*ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
68*ccaff030SJeremy L Thompson                                     user->loadIncrement, NULL, NULL, NULL);
69*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
70*ccaff030SJeremy L Thompson 
71*ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
72*ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
73*ccaff030SJeremy L Thompson 
74*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
75*ccaff030SJeremy L Thompson };
76*ccaff030SJeremy L Thompson 
77*ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES
78*ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
79*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
80*ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
81*ccaff030SJeremy L Thompson 
82*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
83*ccaff030SJeremy L Thompson 
84*ccaff030SJeremy L Thompson   // Use computed BCs
85*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
86*ccaff030SJeremy L Thompson 
87*ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
88*ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
89*ccaff030SJeremy L Thompson 
90*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
91*ccaff030SJeremy L Thompson };
92*ccaff030SJeremy L Thompson 
93*ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian
94*ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) {
95*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
96*ccaff030SJeremy L Thompson   UserMult user;
97*ccaff030SJeremy L Thompson 
98*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
99*ccaff030SJeremy L Thompson 
100*ccaff030SJeremy L Thompson   // Zero boundary values
101*ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
102*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
103*ccaff030SJeremy L Thompson 
104*ccaff030SJeremy L Thompson   // libCEED for local action of Jacobian
105*ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
106*ccaff030SJeremy L Thompson 
107*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
108*ccaff030SJeremy L Thompson };
109*ccaff030SJeremy L Thompson 
110*ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator
111*ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) {
112*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
113*ccaff030SJeremy L Thompson   UserMultProlongRestr user;
114*ccaff030SJeremy L Thompson   PetscScalar *c, *f;
115*ccaff030SJeremy L Thompson 
116*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
117*ccaff030SJeremy L Thompson 
118*ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
119*ccaff030SJeremy L Thompson 
120*ccaff030SJeremy L Thompson   // Global-to-local
121*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
122*ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dmC, X, INSERT_VALUES, user->locVecC);
123*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
124*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
125*ccaff030SJeremy L Thompson 
126*ccaff030SJeremy L Thompson   // Setup CEED vectors
127*ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(user->locVecC, (const PetscScalar **)&c);
128*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
129*ccaff030SJeremy L Thompson   ierr = VecGetArray(user->locVecF, &f); CHKERRQ(ierr);
130*ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecC, CEED_MEM_HOST, CEED_USE_POINTER, c);
131*ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecF, CEED_MEM_HOST, CEED_USE_POINTER, f);
132*ccaff030SJeremy L Thompson 
133*ccaff030SJeremy L Thompson   // Apply CEED operator
134*ccaff030SJeremy L Thompson   CeedOperatorApply(user->opProlong, user->ceedVecC, user->ceedVecF,
135*ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
136*ccaff030SJeremy L Thompson   CeedVectorSyncArray(user->ceedVecF, CEED_MEM_HOST);
137*ccaff030SJeremy L Thompson 
138*ccaff030SJeremy L Thompson   // Restore PETSc vectors
139*ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(user->locVecC, (const PetscScalar **)c);
140*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
141*ccaff030SJeremy L Thompson   ierr = VecRestoreArray(user->locVecF, &f); CHKERRQ(ierr);
142*ccaff030SJeremy L Thompson 
143*ccaff030SJeremy L Thompson   // Multiplicity
144*ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
145*ccaff030SJeremy L Thompson 
146*ccaff030SJeremy L Thompson   // Local-to-global
147*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
148*ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dmF, user->locVecF, ADD_VALUES, Y);
149*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
150*ccaff030SJeremy L Thompson 
151*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
152*ccaff030SJeremy L Thompson }
153*ccaff030SJeremy L Thompson 
154*ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator
155*ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
156*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
157*ccaff030SJeremy L Thompson   UserMultProlongRestr user;
158*ccaff030SJeremy L Thompson   PetscScalar *c, *f;
159*ccaff030SJeremy L Thompson 
160*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
161*ccaff030SJeremy L Thompson 
162*ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
163*ccaff030SJeremy L Thompson 
164*ccaff030SJeremy L Thompson   // Global-to-local
165*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
166*ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dmF, X, INSERT_VALUES, user->locVecF);
167*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
168*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
169*ccaff030SJeremy L Thompson 
170*ccaff030SJeremy L Thompson   // Multiplicity
171*ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
172*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
173*ccaff030SJeremy L Thompson 
174*ccaff030SJeremy L Thompson   // Setup CEED vectors
175*ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(user->locVecF, (const PetscScalar **)&f); CHKERRQ(ierr);
176*ccaff030SJeremy L Thompson   ierr = VecGetArray(user->locVecC, &c); CHKERRQ(ierr);
177*ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecF, CEED_MEM_HOST, CEED_USE_POINTER, f);
178*ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecC, CEED_MEM_HOST, CEED_USE_POINTER, c);
179*ccaff030SJeremy L Thompson 
180*ccaff030SJeremy L Thompson   // Apply CEED operator
181*ccaff030SJeremy L Thompson   CeedOperatorApply(user->opRestrict, user->ceedVecF, user->ceedVecC,
182*ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
183*ccaff030SJeremy L Thompson   CeedVectorSyncArray(user->ceedVecC, CEED_MEM_HOST);
184*ccaff030SJeremy L Thompson 
185*ccaff030SJeremy L Thompson   // Restore PETSc vectors
186*ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(user->locVecF, (const PetscScalar **)&f);
187*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
188*ccaff030SJeremy L Thompson   ierr = VecRestoreArray(user->locVecC, &c); CHKERRQ(ierr);
189*ccaff030SJeremy L Thompson 
190*ccaff030SJeremy L Thompson   // Local-to-global
191*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
192*ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dmC, user->locVecC, ADD_VALUES, Y);
193*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
194*ccaff030SJeremy L Thompson 
195*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
196*ccaff030SJeremy L Thompson };
197*ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator
198*ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
199*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
200*ccaff030SJeremy L Thompson   UserMult user;
201*ccaff030SJeremy L Thompson 
202*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
203*ccaff030SJeremy L Thompson 
204*ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
205*ccaff030SJeremy L Thompson 
206*ccaff030SJeremy L Thompson   // Compute Diagonal via libCEED
207*ccaff030SJeremy L Thompson   CeedVector ceedDiagVec;
208*ccaff030SJeremy L Thompson   const CeedScalar *diagArray;
209*ccaff030SJeremy L Thompson 
210*ccaff030SJeremy L Thompson   // -- Compute Diagonal
211*ccaff030SJeremy L Thompson   CeedOperatorAssembleLinearDiagonal(user->op, &ceedDiagVec,
212*ccaff030SJeremy L Thompson                                      CEED_REQUEST_IMMEDIATE);
213*ccaff030SJeremy L Thompson 
214*ccaff030SJeremy L Thompson   // -- Place in PETSc vector
215*ccaff030SJeremy L Thompson   CeedVectorGetArrayRead(ceedDiagVec, CEED_MEM_HOST, &diagArray);
216*ccaff030SJeremy L Thompson   ierr = VecPlaceArray(user->Xloc, diagArray); CHKERRQ(ierr);
217*ccaff030SJeremy L Thompson 
218*ccaff030SJeremy L Thompson   // -- Local-to-Global
219*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(D); CHKERRQ(ierr);
220*ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, user->Xloc, ADD_VALUES, D); CHKERRQ(ierr);
221*ccaff030SJeremy L Thompson 
222*ccaff030SJeremy L Thompson   // -- Cleanup
223*ccaff030SJeremy L Thompson   ierr = VecResetArray(user->Xloc); CHKERRQ(ierr);
224*ccaff030SJeremy L Thompson   CeedVectorRestoreArrayRead(ceedDiagVec, &diagArray);
225*ccaff030SJeremy L Thompson   CeedVectorDestroy(&ceedDiagVec);
226*ccaff030SJeremy L Thompson 
227*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
228*ccaff030SJeremy L Thompson };
229