v0.9.1
Classes | Functions | Variables
PCMGSetUpViaApproxOrders.hpp File Reference

header of multi-grid solver for p- adaptivity More...

Go to the source code of this file.

Classes

struct  PCMGSubMatrixCtx
 
struct  DMMGViaApproxOrdersCtx
 Structure for DM for multi-grid via approximation orders. More...
 
struct  PCMGSetUpViaApproxOrdersCtx
 Set data structures of MG pre-conditioner via approximation orders. More...
 

Functions

MoFEMErrorCode DMMGViaApproxOrdersGetCtx (DM dm, DMMGViaApproxOrdersCtx **ctx)
 
MoFEMErrorCode DMMGViaApproxOrdersSetAO (DM dm, AO ao)
 Set DM ordering. More...
 
MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize (DM dm, int *size)
 Gets size of coarseningIS in internal data struture DMMGViaApproxOrdersCtx. More...
 
MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS (DM, IS is, Mat A, Mat *subA, bool create_sub_matrix, bool shell_sub_a)
 Push back coarsening level to MG via approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS (DM)
 Pop is form MG via approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS (DM)
 Clear approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS (DM dm, IS *is_vec, int nb_elems, Mat A, int verb=0)
 Replace coarsening IS in DM via approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersGetCtx (DM dm, const DMMGViaApproxOrdersCtx **ctx)
 Get context for DM via approximation orders. More...
 
MoFEMErrorCode DMRegister_MGViaApproxOrders (const char sname[])
 Register DM for Multi-Grid via approximation orders. More...
 
MoFEMErrorCode DMCreate_MGViaApproxOrders (DM dm)
 Create DM data structure for Multi-Grid via approximation orders. More...
 
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders (DM dm, Mat *M)
 Create matrix for Multi-Grid via approximation orders. More...
 
MoFEMErrorCode DMCoarsen_MGViaApproxOrders (DM dm, MPI_Comm comm, DM *dmc)
 Coarsen DM. More...
 
MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders (DM dm1, DM dm2, Mat *mat, Vec *vec)
 Create interpolation matrix between data managers dm1 and dm2. More...
 
MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders (DM dm, Vec *g)
 Create global vector for DMGViaApproxOrders. More...
 
MoFEMErrorCode PCMGSetUpViaApproxOrders (PC pc, PCMGSetUpViaApproxOrdersCtx *ctx, int verb=0)
 Function build MG structure. More...
 

Variables

static const int DMMGVIAAPPROXORDERSCTX_INTERFACE = 1<<2
 
static const MOFEMuuid IDD_DMMGVIAAPPROXORDERSCTX = MOFEMuuid(BitIntefaceId(DMMGVIAAPPROXORDERSCTX_INTERFACE))
 

Detailed Description

header of multi-grid solver for p- adaptivity

MoFEM is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

MoFEM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with MoFEM. If not, see http://www.gnu.org/licenses/

Definition in file PCMGSetUpViaApproxOrders.hpp.

Function Documentation

◆ DMCreate_MGViaApproxOrders()

MoFEMErrorCode DMCreate_MGViaApproxOrders ( DM  dm)

Create DM data structure for Multi-Grid via approximation orders.

It set data structure and operators needed

Parameters
dmDiscrete manager
Returns
Error code

Definition at line 433 of file PCMGSetUpViaApproxOrders.cpp.

