Line data Source code
1 : /**********************************************************************/
2 : /* DO NOT MODIFY THIS HEADER */
3 : /* MAGPIE - Mesoscale Atomistic Glue Program for Integrated Execution */
4 : /* */
5 : /* Copyright 2017 Battelle Energy Alliance, LLC */
6 : /* ALL RIGHTS RESERVED */
7 : /**********************************************************************/
8 :
9 : #pragma once
10 :
11 : #include "GeneralUserObject.h"
12 : #include "MooseVariable.h"
13 : #include "MooseMesh.h"
14 : #include "SystemBase.h"
15 :
16 : #include "spparks/src/library.h"
17 :
18 : class SPPARKSUserObject : public GeneralUserObject
19 : {
20 : public:
21 : static InputParameters validParams();
22 :
23 : SPPARKSUserObject(const InputParameters & parameters);
24 : virtual ~SPPARKSUserObject();
25 :
26 : virtual void initialSetup();
27 0 : virtual void residualSetup() {}
28 9 : virtual void timestepSetup() {}
29 :
30 : virtual void initialize();
31 : virtual void execute();
32 9 : virtual void finalize() {}
33 :
34 : int getIntValue(unsigned int fem_node_id, unsigned int index) const;
35 : Real getDoubleValue(unsigned int fem_node_id, unsigned int index) const;
36 :
37 : protected:
38 : struct FEMID;
39 : typedef FEMID SPPARKSID;
40 :
41 : void initSPPARKS();
42 :
43 : char * runSPPARKSCommand(const std::string & cmd);
44 :
45 9 : Real getSPPARKSTime(Real dt) { return dt / _time_ratio; }
46 :
47 : template <typename T>
48 : Real getValue(const std::map<unsigned int, std::map<unsigned int, T>> & data,
49 : unsigned int fem_node_id,
50 : unsigned int index) const;
51 :
52 : template <typename T>
53 : void getSPPARKSDataPointer(T *& ptr, const std::string & name, unsigned int index);
54 :
55 : template <typename T>
56 : void getSPPARKSPointer(T *& ptr, const std::string & name) const;
57 :
58 : void getSPPARKSData();
59 : void setSPPARKSData();
60 :
61 : template <typename T>
62 : void
63 : getSPPARKSData(std::map<unsigned int, T> & storage, const std::string & name, unsigned int index);
64 :
65 : template <typename T>
66 : void sendRecvSPPARKSData(const T * const data, std::map<unsigned int, T> & storage);
67 :
68 : template <typename T>
69 : void setSPPARKSData(T * data, const std::string & name, unsigned int index, MooseVariable & var);
70 :
71 : // Initiate MOOSE based on data from SPPARKS, added by YF.
72 : void setFEMData();
73 :
74 : template <typename T>
75 : void setFEMData(const std::string & name, unsigned int index, MooseVariable & var);
76 :
77 : template <typename T>
78 : void sendRecvFEMData(const std::map<libMesh::dof_id_type, T> & storage, T * const data);
79 :
80 : /// SPPARKS instance pointer
81 : void * _spparks;
82 :
83 : const std::string & _file;
84 : const bool _spparks_only;
85 : const std::vector<unsigned int> _from_ivar;
86 : const std::vector<unsigned int> _from_dvar;
87 : const std::vector<unsigned int> _to_ivar;
88 : const std::vector<unsigned int> _to_dvar;
89 : const Real _xmin;
90 : const Real _ymin;
91 : const Real _zmin;
92 : const Real _xmax;
93 : const Real _ymax;
94 : const Real _zmax;
95 : bool _init_spparks;
96 : const Real _time_ratio;
97 : bool _initialized;
98 :
99 : // added by YF, do one time SPPARKS run
100 : bool _one_time_run;
101 : int _n_spparks_run;
102 :
103 : int _dim;
104 :
105 : std::vector<MooseVariable *> _int_vars;
106 : std::vector<MooseVariable *> _double_vars;
107 : std::vector<MooseVariable *> _sol_vars; // added by YF
108 :
109 : //
110 : // Communication maps
111 : //
112 :
113 : /// SPPARKSID to vector of procs that need the value
114 : std::map<SPPARKSID, std::vector<unsigned int>> _spparks_to_proc;
115 : /// Processor to list of MOOSE FEM ids
116 : std::map<unsigned int, std::vector<libMesh::dof_id_type>> _sending_proc_to_fem_id;
117 : /// FEMID to vector of procs that need the value
118 : std::map<FEMID, std::vector<unsigned int>> _fem_to_proc;
119 : /// Processor to list of SPPARKS ids
120 : std::map<unsigned int, std::vector<unsigned int>> _sending_proc_to_spparks_id;
121 :
122 : unsigned int _num_local_fem_nodes;
123 : unsigned int _num_local_spparks_nodes;
124 :
125 : /// Local SPPARKSID to local FEMID
126 : std::map<SPPARKSID, FEMID> _spparks_to_fem;
127 :
128 : /// Local FEMID to local SPPARKSID
129 : std::multimap<FEMID, SPPARKSID> _fem_to_spparks;
130 :
131 : /// Maps from variable index to MOOSE FEM node id to value
132 : std::map<unsigned int, std::map<unsigned int, int>> _int_data_for_fem;
133 : std::map<unsigned int, std::map<unsigned int, double>> _double_data_for_fem;
134 :
135 : // Maps from variable index to value (vector index is array index)
136 : std::map<unsigned int, std::vector<int>> _int_data_for_spparks;
137 : std::map<unsigned int, std::vector<double>> _double_data_for_spparks;
138 :
139 : Real _last_time;
140 : };
141 :
142 : struct SPPARKSUserObject::FEMID
143 : {
144 8253 : FEMID(libMesh::dof_id_type ident, const Point & p) : id(ident), coor(p) {}
145 :
146 : libMesh::dof_id_type id;
147 : Point coor;
148 :
149 : // comparison operator for use in maps
150 197976 : bool operator<(const FEMID & rhs) const
151 : {
152 : // return coor < rhs.coor;
153 : const Real tol = 1e-12;
154 :
155 197976 : if (!fuzzyEqual(coor(0), rhs.coor(0), tol))
156 94680 : return (coor(0) + tol) < rhs.coor(0);
157 :
158 103296 : if (!fuzzyEqual(coor(1), rhs.coor(1), tol))
159 56028 : return (coor(1) + tol) < rhs.coor(1);
160 :
161 47268 : return (coor(2) + tol) < rhs.coor(2);
162 : }
163 :
164 : private:
165 197976 : bool fuzzyEqual(Real f, Real s, Real tol) const { return f < (s + tol) && (f > s - tol); }
166 : };
167 :
168 : template <typename T>
169 : Real
170 0 : SPPARKSUserObject::getValue(const std::map<unsigned int, std::map<unsigned int, T>> & data,
171 : unsigned int fem_node_id,
172 : unsigned int index) const
173 : {
174 : // Extract the data
175 : const auto it = data.find(index);
176 0 : if (it == data.end())
177 0 : mooseError("SPPARKSUserObject error: unknown index ", index);
178 :
179 : const auto it2 = it->second.find(fem_node_id);
180 0 : if (it2 == it->second.end())
181 : {
182 0 : mooseWarning("SPPARKSUserObject error: unknown MOOSE FEM node id ", fem_node_id);
183 0 : return 0;
184 : }
185 :
186 0 : return it2->second;
187 : }
188 :
189 : template <typename T>
190 : void
191 36 : SPPARKSUserObject::getSPPARKSDataPointer(T *& ptr, const std::string & name, unsigned int index)
192 : {
193 36 : std::stringstream fullname;
194 : fullname << name;
195 : fullname << index;
196 36 : getSPPARKSPointer(ptr, fullname.str().c_str());
197 36 : }
198 :
199 : template <typename T>
200 : void
201 81 : SPPARKSUserObject::getSPPARKSPointer(T *& ptr, const std::string & name) const
202 : {
203 81 : void * p = spparks_extract(_spparks, name.c_str());
204 81 : if (!p)
205 0 : mooseError("SPPARKS returned NULL pointer for ", name);
206 81 : ptr = static_cast<T *>(p);
207 81 : }
208 :
209 : template <typename T>
210 : void
211 36 : SPPARKSUserObject::getSPPARKSData(std::map<unsigned int, T> & storage,
212 : const std::string & name,
213 : unsigned int index)
214 : {
215 : T * data;
216 72 : getSPPARKSDataPointer(data, name.c_str(), index);
217 :
218 : // Copy data from local SPPARKS node to local MOOSE FEM node
219 : // Index into storage is MOOSE FEM node id.
220 5436 : for (auto & ids : _fem_to_spparks)
221 5400 : storage[ids.first.id] = data[ids.second.id];
222 :
223 : // Copy data across processors
224 36 : if (n_processors() > 1)
225 24 : sendRecvSPPARKSData(data, storage);
226 36 : }
227 :
228 : template <typename T>
229 : void
230 24 : SPPARKSUserObject::sendRecvSPPARKSData(const T * const data, std::map<unsigned int, T> & storage)
231 : {
232 : Parallel::MessageTag comm_tag(101);
233 :
234 : const unsigned int num_recvs = _sending_proc_to_fem_id.size();
235 24 : std::vector<Parallel::Request> recv_request(num_recvs);
236 :
237 : std::map<unsigned int, std::vector<T>>
238 : data_to_me; // sending proc, vector of SPPARKS values (one value per SPPARKS node)
239 : unsigned int offset = 0;
240 48 : for (auto & proc_id : _sending_proc_to_fem_id)
241 : {
242 24 : data_to_me[proc_id.first].resize(proc_id.second.size());
243 24 : _communicator.receive(proc_id.first, data_to_me[proc_id.first], recv_request[offset], comm_tag);
244 24 : ++offset;
245 : }
246 :
247 : std::map<unsigned int, std::vector<T>> data_from_me; // Processor, list of SPPARKS values
248 3096 : for (auto & procs : _spparks_to_proc)
249 7392 : for (auto & proc : procs.second)
250 4320 : data_from_me[proc].push_back(data[procs.first.id]);
251 :
252 48 : for (auto & data : data_from_me)
253 24 : _communicator.send(data.first, data.second, comm_tag);
254 :
255 24 : Parallel::wait(recv_request);
256 :
257 : // Move data into storage
258 48 : for (auto & proc_id : _sending_proc_to_fem_id)
259 : {
260 : const std::vector<libMesh::dof_id_type> & id = proc_id.second;
261 24 : const std::vector<T> & v = data_to_me[proc_id.first];
262 :
263 : // storage is MOOSE FEM node id, value
264 4344 : for (unsigned int j = 0; j < v.size(); ++j)
265 4320 : storage[id[j]] = v[j];
266 : }
267 24 : }
268 :
269 : template <typename T>
270 : void
271 0 : SPPARKSUserObject::setSPPARKSData(T * data,
272 : const std::string & name,
273 : unsigned int index,
274 : MooseVariable & var)
275 : {
276 0 : getSPPARKSDataPointer(data, name, index);
277 :
278 : SystemBase & sys = var.sys();
279 : NumericVector<Number> & solution = sys.solution();
280 :
281 : // Extract MOOSE data
282 : std::map<libMesh::dof_id_type, T> fem_data;
283 0 : for (auto & node : *_fe_problem.mesh().getLocalNodeRange())
284 0 : fem_data[node->id()] = solution(node->dof_number(sys.number(), var.number(), 0));
285 :
286 : // Index into data is SPPARKS node id.
287 0 : for (auto & ids : _fem_to_spparks)
288 0 : data[ids.second.id] = fem_data[ids.first.id];
289 :
290 : // Copy data across processors
291 0 : if (n_processors() > 1)
292 0 : sendRecvFEMData(fem_data, data);
293 0 : }
294 :
295 : template <typename T>
296 : void
297 9 : SPPARKSUserObject::setFEMData(const std::string & name, unsigned int index, MooseVariable & var)
298 : {
299 : // get SPPARKS data
300 : std::map<unsigned int, T> spparks_data;
301 9 : getSPPARKSData(spparks_data, name, index);
302 :
303 : SystemBase & sys = var.sys();
304 : NumericVector<Number> & solution = sys.solution();
305 :
306 : // Extract MOOSE data and set it
307 2439 : for (auto & node : *_fe_problem.mesh().getLocalNodeRange())
308 2430 : solution.set(node->dof_number(sys.number(), var.number(), 0), spparks_data[node->id()]);
309 :
310 9 : solution.close();
311 9 : }
312 :
313 : // template <typename T>
314 : // void
315 : // SPPARKSUserObject::setFEMData(T * data, const std::string & name, unsigned int index,
316 : // MooseVariable & var)
317 : // {
318 : // getSPPARKSDataPointer(data, name, index);
319 : //
320 : // SystemBase & sys = var.sys();
321 : // NumericVector<Number> & solution = sys.solution();
322 : //
323 : // // Extract MOOSE data
324 : // std::map<libMesh::dof_id_type, T> spparks_data;
325 : // ConstNodeRange & node_range = *_fe_problem.mesh().getLocalNodeRange();
326 : //
327 : // // Index into data is SPPARKS node id.
328 : // for (std::multimap<FEMID, SPPARKSID>::const_iterator i = _spparks_to_fem.begin(); i !=
329 : // _spparks_to_fem.end(); ++i)
330 : // spparks_data[i->first.id] = data[i->second.id];
331 : //
332 : // // set data
333 : // for (ConstNodeRange::const_iterator i = node_range.begin(); i < node_range.end(); ++i)
334 : // solution.set((*i)->dof_number(sys.number(), var.number(), 0), spparks_data[(*i)->id()]);
335 : //
336 : // // Copy data across processors
337 : // if (n_processors() > 1)
338 : // sendRecvSPPARKSData(data, spparks_data);
339 : // }
340 :
341 : template <typename T>
342 : void
343 0 : SPPARKSUserObject::sendRecvFEMData(const std::map<libMesh::dof_id_type, T> & storage,
344 : T * const data)
345 : {
346 : Parallel::MessageTag comm_tag(101);
347 :
348 : const unsigned int num_recvs = _sending_proc_to_spparks_id.size();
349 0 : std::vector<Parallel::Request> recv_request(num_recvs);
350 :
351 : // sending proc, vector of MOOSE values (one value per MOOSE FEM node)
352 : std::map<unsigned int, std::vector<T>> data_to_me;
353 : unsigned int offset = 0;
354 0 : for (auto & proc_id : _sending_proc_to_spparks_id)
355 : {
356 0 : data_to_me[proc_id.first].resize(proc_id.second.size());
357 0 : _communicator.receive(proc_id.first, data_to_me[proc_id.first], recv_request[offset], comm_tag);
358 0 : ++offset;
359 : }
360 :
361 : // Processor, list of MOOSE values
362 : std::map<unsigned int, std::vector<T>> data_from_me;
363 0 : for (auto & id_proc : _fem_to_proc)
364 0 : for (unsigned int j = 0; j < id_proc.second.size(); ++j)
365 0 : data_from_me[id_proc.second[j]].push_back(storage.find(id_proc.first.id)->second);
366 :
367 0 : for (auto & data : data_from_me)
368 0 : _communicator.send(data.first, data.second, comm_tag);
369 :
370 0 : Parallel::wait(recv_request);
371 :
372 : // Move data into storage
373 0 : for (auto & proc_id : _sending_proc_to_spparks_id)
374 : {
375 : const std::vector<unsigned int> & id = proc_id.second;
376 0 : const std::vector<T> & v = data_to_me[proc_id.first];
377 0 : if (id.size() != v.size())
378 0 : mooseError("Mismatched communication vectors");
379 :
380 : // data is SPPARKS index, value
381 0 : for (unsigned int j = 0; j < v.size(); ++j)
382 0 : data[id[j]] = v[j];
383 : }
384 0 : }
|