v0.14.0
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
test_cache_on_entities.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  MyStorage
 
struct  OpVolumeSet
 Operator set cache stored data, in this example, global indices, but it can be any structure. More...
 
struct  OpVolumeTest
 Test if cached data can be accessed, and check consistency of data. More...
 
struct  OpVolumeAssemble
 
struct  VolRule
 Set integration rule. More...
 

Typedefs

using VolEle = VolumeElementForcesAndSourcesCore
 
using VolOp = VolEle::UserDataOperator
 
using FaceEle = FaceElementForcesAndSourcesCore
 
using FaceOp = FaceEle::UserDataOperator
 

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 

Typedef Documentation

◆ FaceEle

Definition at line 24 of file test_cache_on_entities.cpp.

◆ FaceOp

Definition at line 25 of file test_cache_on_entities.cpp.

◆ VolEle

Definition at line 22 of file test_cache_on_entities.cpp.

◆ VolOp

Examples
build_large_problem.cpp, and lorentz_force.cpp.

Definition at line 23 of file test_cache_on_entities.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 167 of file test_cache_on_entities.cpp.

167 {
168
169 // initialize petsc
170 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
171
172 try {
173
174 // Create MoAB database
175 moab::Core moab_core;
176 moab::Interface &moab = moab_core;
177
178 // Create MoFEM database and link it to MoAB
179 MoFEM::Core mofem_core(moab);
180 MoFEM::Interface &m_field = mofem_core;
181
182 // Register DM Manager
183 DMType dm_name = "DMMOFEM";
184 CHKERR DMRegister_MoFEM(dm_name);
185
186 // Simple interface
187 Simple *simple_interface;
188 CHKERR m_field.getInterface(simple_interface);
189 {
190 // get options from command line
191 CHKERR simple_interface->getOptions();
192 // load mesh file
193 CHKERR simple_interface->loadFile();
194 // add fields
195 CHKERR simple_interface->addDomainField("FIELD", H1,
197 // set fields order
198 CHKERR simple_interface->setFieldOrder("FIELD", 4);
199 // setup problem
200 CHKERR simple_interface->setUp();
201
202 // get dm
203 auto dm = simple_interface->getDM();
204 // create elements
205 auto domain_fe = boost::make_shared<VolEle>(m_field);
206 // set integration rule
207 domain_fe->getRuleHook = VolRule();
208
209 auto get_skin = [&]() {
210 Range ents;
211 CHKERR m_field.get_moab().get_entities_by_dimension(0, 3, ents, true);
212 Skinner skin(&m_field.get_moab());
213 Range skin_ents;
214 CHKERR skin.find_skin(0, ents, false, skin_ents);
215 // filter not owned entities, those are not on boundary
216 ParallelComm *pcomm =
217 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
218 Range proc_skin;
219 CHKERR pcomm->filter_pstatus(skin_ents,
220 PSTATUS_SHARED | PSTATUS_MULTISHARED,
221 PSTATUS_NOT, -1, &proc_skin);
222 Range adj;
223 for (auto d : {0, 1, 2})
224 CHKERR m_field.get_moab().get_adjacencies(proc_skin, d, false, adj,
225 moab::Interface::UNION);
226 // proc_skin.merge(adj);
227 return proc_skin;
228 };
229
230 auto get_mark_skin_dofs = [&](Range &&skin) {
231 auto problem_manager = m_field.getInterface<ProblemsManager>();
232 auto marker_ptr = boost::make_shared<std::vector<unsigned char>>();
233 problem_manager->markDofs(simple_interface->getProblemName(), ROW,
235 *marker_ptr);
236 return marker_ptr;
237 };
238
241 auto mark_dofs = get_mark_skin_dofs(get_skin());
242
243 // set operator to the volume elements
244 domain_fe->getOpPtrVector().push_back(new OpVolumeSet("FIELD"));
245 domain_fe->getOpPtrVector().push_back(new OpVolumeTest("FIELD"));
246
247 domain_fe->getOpPtrVector().push_back(
248 new OpSetBc("FIELD", true, mark_dofs));
249 domain_fe->getOpPtrVector().push_back(new OpVolumeAssemble("FIELD", v));
250 domain_fe->getOpPtrVector().push_back(new OpUnSetBc("FIELD"));
251
252 // make integration in volume (here real calculations starts)
253 CHKERR DMoFEMLoopFiniteElements(dm, simple_interface->getDomainFEName(),
254 domain_fe);
256
258 CHKERR VecAssemblyBegin(v);
259 CHKERR VecAssemblyEnd(v);
260 CHKERR VecGhostUpdateBegin(v, ADD_VALUES, SCATTER_REVERSE);
261 CHKERR VecGhostUpdateEnd(v, ADD_VALUES, SCATTER_REVERSE);
262 CHKERR VecGhostUpdateBegin(v, INSERT_VALUES, SCATTER_FORWARD);
263 CHKERR VecGhostUpdateEnd(v, INSERT_VALUES, SCATTER_FORWARD);
264
265 CHKERR DMoFEMMeshToLocalVector(simple_interface->getDM(), v,
266 INSERT_VALUES, SCATTER_REVERSE);
267
268 // CHKERR VecView(v,PETSC_VIEWER_STDOUT_WORLD);
269
270 const MoFEM::Problem *problem_ptr;
271 CHKERR DMMoFEMGetProblemPtr(dm, &problem_ptr);
272
273 double *array;
274 CHKERR VecGetArray(v, &array);
275 for (auto dof : *(problem_ptr->getNumeredRowDofsPtr())) {
276 auto loc_idx = dof->getPetscLocalDofIdx();
277
278 if (dof->getFieldData() != array[loc_idx])
279 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
280 "Data inconsistency");
281
282 if ((*mark_dofs)[loc_idx]) {
283 if (array[loc_idx] != 0)
284 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
285 "Dof assembled, but it should be not");
286 } else {
287 if (array[loc_idx] == 0)
288 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
289 "Dof not assembled, but it should be");
290 }
291 }
292
293 CHKERR VecRestoreArray(v, &array);
294 }
295 }
297
298 // finish work cleaning memory, getting statistics, etc.
300
301 return 0;
302}
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
Definition: LogManager.hpp:345
@ ROW
Definition: definitions.h:123
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#define CHKERR
Inline error check.
Definition: definitions.h:535
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:509
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:412
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:572
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1153
const double v
phase velocity of light in medium (cm/ns)
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
Core (interface) class.
Definition: Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Deprecated interface functions.
Set indices on entities on finite element.
keeps basic data about problem
auto & getNumeredRowDofsPtr() const
get access to numeredRowDofsPtr storing DOFs on rows
Problem manager is used to build and partition problems.
MoFEMErrorCode markDofs(const std::string problem_name, RowColData rc, const enum MarkOP op, const Range ents, std::vector< unsigned char > &marker) const
Create vector with marked indices.
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
MoFEMErrorCode addDomainField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:264
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:667
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:362
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition: Simple.hpp:341
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Operator set cache stored data, in this example, global indices, but it can be any structure.
static std::vector< SharedVecInt > entsIndices
This is global static storage.
Test if cached data can be accessed, and check consistency of data.
Set integration rule.
static char help[]

Variable Documentation

◆ help

char help[] = "...\n\n"
static

Definition at line 20 of file test_cache_on_entities.cpp.