433  {
434  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
436  if (!dm->data) {
437  dm->data = new DMMGViaApproxOrdersCtx();
438  } else {
439  ((DMCtx *)(dm->data))->referenceNumber++;
440  }
441  // cerr << "Create " << ((DMCtx*)(dm->data))->referenceNumber << endl;
443  dm->ops->creatematrix = DMCreateMatrix_MGViaApproxOrders;
444  dm->ops->createglobalvector = DMCreateGlobalVector_MGViaApproxOrders;
445  dm->ops->coarsen = DMCoarsen_MGViaApproxOrders;
446  // dm->ops->destroy = DMDestroy_MGViaApproxOrders;
447  dm->ops->createinterpolation = DMCreateInterpolation_MGViaApproxOrders;
448  CHKERR DMKSPSetComputeOperators(dm, ksp_set_operators, NULL);
449  PetscInfo1(dm, "Create DMMGViaApproxOrders reference = %d\n",
450  ((DMCtx *)(dm->data))->referenceNumber);
452 }
MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc)
Coarsen DM.
static MoFEMErrorCode ksp_set_operators(KSP ksp, Mat A, Mat B, void *ctx)
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders(DM dm1, DM dm2, Mat *mat, Vec *vec)
Create interpolation matrix between data managers dm1 and dm2.
MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders(DM dm, Vec *g)
Create global vector for DMGViaApproxOrders.
Structure for DM for multi-grid via approximation orders.
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders(DM dm, Mat *M)
Create matrix for Multi-Grid via approximation orders.
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition: DMMMoFEM.cpp:54
#define CHKERR
Inline error check.
Definition: definitions.h:602
PETSc Discrete Manager data structure.
Definition: DMMoFEM.hpp:877
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ DMCreateGlobalVector_MGViaApproxOrders()

MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders ( DM  dm,
Vec *  g 
)

Create global vector for DMGViaApproxOrders.

Parameters
dmDistributed mesh data structure
greturned pointer to vector
Returns
Error code

Definition at line 599 of file PCMGSetUpViaApproxOrders.cpp.

