xref: /libCEED/examples/petsc/src/matops.c (revision 179e5961c1d96927213c44ef053f7de320144240)
1e83e87a5Sjeremylt #include "../include/matops.h"
22b730f8bSJeremy L Thompson 
3e83e87a5Sjeremylt #include "../include/petscutils.h"
4e83e87a5Sjeremylt 
5e83e87a5Sjeremylt // -----------------------------------------------------------------------------
66c88e6a2Srezgarshakeri // Setup apply operator context data
76c88e6a2Srezgarshakeri // -----------------------------------------------------------------------------
82b730f8bSJeremy L Thompson PetscErrorCode SetupApplyOperatorCtx(MPI_Comm comm, DM dm, Ceed ceed, CeedData ceed_data, Vec X_loc, OperatorApplyContext op_apply_ctx) {
96c88e6a2Srezgarshakeri   PetscFunctionBeginUser;
106c88e6a2Srezgarshakeri 
116c88e6a2Srezgarshakeri   op_apply_ctx->comm  = comm;
126c88e6a2Srezgarshakeri   op_apply_ctx->dm    = dm;
136c88e6a2Srezgarshakeri   op_apply_ctx->X_loc = X_loc;
142b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc, &op_apply_ctx->Y_loc));
156c88e6a2Srezgarshakeri   op_apply_ctx->x_ceed = ceed_data->x_ceed;
166c88e6a2Srezgarshakeri   op_apply_ctx->y_ceed = ceed_data->y_ceed;
176c88e6a2Srezgarshakeri   op_apply_ctx->op     = ceed_data->op_apply;
186c88e6a2Srezgarshakeri   op_apply_ctx->ceed   = ceed;
192b730f8bSJeremy L Thompson 
20ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
216c88e6a2Srezgarshakeri }
226c88e6a2Srezgarshakeri 
236c88e6a2Srezgarshakeri // -----------------------------------------------------------------------------
246c88e6a2Srezgarshakeri // Setup error operator context data
256c88e6a2Srezgarshakeri // -----------------------------------------------------------------------------
262b730f8bSJeremy L Thompson PetscErrorCode SetupErrorOperatorCtx(MPI_Comm comm, DM dm, Ceed ceed, CeedData ceed_data, Vec X_loc, CeedOperator op_error,
276c88e6a2Srezgarshakeri                                      OperatorApplyContext op_error_ctx) {
286c88e6a2Srezgarshakeri   PetscFunctionBeginUser;
296c88e6a2Srezgarshakeri 
306c88e6a2Srezgarshakeri   op_error_ctx->comm  = comm;
316c88e6a2Srezgarshakeri   op_error_ctx->dm    = dm;
326c88e6a2Srezgarshakeri   op_error_ctx->X_loc = X_loc;
332b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc, &op_error_ctx->Y_loc));
346c88e6a2Srezgarshakeri   op_error_ctx->x_ceed = ceed_data->x_ceed;
356c88e6a2Srezgarshakeri   op_error_ctx->y_ceed = ceed_data->y_ceed;
366c88e6a2Srezgarshakeri   op_error_ctx->op     = op_error;
376c88e6a2Srezgarshakeri   op_error_ctx->ceed   = ceed;
382b730f8bSJeremy L Thompson 
39ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
406c88e6a2Srezgarshakeri }
416c88e6a2Srezgarshakeri 
426c88e6a2Srezgarshakeri // -----------------------------------------------------------------------------
43e83e87a5Sjeremylt // This function returns the computed diagonal of the operator
44e83e87a5Sjeremylt // -----------------------------------------------------------------------------
45e83e87a5Sjeremylt PetscErrorCode MatGetDiag(Mat A, Vec D) {
46d4d45553Srezgarshakeri   OperatorApplyContext op_apply_ctx;
47e83e87a5Sjeremylt 
48e83e87a5Sjeremylt   PetscFunctionBeginUser;
49e83e87a5Sjeremylt 
502b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &op_apply_ctx));
51e83e87a5Sjeremylt 
52e83e87a5Sjeremylt   // Compute Diagonal via libCEED
539b072555Sjeremylt   PetscMemType mem_type;
54e83e87a5Sjeremylt 
55*179e5961SZach Atkins   // Place PETSc vector in libCEED vector
56*179e5961SZach Atkins   PetscCall(VecP2C(op_apply_ctx->Y_loc, &mem_type, op_apply_ctx->y_ceed));
57e83e87a5Sjeremylt 
58*179e5961SZach Atkins   // Compute Diagonal
592b730f8bSJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(op_apply_ctx->op, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
60e83e87a5Sjeremylt 
61*179e5961SZach Atkins   // Local-to-Global
62*179e5961SZach Atkins   PetscCall(VecC2P(op_apply_ctx->y_ceed, mem_type, op_apply_ctx->Y_loc));
632b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(D));
642b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(op_apply_ctx->dm, op_apply_ctx->Y_loc, ADD_VALUES, D));
65e83e87a5Sjeremylt 
66ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
67fd8f24faSSebastian Grimberg }
68e83e87a5Sjeremylt 
69e83e87a5Sjeremylt // -----------------------------------------------------------------------------
70ea61e9acSJeremy L Thompson // This function uses libCEED to compute the action of the Laplacian with Dirichlet boundary conditions
71e83e87a5Sjeremylt // -----------------------------------------------------------------------------
722b730f8bSJeremy L Thompson PetscErrorCode ApplyLocal_Ceed(Vec X, Vec Y, OperatorApplyContext op_apply_ctx) {
739b072555Sjeremylt   PetscMemType x_mem_type, y_mem_type;
74e83e87a5Sjeremylt 
75e83e87a5Sjeremylt   PetscFunctionBeginUser;
76e83e87a5Sjeremylt 
77e83e87a5Sjeremylt   // Global-to-local
782b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(op_apply_ctx->dm, X, INSERT_VALUES, op_apply_ctx->X_loc));
79e83e87a5Sjeremylt 
80e83e87a5Sjeremylt   // Setup libCEED vectors
81*179e5961SZach Atkins   PetscCall(VecReadP2C(op_apply_ctx->X_loc, &x_mem_type, op_apply_ctx->x_ceed));
82*179e5961SZach Atkins   PetscCall(VecP2C(op_apply_ctx->Y_loc, &y_mem_type, op_apply_ctx->y_ceed));
83e83e87a5Sjeremylt 
84e83e87a5Sjeremylt   // Apply libCEED operator
852b730f8bSJeremy L Thompson   CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
86e83e87a5Sjeremylt 
87e83e87a5Sjeremylt   // Restore PETSc vectors
88*179e5961SZach Atkins   PetscCall(VecReadC2P(op_apply_ctx->x_ceed, x_mem_type, op_apply_ctx->X_loc));
89*179e5961SZach Atkins   PetscCall(VecC2P(op_apply_ctx->y_ceed, y_mem_type, op_apply_ctx->Y_loc));
90e83e87a5Sjeremylt 
91e83e87a5Sjeremylt   // Local-to-global
922b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Y));
932b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(op_apply_ctx->dm, op_apply_ctx->Y_loc, ADD_VALUES, Y));
94e83e87a5Sjeremylt 
95ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
96fd8f24faSSebastian Grimberg }
97e83e87a5Sjeremylt 
98e83e87a5Sjeremylt // -----------------------------------------------------------------------------
99e83e87a5Sjeremylt // This function wraps the libCEED operator for a MatShell
100e83e87a5Sjeremylt // -----------------------------------------------------------------------------
101e83e87a5Sjeremylt PetscErrorCode MatMult_Ceed(Mat A, Vec X, Vec Y) {
102d4d45553Srezgarshakeri   OperatorApplyContext op_apply_ctx;
103e83e87a5Sjeremylt 
104e83e87a5Sjeremylt   PetscFunctionBeginUser;
105e83e87a5Sjeremylt 
1062b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &op_apply_ctx));
107e83e87a5Sjeremylt 
108e83e87a5Sjeremylt   // libCEED for local action of residual evaluator
1092b730f8bSJeremy L Thompson   PetscCall(ApplyLocal_Ceed(X, Y, op_apply_ctx));
110e83e87a5Sjeremylt 
111ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
112fd8f24faSSebastian Grimberg }
113e83e87a5Sjeremylt 
114e83e87a5Sjeremylt // -----------------------------------------------------------------------------
115e83e87a5Sjeremylt // This function uses libCEED to compute the action of the prolongation operator
116e83e87a5Sjeremylt // -----------------------------------------------------------------------------
117e83e87a5Sjeremylt PetscErrorCode MatMult_Prolong(Mat A, Vec X, Vec Y) {
118d4d45553Srezgarshakeri   ProlongRestrContext pr_restr_ctx;
1199b072555Sjeremylt   PetscMemType        c_mem_type, f_mem_type;
120e83e87a5Sjeremylt 
121e83e87a5Sjeremylt   PetscFunctionBeginUser;
122e83e87a5Sjeremylt 
1232b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &pr_restr_ctx));
124e83e87a5Sjeremylt 
125e83e87a5Sjeremylt   // Global-to-local
1262b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(pr_restr_ctx->loc_vec_c));
1272b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(pr_restr_ctx->dmc, X, INSERT_VALUES, pr_restr_ctx->loc_vec_c));
128e83e87a5Sjeremylt 
129e83e87a5Sjeremylt   // Setup libCEED vectors
130*179e5961SZach Atkins   PetscCall(VecReadP2C(pr_restr_ctx->loc_vec_c, &c_mem_type, pr_restr_ctx->ceed_vec_c));
131*179e5961SZach Atkins   PetscCall(VecP2C(pr_restr_ctx->loc_vec_f, &f_mem_type, pr_restr_ctx->ceed_vec_f));
132e83e87a5Sjeremylt 
133e83e87a5Sjeremylt   // Apply libCEED operator
1342b730f8bSJeremy L Thompson   CeedOperatorApply(pr_restr_ctx->op_prolong, pr_restr_ctx->ceed_vec_c, pr_restr_ctx->ceed_vec_f, CEED_REQUEST_IMMEDIATE);
135e83e87a5Sjeremylt 
136e83e87a5Sjeremylt   // Restore PETSc vectors
137*179e5961SZach Atkins   PetscCall(VecReadC2P(pr_restr_ctx->ceed_vec_c, c_mem_type, pr_restr_ctx->loc_vec_c));
138*179e5961SZach Atkins   PetscCall(VecC2P(pr_restr_ctx->ceed_vec_f, f_mem_type, pr_restr_ctx->loc_vec_f));
139e83e87a5Sjeremylt 
140e83e87a5Sjeremylt   // Multiplicity
1412b730f8bSJeremy L Thompson   PetscCall(VecPointwiseMult(pr_restr_ctx->loc_vec_f, pr_restr_ctx->loc_vec_f, pr_restr_ctx->mult_vec));
142e83e87a5Sjeremylt 
143e83e87a5Sjeremylt   // Local-to-global
1442b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Y));
1452b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(pr_restr_ctx->dmf, pr_restr_ctx->loc_vec_f, ADD_VALUES, Y));
146e83e87a5Sjeremylt 
147ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
148fd8f24faSSebastian Grimberg }
149e83e87a5Sjeremylt 
150e83e87a5Sjeremylt // -----------------------------------------------------------------------------
151e83e87a5Sjeremylt // This function uses libCEED to compute the action of the restriction operator
152e83e87a5Sjeremylt // -----------------------------------------------------------------------------
153e83e87a5Sjeremylt PetscErrorCode MatMult_Restrict(Mat A, Vec X, Vec Y) {
154d4d45553Srezgarshakeri   ProlongRestrContext pr_restr_ctx;
1559b072555Sjeremylt   PetscMemType        c_mem_type, f_mem_type;
156e83e87a5Sjeremylt 
157e83e87a5Sjeremylt   PetscFunctionBeginUser;
158e83e87a5Sjeremylt 
1592b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &pr_restr_ctx));
160e83e87a5Sjeremylt 
161e83e87a5Sjeremylt   // Global-to-local
1622b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(pr_restr_ctx->loc_vec_f));
1632b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(pr_restr_ctx->dmf, X, INSERT_VALUES, pr_restr_ctx->loc_vec_f));
164e83e87a5Sjeremylt 
165e83e87a5Sjeremylt   // Multiplicity
1662b730f8bSJeremy L Thompson   PetscCall(VecPointwiseMult(pr_restr_ctx->loc_vec_f, pr_restr_ctx->loc_vec_f, pr_restr_ctx->mult_vec));
167e83e87a5Sjeremylt 
168e83e87a5Sjeremylt   // Setup libCEED vectors
169*179e5961SZach Atkins   PetscCall(VecReadP2C(pr_restr_ctx->loc_vec_f, &f_mem_type, pr_restr_ctx->ceed_vec_f));
170*179e5961SZach Atkins   PetscCall(VecP2C(pr_restr_ctx->loc_vec_c, &c_mem_type, pr_restr_ctx->ceed_vec_c));
171e83e87a5Sjeremylt 
172e83e87a5Sjeremylt   // Apply CEED operator
1732b730f8bSJeremy L Thompson   CeedOperatorApply(pr_restr_ctx->op_restrict, pr_restr_ctx->ceed_vec_f, pr_restr_ctx->ceed_vec_c, CEED_REQUEST_IMMEDIATE);
174e83e87a5Sjeremylt 
175e83e87a5Sjeremylt   // Restore PETSc vectors
176*179e5961SZach Atkins   PetscCall(VecReadC2P(pr_restr_ctx->ceed_vec_f, f_mem_type, pr_restr_ctx->loc_vec_f));
177*179e5961SZach Atkins   PetscCall(VecC2P(pr_restr_ctx->ceed_vec_c, c_mem_type, pr_restr_ctx->loc_vec_c));
178e83e87a5Sjeremylt 
179e83e87a5Sjeremylt   // Local-to-global
1802b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Y));
1812b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(pr_restr_ctx->dmc, pr_restr_ctx->loc_vec_c, ADD_VALUES, Y));
182e83e87a5Sjeremylt 
183ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
184fd8f24faSSebastian Grimberg }
185e83e87a5Sjeremylt 
186e83e87a5Sjeremylt // -----------------------------------------------------------------------------
187e83e87a5Sjeremylt // This function calculates the error in the final solution
188e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1892b730f8bSJeremy L Thompson PetscErrorCode ComputeL2Error(Vec X, PetscScalar *l2_error, OperatorApplyContext op_error_ctx) {
19038f32c05Srezgarshakeri   Vec E;
191e83e87a5Sjeremylt   PetscFunctionBeginUser;
19238f32c05Srezgarshakeri   PetscCall(VecDuplicate(X, &E));
19338f32c05Srezgarshakeri   PetscCall(ApplyLocal_Ceed(X, E, op_error_ctx));
19438f32c05Srezgarshakeri   PetscScalar error_sq = 1.0;
19538f32c05Srezgarshakeri   PetscCall(VecSum(E, &error_sq));
19638f32c05Srezgarshakeri   *l2_error = sqrt(error_sq);
197e83e87a5Sjeremylt 
198ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
199fd8f24faSSebastian Grimberg }
200e83e87a5Sjeremylt 
201e83e87a5Sjeremylt // -----------------------------------------------------------------------------
202