MPQC 3.0.0-alpha
Loading...
Searching...
No Matches
sr_r12intermediates_util.h
1//
2// sr_r12intermediates_util.h
3//
4// Copyright (C) 2013 Edward Valeev
5//
6// Author: Edward Valeev <evaleev@vt.edu>
7// Maintainer: EV
8//
9// This file is part of the SC Toolkit.
10//
11// The SC Toolkit is free software; you can redistribute it and/or modify
12// it under the terms of the GNU Library General Public License as published by
13// the Free Software Foundation; either version 2, or (at your option)
14// any later version.
15//
16// The SC Toolkit is distributed in the hope that it will be useful,
17// but WITHOUT ANY WARRANTY; without even the implied warranty of
18// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19// GNU Library General Public License for more details.
20//
21// You should have received a copy of the GNU Library General Public License
22// along with the SC Toolkit; see the file COPYING.LIB. If not, write to
23// the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
24//
25// The U.S. Government is granted a limited license as per AL 91-7.
26//
27
28// this should not be included directly, will be included by sr_r12intermediates.h
29
30#ifndef _mpqc_src_lib_chemistry_qc_mbptr12_srr12intermediatesutil_h
31#define _mpqc_src_lib_chemistry_qc_mbptr12_srr12intermediatesutil_h
32
33#include <algorithm>
34#include <numeric>
35#include <TiledArray/bitset.h>
36#include <TiledArray/range1.h>
37
38namespace sc {
39
40 using integer = TA::Range1::index1_type;
41
42 template <typename T>
43 std::shared_ptr<typename SingleReference_R12Intermediates<T>::TArray22>
44 SingleReference_R12Intermediates<T>::ab_O_cd(const std::string& key,
45 const int te_type) {
46
47 typedef typename TArray22::value_type value_type;
48
49 Ref<TwoBodyMOIntsTransform> tform = r12world_->world()->moints_runtime4()->get(key);
50 tform->compute();
51 Ref<DistArray4> darray4 = tform->ints_distarray4();
52 darray4->activate();
53
54 // this only works if all nodes can read integrals
55 std::vector<int> tasks_with_access;
56 const int ntasks_with_access = darray4->tasks_with_access(tasks_with_access);
57 MPQC_ASSERT(ntasks_with_access == world_.nproc());
58 MPQC_ASSERT(ntasks_with_access == r12world_->world()->msg()->n());
59
60 const size_t n1 = darray4->ni();
61 const size_t n2 = darray4->nj();
62 const size_t n3 = darray4->nx();
63 const size_t n4 = darray4->ny();
64
65 // make tiled ranges for occupied indices
66 // N.B. using little trickery using partial sum to create tilesize=1 range
67 std::vector<size_t> i1_hashmarks(n1+1, 1);
68 i1_hashmarks[0] = 0;
69 std::partial_sum(i1_hashmarks.begin(), i1_hashmarks.end(), i1_hashmarks.begin());
70 std::vector<size_t> i2_hashmarks(n2+1, 1);
71 i2_hashmarks[0] = 0;
72 std::partial_sum(i2_hashmarks.begin(), i2_hashmarks.end(), i2_hashmarks.begin());
73
74 std::vector<TA::TiledRange1> hashmarks;
75 hashmarks.push_back(TiledArray::TiledRange1(i1_hashmarks.begin(), i1_hashmarks.end()));
76 hashmarks.push_back(TiledArray::TiledRange1(i2_hashmarks.begin(), i2_hashmarks.end()));
77
78 TiledArray::TiledRange i1i2_trange(hashmarks.begin(), hashmarks.end());
79
80 // make an empty tensor now
81 std::shared_ptr<TArray22> i1i2_g_p1p2(new TArray22(world_, i1i2_trange));
82
83 // make ordinary ranges needed to make p1p2 tiles
84 std::vector<size_t> p1p2_start(2, 0);
85 std::vector<size_t> p1p2_finish(2); p1p2_finish[0] = n3; p1p2_finish[1] = n4;
86 TA::Range p1p2_range(p1p2_start, p1p2_finish);
87
88 { // load all local tiles
89
90 std::vector<double> p1p2_block(n3 * n4);
91
92 std::vector<size_t> i1i2(2, 0);
93 for(i1i2[0] = 0; i1i2[0]<n1; ++i1i2[0]) {
94 for(i1i2[1] = 0; i1i2[1]<n2; ++i1i2[1]) {
95 if (i1i2_g_p1p2->is_local(i1i2)) {
96
97 darray4->retrieve_pair_block(i1i2[0], i1i2[1], te_type, &p1p2_block[0]);
98
99 madness::Future < value_type >
100 tile((value_type(i1i2_g_p1p2->trange().make_tile_range(i1i2),
101 (typename value_type::value_type(p1p2_range,
102 &p1p2_block[0])
103 )
104 )
105 ));
106
107 // Insert the tile into the array
108 i1i2_g_p1p2->set(i1i2, tile);
109
110 darray4->release_pair_block(i1i2[0], i1i2[1], te_type);
111
112 }
113 }
114 }
115
116 }
117
118 return i1i2_g_p1p2;
119 }
120
121 template <typename T>
123 SingleReference_R12Intermediates<T>::_(const std::string& key) {
124 auto v = tarray22_registry_.find(key);
125 if (tarray22_registry_.find(key) == tarray22_registry_.end()) {
126 ParsedTwoBodyFourCenterIntKey pkey(key);
127
129 std::string operset_label;
130 const unsigned int oper_idx = 0;
131 if (pkey.oper() == "g") {
132 operset_label = "ERI";
133 }
134 else if (pkey.oper() == "gr") {
135 operset_label = "R12_m1_G12[0]";
136 }
137 else if (pkey.oper() == "r") {
138 operset_label = "R12_0_G12[0]";
139 }
140 else if (pkey.oper() == "r2") {
141 operset_label = "R12_0_G12[0,0]";
142 }
143 else if (pkey.oper() == "rTr") {
144 operset_label = "G12_T1_G12[0,0]";
145 }
146 else
147 throw ProgrammingError("SingleReference_R12Intermediates<T>::_ : invalid operator key",__FILE__,__LINE__);
148
149 // canonicalize indices, may compute spaces if needed
150 auto bra1 = to_space(pkey.bra1());
151 auto bra2 = to_space(pkey.bra2());
152 auto ket1 = to_space(pkey.ket1());
153 auto ket2 = to_space(pkey.ket2());
154
155 std::string tform_key = ParsedTwoBodyFourCenterIntKey::key(bra1, bra2, ket1, ket2,
156 operset_label, "");
157
158 tarray22_registry_[key] = ab_O_cd(tform_key, oper_idx);
159 v = tarray22_registry_.find(key);
160 }
161
162 return *(v->second);
163 }
164
165 template <typename T>
168
169 auto v = tarray4_registry_.find(key);
170 if (tarray4_registry_.find(key) == tarray4_registry_.end()) {
171
173
175 std::string operset_label;
176 bool rdm2 = false;
177 bool t2 = false;
178 bool l2 = false;
179 const unsigned int oper_idx = 0;
180 if (pkey.oper() == "g") {
181 operset_label = "ERI";
182 }
183 else if (pkey.oper() == "gr") {
184 operset_label = "R12_m1_G12[0]";
185 }
186 else if (pkey.oper() == "r") {
187 operset_label = "R12_0_G12[0]";
188 }
189 else if (pkey.oper() == "r2") {
190 operset_label = "R12_0_G12[0,0]";
191 }
192 else if (pkey.oper() == "rTr") {
193 operset_label = "G12_T1_G12[0,0]";
194 }
195 else if (pkey.oper() == "gamma") {
196 rdm2 = true;
197 }
198 else if (pkey.oper() == "T2") {
199 t2 = true;
200 }
201 else if (pkey.oper() == "L2") {
202 l2 = true;
203 }
204 else
205 throw ProgrammingError("SingleReference_R12Intermediates<T>::ijxy : invalid operator key",__FILE__,__LINE__);
206
207 // canonicalize indices, may compute spaces if needed
208 auto bra1 = to_space(pkey.bra1());
209 auto bra2 = to_space(pkey.bra2());
210 auto ket1 = to_space(pkey.ket1());
211 auto ket2 = to_space(pkey.ket2());
212 auto oreg = r12world_->world()->tfactory()->orbital_registry();
213
214 Ref<DistArray4> darray4;
215 if (rdm2) { // rdm2
216 if (rdm2_.null())
217 throw ProgrammingError("SingleReference_R12Intermediates<T>::ijxy: asked for rdm2, but it had not been given");
218 // if requested spaces don't match exactly, make a new DistArray4
219 auto bra1_space = oreg->value(bra1);
220 auto bra2_space = oreg->value(bra2);
221 auto ket1_space = oreg->value(ket1);
222 auto ket2_space = oreg->value(ket2);
223 const auto& rdm_space = rdm2_->orbs();
224
225 darray4 = rdm2_->da4();
226 if (not (*bra1_space == *rdm_space &&
227 *bra2_space == *rdm_space &&
228 *ket1_space == *rdm_space &&
229 *ket2_space == *rdm_space )) {
230 DistArray4Dimensions dims(1, bra1_space->rank(), bra2_space->rank(), ket1_space->rank(), ket2_space->rank());
231 Ref<DistArray4> darray4_mapped = darray4->clone(dims);
232 std::vector<int> bra1_map = sc::map(*rdm_space, *bra1_space, true);
233 std::vector<int> bra2_map = sc::map(*rdm_space, *bra2_space, true);
234 std::vector<int> ket1_map = sc::map(*rdm_space, *ket1_space, true);
235 std::vector<int> ket2_map = sc::map(*rdm_space, *ket2_space, true);
236
237 std::vector<double> buf(rdm_space->rank() * rdm_space->rank());
238 std::vector<double> map_buf(ket1_space->rank() * ket2_space->rank());
239
240 darray4->activate();
241 darray4_mapped->activate();
242
243 for(int b1=0; b1<bra1_space->rank(); ++b1) {
244 const int bb1 = bra1_map[b1];
245 for(int b2=0; b2<bra2_space->rank(); ++b2) {
246 const int bb2 = bra2_map[b2];
247
248 darray4->retrieve_pair_block(bb1, bb2, 0, &buf[0]);
249
250 const int nk1 = ket1_space->rank();
251 const int nk2 = ket2_space->rank();
252 for(int k1=0, k12=0; k1<nk1; ++k1) {
253 const int kk1 = ket1_map[k1];
254 const int kk1_offset = kk1 * rdm_space->rank();
255 for(int k2=0; k2<nk2; ++k2, ++k12) {
256 const int kk2 = ket2_map[k2];
257 map_buf[k12] = buf[kk1_offset + kk2];
258 }
259 }
260
261 darray4_mapped->store_pair_block(b1, b2, 0, &map_buf[0]);
262 darray4->release_pair_block(bb1, bb2, 0);
263 }
264 }
265 if (darray4->data_persistent()) darray4->deactivate();
266 if (darray4_mapped->data_persistent()) darray4_mapped->deactivate();
267
268 darray4 = darray4_mapped;
269 }
270 }
271 else if (t2) {
272 if (t2_[AlphaBeta].null())
273 throw ProgrammingError("SingleReference_R12Intermediates<T>::ijxy: asked for T2, but it had not been given");
274
275// auto oreg = r12world_->world()->tfactory()->orbital_registry();
276// MPQC_ASSERT(oreg->value(bra1) == r12world_->refwfn()->occ_act(Alpha));
277// MPQC_ASSERT(oreg->value(bra2) == r12world_->refwfn()->occ_act(Beta));
278// MPQC_ASSERT(oreg->value(ket1) == r12world_->refwfn()->uocc_act(Alpha));
279// MPQC_ASSERT(oreg->value(ket2) == r12world_->refwfn()->uocc_act(Beta));
280 darray4 = t2_[AlphaBeta];
281 }
282 else if (l2) {
283 if (l2_[AlphaBeta].null())
284 throw ProgrammingError("SingleReference_R12Intermediates<T>::ijxy: asked for L2, but it had not been given");
285
286 darray4 = l2_[AlphaBeta];
287 }
288 else { // not rdm2 nor t2 nor l2
289 std::string tform_key = ParsedTwoBodyFourCenterIntKey::key(bra1, bra2, ket1, ket2,
290 operset_label, "");
291
292 Ref<TwoBodyMOIntsTransform> tform = r12world_->world()->moints_runtime4()->get(tform_key);
293 tform->compute();
294 darray4 = tform->ints_distarray4();
295 }
296
297 darray4->activate();
298
299 const size_t n1 = darray4->ni();
300 const size_t n2 = darray4->nj();
301 const size_t n3 = darray4->nx();
302 const size_t n4 = darray4->ny();
303
304 // this only works if all nodes can read integrals
305 std::vector<int> tasks_with_access;
306 const int ntasks_with_access = darray4->tasks_with_access(tasks_with_access);
307 MPQC_ASSERT(ntasks_with_access == world_.nproc());
308 MPQC_ASSERT(ntasks_with_access == r12world_->world()->msg()->n());
309
310 // make tiled ranges
311 std::vector<size_t> i_hashmarks = space_hashmarks(bra1);
312 std::vector<size_t> j_hashmarks = space_hashmarks(bra2);
313 std::vector<size_t> x_hashmarks = space_hashmarks(ket1);
314 std::vector<size_t> y_hashmarks = space_hashmarks(ket2);
315
316 std::vector<TA::TiledRange1> hashmarks;
317 hashmarks.push_back(TiledArray::TiledRange1(i_hashmarks.begin(), i_hashmarks.end()));
318 hashmarks.push_back(TiledArray::TiledRange1(j_hashmarks.begin(), j_hashmarks.end()));
319 hashmarks.push_back(TiledArray::TiledRange1(x_hashmarks.begin(), x_hashmarks.end()));
320 hashmarks.push_back(TiledArray::TiledRange1(y_hashmarks.begin(), y_hashmarks.end()));
321 TiledArray::TiledRange ijxy_trange(hashmarks.begin(), hashmarks.end());
322
323 std::shared_ptr<TArray4d> result(new TArray4d(world_, ijxy_trange) );
324 // construct local tiles
325 for(auto t=ijxy_trange.tiles_range().begin();
326 t!=ijxy_trange.tiles_range().end();
327 ++t)
328 if (result->is_local(*t))
329 {
330 std::array<std::size_t, 4> index; std::copy(t->begin(), t->end(), index.begin());
331 madness::Future < typename TArray4d::value_type >
332 tile((DA4_Tile<T>(result.get(), index, darray4, oper_idx)
333 ));
334
335 // Insert the tile into the array
336 result->set(*t, tile);
337 }
338
339 tarray4_registry_[key] = result;
340 v = tarray4_registry_.find(key);
341 }
342
343 return *(v->second);
344 }
345
346 template <typename T>
349
350 auto v = tarray22d_registry_.find(key);
351 if (tarray22d_registry_.find(key) == tarray22d_registry_.end()) {
352
354
356 std::string operset_label;
357 const unsigned int oper_idx = 0;
358 if (pkey.oper() == "g") {
359 operset_label = "ERI";
360 }
361 else if (pkey.oper() == "gr") {
362 operset_label = "R12_m1_G12[0]";
363 }
364 else if (pkey.oper() == "r") {
365 operset_label = "R12_0_G12[0]";
366 }
367 else if (pkey.oper() == "r2") {
368 operset_label = "R12_0_G12[0,0]";
369 }
370 else if (pkey.oper() == "rTr") {
371 operset_label = "G12_T1_G12[0,0]";
372 }
373 else
374 throw ProgrammingError("SingleReference_R12Intermediates<T>::ij_xy : invalid operator key",__FILE__,__LINE__);
375
376 // canonicalize indices, may compute spaces if needed
377 auto bra1 = to_space(pkey.bra1());
378 auto bra2 = to_space(pkey.bra2());
379 auto ket1 = to_space(pkey.ket1());
380 auto ket2 = to_space(pkey.ket2());
381
382 std::string tform_key = ParsedTwoBodyFourCenterIntKey::key(bra1, bra2, ket1, ket2,
383 operset_label, "");
384
385 Ref<TwoBodyMOIntsTransform> tform = r12world_->world()->moints_runtime4()->get(tform_key);
386 tform->compute();
387 Ref<DistArray4> darray4 = tform->ints_distarray4();
388 darray4->activate();
389
390 const size_t n1 = darray4->ni();
391 const size_t n2 = darray4->nj();
392 const size_t n3 = darray4->nx();
393 const size_t n4 = darray4->ny();
394
395 // this only works if all nodes can read integrals
396 std::vector<int> tasks_with_access;
397 const int ntasks_with_access = darray4->tasks_with_access(tasks_with_access);
398 MPQC_ASSERT(ntasks_with_access == world_.nproc());
399 MPQC_ASSERT(ntasks_with_access == r12world_->world()->msg()->n());
400
401 // make tiled ranges
402 // N.B. using little trickery using partial sum to create tilesize=1 range
403 std::vector<size_t> i_hashmarks(n1+1, 1);
404 i_hashmarks[0] = 0;
405 std::partial_sum(i_hashmarks.begin(), i_hashmarks.end(), i_hashmarks.begin());
406 std::vector<size_t> j_hashmarks(n2+1, 1);
407 j_hashmarks[0] = 0;
408 std::partial_sum(j_hashmarks.begin(), j_hashmarks.end(), j_hashmarks.begin());
409
410 std::vector<TA::TiledRange1> hashmarks;
411 hashmarks.push_back(TiledArray::TiledRange1(i_hashmarks.begin(), i_hashmarks.end()));
412 hashmarks.push_back(TiledArray::TiledRange1(j_hashmarks.begin(), j_hashmarks.end()));
413 TiledArray::TiledRange ij_trange(hashmarks.begin(), hashmarks.end());
414
415 std::shared_ptr<TArray22d> result(new TArray22d(world_, ij_trange) );
416 // construct local tiles
417 for(auto t=ij_trange.tiles_range().begin();
418 t!=ij_trange.tiles_range().end();
419 ++t)
420 if (result->is_local(*t))
421 {
422 std::array<std::size_t, 2> index; std::copy(t->begin(), t->end(), index.begin());
423 madness::Future < typename TArray22d::value_type >
424 tile((DA4_Tile34<T>(result.get(), index, darray4, oper_idx)
425 ));
426
427 // Insert the tile into the array
428 result->set(*t, tile);
429 }
430
431 tarray22d_registry_[key] = result;
432 v = tarray22d_registry_.find(key);
433 }
434
435 return *(v->second);
436 }
437
438 template <typename T>
441
442 auto v = tarray2_registry_.find(key);
443 if (v == tarray2_registry_.end()) {
444
445 ParsedOneBodyIntKey pkey(key);
446
447 // canonicalize indices, may compute spaces if needed
448 auto bra_id = to_space(pkey.bra());
449 auto ket_id = to_space(pkey.ket());
450
451 auto oreg = r12world_->world()->tfactory()->orbital_registry();
452 auto freg = r12world_->world()->fockbuild_runtime();
453 auto bra = oreg->value(bra_id);
454 auto ket = oreg->value(ket_id);
455 auto nbra = bra->rank();
456 auto nket = ket->rank();
457
458 // grab the matrix
459 RefSCMatrix operator_matrix;
460 if (pkey.oper() == "J")
461 operator_matrix = freg->get(ParsedOneBodyIntKey::key(bra_id, ket_id, "J"));
462 else if (pkey.oper() == "K")
463 operator_matrix = freg->get(ParsedOneBodyIntKey::key(bra_id, ket_id, "K"));
464 else if (pkey.oper() == "F")
465 operator_matrix = freg->get(ParsedOneBodyIntKey::key(bra_id, ket_id, "F"));
466 else if (pkey.oper() == "h")
467 operator_matrix = freg->get(ParsedOneBodyIntKey::key(bra_id, ket_id, "H"));
468 else if (pkey.oper() == "hJ")
469 operator_matrix = freg->get(ParsedOneBodyIntKey::key(bra_id, ket_id, "H"))
470 + freg->get(ParsedOneBodyIntKey::key(bra_id, ket_id, "J"));
471 else if (pkey.oper().find("mu_") == 0 || pkey.oper().find("q_") == 0 ||
472 pkey.oper().find("dphi_") == 0 || pkey.oper().find("ddphi_") == 0) {
473 operator_matrix = freg->get(ParsedOneBodyIntKey::key(bra_id, ket_id, pkey.oper()));
474 }
475 else if (pkey.oper() == "gamma")
476 operator_matrix = this->rdm1(bra, ket);
477 else if (pkey.oper() == "I") {
478 operator_matrix = bra->coefs()->kit()->matrix(bra->coefs().coldim(), ket->coefs().coldim());
479 if (bra == ket) {
480 operator_matrix.assign(0.0);
481 for(int i=0; i<nbra; ++i)
482 operator_matrix.set_element(i,i,1.0);
483 }
484 else {
485 operator_matrix = sc::overlap(*bra,*ket,operator_matrix->kit());
486 }
487 }
488 else if (pkey.oper() == "T1") {
489 if (ket_id == "a") { // if given t1_ explicitly (CC), make sure its size matches
490 if (t1_) {
491 if (t1_.ncol() != oreg->value(ket_id)->rank())
492 throw ProgrammingError("SingleReference_R12Intermediates::xy() -- T1.ncol() != nvir_act",
493 __FILE__, __LINE__);
494 operator_matrix = t1_;
495 }
496 else if (t1_cabs_) { // if t1_cabs given, this means extract the occ_act x vir_act block
497 if (t1_cabs_.ncol() == oreg->value("a'")->rank()) // doesn't make sense if T1 only include CABS
498 throw ProgrammingError("SingleReference_R12Intermediates::xy() -- asked for <i|T1|a> but T1_cabs does not include conv. virtuals",
499 __FILE__, __LINE__);
500 Ref<OrbitalSpace> vir_act = oreg->value("a");
501 Ref<OrbitalSpace> allvir = oreg->value("A'"); // should exist since t1_cabs_ spans more than just CABS
502 MPQC_ASSERT(t1_cabs_.ncol() == allvir->rank()); // this should be the only other possibility
503 RefSCMatrix t1_subblock = t1_cabs_.kit()->matrix(t1_cabs_.rowdim(),
504 vir_act->coefs()->coldim());
505 std::vector<int> vir_to_allvir = sc::map(*allvir, *vir_act, false);
506 const int nocc_act = t1_subblock.nrow();
507 const int nvir_act = vir_act->rank();
508 for(int i=0; i<nocc_act; ++i) {
509 for(int a=0; a<nvir_act; ++a) {
510 const int Ap = vir_to_allvir[a];
511 t1_subblock.set_element(i, a, t1_cabs_.get_element(i, Ap));
512 }
513 }
514 operator_matrix = t1_subblock;
515 }
516 }
517 else {
518 MPQC_ASSERT(ket_id == "a'" || ket_id == "A'");
519 MPQC_ASSERT(t1_cabs_);
520 MPQC_ASSERT(oreg->value(ket_id)->rank() == t1_cabs_.ncol());
521 operator_matrix = t1_cabs_;
522 }
523 }
524 else if (pkey.oper() == "L1") {
525 MPQC_ASSERT(l1_);
526 MPQC_ASSERT(oreg->value(bra_id)->rank() == l1_.nrow());
527 MPQC_ASSERT(oreg->value(ket_id)->rank() == l1_.ncol());
528 operator_matrix = l1_;
529 }
530 else {
531 std::ostringstream oss;
532 oss << "SingleReference_R12Intermediates<T>::xy -- operator " << pkey.oper() << " is not recognized";
533 throw ProgrammingError(oss.str().c_str(),
534 __FILE__, __LINE__);
535 }
536
537 // make tiled ranges
538 std::vector<size_t> bra_hashmarks = space_hashmarks(bra_id);
539 std::vector<size_t> ket_hashmarks = space_hashmarks(ket_id);
540
541 std::vector<TA::TiledRange1> hashmarks;
542 hashmarks.push_back(TiledArray::TiledRange1(bra_hashmarks.begin(), bra_hashmarks.end()));
543 hashmarks.push_back(TiledArray::TiledRange1(ket_hashmarks.begin(), ket_hashmarks.end()));
544 TiledArray::TiledRange braket_trange(hashmarks.begin(), hashmarks.end());
545
546 std::shared_ptr<TArray2> result(new TArray2(world_, braket_trange) );
547 // construct local tiles
548 for(auto t=result->pmap()->begin();
549 t!=result->pmap()->end();
550 ++t)
551 {
552
553 auto tile_range = result->trange().make_tile_range(*t);
554 std::vector<double> ptr_data(tile_range.volume());
555
556 for(size_t r=tile_range.lobound()[0], rc=0;
557 r != tile_range.upbound()[0];
558 ++r) {
559 for(size_t c=tile_range.lobound()[1];
560 c != tile_range.upbound()[1];
561 ++c, ++rc) {
562 ptr_data[rc] = operator_matrix->get_element(r,c);
563 }
564 }
565
566 madness::Future < typename TArray2::value_type >
567 tile(typename TArray2::value_type(tile_range, &ptr_data[0]));
568
569 // Insert the tile into the array
570 result->set(*t, tile);
571 }
572
573 tarray2_registry_[key] = result;
574 v = tarray2_registry_.find(key);
575 }
576
577 return *(v->second);
578 }
579
580 template <typename T>
581 TA::expressions::TsrExpr<const typename SingleReference_R12Intermediates<T>::TArray4d, true>
583 const TArray4d& tarray4 = ijxy(key);
585
586 std::array<std::string, 4> ind = {{pkey.bra1(), pkey.bra2(), pkey.ket1(), pkey.ket2()}};
587 for(auto v=ind.begin(); v!=ind.end(); ++v) {
588 if (ParsedTransformedOrbitalSpaceKey::valid_key(*v) ) {
590 *v = ptkey.label();
591 }
592 }
593
594 const std::string annotation = ind[0] + "," + ind[1] + "," + ind[2] + "," + ind[3];
595 return tarray4(annotation);
596 }
597
598 template <typename T>
599 TA::expressions::TsrExpr<const typename SingleReference_R12Intermediates<T>::TArray2, true>
601 const TArray2& tarray2 = xy(key);
602 ParsedOneBodyIntKey pkey(key);
603
604 std::array<std::string, 2> ind = {{pkey.bra(), pkey.ket()}};
605 for(auto v=ind.begin(); v!=ind.end(); ++v) {
606 if (ParsedTransformedOrbitalSpaceKey::valid_key(*v) ) {
608 *v = ptkey.label();
609 }
610 }
611
612 const std::string annotation = ind[0] + "," + ind[1];
613 return tarray2(annotation);
614 }
615
616 template <typename T>
618 template <typename T>
620
621 template <typename T>
622 TA::expressions::TsrExpr<const typename SingleReference_R12Intermediates<T>::TArray4Tg, true>
624
626
627 auto v = tarray4tg_registry_.find(key);
628 if (v == tarray4tg_registry_.end()) {
629
630 if (pkey.oper() != "Tg") {
631 std::ostringstream oss;
632 oss << "SingleReference_R12Intermediates<T>::_Tg : " << pkey.oper() << " is an invalid operator key";
633 throw ProgrammingError(oss.str().c_str(), __FILE__, __LINE__);
634 }
635
636 // canonicalize indices, may compute spaces if needed
637 auto bra1 = to_space(pkey.bra1());
638 auto bra2 = to_space(pkey.bra2());
639 auto ket1 = to_space(pkey.ket1());
640 auto ket2 = to_space(pkey.ket2());
641
642 if (not (bra1 == bra2 && bra1 == ket1 && bra1 == ket2)) {
643 throw ProgrammingError("SingleReference_R12Intermediates<T>::_Tg : all spaces must be the same", __FILE__, __LINE__);
644 }
645
646 // make tiled ranges
647 std::vector<size_t> i_hashmarks = space_hashmarks(bra1);
648 std::vector<size_t> j_hashmarks = space_hashmarks(bra2);
649 std::vector<size_t> x_hashmarks = space_hashmarks(ket1);
650 std::vector<size_t> y_hashmarks = space_hashmarks(ket2);
651
652 std::vector<TA::TiledRange1> hashmarks;
653 hashmarks.push_back(TiledArray::TiledRange1(i_hashmarks.begin(), i_hashmarks.end()));
654 hashmarks.push_back(TiledArray::TiledRange1(j_hashmarks.begin(), j_hashmarks.end()));
655 hashmarks.push_back(TiledArray::TiledRange1(x_hashmarks.begin(), x_hashmarks.end()));
656 hashmarks.push_back(TiledArray::TiledRange1(y_hashmarks.begin(), y_hashmarks.end()));
657 TiledArray::TiledRange ijxy_trange(hashmarks.begin(), hashmarks.end());
658
659 std::shared_ptr<TArray4Tg> result(new TArray4Tg(world_, ijxy_trange) );
660 // construct local tiles
661 for(auto t=ijxy_trange.tiles_range().begin();
662 t!=ijxy_trange.tiles_range().end();
663 ++t)
664 if (result->is_local(*t))
665 {
666 std::array<std::size_t, 4> index; std::copy(t->begin(), t->end(), index.begin());
667 madness::Future < typename TArray4Tg::value_type >
668 tile((typename TArray4Tg::value_type(result.get(), index, &tg_s0_gen)
669 ));
670
671 // Insert the tile into the array
672 result->set(*t, tile);
673 }
674
675 tarray4tg_registry_[key] = result;
676 v = tarray4tg_registry_.find(key);
677 }
678 const TArray4Tg& tarray4 = *(v->second);
679
680 std::array<std::string, 4> ind = {{pkey.bra1(), pkey.bra2(), pkey.ket1(), pkey.ket2()}};
681 for(auto v=ind.begin(); v!=ind.end(); ++v) {
682 if (ParsedTransformedOrbitalSpaceKey::valid_key(*v) ) {
684 *v = ptkey.label();
685 }
686 }
687
688 const std::string annotation = ind[0] + "," + ind[1] + "," + ind[2] + "," + ind[3];
689 return tarray4(annotation);
690 }
691
692 namespace expressions {
693
694 template <typename ArgType, bool Transpose = false>
695 struct trace_tensor2_op : public std::binary_function<ArgType, ArgType, TA::Tensor<typename ArgType::numeric_type> > {
696 typedef typename ArgType::numeric_type numeric_type;
697 typedef TA::Tensor<typename ArgType::numeric_type> result_type;
698 typedef result_type value_type;
699
700 static bool is_dense(bool first_dense, bool second_dense) {
701 return true;
702 }
703 static bool is_zero(const bool first_zero, const bool second_zero) {
704 return first_zero || second_zero;
705 }
706
707 template <typename Left, typename Right>
708 static void shape(::TiledArray::detail::Bitset<>& result, const Left& left, const Right& right) {
709 result = ::TiledArray::detail::Bitset<>(0);
710 }
711
712 // applies the scaling factor
713 void scale(numeric_type s) {
714 scale_ = s;
715 }
716
717 result_type
718 operator()(const ArgType& arg1, const ArgType& arg2) const {
719 TA_ASSERT(arg1.size() == 1);
720 TA_ASSERT(arg2.size() == 1);
721 TA_ASSERT(arg1.data()->size() == arg2.data()->size());
722 TA_ASSERT(arg1.range() == arg2.range()); // structure of the arguments must be the same
723
724 numeric_type dot_product;
725 if (not Transpose) // straight dot
726 dot_product = TA::dot(arg1.data()->size(), arg1.data()->data(), arg2.data()->data());
727 else { // transposed dot done as series of row*column products
728 const integer nrow1 = arg1.data()->range().size()[0];
729 const integer ncol1 = arg1.data()->range().size()[1];
730 TA_ASSERT(ncol1 == arg2.data()->range().size()[0]); // ncol1 == nrow2
731 const integer unit_stride = 1;
732 dot_product = 0;
733
734 const numeric_type* ptr1 = arg1.data()->data();
735 const numeric_type* ptr2 = arg2.data()->data();
736 for(integer row1=0; row1<nrow1; ++row1, ptr1+=ncol1, ++ptr2) {
737 dot_product += TA::dot(ncol1, ptr1, unit_stride, ptr2, nrow1);
738 }
739 }
740 //std::cout << dot_product << std::endl;
741
742 result_type result(arg1.range(), &dot_product);
743 return result * scale_;
744 }
745
746 result_type operator()(const TiledArray::ZeroTensor&, const ArgType&) {
747 TA_ASSERT(false); // Should not be used.
748 return result_type();
749 }
750
751 result_type operator()(const ArgType&, const TiledArray::ZeroTensor&) {
752 TA_ASSERT(false); // Should not be used.
753 return result_type();
754 }
755
756 private:
757 numeric_type scale_;
758
759 };
760
761 template <typename ArgType, bool Transpose = false>
762 struct diag_tensor2_op : public std::unary_function<ArgType, TA::Tensor<typename ArgType::numeric_type> > {
763 typedef typename ArgType::numeric_type numeric_type;
764 typedef TA::Tensor<typename ArgType::numeric_type> result_type;
765 typedef result_type value_type;
766
767 static bool is_dense(bool arg_dense) {
768 return true;
769 }
770 static bool is_zero(const bool arg_zero) {
771 return arg_zero;
772 }
773
774 template <typename Arg>
775 static void shape(::TiledArray::detail::Bitset<>& result, const Arg& arg) {
776 result = ::TiledArray::detail::Bitset<>(0);
777 }
778
779 // applies the scaling factor
780 void scale(numeric_type s) {
781 scale_ = s;
782 }
783
784 result_type
785 operator()(const ArgType& arg) const {
786 TA_ASSERT(arg.size() == 1);
787
788 numeric_type result_value = 0;
789 const auto& ij = arg.range().start();
790 const auto& pq_size = arg.data()->range().size();
791 if (not Transpose) { // extract i,j
792 result_value = arg.data()->data()[ij[0] * pq_size[1] + ij[1]];
793 }
794 else { // extract j,i
795 result_value = arg.data()->data()[ij[1] * pq_size[1] + ij[0]];
796 }
797
798 result_type result(arg.range(), &result_value);
799 return result * scale_;
800 }
801
802 result_type operator()(const TiledArray::ZeroTensor&) {
803 TA_ASSERT(false); // Should not be used.
804 return result_type();
805 }
806
807 private:
808 numeric_type scale_;
809 };
810
811 }; // end of namespace sc::expressions
812
813 template <typename T>
814 std::string
815 SingleReference_R12Intermediates<T>::to_space(const std::string& index) {
816
817 bool tformed_space = ParsedTransformedOrbitalSpaceKey::valid_key(index);
818 if (tformed_space) {
820
821 // make sure this space exists
822 const std::string space_label = to_space(pkey.label());
823 std::string space_key = ParsedTransformedOrbitalSpaceKey::key(space_label, pkey.spin(),
824 pkey.support_label(), pkey.support_spin(),
825 pkey.transform_operator());
826 auto oreg = r12world_->world()->tfactory()->orbital_registry();
827 auto freg = r12world_->world()->fockbuild_runtime();
828
829 if (not oreg->key_exists(space_key)) { // compute if needed
830 std::string support_key = ParsedOrbitalSpaceKey::key(to_space(pkey.support_label()), pkey.support_spin());
831 std::string target_key = ParsedOrbitalSpaceKey::key(space_label, pkey.spin());
832
833 Ref<OrbitalSpace> target_space = oreg->value(target_key);
834 Ref<OrbitalSpace> support_space = oreg->value(support_key);
835 RefSCMatrix support_coefs = support_space->coefs();
836
837 RefSCMatrix operator_matrix;
838 switch(pkey.transform_operator()) {
839 case OneBodyOper::J:
840 operator_matrix = freg->get(ParsedOneBodyIntKey::key(support_key, target_key, "J"));
841 break;
842 case OneBodyOper::K:
843 operator_matrix = freg->get(ParsedOneBodyIntKey::key(support_key, target_key, "K"));
844 break;
845 case OneBodyOper::F:
846 operator_matrix = freg->get(ParsedOneBodyIntKey::key(support_key, target_key, "F"));
847 break;
848 case OneBodyOper::hJ:
849 operator_matrix = freg->get(ParsedOneBodyIntKey::key(support_key, target_key, "H"))
850 + freg->get(ParsedOneBodyIntKey::key(support_key, target_key, "J"));
851 break;
853 operator_matrix = this->rdm1(support_space, target_space);
854 break;
855 default:
856 throw ProgrammingError("SingleReference_R12Intermediates<T>::to_space -- transformed spaces only with Fock-like operators are supported now",
857 __FILE__, __LINE__);
858 }
859
860 Ref<OrbitalSpace> tformed_space = new OrbitalSpace(space_key, space_key,
861 target_space, support_coefs * operator_matrix,
862 support_space->basis());
863 operator_matrix = 0;
864 oreg->add(make_keyspace_pair(tformed_space));
865 }
866 return space_key;
867 }
868
869 return to_space_(index);
870 }
871
872 template <typename T>
873 RefSCMatrix SingleReference_R12Intermediates<T>::rdm1(const Ref<OrbitalSpace>& row_space,
874 const Ref<OrbitalSpace>& col_space) const {
875 // get the reference 1-RDM
876 RefSymmSCMatrix gamma_total = r12world_->refwfn()->ordm_occ_sb(Alpha) + r12world_->refwfn()->ordm_occ_sb(Beta);
877 Ref<OrbitalSpace> occ = r12world_->refwfn()->occ_sb(Alpha);
878 // map to the target space and support space
879 std::vector<int> col_to_occ;
880 std::vector<int> row_to_occ;
881 try {
882 const bool require_same_basis = true;
883 col_to_occ = map(*occ, *col_space, require_same_basis);
884 row_to_occ = map(*occ, *row_space, require_same_basis);
885 }
886 catch (sc::CannotConstructMap&) {
887 throw ProgrammingError("requested RDM1 in spaces x and y, but x and y are not supported by the same basis",
888 __FILE__, __LINE__);
889 }
890 RefSCMatrix result = gamma_total.kit()->matrix(row_space->dim(), col_space->dim());
891 result.assign(0.0);
892 const int nr = result.nrow();
893 const int nc = result.ncol();
894 for(int r=0; r<nr; ++r) {
895 const int rr = row_to_occ[r];
896 if (rr >= 0)
897 for(int c=0; c<nc; ++c) {
898 const int cc = col_to_occ[c];
899 if (cc >= 0)
900 result.set_element(r, c, gamma_total.get_element(rr, cc));
901 }
902 }
903
904 return result;
905 }
906
907 template<typename T>
908 std::string SingleReference_R12Intermediates<T>::to_space_(std::string index) {
909
910 index.erase(std::remove_if(index.begin(), index.end(),
911 [](char c){
912 return c >= '0' && c<='9';
913 }
914 ), index.end());
915
916 if (index == "i" || index == "j" || index == "k" || index == "l")
917 return "i";
918 if (index == "m" || index == "n")
919 return "m";
920 if (index == "i'" || index == "j'" || index == "k'" || index == "l'")
921 return "i'";
922 if (index == "a" || index == "b" || index == "c" || index == "d")
923 return "a";
924 if (index == "e" || index == "f")
925 return "e";
926 if (index == "a" || index == "b" || index == "c" || index == "d")
927 return "a";
928 if (index == "p" || index == "q" || index == "r" || index == "s")
929 return "p";
930 if (index == "p'" || index == "q'" || index == "r'" || index == "s'")
931 return "p'";
932 if (index == "a'" || index == "b'" || index == "c'" || index == "d'")
933 return "a'";
934 if (index == "A'" || index == "B'" || index == "C'" || index == "D'")
935 return "A'";
936
937 std::ostringstream oss;
938 oss << "SingleReference_R12Intermediates::to_space(): index label " << index << " not in the current dictionary:";
939 throw ProgrammingError(oss.str().c_str(),
940 __FILE__, __LINE__);
941 }
942
943 template <typename T>
944 std::vector<size_t>
945 SingleReference_R12Intermediates<T>::space_hashmarks(std::string space_label) const {
946
947 const bool tformed_space = ParsedTransformedOrbitalSpaceKey::valid_key(space_label);
948 if (tformed_space) {
949 ParsedTransformedOrbitalSpaceKey ptkey(space_label);
950 space_label = ptkey.label();
951 }
952
953 auto oreg = r12world_->world()->tfactory()->orbital_registry();
954 size_t rank = oreg->value(space_label)->rank();
955 std::vector<size_t> hashmarks(2,0); hashmarks[1] = rank; // all spaces have 1 tile, except
956 if (space_label == "i" || space_label == "m") { // occupied spaces, whose size is determined heuristically for now
957 // TODO create a central registry of space tilings?
958 const int nproc = r12world_->world()->msg()->n();
959 const size_t desired_tilesize = std::max(4 ,
960 std::min( (int)floor((double)rank / sqrt((double)nproc)), 10)
961 );
962 const size_t ntiles = (rank+desired_tilesize-1)/desired_tilesize;
963 const size_t tilesize = (rank + ntiles - 1)/ ntiles;
964 hashmarks.resize(ntiles+1);
965 hashmarks[0] = 0;
966 for(size_t i=1; i<ntiles; ++i)
967 hashmarks[i] = hashmarks[i-1] + tilesize;
968 hashmarks.back() = rank;
969 }
970
971 return hashmarks;
972 }
973
974#if 0
975 namespace {
976 std::vector<std::string>
977 split(const std::string& key,
978 char separator) {
979
980 std::stringstream ss(key);
981 std::string item;
982 std::vector<std::string> tokens;
983 while (std::getline(ss, item, separator)) {
984 tokens.push_back(item);
985 }
986
987 return tokens;
988 }
989
990 }
991
992 template <typename T> template <size_t NDIM>
993 TA::TiledRange
994 SingleReference_R12Intermediates<T>::make_trange(const std::array<std::string, NDIM>& spaces) const {
995
996 std::vector<TA::TiledRange1> hashmarks;
997 for(auto id=spaces.begin(); id!=spaces.end(); ++id) {
998 std::vector<size_t> hashmark = space_hashmarks(*id);
999 hashmarks.push_back(TiledArray::TiledRange1(hashmark.begin(), hashmark.end()));
1000 }
1001 return TA::TiledRange(hashmarks.begin(), hashmarks.end());
1002 }
1003
1004 namespace detail {
1005
1006 template <typename T, typename Index>
1007 void
1008 sieve_tensor(TA::Array<T, 2>* result,
1009 Index t,
1010 madness::Future<typename TA::Array<T, 2>::value_type > source_tile,
1011 TA::Array<T, 2>* source,
1012 Index t_source,
1013 std::vector<int>* row_map,
1014 std::vector<int>* col_map) {
1015
1016 typedef TA::Array<T, 2> TArray2;
1017
1018 auto result_trange = result->trange().make_tile_range(t);
1019 auto source_trange = source->trange().make_tile_range(t_source);
1020
1021 std::vector<T> ptr_result(result_trange.volume(), T(0.0));
1022
1023 for(size_t r=source_trange.start()[0], rc=0;
1024 r != source_trange.finish()[0];
1025 ++r) {
1026 const int rr = row_map->operator[](r);
1027 if (rr != -1) {
1028 for(size_t c=source_trange.start()[1];
1029 c != source_trange.finish()[1];
1030 ++c, ++rc) {
1031
1032 const int cc = col_map->operator[](c);
1033
1034 if (cc != -1) {
1035
1036 MPQC_ASSERT(false);
1037
1038 }
1039 }
1040 }
1041 }
1042
1043 madness::Future < typename TArray2::value_type >
1044 tile(typename TArray2::value_type(result_trange, &ptr_result[0]));
1045
1046 // Insert the tile into the array
1047 result->set(*t, tile);
1048
1049 }
1050 }
1051
1052 template <typename T>
1054 SingleReference_R12Intermediates<T>::sieve(const TArray2& input,
1055 const std::string& output_annotation) {
1056
1057 // find the input key
1058 std::string input_key;
1059 for(auto i = tarray2_registry_.begin();
1060 i!=tarray2_registry_.end(); ++i) {
1061 if (i->second.get() == &input) {
1062 input_key = i->first;
1063 break;
1064 }
1065 }
1066 MPQC_ASSERT(input_key.empty() == false);
1067
1068 std::vector<std::string> output_space_keys = split(output_annotation, ',');
1069 MPQC_ASSERT(output_space_keys.size() == 2);
1070
1071 ParsedOneBodyIntKey input_pkey(input_key);
1072 ParsedOneBodyIntKey output_pkey(ParsedOneBodyIntKey::key(output_space_keys[0],
1073 output_space_keys[1],
1074 input_pkey.oper()));
1075
1076 auto v = tarray2_registry_.find(output_pkey.key());
1077 if (v == tarray2_registry_.end()) {
1078
1079 // canonicalize indices, may compute spaces if needed
1080 auto ibra_id = to_space(input_pkey.bra());
1081 auto iket_id = to_space(input_pkey.ket());
1082 auto obra_id = to_space(output_pkey.bra());
1083 auto oket_id = to_space(output_pkey.ket());
1084
1085 auto oreg = r12world_->world()->tfactory()->orbital_registry();
1086 auto freg = r12world_->world()->fockbuild_runtime();
1087 auto ibra = oreg->value(ibra_id);
1088 auto iket = oreg->value(iket_id);
1089 auto obra = oreg->value(obra_id);
1090 auto oket = oreg->value(oket_id);
1091
1092 auto nbra = obra->rank();
1093 auto nket = oket->rank();
1094
1095 std::vector<int> bra_map = map(*obra, *ibra);
1096 std::vector<int> ket_map = map(*oket, *iket);
1097
1098 // construct the output array
1099 std::array<std::string, 2> ospaces = {{obra_id, oket_id}};
1100 TA::TiledRange output_trange = make_trange(ospaces);
1101 std::shared_ptr<TArray2> result(new TArray2(world_, output_trange) );
1102
1103 // construct tasks to fill in local tiles
1104 for(auto t=result->pmap()->begin();
1105 t!=result->pmap()->end();
1106 ++t)
1107 {
1108
1109 // map t to the index of the tile that contains the data needed for t
1110 auto t_input = *t;
1111 MPQC_ASSERT(false);
1112
1113 world_.taskq.add(& detail::sieve_tensor,
1114 &result, *t, input.find(t_input), &bra_map, &ket_map);
1115 }
1116
1117 tarray2_registry_[output_pkey.key()] = result;
1118 v = tarray2_registry_.find(output_pkey.key());
1119
1120 }
1121
1122 return *(v->second);
1123 }
1124#endif
1125
1126 template <typename T>
1127 std::ostream& operator<<(std::ostream& os, const DA4_Tile<T>& t) {
1128 os << typename DA4_Tile<T>::eval_type(t);
1129 return os;
1130 }
1131
1132 template <typename T>
1133 std::ostream& operator<<(std::ostream& os, const DA4_Tile34<T>& t) {
1134 os << typename DA4_Tile34<T>::eval_type(t);
1135 return os;
1136 }
1137
1138
1139}; // end of namespace sc
1140
1141#endif // end of header guard
1142
1143
1144// Local Variables:
1145// mode: c++
1146// c-file-style: "CLJ-CONDENSED"
1147// End:
Definition orbitalspace.h:532
Tile of a <34> slice of <1234> that's "evaluated" when needed by reading from DistArray4 holding pqrs...
Definition sr_r12intermediates.h:89
Tile of a 4-index tensor that's "evaluated" when needed by reading from DistArray4.
Definition sr_r12intermediates.h:51
Parsed representation of a string key that represents a set of one-body integrals.
Definition fockbuild_runtime.h:219
Parses keys of a "transformed" OrbitalSpace.
Definition orbitalspace.h:635
Parsed representation of a string key that represents a set of 4-center 2-body integrals.
Definition tbint_runtime.h:160
This is thrown when a situations arises that should be impossible.
Definition scexception.h:92
The RefSCMatrix class is a smart pointer to an SCMatrix specialization.
Definition matrix.h:135
A template class that maintains references counts.
Definition ref.h:361
SingleReference_R12Intermediates computes R12/F12 intermediates using MPQC3 runtime.
Definition sr_r12intermediates.h:218
TArray2 rdm1()
returns the 1-particle reduced density matrix
Definition sr_r12intermediates_VXB_diag.h:6783
TA::Array< T, 2 > TArray2
standard 2-index tensor
Definition sr_r12intermediates.h:226
TA::Array< T, 4, DA4_Tile< T > > TArray4d
4-index tensor with lazy tiles
Definition sr_r12intermediates.h:224
TA::Array< T, 4, LazyTensor< T, 4, expressions::TGeminalGenerator< T > > > TArray4Tg
geminal T tensor
Definition sr_r12intermediates.h:237
TArray22d & ij_xy(const std::string &key)
Definition sr_r12intermediates_util.h:348
TA::expressions::TsrExpr< const TArray2, true > _2(const std::string &key)
Given a descriptive key, creates a rank-2 Array of integrals, or other related quantities The syntax ...
Definition sr_r12intermediates_util.h:600
TA::expressions::TsrExpr< const TArray4d, true > _4(const std::string &key)
Given a descriptive key, creates a rank-4 Array of integrals, or other related quantities The syntax ...
Definition sr_r12intermediates_util.h:582
TA::Array< TA::Tensor< T >, 2, DA4_Tile34< T > > TArray22d
2-index tensor of lazy 2-index tensors
Definition sr_r12intermediates.h:234
TArray2 & xy(const std::string &key)
Definition sr_r12intermediates_util.h:440
TArray2 & sieve(const TArray2 &input, const std::string &output_annotation)
sieves x|o1|y -> x'|o1|y' does not throw only if each tile of the result depends only on 1 tile of th...
TA::Array< TA::Tensor< T >, 2 > TArray22
2-index tensor with lazy tiles
Definition sr_r12intermediates.h:232
TA::expressions::TsrExpr< const TArray4Tg, true > _Tg(const std::string &key)
like _4, produces geminal T tensor
Definition sr_r12intermediates_util.h:623
TArray4d & ijxy(const std::string &key)
Definition sr_r12intermediates_util.h:167
std::vector< int > map(const GaussianBasisSet &B, const GaussianBasisSet &A)
same as operator<<, except A does not have to be contained in B, map[a] = -1 if function a is not in ...
std::vector< unsigned int > operator<<(const GaussianBasisSet &B, const GaussianBasisSet &A)
computes a map from basis functions in A to the equivalent basis functions in B.
RefSCMatrix overlap(const OrbitalSpace &space2, const OrbitalSpace &space1, const Ref< SCMatrixKit > &kit=SCMatrixKit::default_matrixkit())
overlap(s2,s1) returns the overlap matrix between s2 and s1.
std::pair< std::string, Ref< OrbitalSpace > > make_keyspace_pair(const Ref< OrbitalSpace > &space, SpinCase1 spin=AnySpinCase1)
helper function to form a key/space pair from a OrbitalSpace
Contains all MPQC code up to version 3.
Definition mpqcin.h:14
Definition distarray4.h:43
@ F
Fock operator.
Definition operator.h:110
@ K
(electronic) exchange
Definition operator.h:109
@ hJ
h+J
Definition operator.h:111
@ gamma
1-body reduced density matrix
Definition operator.h:104
@ J
(electronic) Coulomb
Definition operator.h:108
makes a geminal T tensor
Definition sr_r12intermediates.h:186
Definition sr_r12intermediates_util.h:762
Definition sr_r12intermediates_util.h:695

Generated at Wed Sep 25 2024 02:45:30 for MPQC 3.0.0-alpha using the documentation package Doxygen 1.12.0.