599  {
600 
601  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
603  int leveldown = dm->leveldown;
604  GET_DM_FIELD(dm);
605  if (dm_field->kspOperators.empty()) {
607  } else {
608 #if PETSC_VERSION_GE(3, 5, 3)
609  CHKERR MatCreateVecs(
610  dm_field->kspOperators[dm_field->kspOperators.size() - 1 - leveldown],
611  g, NULL);
612 #else
613  CHKERR MatGetVecs(
614  dm_field->kspOperators[dm_field->kspOperators.size() - 1 - leveldown],
615  g, NULL);
616 #endif
617  }
618  PetscInfo1(dm, "Create global vector DMMGViaApproxOrders leveldown = %d\n",
619  dm->leveldown);
621 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVectorsets the routine to create a global vector associated with the shell DM.
Definition: DMMMoFEM.cpp:1026
#define GET_DM_FIELD(DM)
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ DMCreateInterpolation_MGViaApproxOrders()

MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders ( DM  dm1,
DM  dm2,
Mat *  mat,
Vec *  vec 
)

Create interpolation matrix between data managers dm1 and dm2.

Parameters
dm1Distributed mesh data structure
dm2Distributed mesh data structure
matPointer to returned interpolation matrix
vecPointer to scaling vector here returned NULL
Returns
Error code

Definition at line 500 of file PCMGSetUpViaApproxOrders.cpp.

501  {
502  PetscValidHeaderSpecific(dm1, DM_CLASSID, 1);
503  PetscValidHeaderSpecific(dm2, DM_CLASSID, 1);
505 
506  MPI_Comm comm;
507  CHKERR PetscObjectGetComm((PetscObject)dm1, &comm);
508 
509  int m, n, M, N;
510 
511  DM dm_down = dm1;
512  DM dm_up = dm2;
513 
514  int dm_down_leveldown = dm_down->leveldown;
515  int dm_up_leveldown = dm_up->leveldown;
516 
517  PetscInfo2(dm1,
518  "Create interpolation DMMGViaApproxOrders dm1_leveldown = %d "
519  "dm2_leveldown = %d\n",
520  dm_down_leveldown, dm_up_leveldown);
521 
522  IS is_down, is_up;
523  {
524  // Coarser mesh
525  GET_DM_FIELD(dm_down);
526  if (static_cast<int>(dm_field->coarseningIS.size()) < dm_down_leveldown) {
527  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
528  }
529  is_down = dm_field->coarseningIS[dm_field->coarseningIS.size() - 1 -
530  dm_down_leveldown];
531  CHKERR ISGetSize(is_down, &M);
532  CHKERR ISGetLocalSize(is_down, &m);
533  }
534  {
535  // Finer mesh
536  GET_DM_FIELD(dm_up);
537  if (static_cast<int>(dm_field->coarseningIS.size()) < dm_up_leveldown) {
538  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
539  }
540  is_up =
541  dm_field
542  ->coarseningIS[dm_field->coarseningIS.size() - 1 - dm_up_leveldown];
543  CHKERR ISGetSize(is_up, &N);
544  CHKERR ISGetLocalSize(is_up, &n);
545  }
546 
547  // is_dow rows
548  // is_up columns
549 
550  CHKERR MatCreate(comm, mat);
551  CHKERR MatSetSizes(*mat, m, n, M, N);
552  CHKERR MatSetType(*mat, MATMPIAIJ);
553  CHKERR MatMPIAIJSetPreallocation(*mat, 1, PETSC_NULL, 0, PETSC_NULL);
554 
555  // get matrix layout
556  PetscLayout rmap, cmap;
557  CHKERR MatGetLayouts(*mat, &rmap, &cmap);
558  int rstart, rend, cstart, cend;
559  CHKERR PetscLayoutGetRange(rmap, &rstart, &rend);
560  CHKERR PetscLayoutGetRange(cmap, &cstart, &cend);
561 
562  // if(verb>0) {
563  // PetscSynchronizedPrintf(comm,"level %d row start %d row end
564  // %d\n",kk,rstart,rend); PetscSynchronizedPrintf(comm,"level %d col start
565  // %d col end %d\n",kk,cstart,cend);
566  // }
567 
568  const int *row_indices_ptr, *col_indices_ptr;
569  CHKERR ISGetIndices(is_down, &row_indices_ptr);
570  CHKERR ISGetIndices(is_up, &col_indices_ptr);
571 
572  map<int, int> idx_map;
573  for (int ii = 0; ii < m; ii++) {
574  idx_map[row_indices_ptr[ii]] = rstart + ii;
575  }
576 
577  CHKERR MatZeroEntries(*mat);
578  // FIXME: Use MatCreateMPIAIJWithArrays and set array directly
579  for (int jj = 0; jj < n; jj++) {
580  map<int, int>::iterator mit = idx_map.find(col_indices_ptr[jj]);
581  if (mit != idx_map.end()) {
582  CHKERR MatSetValue(*mat, mit->second, cstart + jj, 1, INSERT_VALUES);
583  }
584  }
585 
586  CHKERR ISRestoreIndices(is_down, &row_indices_ptr);
587  CHKERR ISRestoreIndices(is_up, &col_indices_ptr);
588 
589  CHKERR MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY);
590  CHKERR MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY);
591 
592  if (vec != NULL) {
593  *vec = PETSC_NULL;
594  }
595 
597 }
FTensor::Index< 'n', 2 > n
Definition: PlasticOps.hpp:68
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define GET_DM_FIELD(DM)
FTensor::Index< 'm', 2 > m
Definition: PlasticOps.hpp:67
#define CHKERR
Inline error check.
Definition: definitions.h:602
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:72
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
const int N
Definition: speed_test.cpp:3

◆ DMMGViaApproxOrdersGetCoarseningISSize()

MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize ( DM  dm,
int *  size 
)

Gets size of coarseningIS in internal data struture DMMGViaApproxOrdersCtx.

Parameters
dmDM
sizesize of coarseningIS
Returns
Error code

Definition at line 222 of file PCMGSetUpViaApproxOrders.cpp.

222  {
223  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
225  GET_DM_FIELD(dm);
226  *size = dm_field->coarseningIS.size();
228 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
#define GET_DM_FIELD(DM)

◆ DMMGViaApproxOrdersGetCtx() [1/2]

MoFEMErrorCode DMMGViaApproxOrdersGetCtx ( DM  dm,
DMMGViaApproxOrdersCtx **  ctx 
)

Get DM Ctx

Definition at line 199 of file PCMGSetUpViaApproxOrders.cpp.

199  {
200  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
202  GET_DM_FIELD(dm);
203  *ctx = dm_field;
205 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
#define GET_DM_FIELD(DM)

◆ DMMGViaApproxOrdersGetCtx() [2/2]

MoFEMErrorCode DMMGViaApproxOrdersGetCtx ( DM  dm,
const DMMGViaApproxOrdersCtx **  ctx 
)

Get context for DM via approximation orders.

Parameters
dmthe DM object
ctxdata context
Returns
error code

Definition at line 391 of file PCMGSetUpViaApproxOrders.cpp.

392  {
393  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
395  GET_DM_FIELD(dm);
396  *ctx = dm_field;
398 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
#define GET_DM_FIELD(DM)

◆ DMMGViaApproxOrdersReplaceCoarseningIS()

MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS ( DM  dm,
IS *  is_vec,
int  nb_elems,
Mat  A,
int  verb = 0 
)

Replace coarsening IS in DM via approximation orders.

Parameters
dmdm
is_vecPointer to vector of is
nb_elemsNumber of elements
AFine matrix
Returns
Error code

Definition at line 309 of file PCMGSetUpViaApproxOrders.cpp.

311  {
312  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
314  GET_DM_FIELD(dm);
315  int nb_no_changed = 0;
316  int nb_replaced = 0;
317  int nb_deleted = 0;
318  int nb_added = 0;
319  std::vector<IS>::iterator it;
320  it = dm_field->coarseningIS.begin();
321  int ii = 0;
322  for (; it != dm_field->coarseningIS.end(); it++, ii++) {
323  if (ii < nb_elems) {
324  PetscBool flg;
325  CHKERR ISEqual(*it, is_vec[ii], &flg);
326  if (!flg) {
327  CHKERR ISDestroy(&*it);
328  CHKERR MatDestroy(&dm_field->kspOperators[ii]);
329  *it = is_vec[ii];
330  CHKERR PetscObjectReference((PetscObject)is_vec[ii]);
331  if (ii < nb_elems - 1) {
332  IS is = is_vec[ii];
333  if (dm_field->aO) {
334  CHKERR ISDuplicate(is_vec[ii], &is);
335  CHKERR ISCopy(is_vec[ii], is);
336  CHKERR AOApplicationToPetscIS(dm_field->aO, is);
337  }
338  Mat subA;
339  #if PETSC_VERSION_GE(3, 8, 0)
340  CHKERR MatCreateSubMatrix(A, is, is, MAT_INITIAL_MATRIX, &subA);
341  #else
342  CHKERR MatGetSubMatrix(A, is, is, MAT_INITIAL_MATRIX, &subA);
343  #endif
344  CHKERR PetscObjectReference((PetscObject)subA);
345  dm_field->kspOperators[ii] = subA;
346  CHKERR MatDestroy(&subA);
347  if (dm_field->aO) {
348  CHKERR ISDestroy(&is);
349  }
350  } else {
351  CHKERR PetscObjectReference((PetscObject)A);
352  dm_field->kspOperators[ii] = A;
353  }
354  nb_replaced++;
355  }
356  } else {
357  nb_no_changed++;
358  continue;
359  }
360  }
361  if (static_cast<int>(dm_field->coarseningIS.size()) < nb_elems) {
362  for (; ii < nb_elems - 1; ii++) {
363  Mat subA;
364  CHKERR DMMGViaApproxOrdersPushBackCoarseningIS(dm, is_vec[ii], A, &subA,
365  true, false);
366  CHKERR MatDestroy(&subA);
367  nb_added++;
368  }
369  CHKERR DMMGViaApproxOrdersPushBackCoarseningIS(dm, is_vec[ii], A, &A, false,
370  false);
371  nb_added++;
372  } else {
373  for (; ii < static_cast<int>(dm_field->coarseningIS.size()); ii++) {
375  nb_deleted++;
376  }
377  }
378  MPI_Comm comm;
379  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
380  if (verb > 0) {
381  PetscPrintf(comm,
382  "DMMGViaApproxOrders nb_no_changed = %d, nb_replaced = %d, "
383  "nb_added = %d, nb_deleted = %d, size = %d\n",
384  nb_no_changed, nb_replaced, nb_added, nb_deleted,
385  dm_field->coarseningIS.size());
386  }
387  PetscInfo(dm, "Replace IS to DMMGViaApproxOrders\n");
389 }
MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS(DM dm, IS is, Mat A, Mat *subA, bool create_sub_matrix, bool shell_sub_a)
Push back coarsening level to MG via approximation orders.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define GET_DM_FIELD(DM)
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS(DM dm)
Pop is form MG via approximation orders.
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ DMMGViaApproxOrdersSetAO()

MoFEMErrorCode DMMGViaApproxOrdersSetAO ( DM  dm,
AO  ao 
)

Set DM ordering.

IS can be given is some other ordering, AO will transform indices from coarseningIS ordering to ordering used to construct fine matrix.

Parameters
dm[description]
ao[description]
Returns
[description]

Definition at line 207 of file PCMGSetUpViaApproxOrders.cpp.

207  {
208  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
210  GET_DM_FIELD(dm);
211  if (dm_field->aO) {
212  // std::cerr << dm_field->aO << std::endl;
213  CHKERR AODestroy(&dm_field->aO);
214  // std::cerr << "destroy ao when adding\n";
215  }
216  dm_field->aO = ao;
217  CHKERR PetscObjectReference((PetscObject)ao);
218  // std::cerr << "add ao\n";
220 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define GET_DM_FIELD(DM)
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ PCMGSetUpViaApproxOrders()

MoFEMErrorCode PCMGSetUpViaApproxOrders ( PC  pc,
PCMGSetUpViaApproxOrdersCtx ctx,
int  verb = 0 
)

Function build MG structure.

Parameters
pcMG pre-conditioner http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html
ctxMoFEM data structure for MG
verbverbosity level
Returns
error code

Definition at line 750 of file PCMGSetUpViaApproxOrders.cpp.

751  {
752 
753  PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
755 
756  MPI_Comm comm;
757  CHKERR PetscObjectGetComm((PetscObject)pc, &comm);
758  if (verb > 0) {
759  PetscPrintf(comm, "Start PCMGSetUpViaApproxOrders\n");
760  }
761 
762  CHKERR ctx->getOptions();
763  CHKERR ctx->buildProlongationOperator(true, verb);
764 
765 #if PETSC_VERSION_GE(3, 8, 0)
766  CHKERR PCMGSetGalerkin(pc, PC_MG_GALERKIN_NONE);
767 #else
768  CHKERR PCMGSetGalerkin(pc, PETSC_FALSE);
769 #endif
770 
771  CHKERR PCMGSetLevels(pc, ctx->nbLevels, NULL);
772 
773  if (verb > 0) {
774  PetscPrintf(comm, "End PCMGSetUpViaApproxOrders\n");
775  }
776 
778 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
int nbLevels
number of multi-grid levels
virtual MoFEMErrorCode buildProlongationOperator(bool use_mat_a, int verb=0)
Set up data structures for MG.
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
virtual MoFEMErrorCode getOptions()
get options from line command

Variable Documentation

◆ DMMGVIAAPPROXORDERSCTX_INTERFACE

const int DMMGVIAAPPROXORDERSCTX_INTERFACE = 1<<2
static

Definition at line 21 of file PCMGSetUpViaApproxOrders.hpp.

◆ IDD_DMMGVIAAPPROXORDERSCTX

const MOFEMuuid IDD_DMMGVIAAPPROXORDERSCTX = MOFEMuuid(BitIntefaceId(DMMGVIAAPPROXORDERSCTX_INTERFACE))
static

Definition at line 22 of file PCMGSetUpViaApproxOrders.hpp.